xref: /petsc/src/mat/tests/ex54.c (revision 9566063d113dddea24716c546802770db7481bc0)
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;
16c4762a1bSJed Brown   PetscBool      flg,test_nd0=PETSC_FALSE;
17c4762a1bSJed Brown 
18*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
19*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
20*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
21c4762a1bSJed Brown 
22*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL));
23*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_size",&m,NULL));
24*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL));
25*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL));
26*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-test_nd0",&test_nd0,NULL));
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   /* Create a AIJ matrix A */
29*9566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
30*9566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,m*bs,m*bs,PETSC_DECIDE,PETSC_DECIDE));
31*9566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATAIJ));
32*9566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(A,PETSC_DEFAULT,NULL));
33*9566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(A,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL));
34*9566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
35*9566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE));
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   /* Create a BAIJ matrix B */
38*9566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
39*9566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B,m*bs,m*bs,PETSC_DECIDE,PETSC_DECIDE));
40*9566063dSJacob Faibussowitsch   PetscCall(MatSetType(B,MATBAIJ));
41*9566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetPreallocation(B,bs,PETSC_DEFAULT,NULL));
42*9566063dSJacob Faibussowitsch   PetscCall(MatMPIBAIJSetPreallocation(B,bs,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL));
43*9566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(B));
44*9566063dSJacob Faibussowitsch   PetscCall(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE));
45c4762a1bSJed Brown 
46*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rdm));
47*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
48c4762a1bSJed Brown 
49*9566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A,&rstart,&rend));
50*9566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A,&M,&N));
51c4762a1bSJed Brown   Mbs  = M/bs;
52c4762a1bSJed Brown 
53*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs,&rows));
54*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs,&cols));
55*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs*bs,&vals));
56*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(M,&idx));
57c4762a1bSJed Brown 
58c4762a1bSJed Brown   /* Now set blocks of values */
59c4762a1bSJed Brown   for (i=0; i<40*bs; i++) {
60*9566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rdm,&rval));
61c4762a1bSJed Brown     cols[0] = bs*(int)(PetscRealPart(rval)*Mbs);
62*9566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rdm,&rval));
63c4762a1bSJed Brown     rows[0] = rstart + bs*(int)(PetscRealPart(rval)*m);
64c4762a1bSJed Brown     for (j=1; j<bs; j++) {
65c4762a1bSJed Brown       rows[j] = rows[j-1]+1;
66c4762a1bSJed Brown       cols[j] = cols[j-1]+1;
67c4762a1bSJed Brown     }
68c4762a1bSJed Brown 
69c4762a1bSJed Brown     for (j=0; j<bs*bs; j++) {
70*9566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rdm,&rval));
71c4762a1bSJed Brown       vals[j] = rval;
72c4762a1bSJed Brown     }
73*9566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A,bs,rows,bs,cols,vals,ADD_VALUES));
74*9566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B,bs,rows,bs,cols,vals,ADD_VALUES));
75c4762a1bSJed Brown   }
76c4762a1bSJed Brown 
77*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
78*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
79*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
80*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   /* Test MatIncreaseOverlap() */
83*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nd,&is1));
84*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nd,&is2));
85c4762a1bSJed Brown 
86dd400576SPatrick Sanan   if (rank == 0 && test_nd0) nd = 0; /* test case */
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   for (i=0; i<nd; i++) {
89*9566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rdm,&rval));
90c4762a1bSJed Brown     sz   = (int)(PetscRealPart(rval)*m);
91c4762a1bSJed Brown     for (j=0; j<sz; j++) {
92*9566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rdm,&rval));
93c4762a1bSJed Brown       idx[j*bs] = bs*(int)(PetscRealPart(rval)*Mbs);
94c4762a1bSJed Brown       for (k=1; k<bs; k++) idx[j*bs+k] = idx[j*bs]+k;
95c4762a1bSJed Brown     }
96*9566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF,sz*bs,idx,PETSC_COPY_VALUES,is1+i));
97*9566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF,sz*bs,idx,PETSC_COPY_VALUES,is2+i));
98c4762a1bSJed Brown   }
99*9566063dSJacob Faibussowitsch   PetscCall(MatIncreaseOverlap(A,nd,is1,ov));
100*9566063dSJacob Faibussowitsch   PetscCall(MatIncreaseOverlap(B,nd,is2,ov));
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   for (i=0; i<nd; ++i) {
103*9566063dSJacob Faibussowitsch     PetscCall(ISEqual(is1[i],is2[i],&flg));
104c4762a1bSJed Brown 
105c4762a1bSJed Brown     if (!flg) {
106*9566063dSJacob 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));
107c4762a1bSJed Brown     }
108c4762a1bSJed Brown   }
109c4762a1bSJed Brown 
110c4762a1bSJed Brown   for (i=0; i<nd; ++i) {
111*9566063dSJacob Faibussowitsch     PetscCall(ISSort(is1[i]));
112*9566063dSJacob Faibussowitsch     PetscCall(ISSort(is2[i]));
113c4762a1bSJed Brown   }
114c4762a1bSJed Brown 
115*9566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(B,nd,is2,is2,MAT_INITIAL_MATRIX,&submatB));
116*9566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(A,nd,is1,is1,MAT_INITIAL_MATRIX,&submatA));
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   /* Test MatMult() */
119c4762a1bSJed Brown   for (i=0; i<nd; i++) {
120*9566063dSJacob Faibussowitsch     PetscCall(MatGetSize(submatA[i],&mm,&nn));
121*9566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF,mm,&xx));
122*9566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx,&s1));
123*9566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx,&s2));
124c4762a1bSJed Brown     for (j=0; j<3; j++) {
125*9566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(xx,rdm));
126*9566063dSJacob Faibussowitsch       PetscCall(MatMult(submatA[i],xx,s1));
127*9566063dSJacob Faibussowitsch       PetscCall(MatMult(submatB[i],xx,s2));
128*9566063dSJacob Faibussowitsch       PetscCall(VecNorm(s1,NORM_2,&s1norm));
129*9566063dSJacob Faibussowitsch       PetscCall(VecNorm(s2,NORM_2,&s2norm));
130c4762a1bSJed Brown       rnorm = s2norm-s1norm;
131c4762a1bSJed Brown       if (rnorm<-tol || rnorm>tol) {
132*9566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF,"[%d]Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",rank,(double)s1norm,(double)s2norm));
133c4762a1bSJed Brown       }
134c4762a1bSJed Brown     }
135*9566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&xx));
136*9566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s1));
137*9566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s2));
138c4762a1bSJed Brown   }
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */
141*9566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(A,nd,is1,is1,MAT_REUSE_MATRIX,&submatA));
142*9566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(B,nd,is2,is2,MAT_REUSE_MATRIX,&submatB));
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   /* Test MatMult() */
145c4762a1bSJed Brown   for (i=0; i<nd; i++) {
146*9566063dSJacob Faibussowitsch     PetscCall(MatGetSize(submatA[i],&mm,&nn));
147*9566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF,mm,&xx));
148*9566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx,&s1));
149*9566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx,&s2));
150c4762a1bSJed Brown     for (j=0; j<3; j++) {
151*9566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(xx,rdm));
152*9566063dSJacob Faibussowitsch       PetscCall(MatMult(submatA[i],xx,s1));
153*9566063dSJacob Faibussowitsch       PetscCall(MatMult(submatB[i],xx,s2));
154*9566063dSJacob Faibussowitsch       PetscCall(VecNorm(s1,NORM_2,&s1norm));
155*9566063dSJacob Faibussowitsch       PetscCall(VecNorm(s2,NORM_2,&s2norm));
156c4762a1bSJed Brown       rnorm = s2norm-s1norm;
157c4762a1bSJed Brown       if (rnorm<-tol || rnorm>tol) {
158*9566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF,"[%d]Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",rank,(double)s1norm,(double)s2norm));
159c4762a1bSJed Brown       }
160c4762a1bSJed Brown     }
161*9566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&xx));
162*9566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s1));
163*9566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s2));
164c4762a1bSJed Brown   }
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   /* Free allocated memory */
167c4762a1bSJed Brown   for (i=0; i<nd; ++i) {
168*9566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is1[i]));
169*9566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is2[i]));
170c4762a1bSJed Brown   }
171*9566063dSJacob Faibussowitsch   PetscCall(MatDestroySubMatrices(nd,&submatA));
172*9566063dSJacob Faibussowitsch   PetscCall(MatDestroySubMatrices(nd,&submatB));
173c4762a1bSJed Brown 
174*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(is1));
175*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(is2));
176*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
177*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(rows));
178*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(cols));
179*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(vals));
180*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
181*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
182*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
183*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
184b122ec5aSJacob Faibussowitsch   return 0;
185c4762a1bSJed Brown }
186c4762a1bSJed Brown 
187c4762a1bSJed Brown /*TEST
188c4762a1bSJed Brown 
189c4762a1bSJed Brown    test:
190c4762a1bSJed Brown       nsize: {{1 3}}
191c4762a1bSJed Brown       args: -mat_block_size {{1 3 4 6 8}} -ov {{1 3}} -mat_size {{11 13}} -nd {{7}} ; done
192c4762a1bSJed Brown       output_file: output/ex54.out
193c4762a1bSJed Brown 
194c4762a1bSJed Brown    test:
195c4762a1bSJed Brown       suffix: 2
196c4762a1bSJed Brown       args: -nd 2 -test_nd0
197c4762a1bSJed Brown       output_file: output/ex54.out
198c4762a1bSJed Brown 
199c4762a1bSJed Brown    test:
200c4762a1bSJed Brown       suffix: 3
201c4762a1bSJed Brown       nsize: 3
202c4762a1bSJed Brown       args: -nd 2 -test_nd0
203c4762a1bSJed Brown       output_file: output/ex54.out
204c4762a1bSJed Brown 
205c4762a1bSJed Brown TEST*/
206