xref: /petsc/src/mat/tests/ex54.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatIncreaseOverlap(), MatCreateSubMatrices() for parallel AIJ and BAIJ formats.\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **args)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   Mat            A,B,*submatA,*submatB;
9c4762a1bSJed Brown   PetscInt       bs=1,m=11,ov=1,i,j,k,*rows,*cols,nd=5,*idx,rstart,rend,sz,mm,nn,M,N,Mbs;
10c4762a1bSJed Brown   PetscMPIInt    size,rank;
11c4762a1bSJed Brown   PetscScalar    *vals,rval;
12c4762a1bSJed Brown   IS             *is1,*is2;
13c4762a1bSJed Brown   PetscRandom    rdm;
14c4762a1bSJed Brown   Vec            xx,s1,s2;
15c4762a1bSJed Brown   PetscReal      s1norm,s2norm,rnorm,tol = 100*PETSC_SMALL;
1682b5ce2aSStefano Zampini   PetscBool      flg,test_nd0=PETSC_FALSE, emptynd;
17c4762a1bSJed Brown 
18*327415f7SBarry Smith   PetscFunctionBeginUser;
199566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
22c4762a1bSJed Brown 
239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL));
249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_size",&m,NULL));
259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL));
269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL));
279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-test_nd0",&test_nd0,NULL));
28c4762a1bSJed Brown 
29c4762a1bSJed Brown   /* Create a AIJ matrix A */
309566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
319566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,m*bs,m*bs,PETSC_DECIDE,PETSC_DECIDE));
329566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATAIJ));
339566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(A,PETSC_DEFAULT,NULL));
349566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(A,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL));
359566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
369566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE));
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   /* Create a BAIJ matrix B */
399566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
409566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B,m*bs,m*bs,PETSC_DECIDE,PETSC_DECIDE));
419566063dSJacob Faibussowitsch   PetscCall(MatSetType(B,MATBAIJ));
429566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetPreallocation(B,bs,PETSC_DEFAULT,NULL));
439566063dSJacob Faibussowitsch   PetscCall(MatMPIBAIJSetPreallocation(B,bs,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL));
449566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(B));
459566063dSJacob Faibussowitsch   PetscCall(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE));
46c4762a1bSJed Brown 
479566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rdm));
489566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
49c4762a1bSJed Brown 
509566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A,&rstart,&rend));
519566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A,&M,&N));
52c4762a1bSJed Brown   Mbs  = M/bs;
53c4762a1bSJed Brown 
549566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs,&rows));
559566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs,&cols));
569566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs*bs,&vals));
579566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(M,&idx));
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   /* Now set blocks of values */
60c4762a1bSJed Brown   for (i=0; i<40*bs; i++) {
619566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rdm,&rval));
62c4762a1bSJed Brown     cols[0] = bs*(int)(PetscRealPart(rval)*Mbs);
639566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rdm,&rval));
64c4762a1bSJed Brown     rows[0] = rstart + bs*(int)(PetscRealPart(rval)*m);
65c4762a1bSJed Brown     for (j=1; j<bs; j++) {
66c4762a1bSJed Brown       rows[j] = rows[j-1]+1;
67c4762a1bSJed Brown       cols[j] = cols[j-1]+1;
68c4762a1bSJed Brown     }
69c4762a1bSJed Brown 
70c4762a1bSJed Brown     for (j=0; j<bs*bs; j++) {
719566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rdm,&rval));
72c4762a1bSJed Brown       vals[j] = rval;
73c4762a1bSJed Brown     }
749566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A,bs,rows,bs,cols,vals,ADD_VALUES));
759566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B,bs,rows,bs,cols,vals,ADD_VALUES));
76c4762a1bSJed Brown   }
779566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
789566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
799566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
809566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   /* Test MatIncreaseOverlap() */
839566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nd,&is1));
849566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nd,&is2));
85c4762a1bSJed Brown 
8682b5ce2aSStefano Zampini   emptynd = PETSC_FALSE;
8782b5ce2aSStefano Zampini   if (rank == 0 && test_nd0) emptynd = PETSC_TRUE; /* test case */
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   for (i=0; i<nd; i++) {
909566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rdm,&rval));
91c4762a1bSJed Brown     sz   = (int)(PetscRealPart(rval)*m);
92c4762a1bSJed Brown     for (j=0; j<sz; j++) {
939566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rdm,&rval));
94c4762a1bSJed Brown       idx[j*bs] = bs*(int)(PetscRealPart(rval)*Mbs);
95c4762a1bSJed Brown       for (k=1; k<bs; k++) idx[j*bs+k] = idx[j*bs]+k;
96c4762a1bSJed Brown     }
9782b5ce2aSStefano Zampini     PetscCall(ISCreateGeneral(PETSC_COMM_SELF,emptynd ? 0 : sz*bs,idx,PETSC_COPY_VALUES,is1+i));
9882b5ce2aSStefano Zampini     PetscCall(ISCreateGeneral(PETSC_COMM_SELF,emptynd ? 0 : sz*bs,idx,PETSC_COPY_VALUES,is2+i));
99c4762a1bSJed Brown   }
1009566063dSJacob Faibussowitsch   PetscCall(MatIncreaseOverlap(A,nd,is1,ov));
1019566063dSJacob Faibussowitsch   PetscCall(MatIncreaseOverlap(B,nd,is2,ov));
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   for (i=0; i<nd; ++i) {
1049566063dSJacob Faibussowitsch     PetscCall(ISEqual(is1[i],is2[i],&flg));
105c4762a1bSJed Brown 
106c4762a1bSJed Brown     if (!flg) {
1079566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF,"i=%" PetscInt_FMT ", flg=%d :bs=%" PetscInt_FMT " m=%" PetscInt_FMT " ov=%" PetscInt_FMT " nd=%" PetscInt_FMT " np=%d\n",i,flg,bs,m,ov,nd,size));
108c4762a1bSJed Brown     }
109c4762a1bSJed Brown   }
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   for (i=0; i<nd; ++i) {
1129566063dSJacob Faibussowitsch     PetscCall(ISSort(is1[i]));
1139566063dSJacob Faibussowitsch     PetscCall(ISSort(is2[i]));
114c4762a1bSJed Brown   }
115c4762a1bSJed Brown 
1169566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(B,nd,is2,is2,MAT_INITIAL_MATRIX,&submatB));
1179566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(A,nd,is1,is1,MAT_INITIAL_MATRIX,&submatA));
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   /* Test MatMult() */
120c4762a1bSJed Brown   for (i=0; i<nd; i++) {
1219566063dSJacob Faibussowitsch     PetscCall(MatGetSize(submatA[i],&mm,&nn));
1229566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF,mm,&xx));
1239566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx,&s1));
1249566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx,&s2));
125c4762a1bSJed Brown     for (j=0; j<3; j++) {
1269566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(xx,rdm));
1279566063dSJacob Faibussowitsch       PetscCall(MatMult(submatA[i],xx,s1));
1289566063dSJacob Faibussowitsch       PetscCall(MatMult(submatB[i],xx,s2));
1299566063dSJacob Faibussowitsch       PetscCall(VecNorm(s1,NORM_2,&s1norm));
1309566063dSJacob Faibussowitsch       PetscCall(VecNorm(s2,NORM_2,&s2norm));
131c4762a1bSJed Brown       rnorm = s2norm-s1norm;
132c4762a1bSJed Brown       if (rnorm<-tol || rnorm>tol) {
1339566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF,"[%d]Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",rank,(double)s1norm,(double)s2norm));
134c4762a1bSJed Brown       }
135c4762a1bSJed Brown     }
1369566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&xx));
1379566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s1));
1389566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s2));
139c4762a1bSJed Brown   }
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */
1429566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(A,nd,is1,is1,MAT_REUSE_MATRIX,&submatA));
1439566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(B,nd,is2,is2,MAT_REUSE_MATRIX,&submatB));
144c4762a1bSJed Brown 
145c4762a1bSJed Brown   /* Test MatMult() */
146c4762a1bSJed Brown   for (i=0; i<nd; i++) {
1479566063dSJacob Faibussowitsch     PetscCall(MatGetSize(submatA[i],&mm,&nn));
1489566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF,mm,&xx));
1499566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx,&s1));
1509566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx,&s2));
151c4762a1bSJed Brown     for (j=0; j<3; j++) {
1529566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(xx,rdm));
1539566063dSJacob Faibussowitsch       PetscCall(MatMult(submatA[i],xx,s1));
1549566063dSJacob Faibussowitsch       PetscCall(MatMult(submatB[i],xx,s2));
1559566063dSJacob Faibussowitsch       PetscCall(VecNorm(s1,NORM_2,&s1norm));
1569566063dSJacob Faibussowitsch       PetscCall(VecNorm(s2,NORM_2,&s2norm));
157c4762a1bSJed Brown       rnorm = s2norm-s1norm;
158c4762a1bSJed Brown       if (rnorm<-tol || rnorm>tol) {
1599566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF,"[%d]Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",rank,(double)s1norm,(double)s2norm));
160c4762a1bSJed Brown       }
161c4762a1bSJed Brown     }
1629566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&xx));
1639566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s1));
1649566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s2));
165c4762a1bSJed Brown   }
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   /* Free allocated memory */
168c4762a1bSJed Brown   for (i=0; i<nd; ++i) {
1699566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is1[i]));
1709566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is2[i]));
171c4762a1bSJed Brown   }
1729566063dSJacob Faibussowitsch   PetscCall(MatDestroySubMatrices(nd,&submatA));
1739566063dSJacob Faibussowitsch   PetscCall(MatDestroySubMatrices(nd,&submatB));
174c4762a1bSJed Brown 
1759566063dSJacob Faibussowitsch   PetscCall(PetscFree(is1));
1769566063dSJacob Faibussowitsch   PetscCall(PetscFree(is2));
1779566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
1789566063dSJacob Faibussowitsch   PetscCall(PetscFree(rows));
1799566063dSJacob Faibussowitsch   PetscCall(PetscFree(cols));
1809566063dSJacob Faibussowitsch   PetscCall(PetscFree(vals));
1819566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1829566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
1839566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
1849566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
185b122ec5aSJacob Faibussowitsch   return 0;
186c4762a1bSJed Brown }
187c4762a1bSJed Brown 
188c4762a1bSJed Brown /*TEST
189c4762a1bSJed Brown 
190c4762a1bSJed Brown    test:
191c4762a1bSJed Brown       nsize: {{1 3}}
19282b5ce2aSStefano Zampini       args: -mat_block_size {{1 3 4 6 8}} -ov {{1 3}} -mat_size {{11 13}} -nd 7
193c4762a1bSJed Brown       output_file: output/ex54.out
194c4762a1bSJed Brown 
195c4762a1bSJed Brown    test:
196c4762a1bSJed Brown       suffix: 2
197c4762a1bSJed Brown       args: -nd 2 -test_nd0
198c4762a1bSJed Brown       output_file: output/ex54.out
199c4762a1bSJed Brown 
200c4762a1bSJed Brown    test:
201c4762a1bSJed Brown       suffix: 3
202c4762a1bSJed Brown       nsize: 3
203c4762a1bSJed Brown       args: -nd 2 -test_nd0
204c4762a1bSJed Brown       output_file: output/ex54.out
205c4762a1bSJed Brown 
206c4762a1bSJed Brown TEST*/
207