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; 255f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 262c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-reorder",&reorder,NULL)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-testaij",&TestAIJ,NULL)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-testShift",&TestShift,NULL)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-displ",&displ,NULL)); 33c4762a1bSJed Brown 34c4762a1bSJed Brown n = mbs*bs; 35c4762a1bSJed Brown if (TestAIJ) { /* A is in aij format */ 365f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_WORLD,n,n,nz,NULL,&A)); 37c4762a1bSJed Brown TestBAIJ = PETSC_FALSE; 38c4762a1bSJed Brown } else { /* A is in baij format */ 395f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqBAIJ(PETSC_COMM_WORLD,bs,n,n,nz,NULL,&A)); 40c4762a1bSJed Brown TestAIJ = PETSC_FALSE; 41c4762a1bSJed Brown } 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* Assemble matrix */ 44c4762a1bSJed Brown if (bs == 1) { 455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-test_problem",&prob,NULL)); 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; 505f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 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; 555f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 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; 605f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 61c4762a1bSJed Brown } else if (prob ==2) { /* matrix for the five point stencil */ 62c4762a1bSJed Brown n1 = (PetscInt) (PetscSqrtReal((PetscReal)n) + 0.001); 632c71b3e2SJacob 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; 695f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 70c4762a1bSJed Brown } 71c4762a1bSJed Brown if (i<n1-1) { 72c4762a1bSJed Brown J = Ii + n1; 735f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 74c4762a1bSJed Brown } 75c4762a1bSJed Brown if (j>0) { 76c4762a1bSJed Brown J = Ii - 1; 775f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 78c4762a1bSJed Brown } 79c4762a1bSJed Brown if (j<n1-1) { 80c4762a1bSJed Brown J = Ii + 1; 815f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 82c4762a1bSJed Brown } 835f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&Ii,&four,INSERT_VALUES)); 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; 935f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 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; 985f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES)); 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; 1035f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES)); 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; 1095f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,1,col,value,INSERT_VALUES)); 110c4762a1bSJed Brown col[0]=i; row=i+bs; 1115f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&row,1,col,value,INSERT_VALUES)); 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; 1195f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&row,1,col,value,INSERT_VALUES)); 120c4762a1bSJed Brown } 121c4762a1bSJed Brown } 122c4762a1bSJed Brown 1235f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 125c4762a1bSJed Brown 126c4762a1bSJed Brown /* Test MatConvert */ 1275f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE)); 1285f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&sA)); 1295f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(A,sA,20,&equal)); 130*28b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_USER,"A != sA"); 131c4762a1bSJed Brown 132c4762a1bSJed Brown /* Test MatGetOwnershipRange() */ 1335f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A,&Ii,&J)); 1345f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(sA,&i,&j)); 1352c71b3e2SJacob Faibussowitsch PetscCheckFalse(i-Ii || j-J,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetOwnershipRange() in MatSBAIJ format"); 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* Vectors */ 1385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&rdm)); 1395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rdm)); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&x)); 1415f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&b)); 1425f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&y)); 1435f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(x,rdm)); 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* Test MatReordering() - not work on sbaij matrix */ 146c4762a1bSJed Brown if (reorder) { 1475f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOrdering(A,MATORDERINGRCM,&perm,&cperm)); 148c4762a1bSJed Brown } else { 1495f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOrdering(A,MATORDERINGNATURAL,&perm,&cperm)); 150c4762a1bSJed Brown } 1515f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&cperm)); 152c4762a1bSJed Brown 153c4762a1bSJed Brown /* initialize factinfo */ 1545f80ce2aSJacob Faibussowitsch CHKERRQ(MatFactorInfoInitialize(&factinfo)); 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) { 1675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"AIJ: \n")); 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 1745f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&sC)); 1755f80ce2aSJacob Faibussowitsch CHKERRQ(MatCholeskyFactorSymbolic(sC,A,perm,&factinfo)); 176c4762a1bSJed Brown } else { /* incomplete Cholesky factor */ 177c4762a1bSJed Brown factinfo.fill = 5.0; 178c4762a1bSJed Brown factinfo.levels = lvl; 179c4762a1bSJed Brown 1805f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_ICC,&sC)); 1815f80ce2aSJacob Faibussowitsch CHKERRQ(MatICCFactorSymbolic(sC,A,perm,&factinfo)); 182c4762a1bSJed Brown } 1835f80ce2aSJacob Faibussowitsch CHKERRQ(MatCholeskyFactorNumeric(sC,A,&factinfo)); 184c4762a1bSJed Brown 1855f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,x,b)); 1865f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolve(sC,b,y)); 1875f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&sC)); 188c4762a1bSJed Brown 189c4762a1bSJed Brown /* Check the residual */ 1905f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,neg_one,x)); 1915f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(y,NORM_2,&norm2)); 192c4762a1bSJed Brown 193c4762a1bSJed Brown if (displ) { 1945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," lvl: %" PetscInt_FMT ", residual: %g\n", lvl,(double)norm2)); 195c4762a1bSJed Brown } 196c4762a1bSJed Brown } 197c4762a1bSJed Brown } 198c4762a1bSJed Brown 199c4762a1bSJed Brown /* Test baij matrix A */ 200c4762a1bSJed Brown if (TestBAIJ) { 201c4762a1bSJed Brown if (displ) { 2025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"BAIJ: \n")); 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 2095f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&sC)); 2105f80ce2aSJacob Faibussowitsch CHKERRQ(MatCholeskyFactorSymbolic(sC,A,perm,&factinfo)); 211c4762a1bSJed Brown } else { /* incomplete Cholesky factor */ 212c4762a1bSJed Brown factinfo.fill = 5.0; 213c4762a1bSJed Brown factinfo.levels = lvl; 214c4762a1bSJed Brown 2155f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_ICC,&sC)); 2165f80ce2aSJacob Faibussowitsch CHKERRQ(MatICCFactorSymbolic(sC,A,perm,&factinfo)); 217c4762a1bSJed Brown } 2185f80ce2aSJacob Faibussowitsch CHKERRQ(MatCholeskyFactorNumeric(sC,A,&factinfo)); 219c4762a1bSJed Brown 2205f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,x,b)); 2215f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolve(sC,b,y)); 2225f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&sC)); 223c4762a1bSJed Brown 224c4762a1bSJed Brown /* Check the residual */ 2255f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,neg_one,x)); 2265f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(y,NORM_2,&norm2)); 227c4762a1bSJed Brown if (displ) { 2285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," lvl: %" PetscInt_FMT ", residual: %g\n", lvl,(double)norm2)); 229c4762a1bSJed Brown } 230c4762a1bSJed Brown } 231c4762a1bSJed Brown } 232c4762a1bSJed Brown 233c4762a1bSJed Brown /* Test sbaij matrix sA */ 234c4762a1bSJed Brown if (displ) { 2355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"SBAIJ: \n")); 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 2425f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(sA,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&sC)); 2435f80ce2aSJacob Faibussowitsch CHKERRQ(MatCholeskyFactorSymbolic(sC,sA,perm,&factinfo)); 244c4762a1bSJed Brown } else { /* incomplete Cholesky factor */ 245c4762a1bSJed Brown factinfo.fill = 5.0; 246c4762a1bSJed Brown factinfo.levels = lvl; 247c4762a1bSJed Brown 2485f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(sA,MATSOLVERPETSC,MAT_FACTOR_ICC,&sC)); 2495f80ce2aSJacob Faibussowitsch CHKERRQ(MatICCFactorSymbolic(sC,sA,perm,&factinfo)); 250c4762a1bSJed Brown } 2515f80ce2aSJacob Faibussowitsch CHKERRQ(MatCholeskyFactorNumeric(sC,sA,&factinfo)); 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; 2565f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(sA,MAT_COPY_VALUES,&B)); 2575f80ce2aSJacob Faibussowitsch CHKERRQ(MatICCFactor(B,perm,&factinfo)); 2585f80ce2aSJacob Faibussowitsch CHKERRQ(MatEqual(sC,B,&equal)); 259c4762a1bSJed Brown if (!equal) { 260c4762a1bSJed Brown SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"in-place Cholesky factor != out-place Cholesky factor"); 261c4762a1bSJed Brown } 2625f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 263c4762a1bSJed Brown */ 264c4762a1bSJed Brown } 265c4762a1bSJed Brown 2665f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(sA,x,b)); 2675f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolve(sC,b,y)); 268c4762a1bSJed Brown 269c4762a1bSJed Brown /* Test MatSolves() */ 270c4762a1bSJed Brown if (bs == 1) { 271c4762a1bSJed Brown Vecs xx,bb; 2725f80ce2aSJacob Faibussowitsch CHKERRQ(VecsCreateSeq(PETSC_COMM_SELF,n,4,&xx)); 2735f80ce2aSJacob Faibussowitsch CHKERRQ(VecsDuplicate(xx,&bb)); 2745f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolves(sC,bb,xx)); 2755f80ce2aSJacob Faibussowitsch CHKERRQ(VecsDestroy(xx)); 2765f80ce2aSJacob Faibussowitsch CHKERRQ(VecsDestroy(bb)); 277c4762a1bSJed Brown } 2785f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&sC)); 279c4762a1bSJed Brown 280c4762a1bSJed Brown /* Check the residual */ 2815f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,neg_one,x)); 2825f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(y,NORM_2,&norm2)); 283c4762a1bSJed Brown if (displ) { 2845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," lvl: %" PetscInt_FMT ", residual: %g\n", lvl,(double)norm2)); 285c4762a1bSJed Brown } 286c4762a1bSJed Brown } 287c4762a1bSJed Brown 2885f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&perm)); 2895f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 2905f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&sA)); 2915f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 2925f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); 2935f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b)); 2945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rdm)); 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