Actual source code: ex74.c

  1: /*$Id: ex74.c,v 1.47 2001/08/07 21:30:08 bsmith Exp $*/

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

 5:  #include petscmat.h
  6: extern int MatMissingDiagonal_SeqSBAIJ(Mat);
  7: extern int MatMarkDiagonal_SeqSBAIJ(Mat);

  9: int main(int argc,char **args)
 10: {
 11:   Vec          x,y,b,s1,s2;
 12:   Mat          A;             /* linear system matrix */
 13:   Mat          sA,sC;         /* symmetric part of the matrices */
 14:   int          n,mbs=16,bs=1,nz=3,prob=1;
 15:   int          ierr,i,j,col[3],size,block, row,I,J,n1,*ip_ptr,inc;
 16:   int          lf;           /* level of fill for icc */
 17:   int          *cols1,*cols2;
 18:   PetscReal    norm1,norm2,tol=1.e-10,fill;
 19:   PetscScalar  neg_one = -1.0,four=4.0,value[3],alpha=0.1;
 20:   PetscScalar  *vr1,*vr2,*vr1_wk,*vr2_wk;
 21:   IS           perm, isrow, iscol;
 22:   PetscRandom  rand;
 23:   PetscTruth   getrow=PETSC_FALSE;
 24:   MatInfo      minfo1,minfo2;

 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);

 32:   n = mbs*bs;
 33:   ierr=MatCreateSeqBAIJ(PETSC_COMM_WORLD,bs,n,n,nz,PETSC_NULL, &A);
 34:   ierr=MatCreateSeqSBAIJ(PETSC_COMM_WORLD,bs,n,n,nz,PETSC_NULL, &sA);

 36:   /* Test MatGetOwnershipRange() */
 37:   MatGetOwnershipRange(A,&I,&J);
 38:   MatGetOwnershipRange(sA,&i,&j);
 39:   if (i-I || j-J){
 40:     PetscPrintf(PETSC_COMM_SELF,"Error: MatGetOwnershipRange() in MatSBAIJ formatn");
 41:   }

 43:   /* Assemble matrix */
 44:   if (bs == 1){
 45:     PetscOptionsGetInt(PETSC_NULL,"-test_problem",&prob,PETSC_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:         MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
 52:       }
 53:       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);
 56:       MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);

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

110:       i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs;
111:       value[0]=4.0; value[1] = -1.0;
112:       MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
113:       MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);
114:     }
115:     /* off-diagonal blocks */
116:     value[0]=-1.0;
117:     for (i=0; i<(n/bs-1)*bs; i++){
118:       col[0]=i+bs;
119:       MatSetValues(A,1,&i,1,col,value,INSERT_VALUES);
120:       MatSetValues(sA,1,&i,1,col,value,INSERT_VALUES);
121:       col[0]=i; row=i+bs;
122:       MatSetValues(A,1,&row,1,col,value,INSERT_VALUES);
123:       MatSetValues(sA,1,&row,1,col,value,INSERT_VALUES);
124:     }
125:   }
126:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
127:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
128:   /* PetscPrintf(PETSC_COMM_SELF,"n The Matrix: n");
129:   MatView(A, PETSC_VIEWER_DRAW_WORLD);
130:   MatView(A, PETSC_VIEWER_STDOUT_WORLD); */

132:   MatAssemblyBegin(sA,MAT_FINAL_ASSEMBLY);
133:   MatAssemblyEnd(sA,MAT_FINAL_ASSEMBLY);
134:   /* PetscPrintf(PETSC_COMM_SELF,"n Symmetric Part of Matrix: n");
135:   MatView(sA, PETSC_VIEWER_DRAW_WORLD); 
136:   MatView(sA, PETSC_VIEWER_STDOUT_WORLD); 
137:   */

139:   /* Test MatNorm(), MatDuplicate() */
140:   MatNorm(A,NORM_FROBENIUS,&norm1);
141:   MatDuplicate(sA,MAT_COPY_VALUES,&sC);
142:   MatNorm(sC,NORM_FROBENIUS,&norm2);
143:   MatDestroy(sC);
144:   norm1 -= norm2;
145:   if (norm1<-tol || norm1>tol){
146:     PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm(), fnorm1-fnorm2=%16.14en",norm1);
147:   }
148:   MatNorm(A,NORM_INFINITY,&norm1);
149:   MatNorm(sA,NORM_INFINITY,&norm2);
150:   norm1 -= norm2;
151:   if (norm1<-tol || norm1>tol){
152:     PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm(), inf_norm1-inf_norm2=%16.14en",norm1);
153:   }
154: 
155:   /* Test MatGetInfo(), MatGetSize(), MatGetBlockSize() */
156:   MatGetInfo(A,MAT_LOCAL,&minfo1);
157:   MatGetInfo(sA,MAT_LOCAL,&minfo2);
158:   /*
159:   printf("matrix nonzeros (BAIJ format) = %d, allocated nonzeros= %dn", (int)minfo1.nz_used,(int)minfo1.nz_allocated); 
160:   printf("matrix nonzeros(SBAIJ format) = %d, allocated nonzeros= %dn", (int)minfo2.nz_used,(int)minfo2.nz_allocated); 
161:   */
162:   i = (int) (minfo1.nz_used - minfo2.nz_used);
163:   j = (int) (minfo2.nz_allocated - minfo2.nz_used);
164:   if (i<0 || j<0) {
165:     PetscPrintf(PETSC_COMM_SELF,"Error: MatGetInfo()n");
166:   }

168:   MatGetSize(A,&I,&J);
169:   MatGetSize(sA,&i,&j);
170:   if (i-I || j-J) {
171:     PetscPrintf(PETSC_COMM_SELF,"Error: MatGetSize()n");
172:   }
173: 
174:   MatGetBlockSize(A, &I);
175:   MatGetBlockSize(sA, &i);
176:   if (i-I){
177:     PetscPrintf(PETSC_COMM_SELF,"Error: MatGetBlockSize()n");
178:   }

180:   /* Test MatGetRow() */
181:   if (getrow){
182:     row = n/2;
183:     PetscMalloc(n*sizeof(PetscScalar),&vr1);
184:     vr1_wk = vr1;
185:     PetscMalloc(n*sizeof(PetscScalar),&vr2);
186:     vr2_wk = vr2;
187:     MatGetRow(A,row,&J,&cols1,&vr1);
188:     vr1_wk += J-1;
189:     MatGetRow(sA,row,&j,&cols2,&vr2);
190:     vr2_wk += j-1;
191:     VecCreateSeq(PETSC_COMM_SELF,j,&x);
192: 
193:     for (i=j-1; i>-1; i--){
194:       VecSetValue(x,i,*vr2_wk - *vr1_wk,INSERT_VALUES);
195:       vr2_wk--; vr1_wk--;
196:     }
197:     VecNorm(x,NORM_1,&norm2);
198:     if (norm2<-tol || norm2>tol) {
199:       PetscPrintf(PETSC_COMM_SELF,"Error: MatGetRow()n");
200:     }
201:     VecDestroy(x);
202:     MatRestoreRow(A,row,&J,&cols1,&vr1);
203:     MatRestoreRow(sA,row,&j,&cols2,&vr2);
204:     PetscFree(vr1);
205:     PetscFree(vr2);

207:     /* Test GetSubMatrix() */
208:     /* get a submatrix consisting of every next block row and column of the original matrix */
209:     /* for symm. matrix, iscol=isrow. */
210:     PetscMalloc(n*sizeof(IS),&isrow);
211:     PetscMalloc(n*sizeof(int),&ip_ptr);
212:     j = 0;
213:     for (n1=0; n1<mbs; n1 += 2){ /* n1: block row */
214:       for (i=0; i<bs; i++) ip_ptr[j++] = n1*bs + i;
215:     }
216:     ISCreateGeneral(PETSC_COMM_SELF, j, ip_ptr, &isrow);
217:     /* ISView(isrow, PETSC_VIEWER_STDOUT_SELF); */
218: 
219:     MatGetSubMatrix(sA,isrow,isrow,PETSC_DECIDE,MAT_INITIAL_MATRIX,&sC);
220:     ISDestroy(isrow);
221:     PetscFree(ip_ptr);
222:     printf("sA =n");
223:     MatView(sA,PETSC_VIEWER_STDOUT_WORLD);
224:     printf("submatrix of sA =n");
225:     MatView(sC,PETSC_VIEWER_STDOUT_WORLD);
226:     MatDestroy(sC);
227:   }

229:   /* Test MatDiagonalScale(), MatGetDiagonal(), MatScale() */
230:   PetscRandomCreate(PETSC_COMM_SELF,RANDOM_DEFAULT,&rand);
231:   VecCreateSeq(PETSC_COMM_SELF,n,&x);
232:   VecDuplicate(x,&s1);
233:   VecDuplicate(x,&s2);
234:   VecDuplicate(x,&y);
235:   VecDuplicate(x,&b);
236:   VecSetRandom(rand,x);

238:   MatMarkDiagonal_SeqSBAIJ(sA);
239:   MatMissingDiagonal_SeqSBAIJ(sA);

241:   MatDiagonalScale(A,x,x);
242:   MatDiagonalScale(sA,x,x);

244:   MatGetDiagonal(A,s1);
245:   MatGetDiagonal(sA,s2);
246:   VecNorm(s1,NORM_1,&norm1);
247:   VecNorm(s2,NORM_1,&norm2);
248:   norm1 -= norm2;
249:   if (norm1<-tol || norm1>tol) {
250:     PetscPrintf(PETSC_COMM_SELF,"Error:MatGetDiagonal() n");
251:   }

253:   MatScale(&alpha,A);
254:   MatScale(&alpha,sA);

256:   /* Test MatGetRowMax() */
257:   MatGetRowMax(A,s1);
258:   MatGetRowMax(sA,s2);
259:   VecNorm(s1,NORM_1,&norm1);
260:   VecNorm(s2,NORM_1,&norm2);
261:   norm1 -= norm2;
262:   if (norm1<-tol || norm1>tol) {
263:     PetscPrintf(PETSC_COMM_SELF,"Error:MatGetRowMax() n");
264:   }

266:   /* Test MatMult(), MatMultAdd() */
267:   for (i=0; i<40; i++) {
268:     VecSetRandom(rand,x);
269:     MatMult(A,x,s1);
270:     MatMult(sA,x,s2);
271:     VecNorm(s1,NORM_1,&norm1);
272:     VecNorm(s2,NORM_1,&norm2);
273:     norm1 -= norm2;
274:     if (norm1<-tol || norm1>tol) {
275:       PetscPrintf(PETSC_COMM_SELF,"Error: MatMult(), MatDiagonalScale() or MatScale()n");
276:     }
277:   }

279:   for (i=0; i<40; i++) {
280:     VecSetRandom(rand,x);
281:     VecSetRandom(rand,y);
282:     MatMultAdd(A,x,y,s1);
283:     MatMultAdd(sA,x,y,s2);
284:     VecNorm(s1,NORM_1,&norm1);
285:     VecNorm(s2,NORM_1,&norm2);
286:     norm1 -= norm2;
287:     if (norm1<-tol || norm1>tol) {
288:       PetscPrintf(PETSC_COMM_SELF,"Error:MatMultAdd(), MatDiagonalScale() or MatScale() n");
289:     }
290:   }

292:   /* Test MatCholeskyFactor(), MatICCFactor() with natural ordering */
293:   MatGetOrdering(A,MATORDERING_NATURAL,&perm,&iscol);
294:   ISDestroy(iscol);
295:   norm1 = tol;
296:   inc   = bs;
297:   for (lf=-1; lf<10; lf += inc){
298:     if (lf==-1) {  /* Cholesky factor */
299:       fill = 5.0;
300:       MatCholeskyFactorSymbolic(sA,perm,fill,&sC);
301:     } else {       /* incomplete Cholesky factor */
302:       fill          = 5.0;
303:       MatICCFactorSymbolic(sA,perm,fill,lf,&sC);
304:     }
305:     MatCholeskyFactorNumeric(sA,&sC);
306:     /* MatView(sC, PETSC_VIEWER_DRAW_WORLD); */
307: 
308:     MatMult(sA,x,b);
309:     MatSolve(sC,b,y);
310:     MatDestroy(sC);
311: 
312:     /* Check the error */
313:     VecAXPY(&neg_one,x,y);
314:     VecNorm(y,NORM_2,&norm2);
315:     /* printf("lf: %d, error: %gn", lf,norm2); */
316:     if (10*norm1 < norm2 && lf-inc != -1){
317:       PetscPrintf(PETSC_COMM_SELF,"lf=%d, %d, Norm of error=%g, %gn",lf-inc,lf,norm1,norm2);
318:     }
319:     norm1 = norm2;
320:     if (norm2 < tol && lf != -1) break;
321:   }

323:   ISDestroy(perm);

325:   MatDestroy(A);
326:   MatDestroy(sA);
327:   VecDestroy(x);
328:   VecDestroy(y);
329:   VecDestroy(s1);
330:   VecDestroy(s2);
331:   VecDestroy(b);
332:   PetscRandomDestroy(rand);

334:   PetscFinalize();
335:   return 0;
336: }