xref: /petsc/src/mat/tests/ex74.c (revision 9566063d113dddea24716c546802770db7481bc0)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests the various sequential routines in MATSEQSBAIJ format.\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **args)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   PetscMPIInt    size;
9c4762a1bSJed Brown   Vec            x,y,b,s1,s2;
10c4762a1bSJed Brown   Mat            A;                    /* linear system matrix */
11c4762a1bSJed Brown   Mat            sA,sB,sFactor,B,C;    /* symmetric matrices */
12c4762a1bSJed Brown   PetscInt       n,mbs=16,bs=1,nz=3,prob=1,i,j,k1,k2,col[3],lf,block, row,Ii,J,n1,inc;
13c4762a1bSJed Brown   PetscReal      norm1,norm2,rnorm,tol=10*PETSC_SMALL;
14c4762a1bSJed Brown   PetscScalar    neg_one=-1.0,four=4.0,value[3];
15c4762a1bSJed Brown   IS             perm, iscol;
16c4762a1bSJed Brown   PetscRandom    rdm;
17c4762a1bSJed Brown   PetscBool      doIcc=PETSC_TRUE,equal;
18c4762a1bSJed Brown   MatInfo        minfo1,minfo2;
19c4762a1bSJed Brown   MatFactorInfo  factinfo;
20c4762a1bSJed Brown   MatType        type;
21c4762a1bSJed Brown 
22*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
23*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
242c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
25*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL));
26*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL));
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   n    = mbs*bs;
29*9566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&A));
30*9566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE));
31*9566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATSEQBAIJ));
32*9566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
33*9566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetPreallocation(A,bs,nz,NULL));
34c4762a1bSJed Brown 
35*9566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&sA));
36*9566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(sA,n,n,PETSC_DETERMINE,PETSC_DETERMINE));
37*9566063dSJacob Faibussowitsch   PetscCall(MatSetType(sA,MATSEQSBAIJ));
38*9566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(sA));
39*9566063dSJacob Faibussowitsch   PetscCall(MatGetType(sA,&type));
40*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)sA,MATSEQSBAIJ,&doIcc));
41*9566063dSJacob Faibussowitsch   PetscCall(MatSeqSBAIJSetPreallocation(sA,bs,nz,NULL));
42*9566063dSJacob Faibussowitsch   PetscCall(MatSetOption(sA,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE));
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   /* Test MatGetOwnershipRange() */
45*9566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A,&Ii,&J));
46*9566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(sA,&i,&j));
47c4762a1bSJed Brown   if (i-Ii || j-J) {
48*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatGetOwnershipRange() in MatSBAIJ format\n"));
49c4762a1bSJed Brown   }
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   /* Assemble matrix */
52c4762a1bSJed Brown   if (bs == 1) {
53*9566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetInt(NULL,NULL,"-test_problem",&prob,NULL));
54c4762a1bSJed Brown     if (prob == 1) { /* tridiagonal matrix */
55c4762a1bSJed Brown       value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
56c4762a1bSJed Brown       for (i=1; i<n-1; i++) {
57c4762a1bSJed Brown         col[0] = i-1; col[1] = i; col[2] = i+1;
58*9566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES));
59*9566063dSJacob Faibussowitsch         PetscCall(MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES));
60c4762a1bSJed Brown       }
61c4762a1bSJed Brown       i = n - 1; col[0]=0; col[1] = n - 2; col[2] = n - 1;
62c4762a1bSJed Brown 
63c4762a1bSJed Brown       value[0]= 0.1; value[1]=-1; value[2]=2;
64c4762a1bSJed Brown 
65*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES));
66*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES));
67c4762a1bSJed Brown 
68c4762a1bSJed Brown       i        = 0;
69c4762a1bSJed Brown       col[0]   = n-1;   col[1] = 1;      col[2] = 0;
70c4762a1bSJed Brown       value[0] = 0.1; value[1] = -1.0; value[2] = 2;
71c4762a1bSJed Brown 
72*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES));
73*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES));
74c4762a1bSJed Brown 
75c4762a1bSJed Brown     } else if (prob == 2) { /* matrix for the five point stencil */
76c4762a1bSJed Brown       n1 = (PetscInt) (PetscSqrtReal((PetscReal)n) + 0.001);
772c71b3e2SJacob Faibussowitsch       PetscCheckFalse(n1*n1 - n,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"sqrt(n) must be a positive integer!");
78c4762a1bSJed Brown       for (i=0; i<n1; i++) {
79c4762a1bSJed Brown         for (j=0; j<n1; j++) {
80c4762a1bSJed Brown           Ii = j + n1*i;
81c4762a1bSJed Brown           if (i>0) {
82c4762a1bSJed Brown             J    = Ii - n1;
83*9566063dSJacob Faibussowitsch             PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES));
84*9566063dSJacob Faibussowitsch             PetscCall(MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES));
85c4762a1bSJed Brown           }
86c4762a1bSJed Brown           if (i<n1-1) {
87c4762a1bSJed Brown             J    = Ii + n1;
88*9566063dSJacob Faibussowitsch             PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES));
89*9566063dSJacob Faibussowitsch             PetscCall(MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES));
90c4762a1bSJed Brown           }
91c4762a1bSJed Brown           if (j>0) {
92c4762a1bSJed Brown             J    = Ii - 1;
93*9566063dSJacob Faibussowitsch             PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES));
94*9566063dSJacob Faibussowitsch             PetscCall(MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES));
95c4762a1bSJed Brown           }
96c4762a1bSJed Brown           if (j<n1-1) {
97c4762a1bSJed Brown             J    = Ii + 1;
98*9566063dSJacob Faibussowitsch             PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES));
99*9566063dSJacob Faibussowitsch             PetscCall(MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES));
100c4762a1bSJed Brown           }
101*9566063dSJacob Faibussowitsch           PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&four,INSERT_VALUES));
102*9566063dSJacob Faibussowitsch           PetscCall(MatSetValues(sA,1,&Ii,1,&Ii,&four,INSERT_VALUES));
103c4762a1bSJed Brown         }
104c4762a1bSJed Brown       }
105c4762a1bSJed Brown     }
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   } else { /* bs > 1 */
108c4762a1bSJed Brown     for (block=0; block<n/bs; block++) {
109c4762a1bSJed Brown       /* diagonal blocks */
110c4762a1bSJed Brown       value[0] = -1.0; value[1] = 4.0; value[2] = -1.0;
111c4762a1bSJed Brown       for (i=1+block*bs; i<bs-1+block*bs; i++) {
112c4762a1bSJed Brown         col[0] = i-1; col[1] = i; col[2] = i+1;
113*9566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES));
114*9566063dSJacob Faibussowitsch         PetscCall(MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES));
115c4762a1bSJed Brown       }
116c4762a1bSJed Brown       i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs;
117c4762a1bSJed Brown 
118c4762a1bSJed Brown       value[0]=-1.0; value[1]=4.0;
119c4762a1bSJed Brown 
120*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES));
121*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES));
122c4762a1bSJed Brown 
123c4762a1bSJed Brown       i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs;
124c4762a1bSJed Brown 
125c4762a1bSJed Brown       value[0]=4.0; value[1] = -1.0;
126c4762a1bSJed Brown 
127*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES));
128*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES));
129c4762a1bSJed Brown     }
130c4762a1bSJed Brown     /* off-diagonal blocks */
131c4762a1bSJed Brown     value[0]=-1.0;
132c4762a1bSJed Brown     for (i=0; i<(n/bs-1)*bs; i++) {
133c4762a1bSJed Brown       col[0]=i+bs;
134c4762a1bSJed Brown 
135*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,1,col,value,INSERT_VALUES));
136*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(sA,1,&i,1,col,value,INSERT_VALUES));
137c4762a1bSJed Brown 
138c4762a1bSJed Brown       col[0]=i; row=i+bs;
139c4762a1bSJed Brown 
140*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&row,1,col,value,INSERT_VALUES));
141*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(sA,1,&row,1,col,value,INSERT_VALUES));
142c4762a1bSJed Brown     }
143c4762a1bSJed Brown   }
144*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
145*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
146c4762a1bSJed Brown 
147*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(sA,MAT_FINAL_ASSEMBLY));
148*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(sA,MAT_FINAL_ASSEMBLY));
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   /* Test MatGetInfo() of A and sA */
151*9566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(A,MAT_LOCAL,&minfo1));
152*9566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(sA,MAT_LOCAL,&minfo2));
153c4762a1bSJed Brown   i  = (int) (minfo1.nz_used - minfo2.nz_used);
154c4762a1bSJed Brown   j  = (int) (minfo1.nz_allocated - minfo2.nz_allocated);
155c4762a1bSJed Brown   k1 = (int) (minfo1.nz_allocated - minfo1.nz_used);
156c4762a1bSJed Brown   k2 = (int) (minfo2.nz_allocated - minfo2.nz_used);
157c4762a1bSJed Brown   if (i < 0 || j < 0 || k1 < 0 || k2 < 0) {
158*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error (compare A and sA): MatGetInfo()\n"));
159c4762a1bSJed Brown   }
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   /* Test MatDuplicate() */
162*9566063dSJacob Faibussowitsch   PetscCall(MatNorm(A,NORM_FROBENIUS,&norm1));
163*9566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(sA,MAT_COPY_VALUES,&sB));
164*9566063dSJacob Faibussowitsch   PetscCall(MatEqual(sA,sB,&equal));
16528b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatDuplicate()");
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   /* Test MatNorm() */
168*9566063dSJacob Faibussowitsch   PetscCall(MatNorm(A,NORM_FROBENIUS,&norm1));
169*9566063dSJacob Faibussowitsch   PetscCall(MatNorm(sB,NORM_FROBENIUS,&norm2));
170c4762a1bSJed Brown   rnorm = PetscAbsReal(norm1-norm2)/norm2;
171c4762a1bSJed Brown   if (rnorm > tol) {
172*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_FROBENIUS, NormA=%16.14e NormsB=%16.14e\n",(double)norm1,(double)norm2));
173c4762a1bSJed Brown   }
174*9566063dSJacob Faibussowitsch   PetscCall(MatNorm(A,NORM_INFINITY,&norm1));
175*9566063dSJacob Faibussowitsch   PetscCall(MatNorm(sB,NORM_INFINITY,&norm2));
176c4762a1bSJed Brown   rnorm = PetscAbsReal(norm1-norm2)/norm2;
177c4762a1bSJed Brown   if (rnorm > tol) {
178*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_INFINITY(), NormA=%16.14e NormsB=%16.14e\n",(double)norm1,(double)norm2));
179c4762a1bSJed Brown   }
180*9566063dSJacob Faibussowitsch   PetscCall(MatNorm(A,NORM_1,&norm1));
181*9566063dSJacob Faibussowitsch   PetscCall(MatNorm(sB,NORM_1,&norm2));
182c4762a1bSJed Brown   rnorm = PetscAbsReal(norm1-norm2)/norm2;
183c4762a1bSJed Brown   if (rnorm > tol) {
184*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_INFINITY(), NormA=%16.14e NormsB=%16.14e\n",(double)norm1,(double)norm2));
185c4762a1bSJed Brown   }
186c4762a1bSJed Brown 
187c4762a1bSJed Brown   /* Test MatGetInfo(), MatGetSize(), MatGetBlockSize() */
188*9566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(A,MAT_LOCAL,&minfo1));
189*9566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(sB,MAT_LOCAL,&minfo2));
190c4762a1bSJed Brown   i  = (int) (minfo1.nz_used - minfo2.nz_used);
191c4762a1bSJed Brown   j  = (int) (minfo1.nz_allocated - minfo2.nz_allocated);
192c4762a1bSJed Brown   k1 = (int) (minfo1.nz_allocated - minfo1.nz_used);
193c4762a1bSJed Brown   k2 = (int) (minfo2.nz_allocated - minfo2.nz_used);
194c4762a1bSJed Brown   if (i < 0 || j < 0 || k1 < 0 || k2 < 0) {
195*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error(compare A and sB): MatGetInfo()\n"));
196c4762a1bSJed Brown   }
197c4762a1bSJed Brown 
198*9566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A,&Ii,&J));
199*9566063dSJacob Faibussowitsch   PetscCall(MatGetSize(sB,&i,&j));
200c4762a1bSJed Brown   if (i-Ii || j-J) {
201*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatGetSize()\n"));
202c4762a1bSJed Brown   }
203c4762a1bSJed Brown 
204*9566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A, &Ii));
205*9566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(sB, &i));
206c4762a1bSJed Brown   if (i-Ii) {
207*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatGetBlockSize()\n"));
208c4762a1bSJed Brown   }
209c4762a1bSJed Brown 
210*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&rdm));
211*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
212*9566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&x));
213*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&s1));
214*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&s2));
215*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&y));
216*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&b));
217*9566063dSJacob Faibussowitsch   PetscCall(VecSetRandom(x,rdm));
218c4762a1bSJed Brown 
219c4762a1bSJed Brown   /* Test MatDiagonalScale(), MatGetDiagonal(), MatScale() */
220c4762a1bSJed Brown #if !defined(PETSC_USE_COMPLEX)
221c4762a1bSJed Brown   /* Scaling matrix with complex numbers results non-spd matrix,
222c4762a1bSJed Brown      causing crash of MatForwardSolve() and MatBackwardSolve() */
223*9566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(A,x,x));
224*9566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(sB,x,x));
225*9566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(A,sB,10,&equal));
22628b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatDiagonalScale");
227c4762a1bSJed Brown 
228*9566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(A,s1));
229*9566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(sB,s2));
230*9566063dSJacob Faibussowitsch   PetscCall(VecAXPY(s2,neg_one,s1));
231*9566063dSJacob Faibussowitsch   PetscCall(VecNorm(s2,NORM_1,&norm1));
232c4762a1bSJed Brown   if (norm1>tol) {
233*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatGetDiagonal(), ||s1-s2||=%g\n",(double)norm1));
234c4762a1bSJed Brown   }
235c4762a1bSJed Brown 
236c4762a1bSJed Brown   {
237c4762a1bSJed Brown     PetscScalar alpha=0.1;
238*9566063dSJacob Faibussowitsch     PetscCall(MatScale(A,alpha));
239*9566063dSJacob Faibussowitsch     PetscCall(MatScale(sB,alpha));
240c4762a1bSJed Brown   }
241c4762a1bSJed Brown #endif
242c4762a1bSJed Brown 
243c4762a1bSJed Brown   /* Test MatGetRowMaxAbs() */
244*9566063dSJacob Faibussowitsch   PetscCall(MatGetRowMaxAbs(A,s1,NULL));
245*9566063dSJacob Faibussowitsch   PetscCall(MatGetRowMaxAbs(sB,s2,NULL));
246*9566063dSJacob Faibussowitsch   PetscCall(VecNorm(s1,NORM_1,&norm1));
247*9566063dSJacob Faibussowitsch   PetscCall(VecNorm(s2,NORM_1,&norm2));
248c4762a1bSJed Brown   norm1 -= norm2;
249c4762a1bSJed Brown   if (norm1<-tol || norm1>tol) {
250*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatGetRowMaxAbs() \n"));
251c4762a1bSJed Brown   }
252c4762a1bSJed Brown 
253c4762a1bSJed Brown   /* Test MatMult() */
254c4762a1bSJed Brown   for (i=0; i<40; i++) {
255*9566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(x,rdm));
256*9566063dSJacob Faibussowitsch     PetscCall(MatMult(A,x,s1));
257*9566063dSJacob Faibussowitsch     PetscCall(MatMult(sB,x,s2));
258*9566063dSJacob Faibussowitsch     PetscCall(VecNorm(s1,NORM_1,&norm1));
259*9566063dSJacob Faibussowitsch     PetscCall(VecNorm(s2,NORM_1,&norm2));
260c4762a1bSJed Brown     norm1 -= norm2;
261c4762a1bSJed Brown     if (norm1<-tol || norm1>tol) {
262*9566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatMult(), norm1-norm2: %g\n",(double)norm1));
263c4762a1bSJed Brown     }
264c4762a1bSJed Brown   }
265c4762a1bSJed Brown 
266c4762a1bSJed Brown   /* MatMultAdd() */
267c4762a1bSJed Brown   for (i=0; i<40; i++) {
268*9566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(x,rdm));
269*9566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(y,rdm));
270*9566063dSJacob Faibussowitsch     PetscCall(MatMultAdd(A,x,y,s1));
271*9566063dSJacob Faibussowitsch     PetscCall(MatMultAdd(sB,x,y,s2));
272*9566063dSJacob Faibussowitsch     PetscCall(VecNorm(s1,NORM_1,&norm1));
273*9566063dSJacob Faibussowitsch     PetscCall(VecNorm(s2,NORM_1,&norm2));
274c4762a1bSJed Brown     norm1 -= norm2;
275c4762a1bSJed Brown     if (norm1<-tol || norm1>tol) {
276*9566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatMultAdd(), norm1-norm2: %g\n",(double)norm1));
277c4762a1bSJed Brown     }
278c4762a1bSJed Brown   }
279c4762a1bSJed Brown 
280c20d7725SJed Brown   /* Test MatMatMult() for sbaij and dense matrices */
281*9566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,n,5*n,NULL,&B));
282*9566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(B,rdm));
283*9566063dSJacob Faibussowitsch   PetscCall(MatMatMult(sA,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C));
284*9566063dSJacob Faibussowitsch   PetscCall(MatMatMultEqual(sA,B,C,5*n,&equal));
28528b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error: MatMatMult()");
286*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
287*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
288c4762a1bSJed Brown 
289c4762a1bSJed Brown   /* Test MatCholeskyFactor(), MatICCFactor() with natural ordering */
290*9566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(A,MATORDERINGNATURAL,&perm,&iscol));
291*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscol));
292c4762a1bSJed Brown   norm1 = tol;
293c4762a1bSJed Brown   inc   = bs;
294c4762a1bSJed Brown 
295c4762a1bSJed Brown   /* initialize factinfo */
296*9566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(&factinfo,sizeof(MatFactorInfo)));
297c4762a1bSJed Brown 
298c4762a1bSJed Brown   for (lf=-1; lf<10; lf += inc) {
299c4762a1bSJed Brown     if (lf==-1) {  /* Cholesky factor of sB (duplicate sA) */
300c4762a1bSJed Brown       factinfo.fill = 5.0;
301c4762a1bSJed Brown 
302*9566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(sB,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&sFactor));
303*9566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorSymbolic(sFactor,sB,perm,&factinfo));
304c4762a1bSJed Brown     } else if (!doIcc) break;
305c4762a1bSJed Brown     else {       /* incomplete Cholesky factor */
306c4762a1bSJed Brown       factinfo.fill   = 5.0;
307c4762a1bSJed Brown       factinfo.levels = lf;
308c4762a1bSJed Brown 
309*9566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(sB,MATSOLVERPETSC,MAT_FACTOR_ICC,&sFactor));
310*9566063dSJacob Faibussowitsch       PetscCall(MatICCFactorSymbolic(sFactor,sB,perm,&factinfo));
311c4762a1bSJed Brown     }
312*9566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactorNumeric(sFactor,sB,&factinfo));
313c4762a1bSJed Brown     /* MatView(sFactor, PETSC_VIEWER_DRAW_WORLD); */
314c4762a1bSJed Brown 
315c4762a1bSJed Brown     /* test MatGetDiagonal on numeric factor */
316c4762a1bSJed Brown     /*
317c4762a1bSJed Brown     if (lf == -1) {
318*9566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonal(sFactor,s1));
319c4762a1bSJed Brown       printf(" in ex74.c, diag: \n");
320*9566063dSJacob Faibussowitsch       PetscCall(VecView(s1,PETSC_VIEWER_STDOUT_SELF));
321c4762a1bSJed Brown     }
322c4762a1bSJed Brown     */
323c4762a1bSJed Brown 
324*9566063dSJacob Faibussowitsch     PetscCall(MatMult(sB,x,b));
325c4762a1bSJed Brown 
326c4762a1bSJed Brown     /* test MatForwardSolve() and MatBackwardSolve() */
327c4762a1bSJed Brown     if (lf == -1) {
328*9566063dSJacob Faibussowitsch       PetscCall(MatForwardSolve(sFactor,b,s1));
329*9566063dSJacob Faibussowitsch       PetscCall(MatBackwardSolve(sFactor,s1,s2));
330*9566063dSJacob Faibussowitsch       PetscCall(VecAXPY(s2,neg_one,x));
331*9566063dSJacob Faibussowitsch       PetscCall(VecNorm(s2,NORM_2,&norm2));
332c4762a1bSJed Brown       if (10*norm1 < norm2) {
333*9566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF,"MatForwardSolve and BackwardSolve: Norm of error=%g, bs=%" PetscInt_FMT "\n",(double)norm2,bs));
334c4762a1bSJed Brown       }
335c4762a1bSJed Brown     }
336c4762a1bSJed Brown 
337c4762a1bSJed Brown     /* test MatSolve() */
338*9566063dSJacob Faibussowitsch     PetscCall(MatSolve(sFactor,b,y));
339*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&sFactor));
340c4762a1bSJed Brown     /* Check the error */
341*9566063dSJacob Faibussowitsch     PetscCall(VecAXPY(y,neg_one,x));
342*9566063dSJacob Faibussowitsch     PetscCall(VecNorm(y,NORM_2,&norm2));
343c4762a1bSJed Brown     if (10*norm1 < norm2 && lf-inc != -1) {
344*9566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF,"lf=%" PetscInt_FMT ", %" PetscInt_FMT ", Norm of error=%g, %g\n",lf-inc,lf,(double)norm1,(double)norm2));
345c4762a1bSJed Brown     }
346c4762a1bSJed Brown     norm1 = norm2;
347c4762a1bSJed Brown     if (norm2 < tol && lf != -1) break;
348c4762a1bSJed Brown   }
349c4762a1bSJed Brown 
350c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS)
351*9566063dSJacob Faibussowitsch   PetscCall(MatGetFactor(sA,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&sFactor));
352*9566063dSJacob Faibussowitsch   PetscCall(MatCholeskyFactorSymbolic(sFactor,sA,NULL,NULL));
353*9566063dSJacob Faibussowitsch   PetscCall(MatCholeskyFactorNumeric(sFactor,sA,NULL));
354c4762a1bSJed Brown   for (i=0; i<10; i++) {
355*9566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(b,rdm));
356*9566063dSJacob Faibussowitsch     PetscCall(MatSolve(sFactor,b,y));
357c4762a1bSJed Brown     /* Check the error */
358*9566063dSJacob Faibussowitsch     PetscCall(MatMult(sA,y,x));
359*9566063dSJacob Faibussowitsch     PetscCall(VecAXPY(x,neg_one,b));
360*9566063dSJacob Faibussowitsch     PetscCall(VecNorm(x,NORM_2,&norm2));
361c4762a1bSJed Brown     if (norm2>tol) {
362*9566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatSolve(), norm2: %g\n",(double)norm2));
363c4762a1bSJed Brown     }
364c4762a1bSJed Brown   }
365*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sFactor));
366c4762a1bSJed Brown #endif
367c4762a1bSJed Brown 
368*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
369c4762a1bSJed Brown 
370*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
371*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sB));
372*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sA));
373*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
374*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
375*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s1));
376*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s2));
377*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
378*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
379c4762a1bSJed Brown 
380*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
381b122ec5aSJacob Faibussowitsch   return 0;
382c4762a1bSJed Brown }
383c4762a1bSJed Brown 
384c4762a1bSJed Brown /*TEST
385c4762a1bSJed Brown 
386c4762a1bSJed Brown    test:
387c4762a1bSJed Brown       args: -bs {{1 2 3 4 5 6 7 8}}
388c4762a1bSJed Brown 
389c4762a1bSJed Brown TEST*/
390