Actual source code: ex76.c

petsc-3.10.5 2019-03-28
Report Typos and Errors

  2: static char help[] = "Tests cholesky, icc factorization and solve on sequential aij, baij and sbaij matrices. \n";

  4:  #include <petscmat.h>

  6: int main(int argc,char **args)
  7: {
  8:   Vec            x,y,b;
  9:   Mat            A;           /* linear system matrix */
 10:   Mat            sA,sC;       /* symmetric part of the matrices */
 11:   PetscInt       n,mbs=16,bs=1,nz=3,prob=1,i,j,col[3],block, row,Ii,J,n1,lvl;
 13:   PetscMPIInt    size;
 14:   PetscReal      norm2;
 15:   PetscScalar    neg_one = -1.0,four=4.0,value[3];
 16:   IS             perm,cperm;
 17:   PetscRandom    rdm;
 18:   PetscBool      reorder = PETSC_FALSE,displ = PETSC_FALSE;
 19:   MatFactorInfo  factinfo;
 20:   PetscBool      equal;
 21:   PetscBool      TestAIJ = PETSC_FALSE,TestBAIJ = PETSC_TRUE;
 22:   PetscInt       TestShift=0;

 24:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
 25:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 26:   if (size != 1) SETERRQ(PETSC_COMM_WORLD,1,"This is a uniprocessor example only!");
 27:   PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);
 28:   PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL);
 29:   PetscOptionsGetBool(NULL,NULL,"-reorder",&reorder,NULL);
 30:   PetscOptionsGetBool(NULL,NULL,"-testaij",&TestAIJ,NULL);
 31:   PetscOptionsGetInt(NULL,NULL,"-testShift",&TestShift,NULL);
 32:   PetscOptionsGetBool(NULL,NULL,"-displ",&displ,NULL);

 34:   n = mbs*bs;
 35:   if (TestAIJ) { /* A is in aij format */
 36:     MatCreateSeqAIJ(PETSC_COMM_WORLD,n,n,nz,NULL,&A);
 37:     TestBAIJ = PETSC_FALSE;
 38:   } else { /* A is in baij format */
 39:     ierr    =MatCreateSeqBAIJ(PETSC_COMM_WORLD,bs,n,n,nz,NULL,&A);
 40:     TestAIJ = PETSC_FALSE;
 41:   }

 43:   /* Assemble matrix */
 44:   if (bs == 1) {
 45:     PetscOptionsGetInt(NULL,NULL,"-test_problem",&prob,NULL);
 46:     if (prob == 1) { /* tridiagonal matrix */
 47:       value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
 48:       for (i=1; i<n-1; i++) {
 49:         col[0] = i-1; col[1] = i; col[2] = i+1;
 50:         MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
 51:       }
 52:       i = n - 1; col[0]=0; col[1] = n - 2; col[2] = n - 1;

 54:       value[0]= 0.1; value[1]=-1; value[2]=2;
 55:       MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);

 57:       i = 0; col[0] = 0; col[1] = 1; col[2]=n-1;

 59:       value[0] = 2.0; value[1] = -1.0; value[2]=0.1;
 60:       MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
 61:     } else if (prob ==2) { /* matrix for the five point stencil */
 62:       n1 = (PetscInt) (PetscSqrtReal((PetscReal)n) + 0.001);
 63:       if (n1*n1 - n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"sqrt(n) must be a positive integer!");
 64:       for (i=0; i<n1; i++) {
 65:         for (j=0; j<n1; j++) {
 66:           Ii = j + n1*i;
 67:           if (i>0) {
 68:             J    = Ii - n1;
 69:             MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);
 70:           }
 71:           if (i<n1-1) {
 72:             J    = Ii + n1;
 73:             MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);
 74:           }
 75:           if (j>0) {
 76:             J    = Ii - 1;
 77:             MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);
 78:           }
 79:           if (j<n1-1) {
 80:             J    = Ii + 1;
 81:             MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);
 82:           }
 83:           MatSetValues(A,1,&Ii,1,&Ii,&four,INSERT_VALUES);
 84:         }
 85:       }
 86:     }
 87:   } else { /* bs > 1 */
 88:     for (block=0; block<n/bs; block++) {
 89:       /* diagonal blocks */
 90:       value[0] = -1.0; value[1] = 4.0; value[2] = -1.0;
 91:       for (i=1+block*bs; i<bs-1+block*bs; i++) {
 92:         col[0] = i-1; col[1] = i; col[2] = i+1;
 93:         MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
 94:       }
 95:       i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs;

 97:       value[0]=-1.0; value[1]=4.0;
 98:       MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);

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

102:       value[0]=4.0; value[1] = -1.0;
103:       MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
104:     }
105:     /* off-diagonal blocks */
106:     value[0]=-1.0;
107:     for (i=0; i<(n/bs-1)*bs; i++) {
108:       col[0]=i+bs;
109:       MatSetValues(A,1,&i,1,col,value,INSERT_VALUES);
110:       col[0]=i; row=i+bs;
111:       MatSetValues(A,1,&row,1,col,value,INSERT_VALUES);
112:     }
113:   }

115:   if (TestShift) {
116:     /* set diagonals in the 0-th block as 0 for testing shift numerical factor */
117:     for (i=0; i<bs; i++) {
118:       row  = i; col[0] = i; value[0] = 0.0;
119:       MatSetValues(A,1,&row,1,col,value,INSERT_VALUES);
120:     }
121:   }

123:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
124:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

126:   /* Test MatConvert */
127:   MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE);
128:   MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&sA);
129:   MatMultEqual(A,sA,20,&equal);
130:   if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"A != sA");

132:   /* Test MatGetOwnershipRange() */
133:   MatGetOwnershipRange(A,&Ii,&J);
134:   MatGetOwnershipRange(sA,&i,&j);
135:   if (i-Ii || j-J) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetOwnershipRange() in MatSBAIJ format\n");

137:   /* Vectors */
138:   PetscRandomCreate(PETSC_COMM_SELF,&rdm);
139:   PetscRandomSetFromOptions(rdm);
140:   VecCreateSeq(PETSC_COMM_SELF,n,&x);
141:   VecDuplicate(x,&b);
142:   VecDuplicate(x,&y);
143:   VecSetRandom(x,rdm);

145:   /* Test MatReordering() - not work on sbaij matrix */
146:   if (reorder) {
147:     MatGetOrdering(A,MATORDERINGRCM,&perm,&cperm);
148:   } else {
149:     MatGetOrdering(A,MATORDERINGNATURAL,&perm,&cperm);
150:   }
151:   ISDestroy(&cperm);

153:   /* initialize factinfo */
154:   MatFactorInfoInitialize(&factinfo);
155:   if (TestShift == 1) {
156:     factinfo.shifttype   = (PetscReal)MAT_SHIFT_NONZERO;
157:     factinfo.shiftamount = 0.1;
158:   } else if (TestShift == 2) {
159:     factinfo.shifttype = (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE;
160:   }

162:   /* Test MatCholeskyFactor(), MatICCFactor() */
163:   /*------------------------------------------*/
164:   /* Test aij matrix A */
165:   if (TestAIJ) {
166:     if (displ) {
167:       PetscPrintf(PETSC_COMM_WORLD,"AIJ: \n");
168:     }
169:     i = 0;
170:     for (lvl=-1; lvl<10; lvl++) {
171:       if (lvl==-1) {  /* Cholesky factor */
172:         factinfo.fill = 5.0;

174:         MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&sC);
175:         MatCholeskyFactorSymbolic(sC,A,perm,&factinfo);
176:       } else {       /* incomplete Cholesky factor */
177:         factinfo.fill   = 5.0;
178:         factinfo.levels = lvl;

180:         MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_ICC,&sC);
181:         MatICCFactorSymbolic(sC,A,perm,&factinfo);
182:       }
183:       MatCholeskyFactorNumeric(sC,A,&factinfo);

185:       MatMult(A,x,b);
186:       MatSolve(sC,b,y);
187:       MatDestroy(&sC);

189:       /* Check the residual */
190:       VecAXPY(y,neg_one,x);
191:       VecNorm(y,NORM_2,&norm2);

193:       if (displ) {
194:         PetscPrintf(PETSC_COMM_WORLD,"  lvl: %D, residual: %g\n", lvl,(double)norm2);
195:       }
196:     }
197:   }

199:   /* Test baij matrix A */
200:   if (TestBAIJ) {
201:     if (displ) {
202:       PetscPrintf(PETSC_COMM_WORLD,"BAIJ: \n");
203:     }
204:     i = 0;
205:     for (lvl=-1; lvl<10; lvl++) {
206:       if (lvl==-1) {  /* Cholesky factor */
207:         factinfo.fill = 5.0;

209:         MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&sC);
210:         MatCholeskyFactorSymbolic(sC,A,perm,&factinfo);
211:       } else {       /* incomplete Cholesky factor */
212:         factinfo.fill   = 5.0;
213:         factinfo.levels = lvl;

215:         MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_ICC,&sC);
216:         MatICCFactorSymbolic(sC,A,perm,&factinfo);
217:       }
218:       MatCholeskyFactorNumeric(sC,A,&factinfo);

220:       MatMult(A,x,b);
221:       MatSolve(sC,b,y);
222:       MatDestroy(&sC);

224:       /* Check the residual */
225:       VecAXPY(y,neg_one,x);
226:       VecNorm(y,NORM_2,&norm2);
227:       if (displ) {
228:         PetscPrintf(PETSC_COMM_WORLD,"  lvl: %D, residual: %g\n", lvl,(double)norm2);
229:       }
230:     }
231:   }

233:   /* Test sbaij matrix sA */
234:   if (displ) {
235:     PetscPrintf(PETSC_COMM_WORLD,"SBAIJ: \n");
236:   }
237:   i = 0;
238:   for (lvl=-1; lvl<10; lvl++) {
239:     if (lvl==-1) {  /* Cholesky factor */
240:       factinfo.fill = 5.0;

242:       MatGetFactor(sA,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&sC);
243:       MatCholeskyFactorSymbolic(sC,sA,perm,&factinfo);
244:     } else {       /* incomplete Cholesky factor */
245:       factinfo.fill   = 5.0;
246:       factinfo.levels = lvl;

248:       MatGetFactor(sA,MATSOLVERPETSC,MAT_FACTOR_ICC,&sC);
249:       MatICCFactorSymbolic(sC,sA,perm,&factinfo);
250:     }
251:     MatCholeskyFactorNumeric(sC,sA,&factinfo);

253:     if (lvl==0 && bs==1) { /* Test inplace ICC(0) for sbaij sA - does not work for new datastructure */
254:       /*
255:         Mat B;
256:         MatDuplicate(sA,MAT_COPY_VALUES,&B);
257:         MatICCFactor(B,perm,&factinfo);
258:         MatEqual(sC,B,&equal);
259:         if (!equal) {
260:           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"in-place Cholesky factor != out-place Cholesky factor");
261:         }
262:         MatDestroy(&B);
263:       */
264:     }


267:     MatMult(sA,x,b);
268:     MatSolve(sC,b,y);

270:     /* Test MatSolves() */
271:     if (bs == 1) {
272:       Vecs xx,bb;
273:       VecsCreateSeq(PETSC_COMM_SELF,n,4,&xx);
274:       VecsDuplicate(xx,&bb);
275:       MatSolves(sC,bb,xx);
276:       VecsDestroy(xx);
277:       VecsDestroy(bb);
278:     }
279:     MatDestroy(&sC);

281:     /* Check the residual */
282:     VecAXPY(y,neg_one,x);
283:     VecNorm(y,NORM_2,&norm2);
284:     if (displ) {
285:       PetscPrintf(PETSC_COMM_WORLD,"  lvl: %D, residual: %g\n", lvl,(double)norm2);
286:     }
287:   }

289:   ISDestroy(&perm);
290:   MatDestroy(&A);
291:   MatDestroy(&sA);
292:   VecDestroy(&x);
293:   VecDestroy(&y);
294:   VecDestroy(&b);
295:   PetscRandomDestroy(&rdm);

297:   PetscFinalize();
298:   return ierr;
299: }

301: /*TEST

303:    test:
304:       args: -bs {{1 2 3 4 5 6 7 8}}

306:    test:
307:       suffix: 3
308:       args: -testaij
309:       output_file: output/ex76_1.out

311: TEST*/