Actual source code: ex54.c

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

  3: static char help[] = "Tests MatIncreaseOverlap(), MatGetSubMatrices() for MatBAIJ format.n";

 5:  #include petscmat.h

  7: int main(int argc,char **args)
  8: {
  9:   Mat         A,B,*submatA,*submatB;
 10:   int         bs=1,m=11,ov=1,i,j,k,*rows,*cols,ierr,nd=5,*idx,size;
 11:   int         rank,rstart,rend,sz,mm,nn,M,N,Mbs;
 12:   PetscScalar *vals,rval;
 13:   IS          *is1,*is2;
 14:   PetscRandom rand;
 15:   Vec         xx,s1,s2;
 16:   PetscReal   s1norm,s2norm,rnorm,tol = 1.e-10;
 17:   PetscTruth  flg;

 19:   PetscInitialize(&argc,&args,(char *)0,help);
 20:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 21:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

 23:   PetscOptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);
 24:   PetscOptionsGetInt(PETSC_NULL,"-mat_size",&m,PETSC_NULL);
 25:   PetscOptionsGetInt(PETSC_NULL,"-ov",&ov,PETSC_NULL);
 26:   PetscOptionsGetInt(PETSC_NULL,"-nd",&nd,PETSC_NULL);

 28:   MatCreateMPIBAIJ(PETSC_COMM_WORLD,bs,m*bs,m*bs,PETSC_DECIDE,PETSC_DECIDE,
 29:                           PETSC_DEFAULT,PETSC_NULL,PETSC_DEFAULT,PETSC_NULL,&A);
 30:   MatCreateMPIAIJ(PETSC_COMM_WORLD,m*bs,m*bs,PETSC_DECIDE,PETSC_DECIDE,
 31:                          PETSC_DEFAULT,PETSC_NULL,PETSC_DEFAULT,PETSC_NULL,&B);
 32:   PetscRandomCreate(PETSC_COMM_WORLD,RANDOM_DEFAULT,&rand);

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

 38:   PetscMalloc(bs*sizeof(int),&rows);
 39:   PetscMalloc(bs*sizeof(int),&cols);
 40:   PetscMalloc(bs*bs*sizeof(PetscScalar),&vals);
 41:   PetscMalloc(M*sizeof(PetscScalar),&idx);

 43:   /* Now set blocks of values */
 44:   for (i=0; i<40*bs; i++) {
 45:       PetscRandomGetValue(rand,&rval);
 46:       cols[0] = bs*(int)(PetscRealPart(rval)*Mbs);
 47:       PetscRandomGetValue(rand,&rval);
 48:       rows[0] = rstart + bs*(int)(PetscRealPart(rval)*m);
 49:       for (j=1; j<bs; j++) {
 50:         rows[j] = rows[j-1]+1;
 51:         cols[j] = cols[j-1]+1;
 52:       }

 54:       for (j=0; j<bs*bs; j++) {
 55:         PetscRandomGetValue(rand,&rval);
 56:         vals[j] = rval;
 57:       }
 58:       MatSetValues(A,bs,rows,bs,cols,vals,ADD_VALUES);
 59:       MatSetValues(B,bs,rows,bs,cols,vals,ADD_VALUES);
 60:   }

 62:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 63:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 64:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 65:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

 67:     /* Test MatIncreaseOverlap() */
 68:   PetscMalloc(nd*sizeof(IS **),&is1);
 69:   PetscMalloc(nd*sizeof(IS **),&is2);

 71: 
 72:   for (i=0; i<nd; i++) {
 73:     PetscRandomGetValue(rand,&rval);
 74:     sz = (int)(PetscRealPart(rval)*m);
 75:     for (j=0; j<sz; j++) {
 76:       PetscRandomGetValue(rand,&rval);
 77:       idx[j*bs] = bs*(int)(PetscRealPart(rval)*Mbs);
 78:       for (k=1; k<bs; k++) idx[j*bs+k] = idx[j*bs]+k;
 79:     }
 80:     ISCreateGeneral(PETSC_COMM_SELF,sz*bs,idx,is1+i);
 81:     ISCreateGeneral(PETSC_COMM_SELF,sz*bs,idx,is2+i);
 82:   }
 83:   MatIncreaseOverlap(A,nd,is1,ov);
 84:   MatIncreaseOverlap(B,nd,is2,ov);

 86:   for (i=0; i<nd; ++i) {
 87:     ISEqual(is1[i],is2[i],&flg);

 89:     if (!flg) {
 90:       PetscPrintf(PETSC_COMM_SELF,"i=%d, flg=%d :bs=%d m=%d ov=%d nd=%d np=%dn",i,flg,bs,m,ov,nd,size);
 91:     }
 92:   }

 94:   for (i=0; i<nd; ++i) {
 95:     ISSort(is1[i]);
 96:     ISSort(is2[i]);
 97:   }
 98: 
 99:   MatGetSubMatrices(B,nd,is2,is2,MAT_INITIAL_MATRIX,&submatB);
100:   MatGetSubMatrices(A,nd,is1,is1,MAT_INITIAL_MATRIX,&submatA);


103:   /* Test MatMult() */
104:   for (i=0; i<nd; i++) {
105:     MatGetSize(submatA[i],&mm,&nn);
106:     VecCreateSeq(PETSC_COMM_SELF,mm,&xx);
107:     VecDuplicate(xx,&s1);
108:     VecDuplicate(xx,&s2);
109:     for (j=0; j<3; j++) {
110:       VecSetRandom(rand,xx);
111:       MatMult(submatA[i],xx,s1);
112:       MatMult(submatB[i],xx,s2);
113:       VecNorm(s1,NORM_2,&s1norm);
114:       VecNorm(s2,NORM_2,&s2norm);
115:       rnorm = s2norm-s1norm;
116:       if (rnorm<-tol || rnorm>tol) {
117:         PetscPrintf(PETSC_COMM_SELF,"[%d]Error:MatMult - Norm1=%16.14e Norm2=%16.14en",rank,s1norm,s2norm);
118:       }
119:     }
120:     VecDestroy(xx);
121:     VecDestroy(s1);
122:     VecDestroy(s2);
123:   }

125:   /* Now test MatGetSubmatrices with MAT_REUSE_MATRIX option */
126: 
127:   MatGetSubMatrices(A,nd,is1,is1,MAT_REUSE_MATRIX,&submatA);
128:   MatGetSubMatrices(B,nd,is2,is2,MAT_REUSE_MATRIX,&submatB);

130:   /* Test MatMult() */
131:   for (i=0; i<nd; i++) {
132:     MatGetSize(submatA[i],&mm,&nn);
133:     VecCreateSeq(PETSC_COMM_SELF,mm,&xx);
134:     VecDuplicate(xx,&s1);
135:     VecDuplicate(xx,&s2);
136:     for (j=0; j<3; j++) {
137:       VecSetRandom(rand,xx);
138:       MatMult(submatA[i],xx,s1);
139:       MatMult(submatB[i],xx,s2);
140:       VecNorm(s1,NORM_2,&s1norm);
141:       VecNorm(s2,NORM_2,&s2norm);
142:       rnorm = s2norm-s1norm;
143:       if (rnorm<-tol || rnorm>tol) {
144:         PetscPrintf(PETSC_COMM_SELF,"[%d]Error:MatMult - Norm1=%16.14e Norm2=%16.14en",rank,s1norm,s2norm);
145:       }
146:     }
147:     VecDestroy(xx);
148:     VecDestroy(s1);
149:     VecDestroy(s2);
150:   }
151: 
152:   /* Free allocated memory */
153:   for (i=0; i<nd; ++i) {
154:     ISDestroy(is1[i]);
155:     ISDestroy(is2[i]);
156:     MatDestroy(submatA[i]);
157:     MatDestroy(submatB[i]);
158:  }
159:   PetscFree(is1);
160:   PetscFree(is2);
161:   PetscFree(idx);
162:   PetscFree(rows);
163:   PetscFree(cols);
164:   PetscFree(vals);
165:   MatDestroy(A);
166:   MatDestroy(B);
167:   PetscFree(submatA);
168:   PetscFree(submatB);
169:   PetscRandomDestroy(rand);
170:   PetscFinalize();
171:   return 0;
172: }