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