1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Tests various routines in MatSeqBAIJ format.\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 Mat A,B,C,D,Fact; 9*c4762a1bSJed Brown Vec xx,s1,s2,yy; 10*c4762a1bSJed Brown PetscErrorCode ierr; 11*c4762a1bSJed Brown PetscInt m=45,rows[2],cols[2],bs=1,i,row,col,*idx,M; 12*c4762a1bSJed Brown PetscScalar rval,vals1[4],vals2[4]; 13*c4762a1bSJed Brown PetscRandom rdm; 14*c4762a1bSJed Brown IS is1,is2; 15*c4762a1bSJed Brown PetscReal s1norm,s2norm,rnorm,tol = 1.e-4; 16*c4762a1bSJed Brown PetscBool flg; 17*c4762a1bSJed Brown MatFactorInfo info; 18*c4762a1bSJed Brown 19*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 20*c4762a1bSJed Brown /* Test MatSetValues() and MatGetValues() */ 21*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL);CHKERRQ(ierr); 22*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-mat_size",&m,NULL);CHKERRQ(ierr); 23*c4762a1bSJed Brown M = m*bs; 24*c4762a1bSJed Brown ierr = MatCreateSeqBAIJ(PETSC_COMM_SELF,bs,M,M,1,NULL,&A);CHKERRQ(ierr); 25*c4762a1bSJed Brown ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 26*c4762a1bSJed Brown ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,M,M,15,NULL,&B);CHKERRQ(ierr); 27*c4762a1bSJed Brown ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 28*c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_SELF,&rdm);CHKERRQ(ierr); 29*c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr); 30*c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,M,&xx);CHKERRQ(ierr); 31*c4762a1bSJed Brown ierr = VecDuplicate(xx,&s1);CHKERRQ(ierr); 32*c4762a1bSJed Brown ierr = VecDuplicate(xx,&s2);CHKERRQ(ierr); 33*c4762a1bSJed Brown ierr = VecDuplicate(xx,&yy);CHKERRQ(ierr); 34*c4762a1bSJed Brown 35*c4762a1bSJed Brown /* For each row add atleast 15 elements */ 36*c4762a1bSJed Brown for (row=0; row<M; row++) { 37*c4762a1bSJed Brown for (i=0; i<25*bs; i++) { 38*c4762a1bSJed Brown ierr = PetscRandomGetValue(rdm,&rval);CHKERRQ(ierr); 39*c4762a1bSJed Brown col = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M)); 40*c4762a1bSJed Brown ierr = MatSetValues(A,1,&row,1,&col,&rval,INSERT_VALUES);CHKERRQ(ierr); 41*c4762a1bSJed Brown ierr = MatSetValues(B,1,&row,1,&col,&rval,INSERT_VALUES);CHKERRQ(ierr); 42*c4762a1bSJed Brown } 43*c4762a1bSJed Brown } 44*c4762a1bSJed Brown 45*c4762a1bSJed Brown /* Now set blocks of values */ 46*c4762a1bSJed Brown for (i=0; i<20*bs; i++) { 47*c4762a1bSJed Brown ierr = PetscRandomGetValue(rdm,&rval);CHKERRQ(ierr); 48*c4762a1bSJed Brown cols[0] = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M)); 49*c4762a1bSJed Brown vals1[0] = rval; 50*c4762a1bSJed Brown ierr = PetscRandomGetValue(rdm,&rval);CHKERRQ(ierr); 51*c4762a1bSJed Brown cols[1] = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M)); 52*c4762a1bSJed Brown vals1[1] = rval; 53*c4762a1bSJed Brown ierr = PetscRandomGetValue(rdm,&rval);CHKERRQ(ierr); 54*c4762a1bSJed Brown rows[0] = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M)); 55*c4762a1bSJed Brown vals1[2] = rval; 56*c4762a1bSJed Brown ierr = PetscRandomGetValue(rdm,&rval);CHKERRQ(ierr); 57*c4762a1bSJed Brown rows[1] = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M)); 58*c4762a1bSJed Brown vals1[3] = rval; 59*c4762a1bSJed Brown ierr = MatSetValues(A,2,rows,2,cols,vals1,INSERT_VALUES);CHKERRQ(ierr); 60*c4762a1bSJed Brown ierr = MatSetValues(B,2,rows,2,cols,vals1,INSERT_VALUES);CHKERRQ(ierr); 61*c4762a1bSJed Brown } 62*c4762a1bSJed Brown 63*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 65*c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 66*c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 67*c4762a1bSJed Brown 68*c4762a1bSJed Brown /* Test MatNorm() */ 69*c4762a1bSJed Brown ierr = MatNorm(A,NORM_FROBENIUS,&s1norm);CHKERRQ(ierr); 70*c4762a1bSJed Brown ierr = MatNorm(B,NORM_FROBENIUS,&s2norm);CHKERRQ(ierr); 71*c4762a1bSJed Brown rnorm = PetscAbsReal(s2norm-s1norm)/s2norm; 72*c4762a1bSJed Brown if (rnorm>tol) { 73*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_FROBENIUS()- NormA=%16.14e NormB=%16.14e bs = %D\n",s1norm,s2norm,bs);CHKERRQ(ierr); 74*c4762a1bSJed Brown } 75*c4762a1bSJed Brown ierr = MatNorm(A,NORM_INFINITY,&s1norm);CHKERRQ(ierr); 76*c4762a1bSJed Brown ierr = MatNorm(B,NORM_INFINITY,&s2norm);CHKERRQ(ierr); 77*c4762a1bSJed Brown rnorm = PetscAbsReal(s2norm-s1norm)/s2norm; 78*c4762a1bSJed Brown if (rnorm>tol) { 79*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_INFINITY()- NormA=%16.14e NormB=%16.14e bs = %D\n",s1norm,s2norm,bs);CHKERRQ(ierr); 80*c4762a1bSJed Brown } 81*c4762a1bSJed Brown ierr = MatNorm(A,NORM_1,&s1norm);CHKERRQ(ierr); 82*c4762a1bSJed Brown ierr = MatNorm(B,NORM_1,&s2norm);CHKERRQ(ierr); 83*c4762a1bSJed Brown rnorm = PetscAbsReal(s2norm-s1norm)/s2norm; 84*c4762a1bSJed Brown if (rnorm>tol) { 85*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_NORM_1()- NormA=%16.14e NormB=%16.14e bs = %D\n",s1norm,s2norm,bs);CHKERRQ(ierr); 86*c4762a1bSJed Brown } 87*c4762a1bSJed Brown 88*c4762a1bSJed Brown /* MatShift() */ 89*c4762a1bSJed Brown rval = 10*s1norm; 90*c4762a1bSJed Brown ierr = MatShift(A,rval);CHKERRQ(ierr); 91*c4762a1bSJed Brown ierr = MatShift(B,rval);CHKERRQ(ierr); 92*c4762a1bSJed Brown 93*c4762a1bSJed Brown /* Test MatTranspose() */ 94*c4762a1bSJed Brown ierr = MatTranspose(A,MAT_INPLACE_MATRIX,&A);CHKERRQ(ierr); 95*c4762a1bSJed Brown ierr = MatTranspose(B,MAT_INPLACE_MATRIX,&B);CHKERRQ(ierr); 96*c4762a1bSJed Brown 97*c4762a1bSJed Brown /* Now do MatGetValues() */ 98*c4762a1bSJed Brown for (i=0; i<30; i++) { 99*c4762a1bSJed Brown ierr = PetscRandomGetValue(rdm,&rval);CHKERRQ(ierr); 100*c4762a1bSJed Brown cols[0] = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M)); 101*c4762a1bSJed Brown ierr = PetscRandomGetValue(rdm,&rval);CHKERRQ(ierr); 102*c4762a1bSJed Brown cols[1] = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M)); 103*c4762a1bSJed Brown ierr = PetscRandomGetValue(rdm,&rval);CHKERRQ(ierr); 104*c4762a1bSJed Brown rows[0] = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M)); 105*c4762a1bSJed Brown ierr = PetscRandomGetValue(rdm,&rval);CHKERRQ(ierr); 106*c4762a1bSJed Brown rows[1] = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M)); 107*c4762a1bSJed Brown ierr = MatGetValues(A,2,rows,2,cols,vals1);CHKERRQ(ierr); 108*c4762a1bSJed Brown ierr = MatGetValues(B,2,rows,2,cols,vals2);CHKERRQ(ierr); 109*c4762a1bSJed Brown ierr = PetscArraycmp(vals1,vals2,4,&flg);CHKERRQ(ierr); 110*c4762a1bSJed Brown if (!flg) { 111*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatGetValues bs = %D\n",bs);CHKERRQ(ierr); 112*c4762a1bSJed Brown } 113*c4762a1bSJed Brown } 114*c4762a1bSJed Brown 115*c4762a1bSJed Brown /* Test MatMult(), MatMultAdd() */ 116*c4762a1bSJed Brown for (i=0; i<40; i++) { 117*c4762a1bSJed Brown ierr = VecSetRandom(xx,rdm);CHKERRQ(ierr); 118*c4762a1bSJed Brown ierr = VecSet(s2,0.0);CHKERRQ(ierr); 119*c4762a1bSJed Brown ierr = MatMult(A,xx,s1);CHKERRQ(ierr); 120*c4762a1bSJed Brown ierr = MatMultAdd(A,xx,s2,s2);CHKERRQ(ierr); 121*c4762a1bSJed Brown ierr = VecNorm(s1,NORM_2,&s1norm);CHKERRQ(ierr); 122*c4762a1bSJed Brown ierr = VecNorm(s2,NORM_2,&s2norm);CHKERRQ(ierr); 123*c4762a1bSJed Brown rnorm = s2norm-s1norm; 124*c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 125*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"MatMult not equal to MatMultAdd Norm1=%e Norm2=%e bs = %D\n",s1norm,s2norm,bs);CHKERRQ(ierr); 126*c4762a1bSJed Brown } 127*c4762a1bSJed Brown } 128*c4762a1bSJed Brown 129*c4762a1bSJed Brown /* Test MatMult() */ 130*c4762a1bSJed Brown ierr = MatMultEqual(A,B,10,&flg);CHKERRQ(ierr); 131*c4762a1bSJed Brown if (!flg) { 132*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatMult()\n");CHKERRQ(ierr); 133*c4762a1bSJed Brown } 134*c4762a1bSJed Brown 135*c4762a1bSJed Brown /* Test MatMultAdd() */ 136*c4762a1bSJed Brown ierr = MatMultAddEqual(A,B,10,&flg);CHKERRQ(ierr); 137*c4762a1bSJed Brown if (!flg) { 138*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatMultAdd()\n");CHKERRQ(ierr); 139*c4762a1bSJed Brown } 140*c4762a1bSJed Brown 141*c4762a1bSJed Brown /* Test MatMultTranspose() */ 142*c4762a1bSJed Brown ierr = MatMultTransposeEqual(A,B,10,&flg);CHKERRQ(ierr); 143*c4762a1bSJed Brown if (!flg) { 144*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatMultTranspose()\n");CHKERRQ(ierr); 145*c4762a1bSJed Brown } 146*c4762a1bSJed Brown 147*c4762a1bSJed Brown /* Test MatMultTransposeAdd() */ 148*c4762a1bSJed Brown ierr = MatMultTransposeAddEqual(A,B,10,&flg);CHKERRQ(ierr); 149*c4762a1bSJed Brown if (!flg) { 150*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatMultTransposeAdd()\n");CHKERRQ(ierr); 151*c4762a1bSJed Brown } 152*c4762a1bSJed Brown 153*c4762a1bSJed Brown /* Test MatMatMult() */ 154*c4762a1bSJed Brown ierr = MatCreateSeqDense(PETSC_COMM_SELF,M,40,NULL,&C);CHKERRQ(ierr); 155*c4762a1bSJed Brown ierr = MatSetRandom(C,rdm);CHKERRQ(ierr); 156*c4762a1bSJed Brown ierr = MatMatMult(A,C,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&D);CHKERRQ(ierr); 157*c4762a1bSJed Brown ierr = MatMatMultEqual(A,C,D,40,&flg);CHKERRQ(ierr); 158*c4762a1bSJed Brown if (!flg) { 159*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatMatMult()\n");CHKERRQ(ierr); 160*c4762a1bSJed Brown } 161*c4762a1bSJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 162*c4762a1bSJed Brown ierr = MatCreateSeqDense(PETSC_COMM_SELF,M,40,NULL,&D);CHKERRQ(ierr); 163*c4762a1bSJed Brown ierr = MatMatMult(A,C,MAT_REUSE_MATRIX,PETSC_DEFAULT,&D);CHKERRQ(ierr); 164*c4762a1bSJed Brown ierr = MatMatMultEqual(A,C,D,40,&flg);CHKERRQ(ierr); 165*c4762a1bSJed Brown if (!flg) { 166*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatMatMult()\n");CHKERRQ(ierr); 167*c4762a1bSJed Brown } 168*c4762a1bSJed Brown 169*c4762a1bSJed Brown /* Do LUFactor() on both the matrices */ 170*c4762a1bSJed Brown ierr = PetscMalloc1(M,&idx);CHKERRQ(ierr); 171*c4762a1bSJed Brown for (i=0; i<M; i++) idx[i] = i; 172*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,M,idx,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr); 173*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,M,idx,PETSC_COPY_VALUES,&is2);CHKERRQ(ierr); 174*c4762a1bSJed Brown ierr = PetscFree(idx);CHKERRQ(ierr); 175*c4762a1bSJed Brown ierr = ISSetPermutation(is1);CHKERRQ(ierr); 176*c4762a1bSJed Brown ierr = ISSetPermutation(is2);CHKERRQ(ierr); 177*c4762a1bSJed Brown 178*c4762a1bSJed Brown ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr); 179*c4762a1bSJed Brown 180*c4762a1bSJed Brown info.fill = 2.0; 181*c4762a1bSJed Brown info.dtcol = 0.0; 182*c4762a1bSJed Brown info.zeropivot = 1.e-14; 183*c4762a1bSJed Brown info.pivotinblocks = 1.0; 184*c4762a1bSJed Brown 185*c4762a1bSJed Brown if (bs < 4) { 186*c4762a1bSJed Brown ierr = MatGetFactor(A,"petsc",MAT_FACTOR_LU,&Fact);CHKERRQ(ierr); 187*c4762a1bSJed Brown ierr = MatLUFactorSymbolic(Fact,A,is1,is2,&info);CHKERRQ(ierr); 188*c4762a1bSJed Brown ierr = MatLUFactorNumeric(Fact,A,&info);CHKERRQ(ierr); 189*c4762a1bSJed Brown ierr = VecSetRandom(yy,rdm);CHKERRQ(ierr); 190*c4762a1bSJed Brown ierr = MatForwardSolve(Fact,yy,xx);CHKERRQ(ierr); 191*c4762a1bSJed Brown ierr = MatBackwardSolve(Fact,xx,s1);CHKERRQ(ierr); 192*c4762a1bSJed Brown ierr = MatDestroy(&Fact);CHKERRQ(ierr); 193*c4762a1bSJed Brown ierr = VecScale(s1,-1.0);CHKERRQ(ierr); 194*c4762a1bSJed Brown ierr = MatMultAdd(A,s1,yy,yy);CHKERRQ(ierr); 195*c4762a1bSJed Brown ierr = VecNorm(yy,NORM_2,&rnorm);CHKERRQ(ierr); 196*c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 197*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error:MatForwardSolve/MatBackwardSolve - Norm1=%16.14e bs = %D\n",rnorm,bs);CHKERRQ(ierr); 198*c4762a1bSJed Brown } 199*c4762a1bSJed Brown } 200*c4762a1bSJed Brown 201*c4762a1bSJed Brown ierr = MatLUFactor(B,is1,is2,&info);CHKERRQ(ierr); 202*c4762a1bSJed Brown ierr = MatLUFactor(A,is1,is2,&info);CHKERRQ(ierr); 203*c4762a1bSJed Brown 204*c4762a1bSJed Brown /* Test MatSolveAdd() */ 205*c4762a1bSJed Brown for (i=0; i<10; i++) { 206*c4762a1bSJed Brown ierr = VecSetRandom(xx,rdm);CHKERRQ(ierr); 207*c4762a1bSJed Brown ierr = VecSetRandom(yy,rdm);CHKERRQ(ierr); 208*c4762a1bSJed Brown ierr = MatSolveAdd(B,xx,yy,s2);CHKERRQ(ierr); 209*c4762a1bSJed Brown ierr = MatSolveAdd(A,xx,yy,s1);CHKERRQ(ierr); 210*c4762a1bSJed Brown ierr = VecNorm(s1,NORM_2,&s1norm);CHKERRQ(ierr); 211*c4762a1bSJed Brown ierr = VecNorm(s2,NORM_2,&s2norm);CHKERRQ(ierr); 212*c4762a1bSJed Brown rnorm = s2norm-s1norm; 213*c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 214*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error:MatSolveAdd - Norm1=%16.14e Norm2=%16.14e bs = %D\n",s1norm,s2norm,bs);CHKERRQ(ierr); 215*c4762a1bSJed Brown } 216*c4762a1bSJed Brown } 217*c4762a1bSJed Brown 218*c4762a1bSJed Brown /* Test MatSolveAdd() when x = A'b +x */ 219*c4762a1bSJed Brown for (i=0; i<10; i++) { 220*c4762a1bSJed Brown ierr = VecSetRandom(xx,rdm);CHKERRQ(ierr); 221*c4762a1bSJed Brown ierr = VecSetRandom(s1,rdm);CHKERRQ(ierr); 222*c4762a1bSJed Brown ierr = VecCopy(s2,s1);CHKERRQ(ierr); 223*c4762a1bSJed Brown ierr = MatSolveAdd(B,xx,s2,s2);CHKERRQ(ierr); 224*c4762a1bSJed Brown ierr = MatSolveAdd(A,xx,s1,s1);CHKERRQ(ierr); 225*c4762a1bSJed Brown ierr = VecNorm(s1,NORM_2,&s1norm);CHKERRQ(ierr); 226*c4762a1bSJed Brown ierr = VecNorm(s2,NORM_2,&s2norm);CHKERRQ(ierr); 227*c4762a1bSJed Brown rnorm = s2norm-s1norm; 228*c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 229*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error:MatSolveAdd(same) - Norm1=%16.14e Norm2=%16.14e bs = %D\n",s1norm,s2norm,bs);CHKERRQ(ierr); 230*c4762a1bSJed Brown } 231*c4762a1bSJed Brown } 232*c4762a1bSJed Brown 233*c4762a1bSJed Brown /* Test MatSolve() */ 234*c4762a1bSJed Brown for (i=0; i<10; i++) { 235*c4762a1bSJed Brown ierr = VecSetRandom(xx,rdm);CHKERRQ(ierr); 236*c4762a1bSJed Brown ierr = MatSolve(B,xx,s2);CHKERRQ(ierr); 237*c4762a1bSJed Brown ierr = MatSolve(A,xx,s1);CHKERRQ(ierr); 238*c4762a1bSJed Brown ierr = VecNorm(s1,NORM_2,&s1norm);CHKERRQ(ierr); 239*c4762a1bSJed Brown ierr = VecNorm(s2,NORM_2,&s2norm);CHKERRQ(ierr); 240*c4762a1bSJed Brown rnorm = s2norm-s1norm; 241*c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 242*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error:MatSolve - Norm1=%16.14e Norm2=%16.14e bs = %D\n",s1norm,s2norm,bs);CHKERRQ(ierr); 243*c4762a1bSJed Brown } 244*c4762a1bSJed Brown } 245*c4762a1bSJed Brown 246*c4762a1bSJed Brown /* Test MatSolveTranspose() */ 247*c4762a1bSJed Brown if (bs < 8) { 248*c4762a1bSJed Brown for (i=0; i<10; i++) { 249*c4762a1bSJed Brown ierr = VecSetRandom(xx,rdm);CHKERRQ(ierr); 250*c4762a1bSJed Brown ierr = MatSolveTranspose(B,xx,s2);CHKERRQ(ierr); 251*c4762a1bSJed Brown ierr = MatSolveTranspose(A,xx,s1);CHKERRQ(ierr); 252*c4762a1bSJed Brown ierr = VecNorm(s1,NORM_2,&s1norm);CHKERRQ(ierr); 253*c4762a1bSJed Brown ierr = VecNorm(s2,NORM_2,&s2norm);CHKERRQ(ierr); 254*c4762a1bSJed Brown rnorm = s2norm-s1norm; 255*c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 256*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error:MatSolveTranspose - Norm1=%16.14e Norm2=%16.14e bs = %D\n",s1norm,s2norm,bs);CHKERRQ(ierr); 257*c4762a1bSJed Brown } 258*c4762a1bSJed Brown } 259*c4762a1bSJed Brown } 260*c4762a1bSJed Brown 261*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 262*c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 263*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 264*c4762a1bSJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 265*c4762a1bSJed Brown ierr = VecDestroy(&xx);CHKERRQ(ierr); 266*c4762a1bSJed Brown ierr = VecDestroy(&s1);CHKERRQ(ierr); 267*c4762a1bSJed Brown ierr = VecDestroy(&s2);CHKERRQ(ierr); 268*c4762a1bSJed Brown ierr = VecDestroy(&yy);CHKERRQ(ierr); 269*c4762a1bSJed Brown ierr = ISDestroy(&is1);CHKERRQ(ierr); 270*c4762a1bSJed Brown ierr = ISDestroy(&is2);CHKERRQ(ierr); 271*c4762a1bSJed Brown ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr); 272*c4762a1bSJed Brown ierr = PetscFinalize(); 273*c4762a1bSJed Brown return ierr; 274*c4762a1bSJed Brown } 275*c4762a1bSJed Brown 276*c4762a1bSJed Brown 277*c4762a1bSJed Brown /*TEST 278*c4762a1bSJed Brown 279*c4762a1bSJed Brown test: 280*c4762a1bSJed Brown args: -mat_block_size {{1 2 3 4 5 6 7 8}} 281*c4762a1bSJed Brown 282*c4762a1bSJed Brown TEST*/ 283