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