xref: /petsc/src/mat/tests/ex145.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   PetscInt       m = 5,n,p,i,j,nrows,ncols;
11c4762a1bSJed Brown   PetscScalar    *v,*barray,rval;
12c4762a1bSJed Brown   PetscReal      norm,tol=1.e-11;
13c4762a1bSJed Brown   PetscMPIInt    size,rank;
14c4762a1bSJed Brown   PetscRandom    rand;
15c4762a1bSJed Brown   const PetscInt *rows,*cols;
16c4762a1bSJed Brown   IS             isrows,iscols;
17c4762a1bSJed Brown   PetscBool      mats_view=PETSC_FALSE;
18c4762a1bSJed Brown   MatFactorInfo  finfo;
19c4762a1bSJed Brown 
20*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*) 0,help));
215f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
225f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
23c4762a1bSJed Brown 
245f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rand));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rand));
26c4762a1bSJed Brown 
27c4762a1bSJed Brown   /* Get local dimensions of matrices */
285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
29c4762a1bSJed Brown   n    = m;
305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
31c4762a1bSJed Brown   p    = m/2;
325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-mats_view",&mats_view));
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   /* Create matrix A */
365f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," Create Elemental matrix A\n"));
375f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATELEMENTAL));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
42c4762a1bSJed Brown   /* Set local matrix entries */
435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipIS(A,&isrows,&iscols));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetLocalSize(isrows,&nrows));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(isrows,&rows));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetLocalSize(iscols,&ncols));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(iscols,&cols));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrows*ncols,&v));
49c4762a1bSJed Brown   for (i=0; i<nrows; i++) {
50c4762a1bSJed Brown     for (j=0; j<ncols; j++) {
515f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRandomGetValue(rand,&rval));
52c4762a1bSJed Brown       v[i*ncols+j] = rval;
53c4762a1bSJed Brown     }
54c4762a1bSJed Brown   }
555f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,nrows,rows,ncols,cols,v,INSERT_VALUES));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(isrows,&rows));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(iscols,&cols));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isrows));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&iscols));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(v));
63c4762a1bSJed Brown   if (mats_view) {
645f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "A: nrows %" PetscInt_FMT ", m %" PetscInt_FMT "; ncols %" PetscInt_FMT ", n %" PetscInt_FMT "\n",nrows,m,ncols,n));
655f80ce2aSJacob Faibussowitsch     CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
66c4762a1bSJed Brown   }
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   /* Create rhs matrix B */
695f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," Create rhs matrix B\n"));
705f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B));
715f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(B,m,p,PETSC_DECIDE,PETSC_DECIDE));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(B,MATELEMENTAL));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(B));
745f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(B));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipIS(B,&isrows,&iscols));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetLocalSize(isrows,&nrows));
775f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(isrows,&rows));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetLocalSize(iscols,&ncols));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(iscols,&cols));
805f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrows*ncols,&v));
81c4762a1bSJed Brown   for (i=0; i<nrows; i++) {
82c4762a1bSJed Brown     for (j=0; j<ncols; j++) {
835f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRandomGetValue(rand,&rval));
84c4762a1bSJed Brown       v[i*ncols+j] = rval;
85c4762a1bSJed Brown     }
86c4762a1bSJed Brown   }
875f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(B,nrows,rows,ncols,cols,v,INSERT_VALUES));
885f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(isrows,&rows));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(iscols,&cols));
925f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isrows));
935f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&iscols));
945f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(v));
95c4762a1bSJed Brown   if (mats_view) {
965f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "B: nrows %" PetscInt_FMT ", m %" PetscInt_FMT "; ncols %" PetscInt_FMT ", p %" PetscInt_FMT "\n",nrows,m,ncols,p));
975f80ce2aSJacob Faibussowitsch     CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD));
98c4762a1bSJed Brown   }
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   /* Create rhs vector b and solution x (same size as b) */
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&b));
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(b,m,PETSC_DECIDE));
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(b));
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(b,&barray));
105c4762a1bSJed Brown   for (j=0; j<m; j++) {
1065f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rand,&rval));
107c4762a1bSJed Brown     barray[j] = rval;
108c4762a1bSJed Brown   }
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(b,&barray));
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(b));
1115f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(b));
112c4762a1bSJed Brown   if (mats_view) {
1135f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d] b: m %" PetscInt_FMT "\n",rank,m));
1145f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
1155f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(b,PETSC_VIEWER_STDOUT_WORLD));
116c4762a1bSJed Brown   }
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(b,&x));
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   /* Create matrix X - same size as B */
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," Create solution matrix X\n"));
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&X));
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(X,m,p,PETSC_DECIDE,PETSC_DECIDE));
1235f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(X,MATELEMENTAL));
1245f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(X));
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(X));
1265f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(X,MAT_FINAL_ASSEMBLY));
1275f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(X,MAT_FINAL_ASSEMBLY));
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   /* Cholesky factorization */
130c4762a1bSJed Brown   /*------------------------*/
1315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," Create Elemental matrix Aher\n"));
1325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatHermitianTranspose(A,MAT_INITIAL_MATRIX,&Aher));
1335f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAXPY(Aher,1.0,A,SAME_NONZERO_PATTERN)); /* Aher = A + A^T */
134dd400576SPatrick Sanan   if (rank == 0) { /* add 100.0 to diagonals of Aher to make it spd */
135c4762a1bSJed Brown 
136c4762a1bSJed Brown     /* TODO: Replace this with a call to El::ShiftDiagonal( A, 100.),
137c4762a1bSJed Brown              or at least pre-allocate the right amount of space */
138c4762a1bSJed Brown     PetscInt M,N;
1395f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetSize(Aher,&M,&N));
140c4762a1bSJed Brown     for (i=0; i<M; i++) {
141c4762a1bSJed Brown       rval = 100.0;
1425f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(Aher,1,&i,1,&i,&rval,ADD_VALUES));
143c4762a1bSJed Brown     }
144c4762a1bSJed Brown   }
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(Aher,MAT_FINAL_ASSEMBLY));
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(Aher,MAT_FINAL_ASSEMBLY));
147c4762a1bSJed Brown   if (mats_view) {
1485f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Aher:\n"));
1495f80ce2aSJacob Faibussowitsch     CHKERRQ(MatView(Aher,PETSC_VIEWER_STDOUT_WORLD));
150c4762a1bSJed Brown   }
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   /* Cholesky factorization */
153c4762a1bSJed Brown   /*------------------------*/
1545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," Test Cholesky Solver \n"));
155c4762a1bSJed Brown   /* In-place Cholesky */
156c4762a1bSJed Brown   /* Create matrix factor G, then copy Aher to G */
1575f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&G));
1585f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(G,m,n,PETSC_DECIDE,PETSC_DECIDE));
1595f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(G,MATELEMENTAL));
1605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(G));
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(G));
1625f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(G,MAT_FINAL_ASSEMBLY));
1635f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(G,MAT_FINAL_ASSEMBLY));
1645f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCopy(Aher,G,SAME_NONZERO_PATTERN));
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   /* Only G = U^T * U is implemented for now */
1675f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCholeskyFactor(G,0,0));
168c4762a1bSJed Brown   if (mats_view) {
1695f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Cholesky Factor G:\n"));
1705f80ce2aSJacob Faibussowitsch     CHKERRQ(MatView(G,PETSC_VIEWER_STDOUT_WORLD));
171c4762a1bSJed Brown   }
172c4762a1bSJed Brown 
173c4762a1bSJed Brown   /* Solve U^T * U x = b and U^T * U X = B */
1745f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSolve(G,b,x));
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatSolve(G,B,X));
1765f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&G));
177c4762a1bSJed Brown 
178c4762a1bSJed Brown   /* Out-place Cholesky */
1795f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetFactor(Aher,MATSOLVERELEMENTAL,MAT_FACTOR_CHOLESKY,&G));
1805f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCholeskyFactorSymbolic(G,Aher,0,&finfo));
1815f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCholeskyFactorNumeric(G,Aher,&finfo));
182c4762a1bSJed Brown   if (mats_view) {
1835f80ce2aSJacob Faibussowitsch     CHKERRQ(MatView(G,PETSC_VIEWER_STDOUT_WORLD));
184c4762a1bSJed Brown   }
1855f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSolve(G,b,x));
1865f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatSolve(G,B,X));
1875f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&G));
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   /* Check norm(Aher*x - b) */
1905f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&c));
1915f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(c,m,PETSC_DECIDE));
1925f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(c));
1935f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(Aher,x,c));
1945f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(c,-1.0,b));
1955f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(c,NORM_1,&norm));
196c4762a1bSJed Brown   if (norm > tol) {
1975f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Warning: |Aher*x - b| for Cholesky %g\n",(double)norm));
198c4762a1bSJed Brown   }
199c4762a1bSJed Brown 
200c4762a1bSJed Brown   /* Check norm(Aher*X - B) */
2015f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(Aher,X,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C));
2025f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAXPY(C,-1.0,B,SAME_NONZERO_PATTERN));
2035f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNorm(C,NORM_1,&norm));
204c4762a1bSJed Brown   if (norm > tol) {
2055f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Warning: |Aher*X - B| for Cholesky %g\n",(double)norm));
206c4762a1bSJed Brown   }
207c4762a1bSJed Brown 
208c4762a1bSJed Brown   /* LU factorization */
209c4762a1bSJed Brown   /*------------------*/
2105f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," Test LU Solver \n"));
211c4762a1bSJed Brown   /* In-place LU */
212c4762a1bSJed Brown   /* Create matrix factor F, then copy A to F */
2135f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&F));
2145f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(F,m,n,PETSC_DECIDE,PETSC_DECIDE));
2155f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(F,MATELEMENTAL));
2165f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(F));
2175f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(F));
2185f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(F,MAT_FINAL_ASSEMBLY));
2195f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(F,MAT_FINAL_ASSEMBLY));
2205f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCopy(A,F,SAME_NONZERO_PATTERN));
221c4762a1bSJed Brown   /* Create vector d to test MatSolveAdd() */
2225f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&d));
2235f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(x,d));
224c4762a1bSJed Brown 
225c4762a1bSJed Brown   /* PF=LU or F=LU factorization - perms is ignored by Elemental;
226c4762a1bSJed Brown      set finfo.dtcol !0 or 0 to enable/disable partial pivoting */
227c4762a1bSJed Brown   finfo.dtcol = 0.1;
2285f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLUFactor(F,0,0,&finfo));
229c4762a1bSJed Brown 
230c4762a1bSJed Brown   /* Solve LUX = PB or LUX = B */
2315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSolveAdd(F,b,d,x));
2325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatSolve(F,B,X));
2335f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&F));
234c4762a1bSJed Brown 
235c4762a1bSJed Brown   /* Check norm(A*X - B) */
2365f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&e));
2375f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(e,m,PETSC_DECIDE));
2385f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(e));
2395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(A,x,c));
2405f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(A,d,e));
2415f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(c,-1.0,e));
2425f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(c,-1.0,b));
2435f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(c,NORM_1,&norm));
244c4762a1bSJed Brown   if (norm > tol) {
2455f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Warning: |A*x - b| for LU %g\n",(double)norm));
246c4762a1bSJed Brown   }
247c20d7725SJed Brown   /* Reuse product C; replace Aher with A */
2485f80ce2aSJacob Faibussowitsch   CHKERRQ(MatProductReplaceMats(A,NULL,NULL,C));
2495f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(A,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C));
2505f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAXPY(C,-1.0,B,SAME_NONZERO_PATTERN));
2515f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNorm(C,NORM_1,&norm));
252c4762a1bSJed Brown   if (norm > tol) {
2535f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Warning: |A*X - B| for LU %g\n",(double)norm));
254c4762a1bSJed Brown   }
255c4762a1bSJed Brown 
256c4762a1bSJed Brown   /* Out-place LU */
2575f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetFactor(A,MATSOLVERELEMENTAL,MAT_FACTOR_LU,&F));
2585f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLUFactorSymbolic(F,A,0,0,&finfo));
2595f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLUFactorNumeric(F,A,&finfo));
260c4762a1bSJed Brown   if (mats_view) {
2615f80ce2aSJacob Faibussowitsch     CHKERRQ(MatView(F,PETSC_VIEWER_STDOUT_WORLD));
262c4762a1bSJed Brown   }
2635f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSolve(F,b,x));
2645f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatSolve(F,B,X));
2655f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&F));
266c4762a1bSJed Brown 
267c4762a1bSJed Brown   /* Free space */
2685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
2695f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&Aher));
2705f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
2715f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
2725f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&X));
2735f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&b));
2745f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&c));
2755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&d));
2765f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&e));
2775f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
2785f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rand));
279*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
280*b122ec5aSJacob Faibussowitsch   return 0;
281c4762a1bSJed Brown }
282c4762a1bSJed Brown 
283c4762a1bSJed Brown /*TEST
284c4762a1bSJed Brown 
285c4762a1bSJed Brown    build:
286c4762a1bSJed Brown       requires: elemental
287c4762a1bSJed Brown 
288c4762a1bSJed Brown    test:
289c4762a1bSJed Brown       nsize: 2
290c4762a1bSJed Brown       output_file: output/ex145.out
291c4762a1bSJed Brown 
292c4762a1bSJed Brown    test:
293c4762a1bSJed Brown       suffix: 2
294c4762a1bSJed Brown       nsize: 6
295c4762a1bSJed Brown       output_file: output/ex145.out
296c4762a1bSJed Brown 
297c4762a1bSJed Brown TEST*/
298