Actual source code: ex48.c

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

  3: static char help[] = "Tests the vatious routines in MatBAIJ format.n";

 5:  #include petscmat.h

  7: int main(int argc,char **args)
  8: {
  9:   Mat         A,B;
 10:   Vec         xx,s1,s2,yy;
 11:   int         m=45,ierr,rows[2],cols[2],bs=1,i,row,col,*idx,M;
 12:   PetscScalar rval,vals1[4],vals2[4],zero=0.0;
 13:   PetscRandom rand;
 14:   IS          is1,is2;
 15:   PetscReal   s1norm,s2norm,rnorm,tol = 1.e-10;
 16:   PetscTruth  flg;
 17: 
 18:   PetscInitialize(&argc,&args,(char *)0,help);
 19: 
 20:   /* Test MatSetValues() and MatGetValues() */
 21:   PetscOptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);
 22:   PetscOptionsGetInt(PETSC_NULL,"-mat_size",&m,PETSC_NULL);
 23:   M    = m*bs;
 24:   MatCreateSeqBAIJ(PETSC_COMM_SELF,bs,M,M,1,PETSC_NULL,&A);
 25:   MatCreateSeqAIJ(PETSC_COMM_SELF,M,M,15,PETSC_NULL,&B);
 26:   PetscRandomCreate(PETSC_COMM_SELF,RANDOM_DEFAULT,&rand);
 27:   VecCreateSeq(PETSC_COMM_SELF,M,&xx);
 28:   VecDuplicate(xx,&s1);
 29:   VecDuplicate(xx,&s2);
 30:   VecDuplicate(xx,&yy);
 31: 
 32:   /* For each row add atleast 15 elements */
 33:   for (row=0; row<M; row++) {
 34:     for (i=0; i<25*bs; i++) {
 35:       PetscRandomGetValue(rand,&rval);
 36:       col  = (int)(PetscRealPart(rval)*M);
 37:       MatSetValues(A,1,&row,1,&col,&rval,ADD_VALUES);
 38:       MatSetValues(B,1,&row,1,&col,&rval,ADD_VALUES);
 39:     }
 40:   }
 41: 
 42:   /* Now set blocks of values */
 43:   for (i=0; i<20*bs; i++) {
 44:     PetscRandomGetValue(rand,&rval);
 45:     cols[0] = (int)(PetscRealPart(rval)*M);
 46:     vals1[0] = rval;
 47:     PetscRandomGetValue(rand,&rval);
 48:     cols[1] = (int)(PetscRealPart(rval)*M);
 49:     vals1[1] = rval;
 50:     PetscRandomGetValue(rand,&rval);
 51:     rows[0] = (int)(PetscRealPart(rval)*M);
 52:     vals1[2] = rval;
 53:     PetscRandomGetValue(rand,&rval);
 54:     rows[1] = (int)(PetscRealPart(rval)*M);
 55:     vals1[3] = rval;
 56:     MatSetValues(A,2,rows,2,cols,vals1,ADD_VALUES);
 57:     MatSetValues(B,2,rows,2,cols,vals1,ADD_VALUES);
 58:   }
 59: 
 60:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 61:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 62:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 63:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 64: 
 65:   /* Test MatNorm() */
 66:   MatNorm(A,NORM_FROBENIUS,&s1norm);
 67:   MatNorm(B,NORM_FROBENIUS,&s2norm);
 68:   rnorm = s2norm-s1norm;
 69:   if (rnorm<-tol || rnorm>tol) {
 70:     PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm()- Norm1=%16.14e Norm2=%16.14e bs = %dn",s1norm,s2norm,bs);
 71:   }
 72:   /* MatScale() */
 73:   rval = 10*s1norm;
 74:   MatShift(&rval,A);
 75:   MatShift(&rval,B);
 76: 
 77:   /* Test MatTranspose() */
 78:   MatTranspose(A,PETSC_NULL);
 79:   MatTranspose(B,PETSC_NULL);
 80: 
 81: 
 82:   /* Now do MatGetValues()  */
 83:   for (i=0; i<30; i++) {
 84:     PetscRandomGetValue(rand,&rval);
 85:     cols[0] = (int)(PetscRealPart(rval)*M);
 86:     PetscRandomGetValue(rand,&rval);
 87:     cols[1] = (int)(PetscRealPart(rval)*M);
 88:     PetscRandomGetValue(rand,&rval);
 89:     rows[0] = (int)(PetscRealPart(rval)*M);
 90:     PetscRandomGetValue(rand,&rval);
 91:     rows[1] = (int)(PetscRealPart(rval)*M);
 92:     MatGetValues(A,2,rows,2,cols,vals1);
 93:     MatGetValues(B,2,rows,2,cols,vals2);
 94:     PetscMemcmp(vals1,vals2,4*sizeof(PetscScalar),&flg);
 95:     if (!flg) {
 96:       PetscPrintf(PETSC_COMM_SELF,"Error: MatGetValues bs = %dn",bs);
 97:     }
 98:   }
 99: 
100:   /* Test MatMult(), MatMultAdd() */
101:   for (i=0; i<40; i++) {
102:     VecSetRandom(rand,xx);
103:     VecSet(&zero,s2);
104:     MatMult(A,xx,s1);
105:     MatMultAdd(A,xx,s2,s2);
106:     VecNorm(s1,NORM_2,&s1norm);
107:     VecNorm(s2,NORM_2,&s2norm);
108:     rnorm = s2norm-s1norm;
109:     if (rnorm<-tol || rnorm>tol) {
110:       PetscPrintf(PETSC_COMM_SELF,"MatMult not equal to MatMultAdd Norm1=%e Norm2=%e bs = %dn",s1norm,s2norm,bs);
111:     }
112:   }

114:   /* Test MatMult() */
115:   for (i=0; i<40; i++) {
116:     VecSetRandom(rand,xx);
117:     MatMult(A,xx,s1);
118:     MatMult(B,xx,s2);
119:     VecNorm(s1,NORM_2,&s1norm);
120:     VecNorm(s2,NORM_2,&s2norm);
121:     rnorm = s2norm-s1norm;
122:     if (rnorm<-tol || rnorm>tol) {
123:       PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14e bs = %d n",s1norm,s2norm,bs);
124:     }
125:   }
126: 
127:   /* Test MatMultAdd() */
128:   for (i=0; i<40; i++) {
129:     VecSetRandom(rand,xx);
130:     VecSetRandom(rand,yy);
131:     MatMultAdd(A,xx,yy,s1);
132:     MatMultAdd(B,xx,yy,s2);
133:     VecNorm(s1,NORM_2,&s1norm);
134:     VecNorm(s2,NORM_2,&s2norm);
135:     rnorm = s2norm-s1norm;
136:     if (rnorm<-tol || rnorm>tol) {
137:       PetscPrintf(PETSC_COMM_SELF,"Error:MatMultAdd - Norm1=%16.14e Norm2=%16.14e bs = %dn",s1norm,s2norm,bs);
138:     }
139:   }
140: 
141:   /* Test MatMultTranspose() */
142:   for (i=0; i<40; i++) {
143:     VecSetRandom(rand,xx);
144:     MatMultTranspose(A,xx,s1);
145:     MatMultTranspose(B,xx,s2);
146:     VecNorm(s1,NORM_2,&s1norm);
147:     VecNorm(s2,NORM_2,&s2norm);
148:     rnorm = s2norm-s1norm;
149:     if (rnorm<-tol || rnorm>tol) {
150:       PetscPrintf(PETSC_COMM_SELF,"Error:MatMultTranspose - Norm1=%16.14e Norm2=%16.14e bs = %dn",s1norm,s2norm,bs);
151:     }
152:   }
153:   /* Test MatMultTransposeAdd() */
154:   for (i=0; i<40; i++) {
155:     VecSetRandom(rand,xx);
156:     VecSetRandom(rand,yy);
157:     MatMultTransposeAdd(A,xx,yy,s1);
158:     MatMultTransposeAdd(B,xx,yy,s2);
159:     VecNorm(s1,NORM_2,&s1norm);
160:     VecNorm(s2,NORM_2,&s2norm);
161:     rnorm = s2norm-s1norm;
162:     if (rnorm<-tol || rnorm>tol) {
163:       PetscPrintf(PETSC_COMM_SELF,"Error:MatMultTransposeAdd - Norm1=%16.14e Norm2=%16.14e bs = %dn",s1norm,s2norm,bs);
164:     }
165:   }
166: 
167: 
168:   /* Do LUFactor() on both the matrices */
169:   PetscMalloc(M*sizeof(int),&idx);
170:   for (i=0; i<M; i++) idx[i] = i;
171:   ISCreateGeneral(PETSC_COMM_SELF,M,idx,&is1);
172:   ISCreateGeneral(PETSC_COMM_SELF,M,idx,&is2);
173:   PetscFree(idx);
174:   ISSetPermutation(is1);
175:   ISSetPermutation(is2);
176:   MatLUFactor(B,is1,is2,PETSC_NULL);
177:   MatLUFactor(A,is1,is2,PETSC_NULL);
178: 
179: 
180:   /* Test MatSolveAdd() */
181:   for (i=0; i<40; i++) {
182:     VecSetRandom(rand,xx);
183:     VecSetRandom(rand,yy);
184:     MatSolveAdd(B,xx,yy,s2);
185:     MatSolveAdd(A,xx,yy,s1);
186:     VecNorm(s1,NORM_2,&s1norm);
187:     VecNorm(s2,NORM_2,&s2norm);
188:     rnorm = s2norm-s1norm;
189:     if (rnorm<-tol || rnorm>tol) {
190:       PetscPrintf(PETSC_COMM_SELF,"Error:MatSolveAdd - Norm1=%16.14e Norm2=%16.14e bs = %dn",s1norm,s2norm,bs);
191:     }
192:   }
193: 
194:   /* Test MatSolveAdd() when x = A'b +x */
195:   for (i=0; i<40; i++) {
196:     VecSetRandom(rand,xx);
197:     VecSetRandom(rand,s1);
198:     VecCopy(s2,s1);
199:     MatSolveAdd(B,xx,s2,s2);
200:     MatSolveAdd(A,xx,s1,s1);
201:     VecNorm(s1,NORM_2,&s1norm);
202:     VecNorm(s2,NORM_2,&s2norm);
203:     rnorm = s2norm-s1norm;
204:     if (rnorm<-tol || rnorm>tol) {
205:       PetscPrintf(PETSC_COMM_SELF,"Error:MatSolveAdd(same) - Norm1=%16.14e Norm2=%16.14e bs = %dn",s1norm,s2norm,bs);
206:     }
207:   }
208: 
209:   /* Test MatSolve() */
210:   for (i=0; i<40; i++) {
211:     VecSetRandom(rand,xx);
212:     MatSolve(B,xx,s2);
213:     MatSolve(A,xx,s1);
214:     VecNorm(s1,NORM_2,&s1norm);
215:     VecNorm(s2,NORM_2,&s2norm);
216:     rnorm = s2norm-s1norm;
217:     if (rnorm<-tol || rnorm>tol) {
218:       PetscPrintf(PETSC_COMM_SELF,"Error:MatSolve - Norm1=%16.14e Norm2=%16.14e bs = %dn",s1norm,s2norm,bs);
219:     }
220:   }
221: 
222:   /* Test MatSolveTranspose() */
223:   if (bs < 8) {
224:     for (i=0; i<40; i++) {
225:       VecSetRandom(rand,xx);
226:       MatSolveTranspose(B,xx,s2);
227:       MatSolveTranspose(A,xx,s1);
228:       VecNorm(s1,NORM_2,&s1norm);
229:       VecNorm(s2,NORM_2,&s2norm);
230:       rnorm = s2norm-s1norm;
231:       if (rnorm<-tol || rnorm>tol) {
232:         PetscPrintf(PETSC_COMM_SELF,"Error:MatSolveTranspose - Norm1=%16.14e Norm2=%16.14e bs = %dn",s1norm,s2norm,bs);
233:       }
234:     }
235:   }

237:   MatDestroy(A);
238:   MatDestroy(B);
239:   VecDestroy(xx);
240:   VecDestroy(s1);
241:   VecDestroy(s2);
242:   VecDestroy(yy);
243:   ISDestroy(is1);
244:   ISDestroy(is2);
245:   PetscRandomDestroy(rand);
246:   PetscFinalize();
247:   return 0;
248: }