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