xref: /petsc/src/mat/tests/ex91.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatIncreaseOverlap(), MatCreateSubMatrices() for sequential MatSBAIJ format. Derived from ex51.c\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **args)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   Mat            A,Atrans,sA,*submatA,*submatsA;
9c4762a1bSJed Brown   PetscInt       bs=1,m=43,ov=1,i,j,k,*rows,*cols,M,nd=5,*idx,mm,nn;
10c4762a1bSJed Brown   PetscMPIInt    size;
11c4762a1bSJed Brown   PetscScalar    *vals,rval,one=1.0;
12c4762a1bSJed Brown   IS             *is1,*is2;
13c4762a1bSJed Brown   PetscRandom    rand;
14c4762a1bSJed Brown   Vec            xx,s1,s2;
15c4762a1bSJed Brown   PetscReal      s1norm,s2norm,rnorm,tol = 10*PETSC_SMALL;
16c4762a1bSJed Brown   PetscBool      flg;
17c4762a1bSJed Brown 
18*327415f7SBarry Smith   PetscFunctionBeginUser;
199566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL));
219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_size",&m,NULL));
229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL));
239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL));
24c4762a1bSJed Brown 
25c4762a1bSJed Brown   /* create a SeqBAIJ matrix A */
26c4762a1bSJed Brown   M    = m*bs;
279566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqBAIJ(PETSC_COMM_SELF,bs,M,M,1,NULL,&A));
289566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
299566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&rand));
309566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rand));
31c4762a1bSJed Brown 
329566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs,&rows));
339566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs,&cols));
349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs*bs,&vals));
359566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(M,&idx));
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   /* Now set blocks of random values */
38c4762a1bSJed Brown   /* first, set diagonal blocks as zero */
39c4762a1bSJed Brown   for (j=0; j<bs*bs; j++) vals[j] = 0.0;
40c4762a1bSJed Brown   for (i=0; i<m; i++) {
41c4762a1bSJed Brown     cols[0] = i*bs; rows[0] = i*bs;
42c4762a1bSJed Brown     for (j=1; j<bs; j++) {
43c4762a1bSJed Brown       rows[j] = rows[j-1]+1;
44c4762a1bSJed Brown       cols[j] = cols[j-1]+1;
45c4762a1bSJed Brown     }
469566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A,bs,rows,bs,cols,vals,ADD_VALUES));
47c4762a1bSJed Brown   }
48c4762a1bSJed Brown   /* second, add random blocks */
49c4762a1bSJed Brown   for (i=0; i<20*bs; i++) {
509566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rand,&rval));
51c4762a1bSJed Brown     cols[0] = bs*(int)(PetscRealPart(rval)*m);
529566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rand,&rval));
53c4762a1bSJed Brown     rows[0] = bs*(int)(PetscRealPart(rval)*m);
54c4762a1bSJed Brown     for (j=1; j<bs; j++) {
55c4762a1bSJed Brown       rows[j] = rows[j-1]+1;
56c4762a1bSJed Brown       cols[j] = cols[j-1]+1;
57c4762a1bSJed Brown     }
58c4762a1bSJed Brown 
59c4762a1bSJed Brown     for (j=0; j<bs*bs; j++) {
609566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rand,&rval));
61c4762a1bSJed Brown       vals[j] = rval;
62c4762a1bSJed Brown     }
639566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A,bs,rows,bs,cols,vals,ADD_VALUES));
64c4762a1bSJed Brown   }
65c4762a1bSJed Brown 
669566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
679566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   /* make A a symmetric matrix: A <- A^T + A */
709566063dSJacob Faibussowitsch   PetscCall(MatTranspose(A,MAT_INITIAL_MATRIX, &Atrans));
719566063dSJacob Faibussowitsch   PetscCall(MatAXPY(A,one,Atrans,DIFFERENT_NONZERO_PATTERN));
729566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Atrans));
739566063dSJacob Faibussowitsch   PetscCall(MatTranspose(A,MAT_INITIAL_MATRIX, &Atrans));
749566063dSJacob Faibussowitsch   PetscCall(MatEqual(A, Atrans, &flg));
7528b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"A+A^T is non-symmetric");
769566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Atrans));
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   /* create a SeqSBAIJ matrix sA (= A) */
799566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE));
809566063dSJacob Faibussowitsch   PetscCall(MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&sA));
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   /* Test sA==A through MatMult() */
83c4762a1bSJed Brown   for (i=0; i<nd; i++) {
849566063dSJacob Faibussowitsch     PetscCall(MatGetSize(A,&mm,&nn));
859566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF,mm,&xx));
869566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx,&s1));
879566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx,&s2));
88c4762a1bSJed Brown     for (j=0; j<3; j++) {
899566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(xx,rand));
909566063dSJacob Faibussowitsch       PetscCall(MatMult(A,xx,s1));
919566063dSJacob Faibussowitsch       PetscCall(MatMult(sA,xx,s2));
929566063dSJacob Faibussowitsch       PetscCall(VecNorm(s1,NORM_2,&s1norm));
939566063dSJacob Faibussowitsch       PetscCall(VecNorm(s2,NORM_2,&s2norm));
94c4762a1bSJed Brown       rnorm = s2norm-s1norm;
95c4762a1bSJed Brown       if (rnorm<-tol || rnorm>tol) {
969566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",(double)s1norm,(double)s2norm));
97c4762a1bSJed Brown       }
98c4762a1bSJed Brown     }
999566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&xx));
1009566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s1));
1019566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s2));
102c4762a1bSJed Brown   }
103c4762a1bSJed Brown 
104c4762a1bSJed Brown   /* Test MatIncreaseOverlap() */
1059566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nd,&is1));
1069566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nd,&is2));
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   for (i=0; i<nd; i++) {
1099566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rand,&rval));
110c4762a1bSJed Brown     size = (int)(PetscRealPart(rval)*m);
111c4762a1bSJed Brown     for (j=0; j<size; j++) {
1129566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rand,&rval));
113c4762a1bSJed Brown       idx[j*bs] = bs*(int)(PetscRealPart(rval)*m);
114c4762a1bSJed Brown       for (k=1; k<bs; k++) idx[j*bs+k] = idx[j*bs]+k;
115c4762a1bSJed Brown     }
1169566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF,size*bs,idx,PETSC_COPY_VALUES,is1+i));
1179566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF,size*bs,idx,PETSC_COPY_VALUES,is2+i));
118c4762a1bSJed Brown   }
119c4762a1bSJed Brown   /* for debugging */
120c4762a1bSJed Brown   /*
1219566063dSJacob Faibussowitsch   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_SELF));
1229566063dSJacob Faibussowitsch   PetscCall(MatView(sA,PETSC_VIEWER_STDOUT_SELF));
123c4762a1bSJed Brown   */
124c4762a1bSJed Brown 
1259566063dSJacob Faibussowitsch   PetscCall(MatIncreaseOverlap(A,nd,is1,ov));
1269566063dSJacob Faibussowitsch   PetscCall(MatIncreaseOverlap(sA,nd,is2,ov));
127c4762a1bSJed Brown 
128c4762a1bSJed Brown   for (i=0; i<nd; ++i) {
1299566063dSJacob Faibussowitsch     PetscCall(ISSort(is1[i]));
1309566063dSJacob Faibussowitsch     PetscCall(ISSort(is2[i]));
131c4762a1bSJed Brown   }
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   for (i=0; i<nd; ++i) {
1349566063dSJacob Faibussowitsch     PetscCall(ISEqual(is1[i],is2[i],&flg));
13528b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"i=%" PetscInt_FMT ", is1 != is2",i);
136c4762a1bSJed Brown   }
137c4762a1bSJed Brown 
1389566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(A,nd,is1,is1,MAT_INITIAL_MATRIX,&submatA));
1399566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(sA,nd,is2,is2,MAT_INITIAL_MATRIX,&submatsA));
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   /* Test MatMult() */
142c4762a1bSJed Brown   for (i=0; i<nd; i++) {
1439566063dSJacob Faibussowitsch     PetscCall(MatGetSize(submatA[i],&mm,&nn));
1449566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF,mm,&xx));
1459566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx,&s1));
1469566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx,&s2));
147c4762a1bSJed Brown     for (j=0; j<3; j++) {
1489566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(xx,rand));
1499566063dSJacob Faibussowitsch       PetscCall(MatMult(submatA[i],xx,s1));
1509566063dSJacob Faibussowitsch       PetscCall(MatMult(submatsA[i],xx,s2));
1519566063dSJacob Faibussowitsch       PetscCall(VecNorm(s1,NORM_2,&s1norm));
1529566063dSJacob Faibussowitsch       PetscCall(VecNorm(s2,NORM_2,&s2norm));
153c4762a1bSJed Brown       rnorm = s2norm-s1norm;
154c4762a1bSJed Brown       if (rnorm<-tol || rnorm>tol) {
1559566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",(double)s1norm,(double)s2norm));
156c4762a1bSJed Brown       }
157c4762a1bSJed Brown     }
1589566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&xx));
1599566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s1));
1609566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s2));
161c4762a1bSJed Brown   }
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */
1649566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(A,nd,is1,is1,MAT_REUSE_MATRIX,&submatA));
1659566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(sA,nd,is2,is2,MAT_REUSE_MATRIX,&submatsA));
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   /* Test MatMult() */
168c4762a1bSJed Brown   for (i=0; i<nd; i++) {
1699566063dSJacob Faibussowitsch     PetscCall(MatGetSize(submatA[i],&mm,&nn));
1709566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF,mm,&xx));
1719566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx,&s1));
1729566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx,&s2));
173c4762a1bSJed Brown     for (j=0; j<3; j++) {
1749566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(xx,rand));
1759566063dSJacob Faibussowitsch       PetscCall(MatMult(submatA[i],xx,s1));
1769566063dSJacob Faibussowitsch       PetscCall(MatMult(submatsA[i],xx,s2));
1779566063dSJacob Faibussowitsch       PetscCall(VecNorm(s1,NORM_2,&s1norm));
1789566063dSJacob Faibussowitsch       PetscCall(VecNorm(s2,NORM_2,&s2norm));
179c4762a1bSJed Brown       rnorm = s2norm-s1norm;
180c4762a1bSJed Brown       if (rnorm<-tol || rnorm>tol) {
1819566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",(double)s1norm,(double)s2norm));
182c4762a1bSJed Brown       }
183c4762a1bSJed Brown     }
1849566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&xx));
1859566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s1));
1869566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s2));
187c4762a1bSJed Brown   }
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   /* Free allocated memory */
190c4762a1bSJed Brown   for (i=0; i<nd; ++i) {
1919566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is1[i]));
1929566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is2[i]));
193c4762a1bSJed Brown   }
1949566063dSJacob Faibussowitsch   PetscCall(MatDestroySubMatrices(nd,&submatA));
1959566063dSJacob Faibussowitsch   PetscCall(MatDestroySubMatrices(nd,&submatsA));
196c4762a1bSJed Brown 
1979566063dSJacob Faibussowitsch   PetscCall(PetscFree(is1));
1989566063dSJacob Faibussowitsch   PetscCall(PetscFree(is2));
1999566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
2009566063dSJacob Faibussowitsch   PetscCall(PetscFree(rows));
2019566063dSJacob Faibussowitsch   PetscCall(PetscFree(cols));
2029566063dSJacob Faibussowitsch   PetscCall(PetscFree(vals));
2039566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2049566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sA));
2059566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rand));
2069566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
207b122ec5aSJacob Faibussowitsch   return 0;
208c4762a1bSJed Brown }
209c4762a1bSJed Brown 
210c4762a1bSJed Brown /*TEST
211c4762a1bSJed Brown 
212c4762a1bSJed Brown    test:
213c4762a1bSJed Brown       args: -ov 2
214c4762a1bSJed Brown 
215c4762a1bSJed Brown TEST*/
216