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