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