1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests cholesky, icc factorization and solve on sequential aij, baij and sbaij matrices. \n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char **args) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown Vec x,y,b; 9c4762a1bSJed Brown Mat A; /* linear system matrix */ 10c4762a1bSJed Brown Mat sA,sC; /* symmetric part of the matrices */ 11c4762a1bSJed Brown PetscInt n,mbs=16,bs=1,nz=3,prob=1,i,j,col[3],block, row,Ii,J,n1,lvl; 12c4762a1bSJed Brown PetscErrorCode ierr; 13c4762a1bSJed Brown PetscMPIInt size; 14c4762a1bSJed Brown PetscReal norm2; 15c4762a1bSJed Brown PetscScalar neg_one = -1.0,four=4.0,value[3]; 16c4762a1bSJed Brown IS perm,cperm; 17c4762a1bSJed Brown PetscRandom rdm; 18c4762a1bSJed Brown PetscBool reorder = PETSC_FALSE,displ = PETSC_FALSE; 19c4762a1bSJed Brown MatFactorInfo factinfo; 20c4762a1bSJed Brown PetscBool equal; 21c4762a1bSJed Brown PetscBool TestAIJ = PETSC_FALSE,TestBAIJ = PETSC_TRUE; 22c4762a1bSJed Brown PetscInt TestShift=0; 23c4762a1bSJed Brown 24c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 25ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 26*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 27c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);CHKERRQ(ierr); 28c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL);CHKERRQ(ierr); 29c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-reorder",&reorder,NULL);CHKERRQ(ierr); 30c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-testaij",&TestAIJ,NULL);CHKERRQ(ierr); 31c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-testShift",&TestShift,NULL);CHKERRQ(ierr); 32c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-displ",&displ,NULL);CHKERRQ(ierr); 33c4762a1bSJed Brown 34c4762a1bSJed Brown n = mbs*bs; 35c4762a1bSJed Brown if (TestAIJ) { /* A is in aij format */ 36c4762a1bSJed Brown ierr = MatCreateSeqAIJ(PETSC_COMM_WORLD,n,n,nz,NULL,&A);CHKERRQ(ierr); 37c4762a1bSJed Brown TestBAIJ = PETSC_FALSE; 38c4762a1bSJed Brown } else { /* A is in baij format */ 39c4762a1bSJed Brown ierr = MatCreateSeqBAIJ(PETSC_COMM_WORLD,bs,n,n,nz,NULL,&A);CHKERRQ(ierr); 40c4762a1bSJed Brown TestAIJ = PETSC_FALSE; 41c4762a1bSJed Brown } 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* Assemble matrix */ 44c4762a1bSJed Brown if (bs == 1) { 45c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-test_problem",&prob,NULL);CHKERRQ(ierr); 46c4762a1bSJed Brown if (prob == 1) { /* tridiagonal matrix */ 47c4762a1bSJed Brown value[0] = -1.0; value[1] = 2.0; value[2] = -1.0; 48c4762a1bSJed Brown for (i=1; i<n-1; i++) { 49c4762a1bSJed Brown col[0] = i-1; col[1] = i; col[2] = i+1; 50c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 51c4762a1bSJed Brown } 52c4762a1bSJed Brown i = n - 1; col[0]=0; col[1] = n - 2; col[2] = n - 1; 53c4762a1bSJed Brown 54c4762a1bSJed Brown value[0]= 0.1; value[1]=-1; value[2]=2; 55c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 56c4762a1bSJed Brown 57c4762a1bSJed Brown i = 0; col[0] = 0; col[1] = 1; col[2]=n-1; 58c4762a1bSJed Brown 59c4762a1bSJed Brown value[0] = 2.0; value[1] = -1.0; value[2]=0.1; 60c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 61c4762a1bSJed Brown } else if (prob ==2) { /* matrix for the five point stencil */ 62c4762a1bSJed Brown n1 = (PetscInt) (PetscSqrtReal((PetscReal)n) + 0.001); 63*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(n1*n1 - n,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"sqrt(n) must be a positive integer!"); 64c4762a1bSJed Brown for (i=0; i<n1; i++) { 65c4762a1bSJed Brown for (j=0; j<n1; j++) { 66c4762a1bSJed Brown Ii = j + n1*i; 67c4762a1bSJed Brown if (i>0) { 68c4762a1bSJed Brown J = Ii - n1; 69c4762a1bSJed Brown ierr = MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr); 70c4762a1bSJed Brown } 71c4762a1bSJed Brown if (i<n1-1) { 72c4762a1bSJed Brown J = Ii + n1; 73c4762a1bSJed Brown ierr = MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr); 74c4762a1bSJed Brown } 75c4762a1bSJed Brown if (j>0) { 76c4762a1bSJed Brown J = Ii - 1; 77c4762a1bSJed Brown ierr = MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr); 78c4762a1bSJed Brown } 79c4762a1bSJed Brown if (j<n1-1) { 80c4762a1bSJed Brown J = Ii + 1; 81c4762a1bSJed Brown ierr = MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr); 82c4762a1bSJed Brown } 83c4762a1bSJed Brown ierr = MatSetValues(A,1,&Ii,1,&Ii,&four,INSERT_VALUES);CHKERRQ(ierr); 84c4762a1bSJed Brown } 85c4762a1bSJed Brown } 86c4762a1bSJed Brown } 87c4762a1bSJed Brown } else { /* bs > 1 */ 88c4762a1bSJed Brown for (block=0; block<n/bs; block++) { 89c4762a1bSJed Brown /* diagonal blocks */ 90c4762a1bSJed Brown value[0] = -1.0; value[1] = 4.0; value[2] = -1.0; 91c4762a1bSJed Brown for (i=1+block*bs; i<bs-1+block*bs; i++) { 92c4762a1bSJed Brown col[0] = i-1; col[1] = i; col[2] = i+1; 93c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 94c4762a1bSJed Brown } 95c4762a1bSJed Brown i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs; 96c4762a1bSJed Brown 97c4762a1bSJed Brown value[0]=-1.0; value[1]=4.0; 98c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); 99c4762a1bSJed Brown 100c4762a1bSJed Brown i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs; 101c4762a1bSJed Brown 102c4762a1bSJed Brown value[0]=4.0; value[1] = -1.0; 103c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); 104c4762a1bSJed Brown } 105c4762a1bSJed Brown /* off-diagonal blocks */ 106c4762a1bSJed Brown value[0]=-1.0; 107c4762a1bSJed Brown for (i=0; i<(n/bs-1)*bs; i++) { 108c4762a1bSJed Brown col[0]=i+bs; 109c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,1,col,value,INSERT_VALUES);CHKERRQ(ierr); 110c4762a1bSJed Brown col[0]=i; row=i+bs; 111c4762a1bSJed Brown ierr = MatSetValues(A,1,&row,1,col,value,INSERT_VALUES);CHKERRQ(ierr); 112c4762a1bSJed Brown } 113c4762a1bSJed Brown } 114c4762a1bSJed Brown 115c4762a1bSJed Brown if (TestShift) { 116c4762a1bSJed Brown /* set diagonals in the 0-th block as 0 for testing shift numerical factor */ 117c4762a1bSJed Brown for (i=0; i<bs; i++) { 118c4762a1bSJed Brown row = i; col[0] = i; value[0] = 0.0; 119c4762a1bSJed Brown ierr = MatSetValues(A,1,&row,1,col,value,INSERT_VALUES);CHKERRQ(ierr); 120c4762a1bSJed Brown } 121c4762a1bSJed Brown } 122c4762a1bSJed Brown 123c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 124c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 125c4762a1bSJed Brown 126c4762a1bSJed Brown /* Test MatConvert */ 127c4762a1bSJed Brown ierr = MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr); 128c4762a1bSJed Brown ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&sA);CHKERRQ(ierr); 129c4762a1bSJed Brown ierr = MatMultEqual(A,sA,20,&equal);CHKERRQ(ierr); 130*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!equal,PETSC_COMM_SELF,PETSC_ERR_USER,"A != sA"); 131c4762a1bSJed Brown 132c4762a1bSJed Brown /* Test MatGetOwnershipRange() */ 133c4762a1bSJed Brown ierr = MatGetOwnershipRange(A,&Ii,&J);CHKERRQ(ierr); 134c4762a1bSJed Brown ierr = MatGetOwnershipRange(sA,&i,&j);CHKERRQ(ierr); 135*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(i-Ii || j-J,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetOwnershipRange() in MatSBAIJ format"); 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* Vectors */ 138c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_SELF,&rdm);CHKERRQ(ierr); 139c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr); 140c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,n,&x);CHKERRQ(ierr); 141c4762a1bSJed Brown ierr = VecDuplicate(x,&b);CHKERRQ(ierr); 142c4762a1bSJed Brown ierr = VecDuplicate(x,&y);CHKERRQ(ierr); 143c4762a1bSJed Brown ierr = VecSetRandom(x,rdm);CHKERRQ(ierr); 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* Test MatReordering() - not work on sbaij matrix */ 146c4762a1bSJed Brown if (reorder) { 147c4762a1bSJed Brown ierr = MatGetOrdering(A,MATORDERINGRCM,&perm,&cperm);CHKERRQ(ierr); 148c4762a1bSJed Brown } else { 149c4762a1bSJed Brown ierr = MatGetOrdering(A,MATORDERINGNATURAL,&perm,&cperm);CHKERRQ(ierr); 150c4762a1bSJed Brown } 151c4762a1bSJed Brown ierr = ISDestroy(&cperm);CHKERRQ(ierr); 152c4762a1bSJed Brown 153c4762a1bSJed Brown /* initialize factinfo */ 154c4762a1bSJed Brown ierr = MatFactorInfoInitialize(&factinfo);CHKERRQ(ierr); 155c4762a1bSJed Brown if (TestShift == 1) { 156c4762a1bSJed Brown factinfo.shifttype = (PetscReal)MAT_SHIFT_NONZERO; 157c4762a1bSJed Brown factinfo.shiftamount = 0.1; 158c4762a1bSJed Brown } else if (TestShift == 2) { 159c4762a1bSJed Brown factinfo.shifttype = (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE; 160c4762a1bSJed Brown } 161c4762a1bSJed Brown 162c4762a1bSJed Brown /* Test MatCholeskyFactor(), MatICCFactor() */ 163c4762a1bSJed Brown /*------------------------------------------*/ 164c4762a1bSJed Brown /* Test aij matrix A */ 165c4762a1bSJed Brown if (TestAIJ) { 166c4762a1bSJed Brown if (displ) { 167c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"AIJ: \n");CHKERRQ(ierr); 168c4762a1bSJed Brown } 169c4762a1bSJed Brown i = 0; 170c4762a1bSJed Brown for (lvl=-1; lvl<10; lvl++) { 171c4762a1bSJed Brown if (lvl==-1) { /* Cholesky factor */ 172c4762a1bSJed Brown factinfo.fill = 5.0; 173c4762a1bSJed Brown 174c4762a1bSJed Brown ierr = MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&sC);CHKERRQ(ierr); 175c4762a1bSJed Brown ierr = MatCholeskyFactorSymbolic(sC,A,perm,&factinfo);CHKERRQ(ierr); 176c4762a1bSJed Brown } else { /* incomplete Cholesky factor */ 177c4762a1bSJed Brown factinfo.fill = 5.0; 178c4762a1bSJed Brown factinfo.levels = lvl; 179c4762a1bSJed Brown 180c4762a1bSJed Brown ierr = MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_ICC,&sC);CHKERRQ(ierr); 181c4762a1bSJed Brown ierr = MatICCFactorSymbolic(sC,A,perm,&factinfo);CHKERRQ(ierr); 182c4762a1bSJed Brown } 183c4762a1bSJed Brown ierr = MatCholeskyFactorNumeric(sC,A,&factinfo);CHKERRQ(ierr); 184c4762a1bSJed Brown 185c4762a1bSJed Brown ierr = MatMult(A,x,b);CHKERRQ(ierr); 186c4762a1bSJed Brown ierr = MatSolve(sC,b,y);CHKERRQ(ierr); 187c4762a1bSJed Brown ierr = MatDestroy(&sC);CHKERRQ(ierr); 188c4762a1bSJed Brown 189c4762a1bSJed Brown /* Check the residual */ 190c4762a1bSJed Brown ierr = VecAXPY(y,neg_one,x);CHKERRQ(ierr); 191c4762a1bSJed Brown ierr = VecNorm(y,NORM_2,&norm2);CHKERRQ(ierr); 192c4762a1bSJed Brown 193c4762a1bSJed Brown if (displ) { 1948fffc762SJacob Faibussowitsch ierr = PetscPrintf(PETSC_COMM_WORLD," lvl: %" PetscInt_FMT ", residual: %g\n", lvl,(double)norm2);CHKERRQ(ierr); 195c4762a1bSJed Brown } 196c4762a1bSJed Brown } 197c4762a1bSJed Brown } 198c4762a1bSJed Brown 199c4762a1bSJed Brown /* Test baij matrix A */ 200c4762a1bSJed Brown if (TestBAIJ) { 201c4762a1bSJed Brown if (displ) { 202c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"BAIJ: \n");CHKERRQ(ierr); 203c4762a1bSJed Brown } 204c4762a1bSJed Brown i = 0; 205c4762a1bSJed Brown for (lvl=-1; lvl<10; lvl++) { 206c4762a1bSJed Brown if (lvl==-1) { /* Cholesky factor */ 207c4762a1bSJed Brown factinfo.fill = 5.0; 208c4762a1bSJed Brown 209c4762a1bSJed Brown ierr = MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&sC);CHKERRQ(ierr); 210c4762a1bSJed Brown ierr = MatCholeskyFactorSymbolic(sC,A,perm,&factinfo);CHKERRQ(ierr); 211c4762a1bSJed Brown } else { /* incomplete Cholesky factor */ 212c4762a1bSJed Brown factinfo.fill = 5.0; 213c4762a1bSJed Brown factinfo.levels = lvl; 214c4762a1bSJed Brown 215c4762a1bSJed Brown ierr = MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_ICC,&sC);CHKERRQ(ierr); 216c4762a1bSJed Brown ierr = MatICCFactorSymbolic(sC,A,perm,&factinfo);CHKERRQ(ierr); 217c4762a1bSJed Brown } 218c4762a1bSJed Brown ierr = MatCholeskyFactorNumeric(sC,A,&factinfo);CHKERRQ(ierr); 219c4762a1bSJed Brown 220c4762a1bSJed Brown ierr = MatMult(A,x,b);CHKERRQ(ierr); 221c4762a1bSJed Brown ierr = MatSolve(sC,b,y);CHKERRQ(ierr); 222c4762a1bSJed Brown ierr = MatDestroy(&sC);CHKERRQ(ierr); 223c4762a1bSJed Brown 224c4762a1bSJed Brown /* Check the residual */ 225c4762a1bSJed Brown ierr = VecAXPY(y,neg_one,x);CHKERRQ(ierr); 226c4762a1bSJed Brown ierr = VecNorm(y,NORM_2,&norm2);CHKERRQ(ierr); 227c4762a1bSJed Brown if (displ) { 2288fffc762SJacob Faibussowitsch ierr = PetscPrintf(PETSC_COMM_WORLD," lvl: %" PetscInt_FMT ", residual: %g\n", lvl,(double)norm2);CHKERRQ(ierr); 229c4762a1bSJed Brown } 230c4762a1bSJed Brown } 231c4762a1bSJed Brown } 232c4762a1bSJed Brown 233c4762a1bSJed Brown /* Test sbaij matrix sA */ 234c4762a1bSJed Brown if (displ) { 235c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"SBAIJ: \n");CHKERRQ(ierr); 236c4762a1bSJed Brown } 237c4762a1bSJed Brown i = 0; 238c4762a1bSJed Brown for (lvl=-1; lvl<10; lvl++) { 239c4762a1bSJed Brown if (lvl==-1) { /* Cholesky factor */ 240c4762a1bSJed Brown factinfo.fill = 5.0; 241c4762a1bSJed Brown 242c4762a1bSJed Brown ierr = MatGetFactor(sA,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&sC);CHKERRQ(ierr); 243c4762a1bSJed Brown ierr = MatCholeskyFactorSymbolic(sC,sA,perm,&factinfo);CHKERRQ(ierr); 244c4762a1bSJed Brown } else { /* incomplete Cholesky factor */ 245c4762a1bSJed Brown factinfo.fill = 5.0; 246c4762a1bSJed Brown factinfo.levels = lvl; 247c4762a1bSJed Brown 248c4762a1bSJed Brown ierr = MatGetFactor(sA,MATSOLVERPETSC,MAT_FACTOR_ICC,&sC);CHKERRQ(ierr); 249c4762a1bSJed Brown ierr = MatICCFactorSymbolic(sC,sA,perm,&factinfo);CHKERRQ(ierr); 250c4762a1bSJed Brown } 251c4762a1bSJed Brown ierr = MatCholeskyFactorNumeric(sC,sA,&factinfo);CHKERRQ(ierr); 252c4762a1bSJed Brown 253c4762a1bSJed Brown if (lvl==0 && bs==1) { /* Test inplace ICC(0) for sbaij sA - does not work for new datastructure */ 254c4762a1bSJed Brown /* 255c4762a1bSJed Brown Mat B; 256c4762a1bSJed Brown ierr = MatDuplicate(sA,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 257c4762a1bSJed Brown ierr = MatICCFactor(B,perm,&factinfo);CHKERRQ(ierr); 258c4762a1bSJed Brown ierr = MatEqual(sC,B,&equal);CHKERRQ(ierr); 259c4762a1bSJed Brown if (!equal) { 260c4762a1bSJed Brown SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"in-place Cholesky factor != out-place Cholesky factor"); 261c4762a1bSJed Brown } 262c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 263c4762a1bSJed Brown */ 264c4762a1bSJed Brown } 265c4762a1bSJed Brown 266c4762a1bSJed Brown ierr = MatMult(sA,x,b);CHKERRQ(ierr); 267c4762a1bSJed Brown ierr = MatSolve(sC,b,y);CHKERRQ(ierr); 268c4762a1bSJed Brown 269c4762a1bSJed Brown /* Test MatSolves() */ 270c4762a1bSJed Brown if (bs == 1) { 271c4762a1bSJed Brown Vecs xx,bb; 272c4762a1bSJed Brown ierr = VecsCreateSeq(PETSC_COMM_SELF,n,4,&xx);CHKERRQ(ierr); 273c4762a1bSJed Brown ierr = VecsDuplicate(xx,&bb);CHKERRQ(ierr); 274c4762a1bSJed Brown ierr = MatSolves(sC,bb,xx);CHKERRQ(ierr); 275c4762a1bSJed Brown ierr = VecsDestroy(xx);CHKERRQ(ierr); 276c4762a1bSJed Brown ierr = VecsDestroy(bb);CHKERRQ(ierr); 277c4762a1bSJed Brown } 278c4762a1bSJed Brown ierr = MatDestroy(&sC);CHKERRQ(ierr); 279c4762a1bSJed Brown 280c4762a1bSJed Brown /* Check the residual */ 281c4762a1bSJed Brown ierr = VecAXPY(y,neg_one,x);CHKERRQ(ierr); 282c4762a1bSJed Brown ierr = VecNorm(y,NORM_2,&norm2);CHKERRQ(ierr); 283c4762a1bSJed Brown if (displ) { 2848fffc762SJacob Faibussowitsch ierr = PetscPrintf(PETSC_COMM_WORLD," lvl: %" PetscInt_FMT ", residual: %g\n", lvl,(double)norm2);CHKERRQ(ierr); 285c4762a1bSJed Brown } 286c4762a1bSJed Brown } 287c4762a1bSJed Brown 288c4762a1bSJed Brown ierr = ISDestroy(&perm);CHKERRQ(ierr); 289c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 290c4762a1bSJed Brown ierr = MatDestroy(&sA);CHKERRQ(ierr); 291c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 292c4762a1bSJed Brown ierr = VecDestroy(&y);CHKERRQ(ierr); 293c4762a1bSJed Brown ierr = VecDestroy(&b);CHKERRQ(ierr); 294c4762a1bSJed Brown ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr); 295c4762a1bSJed Brown 296c4762a1bSJed Brown ierr = PetscFinalize(); 297c4762a1bSJed Brown return ierr; 298c4762a1bSJed Brown } 299c4762a1bSJed Brown 300c4762a1bSJed Brown /*TEST 301c4762a1bSJed Brown 302c4762a1bSJed Brown test: 303c4762a1bSJed Brown args: -bs {{1 2 3 4 5 6 7 8}} 304c4762a1bSJed Brown 305c4762a1bSJed Brown test: 306c4762a1bSJed Brown suffix: 3 307c4762a1bSJed Brown args: -testaij 308c4762a1bSJed Brown output_file: output/ex76_1.out 309c4762a1bSJed Brown 310c4762a1bSJed Brown TEST*/ 311