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