xref: /petsc/src/mat/tests/ex91.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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   PetscErrorCode ierr;
11c4762a1bSJed Brown   PetscMPIInt    size;
12c4762a1bSJed Brown   PetscScalar    *vals,rval,one=1.0;
13c4762a1bSJed Brown   IS             *is1,*is2;
14c4762a1bSJed Brown   PetscRandom    rand;
15c4762a1bSJed Brown   Vec            xx,s1,s2;
16c4762a1bSJed Brown   PetscReal      s1norm,s2norm,rnorm,tol = 10*PETSC_SMALL;
17c4762a1bSJed Brown   PetscBool      flg;
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
20*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL));
21*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mat_size",&m,NULL));
22*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL));
23*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL));
24c4762a1bSJed Brown 
25c4762a1bSJed Brown   /* create a SeqBAIJ matrix A */
26c4762a1bSJed Brown   M    = m*bs;
27*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqBAIJ(PETSC_COMM_SELF,bs,M,M,1,NULL,&A));
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&rand));
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rand));
31c4762a1bSJed Brown 
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(bs,&rows));
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(bs,&cols));
34*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(bs*bs,&vals));
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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     }
46*5f80ce2aSJacob Faibussowitsch     CHKERRQ(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++) {
50*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rand,&rval));
51c4762a1bSJed Brown     cols[0] = bs*(int)(PetscRealPart(rval)*m);
52*5f80ce2aSJacob Faibussowitsch     CHKERRQ(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++) {
60*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRandomGetValue(rand,&rval));
61c4762a1bSJed Brown       vals[j] = rval;
62c4762a1bSJed Brown     }
63*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,bs,rows,bs,cols,vals,ADD_VALUES));
64c4762a1bSJed Brown   }
65c4762a1bSJed Brown 
66*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
67*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   /* make A a symmetric matrix: A <- A^T + A */
70*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX, &Atrans));
71*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAXPY(A,one,Atrans,DIFFERENT_NONZERO_PATTERN));
72*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&Atrans));
73*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX, &Atrans));
74*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatEqual(A, Atrans, &flg));
752c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"A+A^T is non-symmetric");
76*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&Atrans));
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   /* create a SeqSBAIJ matrix sA (= A) */
79*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE));
80*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&sA));
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   /* Test sA==A through MatMult() */
83c4762a1bSJed Brown   for (i=0; i<nd; i++) {
84*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetSize(A,&mm,&nn));
85*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,mm,&xx));
86*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(xx,&s1));
87*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(xx,&s2));
88c4762a1bSJed Brown     for (j=0; j<3; j++) {
89*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecSetRandom(xx,rand));
90*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(A,xx,s1));
91*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(sA,xx,s2));
92*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(s1,NORM_2,&s1norm));
93*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(s2,NORM_2,&s2norm));
94c4762a1bSJed Brown       rnorm = s2norm-s1norm;
95c4762a1bSJed Brown       if (rnorm<-tol || rnorm>tol) {
96*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",(double)s1norm,(double)s2norm));
97c4762a1bSJed Brown       }
98c4762a1bSJed Brown     }
99*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&xx));
100*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&s1));
101*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&s2));
102c4762a1bSJed Brown   }
103c4762a1bSJed Brown 
104c4762a1bSJed Brown   /* Test MatIncreaseOverlap() */
105*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nd,&is1));
106*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nd,&is2));
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   for (i=0; i<nd; i++) {
109*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rand,&rval));
110c4762a1bSJed Brown     size = (int)(PetscRealPart(rval)*m);
111c4762a1bSJed Brown     for (j=0; j<size; j++) {
112*5f80ce2aSJacob Faibussowitsch       CHKERRQ(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     }
116*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,size*bs,idx,PETSC_COPY_VALUES,is1+i));
117*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,size*bs,idx,PETSC_COPY_VALUES,is2+i));
118c4762a1bSJed Brown   }
119c4762a1bSJed Brown   /* for debugging */
120c4762a1bSJed Brown   /*
121*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_SELF));
122*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(sA,PETSC_VIEWER_STDOUT_SELF));
123c4762a1bSJed Brown   */
124c4762a1bSJed Brown 
125*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatIncreaseOverlap(A,nd,is1,ov));
126*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatIncreaseOverlap(sA,nd,is2,ov));
127c4762a1bSJed Brown 
128c4762a1bSJed Brown   for (i=0; i<nd; ++i) {
129*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISSort(is1[i]));
130*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISSort(is2[i]));
131c4762a1bSJed Brown   }
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   for (i=0; i<nd; ++i) {
134*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISEqual(is1[i],is2[i],&flg));
1352c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"i=%" PetscInt_FMT ", is1 != is2",i);
136c4762a1bSJed Brown   }
137c4762a1bSJed Brown 
138*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSubMatrices(A,nd,is1,is1,MAT_INITIAL_MATRIX,&submatA));
139*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSubMatrices(sA,nd,is2,is2,MAT_INITIAL_MATRIX,&submatsA));
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   /* Test MatMult() */
142c4762a1bSJed Brown   for (i=0; i<nd; i++) {
143*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetSize(submatA[i],&mm,&nn));
144*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,mm,&xx));
145*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(xx,&s1));
146*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(xx,&s2));
147c4762a1bSJed Brown     for (j=0; j<3; j++) {
148*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecSetRandom(xx,rand));
149*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(submatA[i],xx,s1));
150*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(submatsA[i],xx,s2));
151*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(s1,NORM_2,&s1norm));
152*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(s2,NORM_2,&s2norm));
153c4762a1bSJed Brown       rnorm = s2norm-s1norm;
154c4762a1bSJed Brown       if (rnorm<-tol || rnorm>tol) {
155*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",(double)s1norm,(double)s2norm));
156c4762a1bSJed Brown       }
157c4762a1bSJed Brown     }
158*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&xx));
159*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&s1));
160*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&s2));
161c4762a1bSJed Brown   }
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */
164*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSubMatrices(A,nd,is1,is1,MAT_REUSE_MATRIX,&submatA));
165*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSubMatrices(sA,nd,is2,is2,MAT_REUSE_MATRIX,&submatsA));
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   /* Test MatMult() */
168c4762a1bSJed Brown   for (i=0; i<nd; i++) {
169*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetSize(submatA[i],&mm,&nn));
170*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,mm,&xx));
171*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(xx,&s1));
172*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(xx,&s2));
173c4762a1bSJed Brown     for (j=0; j<3; j++) {
174*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecSetRandom(xx,rand));
175*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(submatA[i],xx,s1));
176*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(submatsA[i],xx,s2));
177*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(s1,NORM_2,&s1norm));
178*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(s2,NORM_2,&s2norm));
179c4762a1bSJed Brown       rnorm = s2norm-s1norm;
180c4762a1bSJed Brown       if (rnorm<-tol || rnorm>tol) {
181*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",(double)s1norm,(double)s2norm));
182c4762a1bSJed Brown       }
183c4762a1bSJed Brown     }
184*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&xx));
185*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&s1));
186*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&s2));
187c4762a1bSJed Brown   }
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   /* Free allocated memory */
190c4762a1bSJed Brown   for (i=0; i<nd; ++i) {
191*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is1[i]));
192*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is2[i]));
193c4762a1bSJed Brown   }
194*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroySubMatrices(nd,&submatA));
195*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroySubMatrices(nd,&submatsA));
196c4762a1bSJed Brown 
197*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(is1));
198*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(is2));
199*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(idx));
200*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(rows));
201*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(cols));
202*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(vals));
203*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
204*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&sA));
205*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rand));
206c4762a1bSJed Brown   ierr = PetscFinalize();
207c4762a1bSJed Brown   return ierr;
208c4762a1bSJed Brown }
209c4762a1bSJed Brown 
210c4762a1bSJed Brown /*TEST
211c4762a1bSJed Brown 
212c4762a1bSJed Brown    test:
213c4762a1bSJed Brown       args: -ov 2
214c4762a1bSJed Brown 
215c4762a1bSJed Brown TEST*/
216