Actual source code: ex92.c

  2: static char help[] = "Tests MatIncreaseOverlap(), MatGetSubMatrices() for parallel MatSBAIJ format.\n";

 4:  #include petscmat.h

  8: int main(int argc,char **args)
  9: {
 10:   Mat            A,Atrans,sA,*submatA,*submatsA;
 12:   PetscMPIInt    size,rank;
 13:   PetscInt       bs=1,mbs=11,ov=1,i,j,k,*rows,*cols,nd=5,*idx,rstart,rend,sz,mm,M,N,Mbs;
 14:   PetscScalar    *vals,rval,one=1.0;
 15:   IS             *is1,*is2;
 16:   PetscRandom    rand;
 17:   Vec            xx,s1,s2;
 18:   PetscReal      s1norm,s2norm,rnorm,tol = 1.e-10;
 19:   PetscTruth     flg;
 20:   int            stages[2];

 22:   PetscInitialize(&argc,&args,(char *)0,help);
 23:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 24:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

 26:   PetscOptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);
 27:   PetscOptionsGetInt(PETSC_NULL,"-mat_size",&mbs,PETSC_NULL);
 28:   PetscOptionsGetInt(PETSC_NULL,"-ov",&ov,PETSC_NULL);
 29:   PetscOptionsGetInt(PETSC_NULL,"-nd",&nd,PETSC_NULL);

 31:   MatCreateMPIBAIJ(PETSC_COMM_WORLD,bs,mbs*bs,mbs*bs,PETSC_DECIDE,PETSC_DECIDE,
 32:                           PETSC_DEFAULT,PETSC_NULL,PETSC_DEFAULT,PETSC_NULL,&A);
 33:   PetscRandomCreate(PETSC_COMM_WORLD,&rand);
 34:   PetscRandomSetFromOptions(rand);

 36:   MatGetOwnershipRange(A,&rstart,&rend);
 37:   MatGetSize(A,&M,&N);
 38:   Mbs  = M/bs;

 40:   PetscMalloc(bs*sizeof(PetscInt),&rows);
 41:   PetscMalloc(bs*sizeof(PetscInt),&cols);
 42:   PetscMalloc(bs*bs*sizeof(PetscScalar),&vals);
 43:   PetscMalloc(M*sizeof(PetscScalar),&idx);

 45:   /* Now set blocks of values */
 46:   for (j=0; j<bs*bs; j++) vals[j] = 0.0;
 47:   for (i=0; i<Mbs; i++){
 48:     cols[0] = i*bs; rows[0] = i*bs;
 49:     for (j=1; j<bs; j++) {
 50:       rows[j] = rows[j-1]+1;
 51:       cols[j] = cols[j-1]+1;
 52:     }
 53:     MatSetValues(A,bs,rows,bs,cols,vals,ADD_VALUES);
 54:   }
 55:   /* second, add random blocks */
 56:   for (i=0; i<20*bs; i++) {
 57:       PetscRandomGetValue(rand,&rval);
 58:       cols[0] = bs*(PetscInt)(PetscRealPart(rval)*Mbs);
 59:       PetscRandomGetValue(rand,&rval);
 60:       rows[0] = rstart + bs*(PetscInt)(PetscRealPart(rval)*mbs);
 61:       for (j=1; j<bs; j++) {
 62:         rows[j] = rows[j-1]+1;
 63:         cols[j] = cols[j-1]+1;
 64:       }

 66:       for (j=0; j<bs*bs; j++) {
 67:         PetscRandomGetValue(rand,&rval);
 68:         vals[j] = rval;
 69:       }
 70:       MatSetValues(A,bs,rows,bs,cols,vals,ADD_VALUES);
 71:   }

 73:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 74:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 75: 
 76:   /* make A a symmetric matrix: A <- A^T + A */
 77:   MatTranspose(A, &Atrans);
 78:   MatAXPY(A,one,Atrans,DIFFERENT_NONZERO_PATTERN);
 79:   MatDestroy(Atrans);
 80:   MatTranspose(A, &Atrans);
 81:   MatEqual(A, Atrans, &flg);
 82:   if (!flg) {
 83:     SETERRQ(1,"A+A^T is non-symmetric");
 84:   }
 85:   MatDestroy(Atrans);
 86:   /* MatView(A,PETSC_VIEWER_STDOUT_WORLD); */

 88:   /* create a SeqSBAIJ matrix sA (= A) */
 89:   MatConvert(A,MATMPISBAIJ,MAT_INITIAL_MATRIX,&sA);
 90:   /* MatView(sA,PETSC_VIEWER_STDOUT_WORLD); */

 92:   /* Test sA==A through MatMult() */
 93:   for (i=0; i<nd; i++) {
 94:     MatGetLocalSize(A, &mm,PETSC_NULL);
 95:     VecCreate(PETSC_COMM_WORLD,&xx);
 96:     VecSetSizes(xx,mm,PETSC_DECIDE);
 97:     VecSetFromOptions(xx);
 98:     VecDuplicate(xx,&s1);
 99:     VecDuplicate(xx,&s2);
100:     for (j=0; j<3; j++) {
101:       VecSetRandom(xx,rand);
102:       MatMult(A,xx,s1);
103:       MatMult(sA,xx,s2);
104:       VecNorm(s1,NORM_2,&s1norm);
105:       VecNorm(s2,NORM_2,&s2norm);
106:       rnorm = s2norm-s1norm;
107:       if (rnorm<-tol || rnorm>tol) {
108:         PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",s1norm,s2norm);
109:       }
110:     }
111:     VecDestroy(xx);
112:     VecDestroy(s1);
113:     VecDestroy(s2);
114:   }
115: 
116:   /* Test MatIncreaseOverlap() */
117:   PetscMalloc(nd*sizeof(IS **),&is1);
118:   PetscMalloc(nd*sizeof(IS **),&is2);

120:   for (i=0; i<nd; i++) {
121:     PetscRandomGetValue(rand,&rval);
122:     sz = (PetscInt)((0.5 + 0.2*PetscRealPart(rval))*mbs); /* 0.5*mbs < sz < 0.7*mbs */
123:     sz /= (size*nd*10);
124:     /*
125:     if (!rank){
126:       PetscPrintf(PETSC_COMM_SELF," [%d] IS sz[%d]: %d\n",rank,i,sz);
127:     } 
128:     */
129:     for (j=0; j<sz; j++) {
130:       PetscRandomGetValue(rand,&rval);
131:       idx[j*bs] = bs*(PetscInt)(PetscRealPart(rval)*Mbs);
132:       for (k=1; k<bs; k++) idx[j*bs+k] = idx[j*bs]+k;
133:     }
134:     ISCreateGeneral(PETSC_COMM_SELF,sz*bs,idx,is1+i);
135:     ISCreateGeneral(PETSC_COMM_SELF,sz*bs,idx,is2+i);
136:   }

138:   PetscLogStageRegister(&stages[0],"MatOv_SBAIJ");
139:   PetscLogStageRegister(&stages[1],"MatOv_ BAIJ");

141:   PetscLogStagePush(stages[0]);
142:   MatIncreaseOverlap(sA,nd,is2,ov);
143:   PetscLogStagePop();

145:   PetscLogStagePush(stages[1]);
146:   MatIncreaseOverlap(A,nd,is1,ov);
147:   PetscLogStagePop();

149:   for (i=0; i<nd; ++i) {
150:     ISEqual(is1[i],is2[i],&flg);
151:     if (!flg ){
152:       if (rank == size){
153:         ISSort(is1[i]);
154:         ISView(is1[i],PETSC_VIEWER_STDOUT_SELF);
155:         ISSort(is2[i]);
156:         ISView(is2[i],PETSC_VIEWER_STDOUT_SELF);
157:       }
158:       SETERRQ1(1,"i=%D, is1 != is2",i);
159:     }
160:   }
161: 
162:   for (i=0; i<nd; ++i) {
163:     ISSort(is1[i]);
164:     ISSort(is2[i]);
165:   }
166:   MatGetSubMatrices(sA,nd,is2,is2,MAT_INITIAL_MATRIX,&submatsA);
167:   MatGetSubMatrices(A,nd,is1,is1,MAT_INITIAL_MATRIX,&submatA);

169:   /* Test MatMult() */
170:   for (i=0; i<nd; i++) {
171:     MatGetSize(submatA[i],&mm,PETSC_NULL);
172:     VecCreateSeq(PETSC_COMM_SELF,mm,&xx);
173:     VecDuplicate(xx,&s1);
174:     VecDuplicate(xx,&s2);
175:     for (j=0; j<3; j++) {
176:       VecSetRandom(xx,rand);
177:       MatMult(submatA[i],xx,s1);
178:       MatMult(submatsA[i],xx,s2);
179:       VecNorm(s1,NORM_2,&s1norm);
180:       VecNorm(s2,NORM_2,&s2norm);
181:       rnorm = s2norm-s1norm;
182:       if (rnorm<-tol || rnorm>tol) {
183:         PetscPrintf(PETSC_COMM_SELF,"[%d]Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",rank,s1norm,s2norm);
184:       }
185:     }
186:     VecDestroy(xx);
187:     VecDestroy(s1);
188:     VecDestroy(s2);
189:   }

191:   /* Now test MatGetSubmatrices with MAT_REUSE_MATRIX option */
192:   MatGetSubMatrices(A,nd,is1,is1,MAT_REUSE_MATRIX,&submatA);
193:   MatGetSubMatrices(sA,nd,is2,is2,MAT_REUSE_MATRIX,&submatsA);

195:   /* Test MatMult() */
196:   for (i=0; i<nd; i++) {
197:     MatGetSize(submatA[i],&mm,PETSC_NULL);
198:     VecCreateSeq(PETSC_COMM_SELF,mm,&xx);
199:     VecDuplicate(xx,&s1);
200:     VecDuplicate(xx,&s2);
201:     for (j=0; j<3; j++) {
202:       VecSetRandom(xx,rand);
203:       MatMult(submatA[i],xx,s1);
204:       MatMult(submatsA[i],xx,s2);
205:       VecNorm(s1,NORM_2,&s1norm);
206:       VecNorm(s2,NORM_2,&s2norm);
207:       rnorm = s2norm-s1norm;
208:       if (rnorm<-tol || rnorm>tol) {
209:         PetscPrintf(PETSC_COMM_SELF,"[%d]Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",rank,s1norm,s2norm);
210:       }
211:     }
212:     VecDestroy(xx);
213:     VecDestroy(s1);
214:     VecDestroy(s2);
215:   }

217:   /* Free allocated memory */
218:   for (i=0; i<nd; ++i) {
219:     ISDestroy(is1[i]);
220:     ISDestroy(is2[i]);
221:     MatDestroy(submatA[i]);
222:     MatDestroy(submatsA[i]);
223:   }
224:   PetscFree(submatA);
225:   PetscFree(submatsA);
226:   PetscFree(is1);
227:   PetscFree(is2);
228:   PetscFree(idx);
229:   PetscFree(rows);
230:   PetscFree(cols);
231:   PetscFree(vals);
232:   MatDestroy(A);
233:   MatDestroy(sA);
234:   PetscRandomDestroy(rand);
235:   PetscFinalize();
236:   return 0;
237: }