Actual source code: ex76.c

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

 4:  #include petscmat.h

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

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

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

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

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

 68:       i = 0; col[0] = 0; col[1] = 1; col[2]=n-1;
 69:       value[0] = 2.0; value[1] = -1.0; value[2]=0.1;
 70:       MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
 71:       MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
 72:     } else if (prob ==2){ /* matrix for the five point stencil */
 73:       n1 = (int) (sqrt((PetscReal)n) + 0.001);
 74:       if (n1*n1 - n) SETERRQ(PETSC_ERR_ARG_WRONG,"sqrt(n) must be a positive interger!");
 75:       for (i=0; i<n1; i++) {
 76:         for (j=0; j<n1; j++) {
 77:           Ii = j + n1*i;
 78:           if (i>0)   {
 79:             J = Ii - n1;
 80:             MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);
 81:             MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES);
 82:           }
 83:           if (i<n1-1) {
 84:             J = Ii + n1;
 85:             MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);
 86:             MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES);
 87:           }
 88:           if (j>0)   {
 89:             J = Ii - 1;
 90:             MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);
 91:             MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES);
 92:           }
 93:           if (j<n1-1) {
 94:             J = Ii + 1;
 95:             MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);
 96:             MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES);
 97:           }
 98:           MatSetValues(A,1,&Ii,1,&Ii,&four,INSERT_VALUES);
 99:           MatSetValues(sA,1,&Ii,1,&Ii,&four,INSERT_VALUES);
100:         }
101:       }
102:     }
103:   } else { /* bs > 1 */
104:     for (block=0; block<n/bs; block++){
105:       /* diagonal blocks */
106:       value[0] = -1.0; value[1] = 4.0; value[2] = -1.0;
107:       for (i=1+block*bs; i<bs-1+block*bs; i++) {
108:         col[0] = i-1; col[1] = i; col[2] = i+1;
109:         MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
110:         MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
111:       }
112:       i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs;
113:       value[0]=-1.0; value[1]=4.0;
114:       MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
115:       MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);

117:       i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs;
118:       value[0]=4.0; value[1] = -1.0;
119:       MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
120:       MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);
121:     }
122:     /* off-diagonal blocks */
123:     value[0]=-1.0;
124:     for (i=0; i<(n/bs-1)*bs; i++){
125:       col[0]=i+bs;
126:       MatSetValues(A,1,&i,1,col,value,INSERT_VALUES);
127:       MatSetValues(sA,1,&i,1,col,value,INSERT_VALUES);
128:       col[0]=i; row=i+bs;
129:       MatSetValues(A,1,&row,1,col,value,INSERT_VALUES);
130:       MatSetValues(sA,1,&row,1,col,value,INSERT_VALUES);
131:     }
132:   }

134:   if (TestShift){
135:      /* set diagonals in the 0-th block as 0 for testing shift numerical factor */
136:      for (i=0; i<bs; i++){
137:        row = i; col[0] = i; value[0] = 0.0;
138:        MatSetValues(A,1,&row,1,col,value,INSERT_VALUES);
139:        MatSetValues(sA,1,&row,1,col,value,INSERT_VALUES);
140:      }
141:    }

143:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
144:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
145: 
146:   MatAssemblyBegin(sA,MAT_FINAL_ASSEMBLY);
147:   MatAssemblyEnd(sA,MAT_FINAL_ASSEMBLY);

149:   MatMultEqual(A,sA,5,&equal);
150:   if (!equal) SETERRQ(PETSC_ERR_USER,"A != sA");

152:   /* Vectors */
153:   PetscRandomCreate(PETSC_COMM_SELF,&rdm);
154:   PetscRandomSetFromOptions(rdm);
155:   VecCreateSeq(PETSC_COMM_SELF,n,&x);
156:   VecDuplicate(x,&b);
157:   VecDuplicate(x,&y);
158:   VecSetRandom(x,rdm);

160:   /* Test MatReordering() on a symmetric ordering */
161:   PetscMalloc(mbs*sizeof(PetscInt),&ip_ptr);
162:   for (i=0; i<mbs; i++) ip_ptr[i] = i;
163:   switch (reorder){
164:   case 0: break;
165:   case 1:
166:     i = ip_ptr[2]; ip_ptr[2] = ip_ptr[mbs-3]; ip_ptr[mbs-3] = i;
167:   case 2:
168:     i = ip_ptr[1]; ip_ptr[1] = ip_ptr[mbs-2]; ip_ptr[mbs-2] = i;
169:   case 3:
170:     i = ip_ptr[0]; ip_ptr[0] = ip_ptr[mbs-1]; ip_ptr[mbs-1] = i;
171:   }

173:   ISCreateGeneral(PETSC_COMM_SELF,mbs,ip_ptr,&perm);
174:   ISSetPermutation(perm);

176:   /* initialize factinfo */
177:   factinfo.shiftnz   = 0.0;
178:   factinfo.shiftpd   = 0.0; /* false */
179:   factinfo.zeropivot = 1.e-12;
180:   if (TestShift == 1){
181:     factinfo.shiftnz = 0.1;
182:   } else if (TestShift == 2){
183:     factinfo.shiftpd = 1.0; /* true*/
184:   }
185: 
186:   /* Test MatCholeskyFactor(), MatICCFactor() */
187:   /*------------------------------------------*/
188:   /* Test aij matrix A */
189:   if (TestAIJ){
190:     if (displ>0){PetscPrintf(PETSC_COMM_SELF,"AIJ: \n");}
191:     i = 0;
192:     for (lvl=-1; lvl<10; lvl++){
193:       if (lvl==-1) {  /* Cholesky factor */
194:         factinfo.fill = 5.0;
195:         MatCholeskyFactorSymbolic(A,perm,&factinfo,&sC);
196:       } else {       /* incomplete Cholesky factor */
197:         factinfo.fill   = 5.0;
198:         factinfo.levels = lvl;
199:         MatICCFactorSymbolic(A,perm,&factinfo,&sC);
200:       }
201:       MatCholeskyFactorNumeric(A,&factinfo,&sC);
202:       MatMult(A,x,b);
203:       MatSolve(sC,b,y);
204:       MatDestroy(sC);

206:       /* Check the error */
207:       VecAXPY(y,neg_one,x);
208:       VecNorm(y,NORM_2,&norm2);
209:       if (displ>0){PetscPrintf(PETSC_COMM_SELF,"  lvl: %d, error: %G\n", lvl,norm2);}
210:       err[i++] = norm2;
211:     }
212:   }
213: 
214:   /* Test baij matrix A */
215:   if (TestBAIJ){
216:     if (displ>0){PetscPrintf(PETSC_COMM_SELF,"BAIJ: \n");}
217:     i = 0;
218:     for (lvl=-1; lvl<10; lvl++){
219:       if (lvl==-1) {  /* Cholesky factor */
220:         factinfo.fill = 5.0;
221:         MatCholeskyFactorSymbolic(A,perm,&factinfo,&sC);
222:       } else {       /* incomplete Cholesky factor */
223:         factinfo.fill   = 5.0;
224:         factinfo.levels = lvl;
225:         MatICCFactorSymbolic(A,perm,&factinfo,&sC);
226:       }
227:       MatCholeskyFactorNumeric(A,&factinfo,&sC);

229:       MatMult(A,x,b);
230:       MatSolve(sC,b,y);
231:       MatDestroy(sC);

233:       /* Check the error */
234:       VecAXPY(y,neg_one,x);
235:       VecNorm(y,NORM_2,&norm2);
236:       if (displ>0){PetscPrintf(PETSC_COMM_SELF,"  lvl: %d, error: %G\n", lvl,norm2);}
237:       err[i++] = norm2;
238:     }
239:   }

241:   /* Test sbaij matrix sA */
242:   if (displ>0){PetscPrintf(PETSC_COMM_SELF,"SBAIJ: \n");}
243:   i = 0;
244:   for (lvl=-1; lvl<10; lvl++){
245:     if (lvl==-1) {  /* Cholesky factor */
246:       factinfo.fill = 5.0;
247:       MatCholeskyFactorSymbolic(sA,perm,&factinfo,&sC);
248:     } else {       /* incomplete Cholesky factor */
249:       factinfo.fill   = 5.0;
250:       factinfo.levels = lvl;
251:       MatICCFactorSymbolic(sA,perm,&factinfo,&sC);
252:     }
253:     MatCholeskyFactorNumeric(sA,&factinfo,&sC);
254:     MatMult(sA,x,b);
255:     MatSolve(sC,b,y);

257:     /* Test MatSolves() */
258:     if (bs == 1) {
259:       Vecs xx,bb;
260:       VecsCreateSeq(PETSC_COMM_SELF,n,4,&xx);
261:       VecsDuplicate(xx,&bb);
262:       MatSolves(sC,bb,xx);
263:       VecsDestroy(xx);
264:       VecsDestroy(bb);
265:     }
266:     MatDestroy(sC);

268:     /* Check the error */
269:     VecAXPY(y,neg_one,x);
270:     VecNorm(y,NORM_2,&norm2);
271:     if (displ>0){PetscPrintf(PETSC_COMM_SELF,"  lvl: %d, error: %G\n", lvl,norm2); }
272:     err[i] -= norm2;
273:     if (err[i] > tol) SETERRQ2(PETSC_ERR_USER," level: %d, err: %G\n", lvl,err[i]);
274:   }

276:   ISDestroy(perm);
277:   PetscFree(ip_ptr);
278:   MatDestroy(A);
279:   MatDestroy(sA);
280:   VecDestroy(x);
281:   VecDestroy(y);
282:   VecDestroy(b);
283:   PetscRandomDestroy(rdm);

285:   PetscFinalize();
286:   return 0;
287: }