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 209566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*) 0,help)); 219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 23c4762a1bSJed Brown 249566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 259566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 26c4762a1bSJed Brown 27c4762a1bSJed Brown /* Get local dimensions of matrices */ 289566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 29c4762a1bSJed Brown n = m; 309566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 31c4762a1bSJed Brown p = m/2; 329566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL)); 339566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-mats_view",&mats_view)); 34c4762a1bSJed Brown 35c4762a1bSJed Brown /* Create matrix A */ 369566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," Create Elemental matrix A\n")); 379566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 389566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE)); 399566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATELEMENTAL)); 409566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 419566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 42c4762a1bSJed Brown /* Set local matrix entries */ 439566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipIS(A,&isrows,&iscols)); 449566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrows,&nrows)); 459566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrows,&rows)); 469566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iscols,&ncols)); 479566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iscols,&cols)); 489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrows*ncols,&v)); 49c4762a1bSJed Brown for (i=0; i<nrows; i++) { 50c4762a1bSJed Brown for (j=0; j<ncols; j++) { 519566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&rval)); 52c4762a1bSJed Brown v[i*ncols+j] = rval; 53c4762a1bSJed Brown } 54c4762a1bSJed Brown } 559566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,nrows,rows,ncols,cols,v,INSERT_VALUES)); 569566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 579566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 589566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrows,&rows)); 599566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscols,&cols)); 609566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrows)); 619566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscols)); 629566063dSJacob Faibussowitsch PetscCall(PetscFree(v)); 63c4762a1bSJed Brown if (mats_view) { 649566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "A: nrows %" PetscInt_FMT ", m %" PetscInt_FMT "; ncols %" PetscInt_FMT ", n %" PetscInt_FMT "\n",nrows,m,ncols,n)); 659566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 66c4762a1bSJed Brown } 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* Create rhs matrix B */ 699566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," Create rhs matrix B\n")); 709566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&B)); 719566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B,m,p,PETSC_DECIDE,PETSC_DECIDE)); 729566063dSJacob Faibussowitsch PetscCall(MatSetType(B,MATELEMENTAL)); 739566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 749566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 759566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipIS(B,&isrows,&iscols)); 769566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrows,&nrows)); 779566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrows,&rows)); 789566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iscols,&ncols)); 799566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iscols,&cols)); 809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrows*ncols,&v)); 81c4762a1bSJed Brown for (i=0; i<nrows; i++) { 82c4762a1bSJed Brown for (j=0; j<ncols; j++) { 839566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&rval)); 84c4762a1bSJed Brown v[i*ncols+j] = rval; 85c4762a1bSJed Brown } 86c4762a1bSJed Brown } 879566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,nrows,rows,ncols,cols,v,INSERT_VALUES)); 889566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 899566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 909566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrows,&rows)); 919566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscols,&cols)); 929566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrows)); 939566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscols)); 949566063dSJacob Faibussowitsch PetscCall(PetscFree(v)); 95c4762a1bSJed Brown if (mats_view) { 969566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "B: nrows %" PetscInt_FMT ", m %" PetscInt_FMT "; ncols %" PetscInt_FMT ", p %" PetscInt_FMT "\n",nrows,m,ncols,p)); 979566063dSJacob 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) */ 1019566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&b)); 1029566063dSJacob Faibussowitsch PetscCall(VecSetSizes(b,m,PETSC_DECIDE)); 1039566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(b)); 1049566063dSJacob Faibussowitsch PetscCall(VecGetArray(b,&barray)); 105c4762a1bSJed Brown for (j=0; j<m; j++) { 1069566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&rval)); 107c4762a1bSJed Brown barray[j] = rval; 108c4762a1bSJed Brown } 1099566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(b,&barray)); 1109566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(b)); 1119566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(b)); 112c4762a1bSJed Brown if (mats_view) { 1139566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d] b: m %" PetscInt_FMT "\n",rank,m)); 1149566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 1159566063dSJacob Faibussowitsch PetscCall(VecView(b,PETSC_VIEWER_STDOUT_WORLD)); 116c4762a1bSJed Brown } 1179566063dSJacob Faibussowitsch PetscCall(VecDuplicate(b,&x)); 118c4762a1bSJed Brown 119c4762a1bSJed Brown /* Create matrix X - same size as B */ 1209566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," Create solution matrix X\n")); 1219566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&X)); 1229566063dSJacob Faibussowitsch PetscCall(MatSetSizes(X,m,p,PETSC_DECIDE,PETSC_DECIDE)); 1239566063dSJacob Faibussowitsch PetscCall(MatSetType(X,MATELEMENTAL)); 1249566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(X)); 1259566063dSJacob Faibussowitsch PetscCall(MatSetUp(X)); 1269566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(X,MAT_FINAL_ASSEMBLY)); 1279566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(X,MAT_FINAL_ASSEMBLY)); 128c4762a1bSJed Brown 129c4762a1bSJed Brown /* Cholesky factorization */ 130c4762a1bSJed Brown /*------------------------*/ 1319566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," Create Elemental matrix Aher\n")); 1329566063dSJacob Faibussowitsch PetscCall(MatHermitianTranspose(A,MAT_INITIAL_MATRIX,&Aher)); 1339566063dSJacob 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; 1399566063dSJacob Faibussowitsch PetscCall(MatGetSize(Aher,&M,&N)); 140c4762a1bSJed Brown for (i=0; i<M; i++) { 141c4762a1bSJed Brown rval = 100.0; 1429566063dSJacob Faibussowitsch PetscCall(MatSetValues(Aher,1,&i,1,&i,&rval,ADD_VALUES)); 143c4762a1bSJed Brown } 144c4762a1bSJed Brown } 1459566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Aher,MAT_FINAL_ASSEMBLY)); 1469566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Aher,MAT_FINAL_ASSEMBLY)); 147c4762a1bSJed Brown if (mats_view) { 1489566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Aher:\n")); 1499566063dSJacob Faibussowitsch PetscCall(MatView(Aher,PETSC_VIEWER_STDOUT_WORLD)); 150c4762a1bSJed Brown } 151c4762a1bSJed Brown 152c4762a1bSJed Brown /* Cholesky factorization */ 153c4762a1bSJed Brown /*------------------------*/ 1549566063dSJacob 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 */ 1579566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&G)); 1589566063dSJacob Faibussowitsch PetscCall(MatSetSizes(G,m,n,PETSC_DECIDE,PETSC_DECIDE)); 1599566063dSJacob Faibussowitsch PetscCall(MatSetType(G,MATELEMENTAL)); 1609566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(G)); 1619566063dSJacob Faibussowitsch PetscCall(MatSetUp(G)); 1629566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(G,MAT_FINAL_ASSEMBLY)); 1639566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(G,MAT_FINAL_ASSEMBLY)); 1649566063dSJacob Faibussowitsch PetscCall(MatCopy(Aher,G,SAME_NONZERO_PATTERN)); 165c4762a1bSJed Brown 166c4762a1bSJed Brown /* Only G = U^T * U is implemented for now */ 1679566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactor(G,0,0)); 168c4762a1bSJed Brown if (mats_view) { 1699566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Cholesky Factor G:\n")); 1709566063dSJacob 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 */ 1749566063dSJacob Faibussowitsch PetscCall(MatSolve(G,b,x)); 1759566063dSJacob Faibussowitsch PetscCall(MatMatSolve(G,B,X)); 1769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&G)); 177c4762a1bSJed Brown 178c4762a1bSJed Brown /* Out-place Cholesky */ 1799566063dSJacob Faibussowitsch PetscCall(MatGetFactor(Aher,MATSOLVERELEMENTAL,MAT_FACTOR_CHOLESKY,&G)); 1809566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(G,Aher,0,&finfo)); 1819566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(G,Aher,&finfo)); 182*1baa6e33SBarry Smith if (mats_view) PetscCall(MatView(G,PETSC_VIEWER_STDOUT_WORLD)); 1839566063dSJacob Faibussowitsch PetscCall(MatSolve(G,b,x)); 1849566063dSJacob Faibussowitsch PetscCall(MatMatSolve(G,B,X)); 1859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&G)); 186c4762a1bSJed Brown 187c4762a1bSJed Brown /* Check norm(Aher*x - b) */ 1889566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&c)); 1899566063dSJacob Faibussowitsch PetscCall(VecSetSizes(c,m,PETSC_DECIDE)); 1909566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(c)); 1919566063dSJacob Faibussowitsch PetscCall(MatMult(Aher,x,c)); 1929566063dSJacob Faibussowitsch PetscCall(VecAXPY(c,-1.0,b)); 1939566063dSJacob Faibussowitsch PetscCall(VecNorm(c,NORM_1,&norm)); 194c4762a1bSJed Brown if (norm > tol) { 1959566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: |Aher*x - b| for Cholesky %g\n",(double)norm)); 196c4762a1bSJed Brown } 197c4762a1bSJed Brown 198c4762a1bSJed Brown /* Check norm(Aher*X - B) */ 1999566063dSJacob Faibussowitsch PetscCall(MatMatMult(Aher,X,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C)); 2009566063dSJacob Faibussowitsch PetscCall(MatAXPY(C,-1.0,B,SAME_NONZERO_PATTERN)); 2019566063dSJacob Faibussowitsch PetscCall(MatNorm(C,NORM_1,&norm)); 202c4762a1bSJed Brown if (norm > tol) { 2039566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: |Aher*X - B| for Cholesky %g\n",(double)norm)); 204c4762a1bSJed Brown } 205c4762a1bSJed Brown 206c4762a1bSJed Brown /* LU factorization */ 207c4762a1bSJed Brown /*------------------*/ 2089566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," Test LU Solver \n")); 209c4762a1bSJed Brown /* In-place LU */ 210c4762a1bSJed Brown /* Create matrix factor F, then copy A to F */ 2119566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&F)); 2129566063dSJacob Faibussowitsch PetscCall(MatSetSizes(F,m,n,PETSC_DECIDE,PETSC_DECIDE)); 2139566063dSJacob Faibussowitsch PetscCall(MatSetType(F,MATELEMENTAL)); 2149566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(F)); 2159566063dSJacob Faibussowitsch PetscCall(MatSetUp(F)); 2169566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(F,MAT_FINAL_ASSEMBLY)); 2179566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(F,MAT_FINAL_ASSEMBLY)); 2189566063dSJacob Faibussowitsch PetscCall(MatCopy(A,F,SAME_NONZERO_PATTERN)); 219c4762a1bSJed Brown /* Create vector d to test MatSolveAdd() */ 2209566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&d)); 2219566063dSJacob Faibussowitsch PetscCall(VecCopy(x,d)); 222c4762a1bSJed Brown 223c4762a1bSJed Brown /* PF=LU or F=LU factorization - perms is ignored by Elemental; 224c4762a1bSJed Brown set finfo.dtcol !0 or 0 to enable/disable partial pivoting */ 225c4762a1bSJed Brown finfo.dtcol = 0.1; 2269566063dSJacob Faibussowitsch PetscCall(MatLUFactor(F,0,0,&finfo)); 227c4762a1bSJed Brown 228c4762a1bSJed Brown /* Solve LUX = PB or LUX = B */ 2299566063dSJacob Faibussowitsch PetscCall(MatSolveAdd(F,b,d,x)); 2309566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F,B,X)); 2319566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 232c4762a1bSJed Brown 233c4762a1bSJed Brown /* Check norm(A*X - B) */ 2349566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&e)); 2359566063dSJacob Faibussowitsch PetscCall(VecSetSizes(e,m,PETSC_DECIDE)); 2369566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(e)); 2379566063dSJacob Faibussowitsch PetscCall(MatMult(A,x,c)); 2389566063dSJacob Faibussowitsch PetscCall(MatMult(A,d,e)); 2399566063dSJacob Faibussowitsch PetscCall(VecAXPY(c,-1.0,e)); 2409566063dSJacob Faibussowitsch PetscCall(VecAXPY(c,-1.0,b)); 2419566063dSJacob Faibussowitsch PetscCall(VecNorm(c,NORM_1,&norm)); 242c4762a1bSJed Brown if (norm > tol) { 2439566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: |A*x - b| for LU %g\n",(double)norm)); 244c4762a1bSJed Brown } 245c20d7725SJed Brown /* Reuse product C; replace Aher with A */ 2469566063dSJacob Faibussowitsch PetscCall(MatProductReplaceMats(A,NULL,NULL,C)); 2479566063dSJacob Faibussowitsch PetscCall(MatMatMult(A,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C)); 2489566063dSJacob Faibussowitsch PetscCall(MatAXPY(C,-1.0,B,SAME_NONZERO_PATTERN)); 2499566063dSJacob Faibussowitsch PetscCall(MatNorm(C,NORM_1,&norm)); 250c4762a1bSJed Brown if (norm > tol) { 2519566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: |A*X - B| for LU %g\n",(double)norm)); 252c4762a1bSJed Brown } 253c4762a1bSJed Brown 254c4762a1bSJed Brown /* Out-place LU */ 2559566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A,MATSOLVERELEMENTAL,MAT_FACTOR_LU,&F)); 2569566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(F,A,0,0,&finfo)); 2579566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(F,A,&finfo)); 258*1baa6e33SBarry Smith if (mats_view) PetscCall(MatView(F,PETSC_VIEWER_STDOUT_WORLD)); 2599566063dSJacob Faibussowitsch PetscCall(MatSolve(F,b,x)); 2609566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F,B,X)); 2619566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 262c4762a1bSJed Brown 263c4762a1bSJed Brown /* Free space */ 2649566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2659566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Aher)); 2669566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 2679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 2689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&X)); 2699566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 2709566063dSJacob Faibussowitsch PetscCall(VecDestroy(&c)); 2719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&d)); 2729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&e)); 2739566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2749566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 2759566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 276b122ec5aSJacob Faibussowitsch return 0; 277c4762a1bSJed Brown } 278c4762a1bSJed Brown 279c4762a1bSJed Brown /*TEST 280c4762a1bSJed Brown 281c4762a1bSJed Brown build: 282c4762a1bSJed Brown requires: elemental 283c4762a1bSJed Brown 284c4762a1bSJed Brown test: 285c4762a1bSJed Brown nsize: 2 286c4762a1bSJed Brown output_file: output/ex145.out 287c4762a1bSJed Brown 288c4762a1bSJed Brown test: 289c4762a1bSJed Brown suffix: 2 290c4762a1bSJed Brown nsize: 6 291c4762a1bSJed Brown output_file: output/ex145.out 292c4762a1bSJed Brown 293c4762a1bSJed Brown TEST*/ 294