Actual source code: test4.c

slepc-3.11.2 2019-07-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */

 11: static char help[] = "Test BV operations, changing the number of active columns.\n\n";

 13: #include <slepcbv.h>

 15: int main(int argc,char **argv)
 16: {
 18:   Vec            t,v;
 19:   Mat            Q,M;
 20:   BV             X,Y;
 21:   PetscInt       i,j,n=10,kx=6,lx=3,ky=5,ly=2;
 22:   PetscScalar    *q,*z;
 23:   PetscReal      nrm;
 24:   PetscViewer    view;
 25:   PetscBool      verbose,trans;

 27:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 28:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 29:   PetscOptionsGetInt(NULL,NULL,"-kx",&kx,NULL);
 30:   PetscOptionsGetInt(NULL,NULL,"-lx",&lx,NULL);
 31:   PetscOptionsGetInt(NULL,NULL,"-ky",&ky,NULL);
 32:   PetscOptionsGetInt(NULL,NULL,"-ly",&ly,NULL);
 33:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 34:   PetscPrintf(PETSC_COMM_WORLD,"First BV with %D active columns (%D leading columns) of dimension %D.\n",kx,lx,n);
 35:   PetscPrintf(PETSC_COMM_WORLD,"Second BV with %D active columns (%D leading columns) of dimension %D.\n",ky,ly,n);

 37:   /* Create template vector */
 38:   VecCreate(PETSC_COMM_WORLD,&t);
 39:   VecSetSizes(t,PETSC_DECIDE,n);
 40:   VecSetFromOptions(t);

 42:   /* Create BV object X */
 43:   BVCreate(PETSC_COMM_WORLD,&X);
 44:   PetscObjectSetName((PetscObject)X,"X");
 45:   BVSetSizesFromVec(X,t,kx+2);  /* two extra columns to test active columns */
 46:   BVSetFromOptions(X);
 47:   BVSetActiveColumns(X,lx,kx);

 49:   /* Set up viewer */
 50:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
 51:   if (verbose) {
 52:     PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
 53:   }

 55:   /* Fill X entries */
 56:   for (j=0;j<kx+2;j++) {
 57:     BVGetColumn(X,j,&v);
 58:     VecSet(v,0.0);
 59:     for (i=0;i<4;i++) {
 60:       if (i+j<n) {
 61:         VecSetValue(v,i+j,(PetscScalar)(3*i+j-2),INSERT_VALUES);
 62:       }
 63:     }
 64:     VecAssemblyBegin(v);
 65:     VecAssemblyEnd(v);
 66:     BVRestoreColumn(X,j,&v);
 67:   }
 68:   if (verbose) {
 69:     BVView(X,view);
 70:   }

 72:   /* Create BV object Y */
 73:   BVCreate(PETSC_COMM_WORLD,&Y);
 74:   PetscObjectSetName((PetscObject)Y,"Y");
 75:   BVSetSizesFromVec(Y,t,ky+1);
 76:   BVSetFromOptions(Y);
 77:   BVSetActiveColumns(Y,ly,ky);

 79:   /* Fill Y entries */
 80:   for (j=0;j<ky+1;j++) {
 81:     BVGetColumn(Y,j,&v);
 82:     VecSet(v,(PetscScalar)(j+1)/4.0);
 83:     BVRestoreColumn(Y,j,&v);
 84:   }
 85:   if (verbose) {
 86:     BVView(Y,view);
 87:   }

 89:   /* Create Mat */
 90:   MatCreateSeqDense(PETSC_COMM_SELF,kx,ky,NULL,&Q);
 91:   PetscObjectSetName((PetscObject)Q,"Q");
 92:   MatDenseGetArray(Q,&q);
 93:   for (i=0;i<kx;i++)
 94:     for (j=0;j<ky;j++)
 95:       q[i+j*kx] = (i<j)? 2.0: -0.5;
 96:   MatDenseRestoreArray(Q,&q);
 97:   if (verbose) {
 98:     MatView(Q,NULL);
 99:   }

101:   /* Test BVResize */
102:   BVResize(X,kx+4,PETSC_TRUE);

104:   /* Test BVMult */
105:   BVMult(Y,2.0,0.5,X,Q);
106:   if (verbose) {
107:     PetscPrintf(PETSC_COMM_WORLD,"After BVMult - - - - - - - - -\n");
108:     BVView(Y,view);
109:   }

111:   /* Test BVMultVec */
112:   BVGetColumn(Y,0,&v);
113:   PetscMalloc1(kx-lx,&z);
114:   z[0] = 2.0;
115:   for (i=1;i<kx-lx;i++) z[i] = -0.5*z[i-1];
116:   BVMultVec(X,-1.0,1.0,v,z);
117:   PetscFree(z);
118:   BVRestoreColumn(Y,0,&v);
119:   if (verbose) {
120:     PetscPrintf(PETSC_COMM_WORLD,"After BVMultVec - - - - - - -\n");
121:     BVView(Y,view);
122:   }

124:   /* Test BVDot */
125:   MatCreateSeqDense(PETSC_COMM_SELF,ky,kx,NULL,&M);
126:   PetscObjectSetName((PetscObject)M,"M");
127:   BVDot(X,Y,M);
128:   if (verbose) {
129:     PetscPrintf(PETSC_COMM_WORLD,"After BVDot - - - - - - - - -\n");
130:     MatView(M,NULL);
131:   }

133:   /* Test BVDotVec */
134:   BVGetColumn(Y,0,&v);
135:   PetscMalloc1(kx-lx,&z);
136:   BVDotVec(X,v,z);
137:   BVRestoreColumn(Y,0,&v);
138:   if (verbose) {
139:     PetscPrintf(PETSC_COMM_WORLD,"After BVDotVec - - - - - - -\n");
140:     VecCreateSeqWithArray(PETSC_COMM_SELF,1,kx-lx,z,&v);
141:     PetscObjectSetName((PetscObject)v,"z");
142:     VecView(v,view);
143:     VecDestroy(&v);
144:   }
145:   PetscFree(z);

147:   /* Test BVMultInPlace and BVScale */
148:   PetscOptionsHasName(NULL,NULL,"-trans",&trans);
149:   if (trans) {
150:     Mat Qt;
151:     MatTranspose(Q,MAT_INITIAL_MATRIX,&Qt);
152:     BVMultInPlaceTranspose(X,Qt,lx+1,ky);
153:     MatDestroy(&Qt);
154:   } else {
155:     BVMultInPlace(X,Q,lx+1,ky);
156:   }
157:   BVScale(X,2.0);
158:   if (verbose) {
159:     PetscPrintf(PETSC_COMM_WORLD,"After BVMultInPlace - - - - -\n");
160:     BVView(X,view);
161:   }

163:   /* Test BVNorm */
164:   BVNormColumn(X,lx,NORM_2,&nrm);
165:   PetscPrintf(PETSC_COMM_WORLD,"2-Norm of X[%D] = %g\n",lx,(double)nrm);
166:   BVNorm(X,NORM_FROBENIUS,&nrm);
167:   PetscPrintf(PETSC_COMM_WORLD,"Frobenius Norm of X = %g\n",(double)nrm);

169:   BVDestroy(&X);
170:   BVDestroy(&Y);
171:   MatDestroy(&Q);
172:   MatDestroy(&M);
173:   VecDestroy(&t);
174:   SlepcFinalize();
175:   return ierr;
176: }

178: /*TEST

180:    test:
181:       suffix: 1
182:       nsize: 1
183:       args: -n 18 -kx 12 -ky 8 -bv_type {{vecs contiguous svec mat}shared output}
184:       output_file: output/test4_1.out

186:    test:
187:       suffix: 1_vecs_vmip
188:       nsize: 1
189:       args: -n 18 -kx 12 -ky 8 -bv_type vecs -bv_vecs_vmip 1
190:       output_file: output/test4_1.out

192:    test:
193:       suffix: 1_cuda
194:       nsize: 1
195:       args: -n 18 -kx 12 -ky 8 -bv_type svec -vec_type cuda
196:       requires: cuda
197:       output_file: output/test4_1.out

199:    test:
200:       suffix: 2
201:       nsize: 1
202:       args: -n 18 -kx 12 -ky 8 -bv_type {{vecs contiguous svec mat}shared output} -trans
203:       output_file: output/test4_1.out

205: TEST*/