xref: /petsc/src/mat/tests/ex145.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Tests LU, Cholesky factorization and MatMatSolve() for an Elemental dense matrix.\n\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscmat.h>
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown int main(int argc,char **argv)
7*c4762a1bSJed Brown {
8*c4762a1bSJed Brown   Mat            A,F,B,X,C,Aher,G;
9*c4762a1bSJed Brown   Vec            b,x,c,d,e;
10*c4762a1bSJed Brown   PetscErrorCode ierr;
11*c4762a1bSJed Brown   PetscInt       m = 5,n,p,i,j,nrows,ncols;
12*c4762a1bSJed Brown   PetscScalar    *v,*barray,rval;
13*c4762a1bSJed Brown   PetscReal      norm,tol=1.e-11;
14*c4762a1bSJed Brown   PetscMPIInt    size,rank;
15*c4762a1bSJed Brown   PetscRandom    rand;
16*c4762a1bSJed Brown   const PetscInt *rows,*cols;
17*c4762a1bSJed Brown   IS             isrows,iscols;
18*c4762a1bSJed Brown   PetscBool      mats_view=PETSC_FALSE;
19*c4762a1bSJed Brown   MatFactorInfo  finfo;
20*c4762a1bSJed Brown 
21*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*) 0,help);if (ierr) return ierr;
22*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
23*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
24*c4762a1bSJed Brown 
25*c4762a1bSJed Brown   ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rand);CHKERRQ(ierr);
26*c4762a1bSJed Brown   ierr = PetscRandomSetFromOptions(rand);CHKERRQ(ierr);
27*c4762a1bSJed Brown 
28*c4762a1bSJed Brown   /* Get local dimensions of matrices */
29*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr);
30*c4762a1bSJed Brown   n    = m;
31*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr);
32*c4762a1bSJed Brown   p    = m/2;
33*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL);CHKERRQ(ierr);
34*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-mats_view",&mats_view);CHKERRQ(ierr);
35*c4762a1bSJed Brown 
36*c4762a1bSJed Brown   /* Create matrix A */
37*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD," Create Elemental matrix A\n");CHKERRQ(ierr);
38*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
39*c4762a1bSJed Brown   ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
40*c4762a1bSJed Brown   ierr = MatSetType(A,MATELEMENTAL);CHKERRQ(ierr);
41*c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
42*c4762a1bSJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
43*c4762a1bSJed Brown   /* Set local matrix entries */
44*c4762a1bSJed Brown   ierr = MatGetOwnershipIS(A,&isrows,&iscols);CHKERRQ(ierr);
45*c4762a1bSJed Brown   ierr = ISGetLocalSize(isrows,&nrows);CHKERRQ(ierr);
46*c4762a1bSJed Brown   ierr = ISGetIndices(isrows,&rows);CHKERRQ(ierr);
47*c4762a1bSJed Brown   ierr = ISGetLocalSize(iscols,&ncols);CHKERRQ(ierr);
48*c4762a1bSJed Brown   ierr = ISGetIndices(iscols,&cols);CHKERRQ(ierr);
49*c4762a1bSJed Brown   ierr = PetscMalloc1(nrows*ncols,&v);CHKERRQ(ierr);
50*c4762a1bSJed Brown   for (i=0; i<nrows; i++) {
51*c4762a1bSJed Brown     for (j=0; j<ncols; j++) {
52*c4762a1bSJed Brown       ierr         = PetscRandomGetValue(rand,&rval);CHKERRQ(ierr);
53*c4762a1bSJed Brown       v[i*ncols+j] = rval;
54*c4762a1bSJed Brown     }
55*c4762a1bSJed Brown   }
56*c4762a1bSJed Brown   ierr = MatSetValues(A,nrows,rows,ncols,cols,v,INSERT_VALUES);CHKERRQ(ierr);
57*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
58*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
59*c4762a1bSJed Brown   ierr = ISRestoreIndices(isrows,&rows);CHKERRQ(ierr);
60*c4762a1bSJed Brown   ierr = ISRestoreIndices(iscols,&cols);CHKERRQ(ierr);
61*c4762a1bSJed Brown   ierr = ISDestroy(&isrows);CHKERRQ(ierr);
62*c4762a1bSJed Brown   ierr = ISDestroy(&iscols);CHKERRQ(ierr);
63*c4762a1bSJed Brown   ierr = PetscFree(v);CHKERRQ(ierr);
64*c4762a1bSJed Brown   if (mats_view) {
65*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD, "A: nrows %d, m %d; ncols %d, n %d\n",nrows,m,ncols,n);CHKERRQ(ierr);
66*c4762a1bSJed Brown     ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
67*c4762a1bSJed Brown   }
68*c4762a1bSJed Brown 
69*c4762a1bSJed Brown   /* Create rhs matrix B */
70*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD," Create rhs matrix B\n");CHKERRQ(ierr);
71*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr);
72*c4762a1bSJed Brown   ierr = MatSetSizes(B,m,p,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
73*c4762a1bSJed Brown   ierr = MatSetType(B,MATELEMENTAL);CHKERRQ(ierr);
74*c4762a1bSJed Brown   ierr = MatSetFromOptions(B);CHKERRQ(ierr);
75*c4762a1bSJed Brown   ierr = MatSetUp(B);CHKERRQ(ierr);
76*c4762a1bSJed Brown   ierr = MatGetOwnershipIS(B,&isrows,&iscols);CHKERRQ(ierr);
77*c4762a1bSJed Brown   ierr = ISGetLocalSize(isrows,&nrows);CHKERRQ(ierr);
78*c4762a1bSJed Brown   ierr = ISGetIndices(isrows,&rows);CHKERRQ(ierr);
79*c4762a1bSJed Brown   ierr = ISGetLocalSize(iscols,&ncols);CHKERRQ(ierr);
80*c4762a1bSJed Brown   ierr = ISGetIndices(iscols,&cols);CHKERRQ(ierr);
81*c4762a1bSJed Brown   ierr = PetscMalloc1(nrows*ncols,&v);CHKERRQ(ierr);
82*c4762a1bSJed Brown   for (i=0; i<nrows; i++) {
83*c4762a1bSJed Brown     for (j=0; j<ncols; j++) {
84*c4762a1bSJed Brown       ierr         = PetscRandomGetValue(rand,&rval);CHKERRQ(ierr);
85*c4762a1bSJed Brown       v[i*ncols+j] = rval;
86*c4762a1bSJed Brown     }
87*c4762a1bSJed Brown   }
88*c4762a1bSJed Brown   ierr = MatSetValues(B,nrows,rows,ncols,cols,v,INSERT_VALUES);CHKERRQ(ierr);
89*c4762a1bSJed Brown   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
90*c4762a1bSJed Brown   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
91*c4762a1bSJed Brown   ierr = ISRestoreIndices(isrows,&rows);CHKERRQ(ierr);
92*c4762a1bSJed Brown   ierr = ISRestoreIndices(iscols,&cols);CHKERRQ(ierr);
93*c4762a1bSJed Brown   ierr = ISDestroy(&isrows);CHKERRQ(ierr);
94*c4762a1bSJed Brown   ierr = ISDestroy(&iscols);CHKERRQ(ierr);
95*c4762a1bSJed Brown   ierr = PetscFree(v);CHKERRQ(ierr);
96*c4762a1bSJed Brown   if (mats_view) {
97*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD, "B: nrows %d, m %d; ncols %d, p %d\n",nrows,m,ncols,p);CHKERRQ(ierr);
98*c4762a1bSJed Brown     ierr = MatView(B,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
99*c4762a1bSJed Brown   }
100*c4762a1bSJed Brown 
101*c4762a1bSJed Brown   /* Create rhs vector b and solution x (same size as b) */
102*c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&b);CHKERRQ(ierr);
103*c4762a1bSJed Brown   ierr = VecSetSizes(b,m,PETSC_DECIDE);CHKERRQ(ierr);
104*c4762a1bSJed Brown   ierr = VecSetFromOptions(b);CHKERRQ(ierr);
105*c4762a1bSJed Brown   ierr = VecGetArray(b,&barray);CHKERRQ(ierr);
106*c4762a1bSJed Brown   for (j=0; j<m; j++) {
107*c4762a1bSJed Brown     ierr      = PetscRandomGetValue(rand,&rval);CHKERRQ(ierr);
108*c4762a1bSJed Brown     barray[j] = rval;
109*c4762a1bSJed Brown   }
110*c4762a1bSJed Brown   ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr);
111*c4762a1bSJed Brown   ierr = VecAssemblyBegin(b);CHKERRQ(ierr);
112*c4762a1bSJed Brown   ierr = VecAssemblyEnd(b);CHKERRQ(ierr);
113*c4762a1bSJed Brown   if (mats_view) {
114*c4762a1bSJed Brown     ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d] b: m %d\n",rank,m);CHKERRQ(ierr);
115*c4762a1bSJed Brown     ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr);
116*c4762a1bSJed Brown     ierr = VecView(b,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
117*c4762a1bSJed Brown   }
118*c4762a1bSJed Brown   ierr = VecDuplicate(b,&x);CHKERRQ(ierr);
119*c4762a1bSJed Brown 
120*c4762a1bSJed Brown   /* Create matrix X - same size as B */
121*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD," Create solution matrix X\n");CHKERRQ(ierr);
122*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&X);CHKERRQ(ierr);
123*c4762a1bSJed Brown   ierr = MatSetSizes(X,m,p,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
124*c4762a1bSJed Brown   ierr = MatSetType(X,MATELEMENTAL);CHKERRQ(ierr);
125*c4762a1bSJed Brown   ierr = MatSetFromOptions(X);CHKERRQ(ierr);
126*c4762a1bSJed Brown   ierr = MatSetUp(X);CHKERRQ(ierr);
127*c4762a1bSJed Brown   ierr = MatAssemblyBegin(X,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
128*c4762a1bSJed Brown   ierr = MatAssemblyEnd(X,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
129*c4762a1bSJed Brown 
130*c4762a1bSJed Brown   /* Cholesky factorization */
131*c4762a1bSJed Brown   /*------------------------*/
132*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD," Create Elemental matrix Aher\n");CHKERRQ(ierr);
133*c4762a1bSJed Brown   ierr = MatHermitianTranspose(A,MAT_INITIAL_MATRIX,&Aher);CHKERRQ(ierr);
134*c4762a1bSJed Brown   ierr = MatAXPY(Aher,1.0,A,SAME_NONZERO_PATTERN);CHKERRQ(ierr); /* Aher = A + A^T */
135*c4762a1bSJed Brown   if(!rank) { /* add 100.0 to diagonals of Aher to make it spd */
136*c4762a1bSJed Brown 
137*c4762a1bSJed Brown     /* TODO: Replace this with a call to El::ShiftDiagonal( A, 100. ),
138*c4762a1bSJed Brown              or at least pre-allocate the right amount of space */
139*c4762a1bSJed Brown     PetscInt M,N;
140*c4762a1bSJed Brown     ierr = MatGetSize(Aher,&M,&N);CHKERRQ(ierr);
141*c4762a1bSJed Brown     for (i=0; i<M; i++) {
142*c4762a1bSJed Brown       rval = 100.0;
143*c4762a1bSJed Brown       ierr = MatSetValues(Aher,1,&i,1,&i,&rval,ADD_VALUES);CHKERRQ(ierr);
144*c4762a1bSJed Brown     }
145*c4762a1bSJed Brown   }
146*c4762a1bSJed Brown   ierr = MatAssemblyBegin(Aher,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
147*c4762a1bSJed Brown   ierr = MatAssemblyEnd(Aher,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
148*c4762a1bSJed Brown   if (mats_view) {
149*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD, "Aher:\n");CHKERRQ(ierr);
150*c4762a1bSJed Brown     ierr = MatView(Aher,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
151*c4762a1bSJed Brown   }
152*c4762a1bSJed Brown 
153*c4762a1bSJed Brown   /* Cholesky factorization */
154*c4762a1bSJed Brown   /*------------------------*/
155*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD," Test Cholesky Solver \n");CHKERRQ(ierr);
156*c4762a1bSJed Brown   /* In-place Cholesky */
157*c4762a1bSJed Brown   /* Create matrix factor G, then copy Aher to G */
158*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&G);CHKERRQ(ierr);
159*c4762a1bSJed Brown   ierr = MatSetSizes(G,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
160*c4762a1bSJed Brown   ierr = MatSetType(G,MATELEMENTAL);CHKERRQ(ierr);
161*c4762a1bSJed Brown   ierr = MatSetFromOptions(G);CHKERRQ(ierr);
162*c4762a1bSJed Brown   ierr = MatSetUp(G);CHKERRQ(ierr);
163*c4762a1bSJed Brown   ierr = MatAssemblyBegin(G,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
164*c4762a1bSJed Brown   ierr = MatAssemblyEnd(G,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
165*c4762a1bSJed Brown   ierr = MatCopy(Aher,G,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
166*c4762a1bSJed Brown 
167*c4762a1bSJed Brown   /* Only G = U^T * U is implemented for now */
168*c4762a1bSJed Brown   ierr = MatCholeskyFactor(G,0,0);CHKERRQ(ierr);
169*c4762a1bSJed Brown   if (mats_view) {
170*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD, "Cholesky Factor G:\n");CHKERRQ(ierr);
171*c4762a1bSJed Brown     ierr = MatView(G,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
172*c4762a1bSJed Brown   }
173*c4762a1bSJed Brown 
174*c4762a1bSJed Brown   /* Solve U^T * U x = b and U^T * U X = B */
175*c4762a1bSJed Brown   ierr = MatSolve(G,b,x);CHKERRQ(ierr);
176*c4762a1bSJed Brown   ierr = MatMatSolve(G,B,X);CHKERRQ(ierr);
177*c4762a1bSJed Brown   ierr = MatDestroy(&G);CHKERRQ(ierr);
178*c4762a1bSJed Brown 
179*c4762a1bSJed Brown   /* Out-place Cholesky */
180*c4762a1bSJed Brown   ierr = MatGetFactor(Aher,MATSOLVERELEMENTAL,MAT_FACTOR_CHOLESKY,&G);CHKERRQ(ierr);
181*c4762a1bSJed Brown   ierr = MatCholeskyFactorSymbolic(G,Aher,0,&finfo);CHKERRQ(ierr);
182*c4762a1bSJed Brown   ierr = MatCholeskyFactorNumeric(G,Aher,&finfo);CHKERRQ(ierr);
183*c4762a1bSJed Brown   if (mats_view) {
184*c4762a1bSJed Brown     ierr = MatView(G,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
185*c4762a1bSJed Brown   }
186*c4762a1bSJed Brown   ierr = MatSolve(G,b,x);CHKERRQ(ierr);
187*c4762a1bSJed Brown   ierr = MatMatSolve(G,B,X);CHKERRQ(ierr);
188*c4762a1bSJed Brown   ierr = MatDestroy(&G);CHKERRQ(ierr);
189*c4762a1bSJed Brown 
190*c4762a1bSJed Brown   /* Check norm(Aher*x - b) */
191*c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&c);CHKERRQ(ierr);
192*c4762a1bSJed Brown   ierr = VecSetSizes(c,m,PETSC_DECIDE);CHKERRQ(ierr);
193*c4762a1bSJed Brown   ierr = VecSetFromOptions(c);CHKERRQ(ierr);
194*c4762a1bSJed Brown   ierr = MatMult(Aher,x,c);CHKERRQ(ierr);
195*c4762a1bSJed Brown   ierr = VecAXPY(c,-1.0,b);CHKERRQ(ierr);
196*c4762a1bSJed Brown   ierr = VecNorm(c,NORM_1,&norm);CHKERRQ(ierr);
197*c4762a1bSJed Brown   if (norm > tol) {
198*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Warning: |Aher*x - b| for Cholesky %g\n",(double)norm);CHKERRQ(ierr);
199*c4762a1bSJed Brown   }
200*c4762a1bSJed Brown 
201*c4762a1bSJed Brown   /* Check norm(Aher*X - B) */
202*c4762a1bSJed Brown   ierr = MatMatMult(Aher,X,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C);CHKERRQ(ierr);
203*c4762a1bSJed Brown   ierr = MatAXPY(C,-1.0,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
204*c4762a1bSJed Brown   ierr = MatNorm(C,NORM_1,&norm);CHKERRQ(ierr);
205*c4762a1bSJed Brown   if (norm > tol) {
206*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Warning: |Aher*X - B| for Cholesky %g\n",(double)norm);CHKERRQ(ierr);
207*c4762a1bSJed Brown   }
208*c4762a1bSJed Brown 
209*c4762a1bSJed Brown   /* LU factorization */
210*c4762a1bSJed Brown   /*------------------*/
211*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD," Test LU Solver \n");CHKERRQ(ierr);
212*c4762a1bSJed Brown   /* In-place LU */
213*c4762a1bSJed Brown   /* Create matrix factor F, then copy A to F */
214*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&F);CHKERRQ(ierr);
215*c4762a1bSJed Brown   ierr = MatSetSizes(F,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
216*c4762a1bSJed Brown   ierr = MatSetType(F,MATELEMENTAL);CHKERRQ(ierr);
217*c4762a1bSJed Brown   ierr = MatSetFromOptions(F);CHKERRQ(ierr);
218*c4762a1bSJed Brown   ierr = MatSetUp(F);CHKERRQ(ierr);
219*c4762a1bSJed Brown   ierr = MatAssemblyBegin(F,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
220*c4762a1bSJed Brown   ierr = MatAssemblyEnd(F,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
221*c4762a1bSJed Brown   ierr = MatCopy(A,F,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
222*c4762a1bSJed Brown   /* Create vector d to test MatSolveAdd() */
223*c4762a1bSJed Brown   ierr = VecDuplicate(x,&d);CHKERRQ(ierr);
224*c4762a1bSJed Brown   ierr = VecCopy(x,d);CHKERRQ(ierr);
225*c4762a1bSJed Brown 
226*c4762a1bSJed Brown   /* PF=LU or F=LU factorization - perms is ignored by Elemental;
227*c4762a1bSJed Brown      set finfo.dtcol !0 or 0 to enable/disable partial pivoting */
228*c4762a1bSJed Brown   finfo.dtcol = 0.1;
229*c4762a1bSJed Brown   ierr        = MatLUFactor(F,0,0,&finfo);CHKERRQ(ierr);
230*c4762a1bSJed Brown 
231*c4762a1bSJed Brown   /* Solve LUX = PB or LUX = B */
232*c4762a1bSJed Brown   ierr = MatSolveAdd(F,b,d,x);CHKERRQ(ierr);
233*c4762a1bSJed Brown   ierr = MatMatSolve(F,B,X);CHKERRQ(ierr);
234*c4762a1bSJed Brown   ierr = MatDestroy(&F);CHKERRQ(ierr);
235*c4762a1bSJed Brown 
236*c4762a1bSJed Brown   /* Check norm(A*X - B) */
237*c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&e);CHKERRQ(ierr);
238*c4762a1bSJed Brown   ierr = VecSetSizes(e,m,PETSC_DECIDE);CHKERRQ(ierr);
239*c4762a1bSJed Brown   ierr = VecSetFromOptions(e);CHKERRQ(ierr);
240*c4762a1bSJed Brown   ierr = MatMult(A,x,c);CHKERRQ(ierr);
241*c4762a1bSJed Brown   ierr = MatMult(A,d,e);CHKERRQ(ierr);
242*c4762a1bSJed Brown   ierr = VecAXPY(c,-1.0,e);CHKERRQ(ierr);
243*c4762a1bSJed Brown   ierr = VecAXPY(c,-1.0,b);CHKERRQ(ierr);
244*c4762a1bSJed Brown   ierr = VecNorm(c,NORM_1,&norm);CHKERRQ(ierr);
245*c4762a1bSJed Brown   if (norm > tol) {
246*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Warning: |A*x - b| for LU %g\n",(double)norm);CHKERRQ(ierr);
247*c4762a1bSJed Brown   }
248*c4762a1bSJed Brown   ierr = MatMatMult(A,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C);CHKERRQ(ierr);
249*c4762a1bSJed Brown   ierr = MatAXPY(C,-1.0,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
250*c4762a1bSJed Brown   ierr = MatNorm(C,NORM_1,&norm);CHKERRQ(ierr);
251*c4762a1bSJed Brown   if (norm > tol) {
252*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Warning: |A*X - B| for LU %g\n",(double)norm);CHKERRQ(ierr);
253*c4762a1bSJed Brown   }
254*c4762a1bSJed Brown 
255*c4762a1bSJed Brown   /* Out-place LU */
256*c4762a1bSJed Brown   ierr = MatGetFactor(A,MATSOLVERELEMENTAL,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
257*c4762a1bSJed Brown   ierr = MatLUFactorSymbolic(F,A,0,0,&finfo);CHKERRQ(ierr);
258*c4762a1bSJed Brown   ierr = MatLUFactorNumeric(F,A,&finfo);CHKERRQ(ierr);
259*c4762a1bSJed Brown   if (mats_view) {
260*c4762a1bSJed Brown     ierr = MatView(F,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
261*c4762a1bSJed Brown   }
262*c4762a1bSJed Brown   ierr = MatSolve(F,b,x);CHKERRQ(ierr);
263*c4762a1bSJed Brown   ierr = MatMatSolve(F,B,X);CHKERRQ(ierr);
264*c4762a1bSJed Brown   ierr = MatDestroy(&F);CHKERRQ(ierr);
265*c4762a1bSJed Brown 
266*c4762a1bSJed Brown   /* Free space */
267*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
268*c4762a1bSJed Brown   ierr = MatDestroy(&Aher);CHKERRQ(ierr);
269*c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
270*c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
271*c4762a1bSJed Brown   ierr = MatDestroy(&X);CHKERRQ(ierr);
272*c4762a1bSJed Brown   ierr = VecDestroy(&b);CHKERRQ(ierr);
273*c4762a1bSJed Brown   ierr = VecDestroy(&c);CHKERRQ(ierr);
274*c4762a1bSJed Brown   ierr = VecDestroy(&d);CHKERRQ(ierr);
275*c4762a1bSJed Brown   ierr = VecDestroy(&e);CHKERRQ(ierr);
276*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
277*c4762a1bSJed Brown   ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr);
278*c4762a1bSJed Brown   ierr = PetscFinalize();
279*c4762a1bSJed Brown   return ierr;
280*c4762a1bSJed Brown }
281*c4762a1bSJed Brown 
282*c4762a1bSJed Brown /*TEST
283*c4762a1bSJed Brown 
284*c4762a1bSJed Brown    build:
285*c4762a1bSJed Brown       requires: elemental
286*c4762a1bSJed Brown 
287*c4762a1bSJed Brown    test:
288*c4762a1bSJed Brown       nsize: 2
289*c4762a1bSJed Brown       output_file: output/ex145.out
290*c4762a1bSJed Brown 
291*c4762a1bSJed Brown    test:
292*c4762a1bSJed Brown       suffix: 2
293*c4762a1bSJed Brown       nsize: 6
294*c4762a1bSJed Brown       output_file: output/ex145.out
295*c4762a1bSJed Brown 
296*c4762a1bSJed Brown TEST*/
297