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; 22c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 23c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 24c4762a1bSJed Brown 25c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rand);CHKERRQ(ierr); 26c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rand);CHKERRQ(ierr); 27c4762a1bSJed Brown 28c4762a1bSJed Brown /* Get local dimensions of matrices */ 29c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr); 30c4762a1bSJed Brown n = m; 31c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr); 32c4762a1bSJed Brown p = m/2; 33c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL);CHKERRQ(ierr); 34c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-mats_view",&mats_view);CHKERRQ(ierr); 35c4762a1bSJed Brown 36c4762a1bSJed Brown /* Create matrix A */ 37c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," Create Elemental matrix A\n");CHKERRQ(ierr); 38c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 39c4762a1bSJed Brown ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 40c4762a1bSJed Brown ierr = MatSetType(A,MATELEMENTAL);CHKERRQ(ierr); 41c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 42c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 43c4762a1bSJed Brown /* Set local matrix entries */ 44c4762a1bSJed Brown ierr = MatGetOwnershipIS(A,&isrows,&iscols);CHKERRQ(ierr); 45c4762a1bSJed Brown ierr = ISGetLocalSize(isrows,&nrows);CHKERRQ(ierr); 46c4762a1bSJed Brown ierr = ISGetIndices(isrows,&rows);CHKERRQ(ierr); 47c4762a1bSJed Brown ierr = ISGetLocalSize(iscols,&ncols);CHKERRQ(ierr); 48c4762a1bSJed Brown ierr = ISGetIndices(iscols,&cols);CHKERRQ(ierr); 49c4762a1bSJed Brown ierr = PetscMalloc1(nrows*ncols,&v);CHKERRQ(ierr); 50c4762a1bSJed Brown for (i=0; i<nrows; i++) { 51c4762a1bSJed Brown for (j=0; j<ncols; j++) { 52c4762a1bSJed Brown ierr = PetscRandomGetValue(rand,&rval);CHKERRQ(ierr); 53c4762a1bSJed Brown v[i*ncols+j] = rval; 54c4762a1bSJed Brown } 55c4762a1bSJed Brown } 56c4762a1bSJed Brown ierr = MatSetValues(A,nrows,rows,ncols,cols,v,INSERT_VALUES);CHKERRQ(ierr); 57c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 58c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 59c4762a1bSJed Brown ierr = ISRestoreIndices(isrows,&rows);CHKERRQ(ierr); 60c4762a1bSJed Brown ierr = ISRestoreIndices(iscols,&cols);CHKERRQ(ierr); 61c4762a1bSJed Brown ierr = ISDestroy(&isrows);CHKERRQ(ierr); 62c4762a1bSJed Brown ierr = ISDestroy(&iscols);CHKERRQ(ierr); 63c4762a1bSJed Brown ierr = PetscFree(v);CHKERRQ(ierr); 64c4762a1bSJed Brown if (mats_view) { 65c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "A: nrows %d, m %d; ncols %d, n %d\n",nrows,m,ncols,n);CHKERRQ(ierr); 66c4762a1bSJed Brown ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 67c4762a1bSJed Brown } 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* Create rhs matrix B */ 70c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," Create rhs matrix B\n");CHKERRQ(ierr); 71c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr); 72c4762a1bSJed Brown ierr = MatSetSizes(B,m,p,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 73c4762a1bSJed Brown ierr = MatSetType(B,MATELEMENTAL);CHKERRQ(ierr); 74c4762a1bSJed Brown ierr = MatSetFromOptions(B);CHKERRQ(ierr); 75c4762a1bSJed Brown ierr = MatSetUp(B);CHKERRQ(ierr); 76c4762a1bSJed Brown ierr = MatGetOwnershipIS(B,&isrows,&iscols);CHKERRQ(ierr); 77c4762a1bSJed Brown ierr = ISGetLocalSize(isrows,&nrows);CHKERRQ(ierr); 78c4762a1bSJed Brown ierr = ISGetIndices(isrows,&rows);CHKERRQ(ierr); 79c4762a1bSJed Brown ierr = ISGetLocalSize(iscols,&ncols);CHKERRQ(ierr); 80c4762a1bSJed Brown ierr = ISGetIndices(iscols,&cols);CHKERRQ(ierr); 81c4762a1bSJed Brown ierr = PetscMalloc1(nrows*ncols,&v);CHKERRQ(ierr); 82c4762a1bSJed Brown for (i=0; i<nrows; i++) { 83c4762a1bSJed Brown for (j=0; j<ncols; j++) { 84c4762a1bSJed Brown ierr = PetscRandomGetValue(rand,&rval);CHKERRQ(ierr); 85c4762a1bSJed Brown v[i*ncols+j] = rval; 86c4762a1bSJed Brown } 87c4762a1bSJed Brown } 88c4762a1bSJed Brown ierr = MatSetValues(B,nrows,rows,ncols,cols,v,INSERT_VALUES);CHKERRQ(ierr); 89c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 90c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 91c4762a1bSJed Brown ierr = ISRestoreIndices(isrows,&rows);CHKERRQ(ierr); 92c4762a1bSJed Brown ierr = ISRestoreIndices(iscols,&cols);CHKERRQ(ierr); 93c4762a1bSJed Brown ierr = ISDestroy(&isrows);CHKERRQ(ierr); 94c4762a1bSJed Brown ierr = ISDestroy(&iscols);CHKERRQ(ierr); 95c4762a1bSJed Brown ierr = PetscFree(v);CHKERRQ(ierr); 96c4762a1bSJed Brown if (mats_view) { 97c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "B: nrows %d, m %d; ncols %d, p %d\n",nrows,m,ncols,p);CHKERRQ(ierr); 98c4762a1bSJed Brown ierr = MatView(B,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 99c4762a1bSJed Brown } 100c4762a1bSJed Brown 101c4762a1bSJed Brown /* Create rhs vector b and solution x (same size as b) */ 102c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&b);CHKERRQ(ierr); 103c4762a1bSJed Brown ierr = VecSetSizes(b,m,PETSC_DECIDE);CHKERRQ(ierr); 104c4762a1bSJed Brown ierr = VecSetFromOptions(b);CHKERRQ(ierr); 105c4762a1bSJed Brown ierr = VecGetArray(b,&barray);CHKERRQ(ierr); 106c4762a1bSJed Brown for (j=0; j<m; j++) { 107c4762a1bSJed Brown ierr = PetscRandomGetValue(rand,&rval);CHKERRQ(ierr); 108c4762a1bSJed Brown barray[j] = rval; 109c4762a1bSJed Brown } 110c4762a1bSJed Brown ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr); 111c4762a1bSJed Brown ierr = VecAssemblyBegin(b);CHKERRQ(ierr); 112c4762a1bSJed Brown ierr = VecAssemblyEnd(b);CHKERRQ(ierr); 113c4762a1bSJed Brown if (mats_view) { 114c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d] b: m %d\n",rank,m);CHKERRQ(ierr); 115c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr); 116c4762a1bSJed Brown ierr = VecView(b,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 117c4762a1bSJed Brown } 118c4762a1bSJed Brown ierr = VecDuplicate(b,&x);CHKERRQ(ierr); 119c4762a1bSJed Brown 120c4762a1bSJed Brown /* Create matrix X - same size as B */ 121c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," Create solution matrix X\n");CHKERRQ(ierr); 122c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&X);CHKERRQ(ierr); 123c4762a1bSJed Brown ierr = MatSetSizes(X,m,p,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 124c4762a1bSJed Brown ierr = MatSetType(X,MATELEMENTAL);CHKERRQ(ierr); 125c4762a1bSJed Brown ierr = MatSetFromOptions(X);CHKERRQ(ierr); 126c4762a1bSJed Brown ierr = MatSetUp(X);CHKERRQ(ierr); 127c4762a1bSJed Brown ierr = MatAssemblyBegin(X,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 128c4762a1bSJed Brown ierr = MatAssemblyEnd(X,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 129c4762a1bSJed Brown 130c4762a1bSJed Brown /* Cholesky factorization */ 131c4762a1bSJed Brown /*------------------------*/ 132c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," Create Elemental matrix Aher\n");CHKERRQ(ierr); 133c4762a1bSJed Brown ierr = MatHermitianTranspose(A,MAT_INITIAL_MATRIX,&Aher);CHKERRQ(ierr); 134c4762a1bSJed Brown ierr = MatAXPY(Aher,1.0,A,SAME_NONZERO_PATTERN);CHKERRQ(ierr); /* Aher = A + A^T */ 135c4762a1bSJed Brown if(!rank) { /* 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; 140c4762a1bSJed Brown ierr = MatGetSize(Aher,&M,&N);CHKERRQ(ierr); 141c4762a1bSJed Brown for (i=0; i<M; i++) { 142c4762a1bSJed Brown rval = 100.0; 143c4762a1bSJed Brown ierr = MatSetValues(Aher,1,&i,1,&i,&rval,ADD_VALUES);CHKERRQ(ierr); 144c4762a1bSJed Brown } 145c4762a1bSJed Brown } 146c4762a1bSJed Brown ierr = MatAssemblyBegin(Aher,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 147c4762a1bSJed Brown ierr = MatAssemblyEnd(Aher,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 148c4762a1bSJed Brown if (mats_view) { 149c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Aher:\n");CHKERRQ(ierr); 150c4762a1bSJed Brown ierr = MatView(Aher,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 151c4762a1bSJed Brown } 152c4762a1bSJed Brown 153c4762a1bSJed Brown /* Cholesky factorization */ 154c4762a1bSJed Brown /*------------------------*/ 155c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," Test Cholesky Solver \n");CHKERRQ(ierr); 156c4762a1bSJed Brown /* In-place Cholesky */ 157c4762a1bSJed Brown /* Create matrix factor G, then copy Aher to G */ 158c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&G);CHKERRQ(ierr); 159c4762a1bSJed Brown ierr = MatSetSizes(G,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 160c4762a1bSJed Brown ierr = MatSetType(G,MATELEMENTAL);CHKERRQ(ierr); 161c4762a1bSJed Brown ierr = MatSetFromOptions(G);CHKERRQ(ierr); 162c4762a1bSJed Brown ierr = MatSetUp(G);CHKERRQ(ierr); 163c4762a1bSJed Brown ierr = MatAssemblyBegin(G,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 164c4762a1bSJed Brown ierr = MatAssemblyEnd(G,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 165c4762a1bSJed Brown ierr = MatCopy(Aher,G,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 166c4762a1bSJed Brown 167c4762a1bSJed Brown /* Only G = U^T * U is implemented for now */ 168c4762a1bSJed Brown ierr = MatCholeskyFactor(G,0,0);CHKERRQ(ierr); 169c4762a1bSJed Brown if (mats_view) { 170c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Cholesky Factor G:\n");CHKERRQ(ierr); 171c4762a1bSJed Brown ierr = MatView(G,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 172c4762a1bSJed Brown } 173c4762a1bSJed Brown 174c4762a1bSJed Brown /* Solve U^T * U x = b and U^T * U X = B */ 175c4762a1bSJed Brown ierr = MatSolve(G,b,x);CHKERRQ(ierr); 176c4762a1bSJed Brown ierr = MatMatSolve(G,B,X);CHKERRQ(ierr); 177c4762a1bSJed Brown ierr = MatDestroy(&G);CHKERRQ(ierr); 178c4762a1bSJed Brown 179c4762a1bSJed Brown /* Out-place Cholesky */ 180c4762a1bSJed Brown ierr = MatGetFactor(Aher,MATSOLVERELEMENTAL,MAT_FACTOR_CHOLESKY,&G);CHKERRQ(ierr); 181c4762a1bSJed Brown ierr = MatCholeskyFactorSymbolic(G,Aher,0,&finfo);CHKERRQ(ierr); 182c4762a1bSJed Brown ierr = MatCholeskyFactorNumeric(G,Aher,&finfo);CHKERRQ(ierr); 183c4762a1bSJed Brown if (mats_view) { 184c4762a1bSJed Brown ierr = MatView(G,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 185c4762a1bSJed Brown } 186c4762a1bSJed Brown ierr = MatSolve(G,b,x);CHKERRQ(ierr); 187c4762a1bSJed Brown ierr = MatMatSolve(G,B,X);CHKERRQ(ierr); 188c4762a1bSJed Brown ierr = MatDestroy(&G);CHKERRQ(ierr); 189c4762a1bSJed Brown 190c4762a1bSJed Brown /* Check norm(Aher*x - b) */ 191c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&c);CHKERRQ(ierr); 192c4762a1bSJed Brown ierr = VecSetSizes(c,m,PETSC_DECIDE);CHKERRQ(ierr); 193c4762a1bSJed Brown ierr = VecSetFromOptions(c);CHKERRQ(ierr); 194c4762a1bSJed Brown ierr = MatMult(Aher,x,c);CHKERRQ(ierr); 195c4762a1bSJed Brown ierr = VecAXPY(c,-1.0,b);CHKERRQ(ierr); 196c4762a1bSJed Brown ierr = VecNorm(c,NORM_1,&norm);CHKERRQ(ierr); 197c4762a1bSJed Brown if (norm > tol) { 198c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Warning: |Aher*x - b| for Cholesky %g\n",(double)norm);CHKERRQ(ierr); 199c4762a1bSJed Brown } 200c4762a1bSJed Brown 201c4762a1bSJed Brown /* Check norm(Aher*X - B) */ 202c4762a1bSJed Brown ierr = MatMatMult(Aher,X,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C);CHKERRQ(ierr); 203c4762a1bSJed Brown ierr = MatAXPY(C,-1.0,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 204c4762a1bSJed Brown ierr = MatNorm(C,NORM_1,&norm);CHKERRQ(ierr); 205c4762a1bSJed Brown if (norm > tol) { 206c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Warning: |Aher*X - B| for Cholesky %g\n",(double)norm);CHKERRQ(ierr); 207c4762a1bSJed Brown } 208c4762a1bSJed Brown 209c4762a1bSJed Brown /* LU factorization */ 210c4762a1bSJed Brown /*------------------*/ 211c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," Test LU Solver \n");CHKERRQ(ierr); 212c4762a1bSJed Brown /* In-place LU */ 213c4762a1bSJed Brown /* Create matrix factor F, then copy A to F */ 214c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&F);CHKERRQ(ierr); 215c4762a1bSJed Brown ierr = MatSetSizes(F,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 216c4762a1bSJed Brown ierr = MatSetType(F,MATELEMENTAL);CHKERRQ(ierr); 217c4762a1bSJed Brown ierr = MatSetFromOptions(F);CHKERRQ(ierr); 218c4762a1bSJed Brown ierr = MatSetUp(F);CHKERRQ(ierr); 219c4762a1bSJed Brown ierr = MatAssemblyBegin(F,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 220c4762a1bSJed Brown ierr = MatAssemblyEnd(F,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 221c4762a1bSJed Brown ierr = MatCopy(A,F,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 222c4762a1bSJed Brown /* Create vector d to test MatSolveAdd() */ 223c4762a1bSJed Brown ierr = VecDuplicate(x,&d);CHKERRQ(ierr); 224c4762a1bSJed Brown ierr = VecCopy(x,d);CHKERRQ(ierr); 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; 229c4762a1bSJed Brown ierr = MatLUFactor(F,0,0,&finfo);CHKERRQ(ierr); 230c4762a1bSJed Brown 231c4762a1bSJed Brown /* Solve LUX = PB or LUX = B */ 232c4762a1bSJed Brown ierr = MatSolveAdd(F,b,d,x);CHKERRQ(ierr); 233c4762a1bSJed Brown ierr = MatMatSolve(F,B,X);CHKERRQ(ierr); 234c4762a1bSJed Brown ierr = MatDestroy(&F);CHKERRQ(ierr); 235c4762a1bSJed Brown 236c4762a1bSJed Brown /* Check norm(A*X - B) */ 237c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&e);CHKERRQ(ierr); 238c4762a1bSJed Brown ierr = VecSetSizes(e,m,PETSC_DECIDE);CHKERRQ(ierr); 239c4762a1bSJed Brown ierr = VecSetFromOptions(e);CHKERRQ(ierr); 240c4762a1bSJed Brown ierr = MatMult(A,x,c);CHKERRQ(ierr); 241c4762a1bSJed Brown ierr = MatMult(A,d,e);CHKERRQ(ierr); 242c4762a1bSJed Brown ierr = VecAXPY(c,-1.0,e);CHKERRQ(ierr); 243c4762a1bSJed Brown ierr = VecAXPY(c,-1.0,b);CHKERRQ(ierr); 244c4762a1bSJed Brown ierr = VecNorm(c,NORM_1,&norm);CHKERRQ(ierr); 245c4762a1bSJed Brown if (norm > tol) { 246c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Warning: |A*x - b| for LU %g\n",(double)norm);CHKERRQ(ierr); 247c4762a1bSJed Brown } 248*c20d7725SJed Brown /* Reuse product C; replace Aher with A */ 249*c20d7725SJed Brown ierr = MatProductReplaceMats(A,NULL,NULL,C);CHKERRQ(ierr); 250c4762a1bSJed Brown ierr = MatMatMult(A,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C);CHKERRQ(ierr); 251c4762a1bSJed Brown ierr = MatAXPY(C,-1.0,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 252c4762a1bSJed Brown ierr = MatNorm(C,NORM_1,&norm);CHKERRQ(ierr); 253c4762a1bSJed Brown if (norm > tol) { 254c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Warning: |A*X - B| for LU %g\n",(double)norm);CHKERRQ(ierr); 255c4762a1bSJed Brown } 256c4762a1bSJed Brown 257c4762a1bSJed Brown /* Out-place LU */ 258c4762a1bSJed Brown ierr = MatGetFactor(A,MATSOLVERELEMENTAL,MAT_FACTOR_LU,&F);CHKERRQ(ierr); 259c4762a1bSJed Brown ierr = MatLUFactorSymbolic(F,A,0,0,&finfo);CHKERRQ(ierr); 260c4762a1bSJed Brown ierr = MatLUFactorNumeric(F,A,&finfo);CHKERRQ(ierr); 261c4762a1bSJed Brown if (mats_view) { 262c4762a1bSJed Brown ierr = MatView(F,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 263c4762a1bSJed Brown } 264c4762a1bSJed Brown ierr = MatSolve(F,b,x);CHKERRQ(ierr); 265c4762a1bSJed Brown ierr = MatMatSolve(F,B,X);CHKERRQ(ierr); 266c4762a1bSJed Brown ierr = MatDestroy(&F);CHKERRQ(ierr); 267c4762a1bSJed Brown 268c4762a1bSJed Brown /* Free space */ 269c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 270c4762a1bSJed Brown ierr = MatDestroy(&Aher);CHKERRQ(ierr); 271c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 272c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 273c4762a1bSJed Brown ierr = MatDestroy(&X);CHKERRQ(ierr); 274c4762a1bSJed Brown ierr = VecDestroy(&b);CHKERRQ(ierr); 275c4762a1bSJed Brown ierr = VecDestroy(&c);CHKERRQ(ierr); 276c4762a1bSJed Brown ierr = VecDestroy(&d);CHKERRQ(ierr); 277c4762a1bSJed Brown ierr = VecDestroy(&e);CHKERRQ(ierr); 278c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 279c4762a1bSJed Brown ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr); 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