Actual source code: ex74.c

petsc-3.8.4 2018-03-24
Report Typos and Errors

  2: static char help[] = "Tests the various sequential routines in MatSBAIJ format.\n";

  4:  #include <petscmat.h>

  6: int main(int argc,char **args)
  7: {
  8:   PetscMPIInt    size;
 10:   Vec            x,y,b,s1,s2;
 11:   Mat            A;                    /* linear system matrix */
 12:   Mat            sA,sB,sFactor;        /* symmetric matrices */
 13:   PetscInt       n,mbs=16,bs=1,nz=3,prob=1,i,j,k1,k2,col[3],lf,block, row,Ii,J,n1,inc;
 14:   PetscReal      norm1,norm2,rnorm,tol=PETSC_SMALL;
 15:   PetscScalar    neg_one = -1.0,four=4.0,value[3];
 16:   IS             perm, iscol;
 17:   PetscRandom    rdm;
 18:   PetscBool      doIcc=PETSC_TRUE,equal;
 19:   MatInfo        minfo1,minfo2;
 20:   MatFactorInfo  factinfo;
 21:   MatType        type;

 23:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
 24:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 25:   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
 26:   PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);
 27:   PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL);

 29:   n    = mbs*bs;
 30:   MatCreate(PETSC_COMM_SELF,&A);
 31:   MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
 32:   MatSetType(A,MATSEQBAIJ);
 33:   MatSetFromOptions(A);
 34:   MatSeqBAIJSetPreallocation(A,bs,nz,NULL);

 36:   MatCreate(PETSC_COMM_SELF,&sA);
 37:   MatSetSizes(sA,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
 38:   MatSetType(sA,MATSEQSBAIJ);
 39:   MatSetFromOptions(sA);
 40:   MatGetType(sA,&type);
 41:   PetscObjectTypeCompare((PetscObject)sA,MATSEQSBAIJ,&doIcc);
 42:   MatSeqSBAIJSetPreallocation(sA,bs,nz,NULL);
 43:   MatSetOption(sA,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);

 45:   /* Test MatGetOwnershipRange() */
 46:   MatGetOwnershipRange(A,&Ii,&J);
 47:   MatGetOwnershipRange(sA,&i,&j);
 48:   if (i-Ii || j-J) {
 49:     PetscPrintf(PETSC_COMM_SELF,"Error: MatGetOwnershipRange() in MatSBAIJ format\n");
 50:   }

 52:   /* Assemble matrix */
 53:   if (bs == 1) {
 54:     PetscOptionsGetInt(NULL,NULL,"-test_problem",&prob,NULL);
 55:     if (prob == 1) { /* tridiagonal matrix */
 56:       value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
 57:       for (i=1; i<n-1; i++) {
 58:         col[0] = i-1; col[1] = i; col[2] = i+1;
 59:         MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
 60:         MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
 61:       }
 62:       i = n - 1; col[0]=0; col[1] = n - 2; col[2] = n - 1;

 64:       value[0]= 0.1; value[1]=-1; value[2]=2;

 66:       MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
 67:       MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);

 69:       i        = 0;
 70:       col[0]   = n-1;   col[1] = 1;      col[2] = 0;
 71:       value[0] = 0.1; value[1] = -1.0; value[2] = 2;

 73:       MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
 74:       MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);

 76:     } else if (prob ==2) { /* matrix for the five point stencil */
 77:       n1 = (PetscInt) (PetscSqrtReal((PetscReal)n) + 0.001);
 78:       if (n1*n1 - n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"sqrt(n) must be a positive integer!");
 79:       for (i=0; i<n1; i++) {
 80:         for (j=0; j<n1; j++) {
 81:           Ii = j + n1*i;
 82:           if (i>0) {
 83:             J    = Ii - n1;
 84:             MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);
 85:             MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES);
 86:           }
 87:           if (i<n1-1) {
 88:             J    = Ii + n1;
 89:             MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);
 90:             MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES);
 91:           }
 92:           if (j>0) {
 93:             J    = Ii - 1;
 94:             MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);
 95:             MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES);
 96:           }
 97:           if (j<n1-1) {
 98:             J    = Ii + 1;
 99:             MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);
100:             MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES);
101:           }
102:           MatSetValues(A,1,&Ii,1,&Ii,&four,INSERT_VALUES);
103:           MatSetValues(sA,1,&Ii,1,&Ii,&four,INSERT_VALUES);
104:         }
105:       }
106:     }

108:   } else { /* bs > 1 */
109:     for (block=0; block<n/bs; block++) {
110:       /* diagonal blocks */
111:       value[0] = -1.0; value[1] = 4.0; value[2] = -1.0;
112:       for (i=1+block*bs; i<bs-1+block*bs; i++) {
113:         col[0] = i-1; col[1] = i; col[2] = i+1;
114:         MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
115:         MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
116:       }
117:       i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs;

119:       value[0]=-1.0; value[1]=4.0;

121:       MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
122:       MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);

124:       i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs;

126:       value[0]=4.0; value[1] = -1.0;

128:       MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
129:       MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);
130:     }
131:     /* off-diagonal blocks */
132:     value[0]=-1.0;
133:     for (i=0; i<(n/bs-1)*bs; i++) {
134:       col[0]=i+bs;

136:       MatSetValues(A,1,&i,1,col,value,INSERT_VALUES);
137:       MatSetValues(sA,1,&i,1,col,value,INSERT_VALUES);

139:       col[0]=i; row=i+bs;

141:       MatSetValues(A,1,&row,1,col,value,INSERT_VALUES);
142:       MatSetValues(sA,1,&row,1,col,value,INSERT_VALUES);
143:     }
144:   }
145:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
146:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

148:   MatAssemblyBegin(sA,MAT_FINAL_ASSEMBLY);
149:   MatAssemblyEnd(sA,MAT_FINAL_ASSEMBLY);

151:   /* Test MatGetInfo() of A and sA */
152:   MatGetInfo(A,MAT_LOCAL,&minfo1);
153:   MatGetInfo(sA,MAT_LOCAL,&minfo2);
154:   /*
155:   printf("A matrix nonzeros (BAIJ format) = %d, allocated nonzeros= %d\n", (int)minfo1.nz_used,(int)minfo1.nz_allocated);
156:   printf("sA matrix nonzeros(SBAIJ format) = %d, allocated nonzeros= %d\n", (int)minfo2.nz_used,(int)minfo2.nz_allocated);
157:   */
158:   i  = (int) (minfo1.nz_used - minfo2.nz_used);
159:   j  = (int) (minfo1.nz_allocated - minfo2.nz_allocated);
160:   k1 = (int) (minfo1.nz_allocated - minfo1.nz_used);
161:   k2 = (int) (minfo2.nz_allocated - minfo2.nz_used);
162:   if (i < 0 || j < 0 || k1 < 0 || k2 < 0) {
163:     PetscPrintf(PETSC_COMM_SELF,"Error (compare A and sA): MatGetInfo()\n");
164:   }

166:   /* Test MatDuplicate() */
167:   MatNorm(A,NORM_FROBENIUS,&norm1);
168:   MatDuplicate(sA,MAT_COPY_VALUES,&sB);
169:   MatEqual(sA,sB,&equal);
170:   if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatDuplicate()");

172:   /* Test MatNorm() */
173:   MatNorm(A,NORM_FROBENIUS,&norm1);
174:   MatNorm(sB,NORM_FROBENIUS,&norm2);
175:   rnorm = PetscAbsReal(norm1-norm2)/norm2;
176:   if (rnorm > tol) {
177:     PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_FROBENIUS, NormA=%16.14e NormsB=%16.14e\n",norm1,norm2);
178:   }
179:   MatNorm(A,NORM_INFINITY,&norm1);
180:   MatNorm(sB,NORM_INFINITY,&norm2);
181:   rnorm = PetscAbsReal(norm1-norm2)/norm2;
182:   if (rnorm > tol) {
183:     PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_INFINITY(), NormA=%16.14e NormsB=%16.14e\n",norm1,norm2);
184:   }
185:   MatNorm(A,NORM_1,&norm1);
186:   MatNorm(sB,NORM_1,&norm2);
187:   rnorm = PetscAbsReal(norm1-norm2)/norm2;
188:   if (rnorm > tol) {
189:     PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_INFINITY(), NormA=%16.14e NormsB=%16.14e\n",norm1,norm2);
190:   }

192:   /* Test MatGetInfo(), MatGetSize(), MatGetBlockSize() */
193:   MatGetInfo(A,MAT_LOCAL,&minfo1);
194:   MatGetInfo(sB,MAT_LOCAL,&minfo2);
195:   /*
196:   printf("matrix nonzeros (BAIJ format) = %d, allocated nonzeros= %d\n", (int)minfo1.nz_used,(int)minfo1.nz_allocated);
197:   printf("matrix nonzeros(SBAIJ format) = %d, allocated nonzeros= %d\n", (int)minfo2.nz_used,(int)minfo2.nz_allocated);
198:   */
199:   i  = (int) (minfo1.nz_used - minfo2.nz_used);
200:   j  = (int) (minfo1.nz_allocated - minfo2.nz_allocated);
201:   k1 = (int) (minfo1.nz_allocated - minfo1.nz_used);
202:   k2 = (int) (minfo2.nz_allocated - minfo2.nz_used);
203:   if (i < 0 || j < 0 || k1 < 0 || k2 < 0) {
204:     PetscPrintf(PETSC_COMM_SELF,"Error(compare A and sB): MatGetInfo()\n");
205:   }

207:   MatGetSize(A,&Ii,&J);
208:   MatGetSize(sB,&i,&j);
209:   if (i-Ii || j-J) {
210:     PetscPrintf(PETSC_COMM_SELF,"Error: MatGetSize()\n");
211:   }

213:   MatGetBlockSize(A, &Ii);
214:   MatGetBlockSize(sB, &i);
215:   if (i-Ii) {
216:     PetscPrintf(PETSC_COMM_SELF,"Error: MatGetBlockSize()\n");
217:   }

219:   PetscRandomCreate(PETSC_COMM_SELF,&rdm);
220:   PetscRandomSetFromOptions(rdm);
221:   VecCreateSeq(PETSC_COMM_SELF,n,&x);
222:   VecDuplicate(x,&s1);
223:   VecDuplicate(x,&s2);
224:   VecDuplicate(x,&y);
225:   VecDuplicate(x,&b);
226:   VecSetRandom(x,rdm);

228:   /* Test MatDiagonalScale(), MatGetDiagonal(), MatScale() */
229: #if !defined(PETSC_USE_COMPLEX)
230:   /* Scaling matrix with complex numbers results non-spd matrix,
231:      causing crash of MatForwardSolve() and MatBackwardSolve() */
232:   MatDiagonalScale(A,x,x);
233:   MatDiagonalScale(sB,x,x);
234:   MatMultEqual(A,sB,10,&equal);
235:   if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatDiagonalScale");

237:   MatGetDiagonal(A,s1);
238:   MatGetDiagonal(sB,s2);
239:   VecAXPY(s2,neg_one,s1);
240:   VecNorm(s2,NORM_1,&norm1);
241:   if (norm1>tol) {
242:     PetscPrintf(PETSC_COMM_SELF,"Error:MatGetDiagonal(), ||s1-s2||=%g\n",(double)norm1);
243:   }

245:   {
246:     PetscScalar alpha=0.1;
247:     MatScale(A,alpha);
248:     MatScale(sB,alpha);
249:   }
250: #endif

252:   /* Test MatGetRowMaxAbs() */
253:   MatGetRowMaxAbs(A,s1,NULL);
254:   MatGetRowMaxAbs(sB,s2,NULL);
255:   VecNorm(s1,NORM_1,&norm1);
256:   VecNorm(s2,NORM_1,&norm2);
257:   norm1 -= norm2;
258:   if (norm1<-tol || norm1>tol) {
259:     PetscPrintf(PETSC_COMM_SELF,"Error:MatGetRowMaxAbs() \n");
260:   }

262:   /* Test MatMult() */
263:   for (i=0; i<40; i++) {
264:     VecSetRandom(x,rdm);
265:     MatMult(A,x,s1);
266:     MatMult(sB,x,s2);
267:     VecNorm(s1,NORM_1,&norm1);
268:     VecNorm(s2,NORM_1,&norm2);
269:     norm1 -= norm2;
270:     if (norm1<-tol || norm1>tol) {
271:       PetscPrintf(PETSC_COMM_SELF,"Error: MatMult(), norm1-norm2: %g\n",(double)norm1);
272:     }
273:   }

275:   /* MatMultAdd() */
276:   for (i=0; i<40; i++) {
277:     VecSetRandom(x,rdm);
278:     VecSetRandom(y,rdm);
279:     MatMultAdd(A,x,y,s1);
280:     MatMultAdd(sB,x,y,s2);
281:     VecNorm(s1,NORM_1,&norm1);
282:     VecNorm(s2,NORM_1,&norm2);
283:     norm1 -= norm2;
284:     if (norm1<-tol || norm1>tol) {
285:       PetscPrintf(PETSC_COMM_SELF,"Error:MatMultAdd(),  norm1-norm2: %g\n",(double)norm1);
286:     }
287:   }

289:   /* Test MatCholeskyFactor(), MatICCFactor() with natural ordering */
290:   MatGetOrdering(A,MATORDERINGNATURAL,&perm,&iscol);
291:   ISDestroy(&iscol);
292:   norm1 = tol;
293:   inc   = bs;

295:   /* initialize factinfo */
296:   PetscMemzero(&factinfo,sizeof(MatFactorInfo));

298:   for (lf=-1; lf<10; lf += inc) {
299:     if (lf==-1) {  /* Cholesky factor of sB (duplicate sA) */
300:       factinfo.fill = 5.0;

302:       MatGetFactor(sB,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&sFactor);
303:       MatCholeskyFactorSymbolic(sFactor,sB,perm,&factinfo);
304:     } else if (!doIcc) break;
305:     else {       /* incomplete Cholesky factor */
306:       factinfo.fill   = 5.0;
307:       factinfo.levels = lf;

309:       MatGetFactor(sB,MATSOLVERPETSC,MAT_FACTOR_ICC,&sFactor);
310:       MatICCFactorSymbolic(sFactor,sB,perm,&factinfo);
311:     }
312:     MatCholeskyFactorNumeric(sFactor,sB,&factinfo);
313:     /* MatView(sFactor, PETSC_VIEWER_DRAW_WORLD); */

315:     /* test MatGetDiagonal on numeric factor */
316:     /*
317:     if (lf == -1) {
318:       MatGetDiagonal(sFactor,s1);
319:       printf(" in ex74.c, diag: \n");
320:       VecView(s1,PETSC_VIEWER_STDOUT_SELF);
321:     }
322:     */

324:     MatMult(sB,x,b);

326:     /* test MatForwardSolve() and MatBackwardSolve() */
327:     if (lf == -1) {
328:       MatForwardSolve(sFactor,b,s1);
329:       MatBackwardSolve(sFactor,s1,s2);
330:       VecAXPY(s2,neg_one,x);
331:       VecNorm(s2,NORM_2,&norm2);
332:       if (10*norm1 < norm2) {
333:         PetscPrintf(PETSC_COMM_SELF,"MatForwardSolve and BackwardSolve: Norm of error=%g, bs=%D\n",(double)norm2,bs);
334:       }
335:     }

337:     /* test MatSolve() */
338:     MatSolve(sFactor,b,y);
339:     MatDestroy(&sFactor);
340:     /* Check the error */
341:     VecAXPY(y,neg_one,x);
342:     VecNorm(y,NORM_2,&norm2);
343:     if (10*norm1 < norm2 && lf-inc != -1) {
344:       PetscPrintf(PETSC_COMM_SELF,"lf=%D, %D, Norm of error=%g, %g\n",lf-inc,lf,(double)norm1,(double)norm2);
345:     }
346:     norm1 = norm2;
347:     if (norm2 < tol && lf != -1) break;
348:   }

350:   ISDestroy(&perm);

352:   MatDestroy(&A);
353:   MatDestroy(&sB);
354:   MatDestroy(&sA);
355:   VecDestroy(&x);
356:   VecDestroy(&y);
357:   VecDestroy(&s1);
358:   VecDestroy(&s2);
359:   VecDestroy(&b);
360:   PetscRandomDestroy(&rdm);

362:   PetscFinalize();
363:   return ierr;
364: }