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