Actual source code: test11.c

slepc-3.18.3 2023-03-24
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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 block orthogonalization.\n\n";

 13: #include <slepcbv.h>

 15: /*
 16:    Compute the Frobenius norm ||A(l:k,l:k)-diag||_F
 17:  */
 18: PetscErrorCode MyMatNorm(Mat A,PetscInt lda,PetscInt l,PetscInt k,PetscScalar diag,PetscReal *norm)
 19: {
 20:   PetscInt          i,j;
 21:   const PetscScalar *pA;
 22:   PetscReal         s,val;

 25:   MatDenseGetArrayRead(A,&pA);
 26:   s = 0.0;
 27:   for (i=l;i<k;i++) {
 28:     for (j=l;j<k;j++) {
 29:       val = (i==j)? PetscAbsScalar(pA[i+j*lda]-diag): PetscAbsScalar(pA[i+j*lda]);
 30:       s += val*val;
 31:     }
 32:   }
 33:   *norm = PetscSqrtReal(s);
 34:   MatDenseRestoreArrayRead(A,&pA);
 35:   return 0;
 36: }

 38: int main(int argc,char **argv)
 39: {
 40:   BV             X,Y,Z,cached;
 41:   Mat            B=NULL,M,R=NULL;
 42:   Vec            v,t;
 43:   PetscInt       i,j,n=20,l=2,k=8,Istart,Iend;
 44:   PetscViewer    view;
 45:   PetscBool      withb,resid,rand,verbose;
 46:   PetscReal      norm;
 47:   PetscScalar    alpha;

 50:   SlepcInitialize(&argc,&argv,(char*)0,help);
 51:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 52:   PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL);
 53:   PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
 54:   PetscOptionsHasName(NULL,NULL,"-withb",&withb);
 55:   PetscOptionsHasName(NULL,NULL,"-resid",&resid);
 56:   PetscOptionsHasName(NULL,NULL,"-rand",&rand);
 57:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 58:   PetscPrintf(PETSC_COMM_WORLD,"Test BV block orthogonalization (length %" PetscInt_FMT ", l=%" PetscInt_FMT ", k=%" PetscInt_FMT ")%s.\n",n,l,k,withb?" with non-standard inner product":"");

 60:   /* Create template vector */
 61:   VecCreate(PETSC_COMM_WORLD,&t);
 62:   VecSetSizes(t,PETSC_DECIDE,n);
 63:   VecSetFromOptions(t);

 65:   /* Create BV object X */
 66:   BVCreate(PETSC_COMM_WORLD,&X);
 67:   PetscObjectSetName((PetscObject)X,"X");
 68:   BVSetSizesFromVec(X,t,k);
 69:   BVSetFromOptions(X);

 71:   /* Set up viewer */
 72:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
 73:   if (verbose) PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);

 75:   /* Fill X entries */
 76:   if (rand) BVSetRandom(X);
 77:   else {
 78:     for (j=0;j<k;j++) {
 79:       BVGetColumn(X,j,&v);
 80:       VecSet(v,0.0);
 81:       for (i=0;i<=n/2;i++) {
 82:         if (i+j<n) {
 83:           alpha = (3.0*i+j-2)/(2*(i+j+1));
 84:           VecSetValue(v,i+j,alpha,INSERT_VALUES);
 85:         }
 86:       }
 87:       VecAssemblyBegin(v);
 88:       VecAssemblyEnd(v);
 89:       BVRestoreColumn(X,j,&v);
 90:     }
 91:   }
 92:   if (verbose) BVView(X,view);

 94:   if (withb) {
 95:     /* Create inner product matrix */
 96:     MatCreate(PETSC_COMM_WORLD,&B);
 97:     MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
 98:     MatSetFromOptions(B);
 99:     MatSetUp(B);
100:     PetscObjectSetName((PetscObject)B,"B");

102:     MatGetOwnershipRange(B,&Istart,&Iend);
103:     for (i=Istart;i<Iend;i++) {
104:       if (i>0) MatSetValue(B,i,i-1,-1.0,INSERT_VALUES);
105:       if (i<n-1) MatSetValue(B,i,i+1,-1.0,INSERT_VALUES);
106:       MatSetValue(B,i,i,2.0,INSERT_VALUES);
107:     }
108:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
109:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
110:     if (verbose) MatView(B,view);
111:     BVSetMatrix(X,B,PETSC_FALSE);
112:   }

114:   /* Create copy on Y */
115:   BVDuplicate(X,&Y);
116:   PetscObjectSetName((PetscObject)Y,"Y");
117:   BVCopy(X,Y);
118:   MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M);

120:   if (resid) {
121:     /* Create matrix R to store triangular factor */
122:     MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&R);
123:     PetscObjectSetName((PetscObject)R,"R");
124:   }

126:   if (l>0) {
127:     /* First orthogonalize leading columns */
128:     PetscPrintf(PETSC_COMM_WORLD,"Orthogonalizing leading columns\n");
129:     BVSetActiveColumns(Y,0,l);
130:     BVSetActiveColumns(X,0,l);
131:     BVOrthogonalize(Y,R);
132:     if (verbose) {
133:       BVView(Y,view);
134:       if (resid) MatView(R,view);
135:     }

137:     if (withb) {
138:       /* Extract cached BV and check it is equal to B*X */
139:       BVGetCachedBV(Y,&cached);
140:       BVDuplicate(X,&Z);
141:       BVSetMatrix(Z,NULL,PETSC_FALSE);
142:       BVSetActiveColumns(Z,0,l);
143:       BVCopy(X,Z);
144:       BVMatMult(X,B,Z);
145:       BVMult(Z,-1.0,1.0,cached,NULL);
146:       BVNorm(Z,NORM_FROBENIUS,&norm);
147:       if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"  Difference ||cached-BX|| < 100*eps\n");
148:       else PetscPrintf(PETSC_COMM_WORLD,"  Difference ||cached-BX||: %g\n",(double)norm);
149:       BVDestroy(&Z);
150:     }

152:     /* Check orthogonality */
153:     BVDot(Y,Y,M);
154:     MyMatNorm(M,k,0,l,1.0,&norm);
155:     if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"  Level of orthogonality of Q1 < 100*eps\n");
156:     else PetscPrintf(PETSC_COMM_WORLD,"  Level of orthogonality of Q1: %g\n",(double)norm);

158:     if (resid) {
159:       /* Check residual */
160:       BVDuplicate(X,&Z);
161:       BVSetMatrix(Z,NULL,PETSC_FALSE);
162:       BVSetActiveColumns(Z,0,l);
163:       BVCopy(X,Z);
164:       BVMult(Z,-1.0,1.0,Y,R);
165:       BVNorm(Z,NORM_FROBENIUS,&norm);
166:       if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"  Residual ||X1-Q1*R11|| < 100*eps\n");
167:       else PetscPrintf(PETSC_COMM_WORLD,"  Residual ||X1-Q1*R11||: %g\n",(double)norm);
168:       BVDestroy(&Z);
169:     }

171:   }

173:   /* Now orthogonalize the rest of columns */
174:   PetscPrintf(PETSC_COMM_WORLD,"Orthogonalizing active columns\n");
175:   BVSetActiveColumns(Y,l,k);
176:   BVSetActiveColumns(X,l,k);
177:   BVOrthogonalize(Y,R);
178:   if (verbose) {
179:     BVView(Y,view);
180:     if (resid) MatView(R,view);
181:   }

183:   if (l>0) {
184:     /* Check orthogonality */
185:     BVDot(Y,Y,M);
186:     MyMatNorm(M,k,l,k,1.0,&norm);
187:     if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"  Level of orthogonality of Q2 < 100*eps\n");
188:     else PetscPrintf(PETSC_COMM_WORLD,"  Level of orthogonality of Q2: %g\n",(double)norm);
189:   }

191:   /* Check the complete decomposition */
192:   PetscPrintf(PETSC_COMM_WORLD,"Overall decomposition\n");
193:   BVSetActiveColumns(Y,0,k);
194:   BVSetActiveColumns(X,0,k);

196:   /* Check orthogonality */
197:   BVDot(Y,Y,M);
198:   MyMatNorm(M,k,0,k,1.0,&norm);
199:   if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"  Level of orthogonality of Q < 100*eps\n");
200:   else PetscPrintf(PETSC_COMM_WORLD,"  Level of orthogonality of Q: %g\n",(double)norm);

202:   if (resid) {
203:     /* Check residual */
204:     BVMult(X,-1.0,1.0,Y,R);
205:     BVSetMatrix(X,NULL,PETSC_FALSE);
206:     BVNorm(X,NORM_FROBENIUS,&norm);
207:     if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"  Residual ||X-Q*R|| < 100*eps\n");
208:     else PetscPrintf(PETSC_COMM_WORLD,"  Residual ||X-Q*R||: %g\n",(double)norm);
209:     MatDestroy(&R);
210:   }

212:   if (B) MatDestroy(&B);
213:   MatDestroy(&M);
214:   BVDestroy(&X);
215:   BVDestroy(&Y);
216:   VecDestroy(&t);
217:   SlepcFinalize();
218:   return 0;
219: }

221: /*TEST

223:    testset:
224:       args: -bv_orthog_block {{gs chol tsqr tsqrchol svqb}}
225:       nsize: 2
226:       output_file: output/test11_1.out
227:       test:
228:          suffix: 1
229:          args: -bv_type {{vecs contiguous svec mat}shared output}
230:       test:
231:          suffix: 1_cuda
232:          args: -bv_type svec -vec_type cuda
233:          requires: cuda

235:    testset:
236:       args: -withb -bv_orthog_block {{gs chol svqb}}
237:       nsize: 2
238:       output_file: output/test11_4.out
239:       test:
240:          suffix: 4
241:          args: -bv_type {{vecs contiguous svec mat}shared output}
242:       test:
243:          suffix: 4_cuda
244:          args: -bv_type svec -vec_type cuda -mat_type aijcusparse
245:          requires: cuda

247:    testset:
248:       args: -resid -bv_orthog_block {{gs chol tsqr tsqrchol svqb}}
249:       nsize: 2
250:       output_file: output/test11_6.out
251:       test:
252:          suffix: 6
253:          args: -bv_type {{vecs contiguous svec mat}shared output}
254:       test:
255:          suffix: 6_cuda
256:          args: -bv_type svec -vec_type cuda
257:          requires: cuda

259:    testset:
260:       args: -resid -withb -bv_orthog_block {{gs chol svqb}}
261:       nsize: 2
262:       output_file: output/test11_9.out
263:       test:
264:          suffix: 9
265:          args: -bv_type {{vecs contiguous svec mat}shared output}
266:       test:
267:          suffix: 9_cuda
268:          args: -bv_type svec -vec_type cuda -mat_type aijcusparse
269:          requires: cuda

271:    testset:
272:       args: -bv_orthog_block tsqr
273:       nsize: 7
274:       output_file: output/test11_1.out
275:       test:
276:          suffix: 11
277:          args: -bv_type {{vecs contiguous svec mat}shared output}
278:          requires: !defined(PETSCTEST_VALGRIND)
279:       test:
280:          suffix: 11_cuda
281:          TODO: too many processes accessing the GPU
282:          args: -bv_type svec -vec_type cuda
283:          requires: cuda !defined(PETSCTEST_VALGRIND)

285:    testset:
286:       args: -resid -n 180 -l 0 -k 7 -bv_orthog_block tsqr
287:       nsize: 9
288:       output_file: output/test11_12.out
289:       test:
290:          suffix: 12
291:          args: -bv_type {{vecs contiguous svec mat}shared output}
292:          requires: !single !defined(PETSCTEST_VALGRIND)
293:       test:
294:          suffix: 12_cuda
295:          TODO: too many processes accessing the GPU
296:          args: -bv_type svec -vec_type cuda
297:          requires: cuda !single !defined(PETSCTEST_VALGRIND)

299: TEST*/