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*327415f7SBarry Smith PetscFunctionBeginUser; 199566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 20c4762a1bSJed Brown /* Test MatSetValues() and MatGetValues() */ 219566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL)); 229566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_size",&m,NULL)); 23c4762a1bSJed Brown M = m*bs; 249566063dSJacob Faibussowitsch PetscCall(MatCreateSeqBAIJ(PETSC_COMM_SELF,bs,M,M,1,NULL,&A)); 259566063dSJacob Faibussowitsch PetscCall(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 269566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,M,M,15,NULL,&B)); 279566063dSJacob Faibussowitsch PetscCall(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 289566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&rdm)); 299566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 309566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,M,&xx)); 319566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx,&s1)); 329566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx,&s2)); 339566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx,&yy)); 34c4762a1bSJed Brown 35c4762a1bSJed Brown /* For each row add at least 15 elements */ 36c4762a1bSJed Brown for (row=0; row<M; row++) { 37c4762a1bSJed Brown for (i=0; i<25*bs; i++) { 389566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm,&rval)); 39c4762a1bSJed Brown col = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M)); 409566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&row,1,&col,&rval,INSERT_VALUES)); 419566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,1,&row,1,&col,&rval,INSERT_VALUES)); 42c4762a1bSJed Brown } 43c4762a1bSJed Brown } 44c4762a1bSJed Brown 45c4762a1bSJed Brown /* Now set blocks of values */ 46c4762a1bSJed Brown for (i=0; i<20*bs; i++) { 479566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm,&rval)); 48c4762a1bSJed Brown cols[0] = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M)); 49c4762a1bSJed Brown vals1[0] = rval; 509566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm,&rval)); 51c4762a1bSJed Brown cols[1] = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M)); 52c4762a1bSJed Brown vals1[1] = rval; 539566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm,&rval)); 54c4762a1bSJed Brown rows[0] = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M)); 55c4762a1bSJed Brown vals1[2] = rval; 569566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm,&rval)); 57c4762a1bSJed Brown rows[1] = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M)); 58c4762a1bSJed Brown vals1[3] = rval; 599566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,2,rows,2,cols,vals1,INSERT_VALUES)); 609566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,2,rows,2,cols,vals1,INSERT_VALUES)); 61c4762a1bSJed Brown } 62c4762a1bSJed Brown 639566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 649566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 659566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 669566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* Test MatNorm() */ 699566063dSJacob Faibussowitsch PetscCall(MatNorm(A,NORM_FROBENIUS,&s1norm)); 709566063dSJacob Faibussowitsch PetscCall(MatNorm(B,NORM_FROBENIUS,&s2norm)); 71c4762a1bSJed Brown rnorm = PetscAbsReal(s2norm-s1norm)/s2norm; 72c4762a1bSJed Brown if (rnorm>tol) { 739566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_FROBENIUS()- NormA=%16.14e NormB=%16.14e bs = %" PetscInt_FMT "\n",(double)s1norm,(double)s2norm,bs)); 74c4762a1bSJed Brown } 759566063dSJacob Faibussowitsch PetscCall(MatNorm(A,NORM_INFINITY,&s1norm)); 769566063dSJacob Faibussowitsch PetscCall(MatNorm(B,NORM_INFINITY,&s2norm)); 77c4762a1bSJed Brown rnorm = PetscAbsReal(s2norm-s1norm)/s2norm; 78c4762a1bSJed Brown if (rnorm>tol) { 799566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_INFINITY()- NormA=%16.14e NormB=%16.14e bs = %" PetscInt_FMT "\n",(double)s1norm,(double)s2norm,bs)); 80c4762a1bSJed Brown } 819566063dSJacob Faibussowitsch PetscCall(MatNorm(A,NORM_1,&s1norm)); 829566063dSJacob Faibussowitsch PetscCall(MatNorm(B,NORM_1,&s2norm)); 83c4762a1bSJed Brown rnorm = PetscAbsReal(s2norm-s1norm)/s2norm; 84c4762a1bSJed Brown if (rnorm>tol) { 859566063dSJacob 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)); 86c4762a1bSJed Brown } 87c4762a1bSJed Brown 88c4762a1bSJed Brown /* MatShift() */ 89c4762a1bSJed Brown rval = 10*s1norm; 909566063dSJacob Faibussowitsch PetscCall(MatShift(A,rval)); 919566063dSJacob Faibussowitsch PetscCall(MatShift(B,rval)); 92c4762a1bSJed Brown 93c4762a1bSJed Brown /* Test MatTranspose() */ 949566063dSJacob Faibussowitsch PetscCall(MatTranspose(A,MAT_INPLACE_MATRIX,&A)); 959566063dSJacob Faibussowitsch PetscCall(MatTranspose(B,MAT_INPLACE_MATRIX,&B)); 96c4762a1bSJed Brown 97c4762a1bSJed Brown /* Now do MatGetValues() */ 98c4762a1bSJed Brown for (i=0; i<30; i++) { 999566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm,&rval)); 100c4762a1bSJed Brown cols[0] = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M)); 1019566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm,&rval)); 102c4762a1bSJed Brown cols[1] = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M)); 1039566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm,&rval)); 104c4762a1bSJed Brown rows[0] = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M)); 1059566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm,&rval)); 106c4762a1bSJed Brown rows[1] = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M)); 1079566063dSJacob Faibussowitsch PetscCall(MatGetValues(A,2,rows,2,cols,vals1)); 1089566063dSJacob Faibussowitsch PetscCall(MatGetValues(B,2,rows,2,cols,vals2)); 1099566063dSJacob Faibussowitsch PetscCall(PetscArraycmp(vals1,vals2,4,&flg)); 110c4762a1bSJed Brown if (!flg) { 1119566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatGetValues bs = %" PetscInt_FMT "\n",bs)); 112c4762a1bSJed Brown } 113c4762a1bSJed Brown } 114c4762a1bSJed Brown 115c4762a1bSJed Brown /* Test MatMult(), MatMultAdd() */ 116c4762a1bSJed Brown for (i=0; i<40; i++) { 1179566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx,rdm)); 1189566063dSJacob Faibussowitsch PetscCall(VecSet(s2,0.0)); 1199566063dSJacob Faibussowitsch PetscCall(MatMult(A,xx,s1)); 1209566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A,xx,s2,s2)); 1219566063dSJacob Faibussowitsch PetscCall(VecNorm(s1,NORM_2,&s1norm)); 1229566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_2,&s2norm)); 123c4762a1bSJed Brown rnorm = s2norm-s1norm; 124c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 1259566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"MatMult not equal to MatMultAdd Norm1=%e Norm2=%e bs = %" PetscInt_FMT "\n",(double)s1norm,(double)s2norm,bs)); 126c4762a1bSJed Brown } 127c4762a1bSJed Brown } 128c4762a1bSJed Brown 129c4762a1bSJed Brown /* Test MatMult() */ 1309566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A,B,10,&flg)); 131c4762a1bSJed Brown if (!flg) { 1329566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatMult()\n")); 133c4762a1bSJed Brown } 134c4762a1bSJed Brown 135c4762a1bSJed Brown /* Test MatMultAdd() */ 1369566063dSJacob Faibussowitsch PetscCall(MatMultAddEqual(A,B,10,&flg)); 137c4762a1bSJed Brown if (!flg) { 1389566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatMultAdd()\n")); 139c4762a1bSJed Brown } 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* Test MatMultTranspose() */ 1429566063dSJacob Faibussowitsch PetscCall(MatMultTransposeEqual(A,B,10,&flg)); 143c4762a1bSJed Brown if (!flg) { 1449566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatMultTranspose()\n")); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* Test MatMultTransposeAdd() */ 1489566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAddEqual(A,B,10,&flg)); 149c4762a1bSJed Brown if (!flg) { 1509566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatMultTransposeAdd()\n")); 151c4762a1bSJed Brown } 152c4762a1bSJed Brown 153c4762a1bSJed Brown /* Test MatMatMult() */ 1549566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,M,40,NULL,&C)); 1559566063dSJacob Faibussowitsch PetscCall(MatSetRandom(C,rdm)); 1569566063dSJacob Faibussowitsch PetscCall(MatMatMult(A,C,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&D)); 1579566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(A,C,D,40,&flg)); 158c4762a1bSJed Brown if (!flg) { 1599566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatMatMult()\n")); 160c4762a1bSJed Brown } 1619566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 1629566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,M,40,NULL,&D)); 1639566063dSJacob Faibussowitsch PetscCall(MatMatMult(A,C,MAT_REUSE_MATRIX,PETSC_DEFAULT,&D)); 1649566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(A,C,D,40,&flg)); 165c4762a1bSJed Brown if (!flg) { 1669566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatMatMult()\n")); 167c4762a1bSJed Brown } 168c4762a1bSJed Brown 169c4762a1bSJed Brown /* Do LUFactor() on both the matrices */ 1709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M,&idx)); 171c4762a1bSJed Brown for (i=0; i<M; i++) idx[i] = i; 1729566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,M,idx,PETSC_COPY_VALUES,&is1)); 1739566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,M,idx,PETSC_COPY_VALUES,&is2)); 1749566063dSJacob Faibussowitsch PetscCall(PetscFree(idx)); 1759566063dSJacob Faibussowitsch PetscCall(ISSetPermutation(is1)); 1769566063dSJacob Faibussowitsch PetscCall(ISSetPermutation(is2)); 177c4762a1bSJed Brown 1789566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&info)); 179c4762a1bSJed Brown 180c4762a1bSJed Brown info.fill = 2.0; 181c4762a1bSJed Brown info.dtcol = 0.0; 182c4762a1bSJed Brown info.zeropivot = 1.e-14; 183c4762a1bSJed Brown info.pivotinblocks = 1.0; 184c4762a1bSJed Brown 185c4762a1bSJed Brown if (bs < 4) { 1869566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A,"petsc",MAT_FACTOR_LU,&Fact)); 1879566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(Fact,A,is1,is2,&info)); 1889566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(Fact,A,&info)); 1899566063dSJacob Faibussowitsch PetscCall(VecSetRandom(yy,rdm)); 1909566063dSJacob Faibussowitsch PetscCall(MatForwardSolve(Fact,yy,xx)); 1919566063dSJacob Faibussowitsch PetscCall(MatBackwardSolve(Fact,xx,s1)); 1929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Fact)); 1939566063dSJacob Faibussowitsch PetscCall(VecScale(s1,-1.0)); 1949566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A,s1,yy,yy)); 1959566063dSJacob Faibussowitsch PetscCall(VecNorm(yy,NORM_2,&rnorm)); 196c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 1979566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatForwardSolve/MatBackwardSolve - Norm1=%16.14e bs = %" PetscInt_FMT "\n",(double)rnorm,bs)); 198c4762a1bSJed Brown } 199c4762a1bSJed Brown } 200c4762a1bSJed Brown 2019566063dSJacob Faibussowitsch PetscCall(MatLUFactor(B,is1,is2,&info)); 2029566063dSJacob Faibussowitsch PetscCall(MatLUFactor(A,is1,is2,&info)); 203c4762a1bSJed Brown 204c4762a1bSJed Brown /* Test MatSolveAdd() */ 205c4762a1bSJed Brown for (i=0; i<10; i++) { 2069566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx,rdm)); 2079566063dSJacob Faibussowitsch PetscCall(VecSetRandom(yy,rdm)); 2089566063dSJacob Faibussowitsch PetscCall(MatSolveAdd(B,xx,yy,s2)); 2099566063dSJacob Faibussowitsch PetscCall(MatSolveAdd(A,xx,yy,s1)); 2109566063dSJacob Faibussowitsch PetscCall(VecNorm(s1,NORM_2,&s1norm)); 2119566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_2,&s2norm)); 212c4762a1bSJed Brown rnorm = s2norm-s1norm; 213c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 2149566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatSolveAdd - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n",(double)s1norm,(double)s2norm,bs)); 215c4762a1bSJed Brown } 216c4762a1bSJed Brown } 217c4762a1bSJed Brown 218c4762a1bSJed Brown /* Test MatSolveAdd() when x = A'b +x */ 219c4762a1bSJed Brown for (i=0; i<10; i++) { 2209566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx,rdm)); 2219566063dSJacob Faibussowitsch PetscCall(VecSetRandom(s1,rdm)); 2229566063dSJacob Faibussowitsch PetscCall(VecCopy(s2,s1)); 2239566063dSJacob Faibussowitsch PetscCall(MatSolveAdd(B,xx,s2,s2)); 2249566063dSJacob Faibussowitsch PetscCall(MatSolveAdd(A,xx,s1,s1)); 2259566063dSJacob Faibussowitsch PetscCall(VecNorm(s1,NORM_2,&s1norm)); 2269566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_2,&s2norm)); 227c4762a1bSJed Brown rnorm = s2norm-s1norm; 228c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 2299566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatSolveAdd(same) - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n",(double)s1norm,(double)s2norm,bs)); 230c4762a1bSJed Brown } 231c4762a1bSJed Brown } 232c4762a1bSJed Brown 233c4762a1bSJed Brown /* Test MatSolve() */ 234c4762a1bSJed Brown for (i=0; i<10; i++) { 2359566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx,rdm)); 2369566063dSJacob Faibussowitsch PetscCall(MatSolve(B,xx,s2)); 2379566063dSJacob Faibussowitsch PetscCall(MatSolve(A,xx,s1)); 2389566063dSJacob Faibussowitsch PetscCall(VecNorm(s1,NORM_2,&s1norm)); 2399566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_2,&s2norm)); 240c4762a1bSJed Brown rnorm = s2norm-s1norm; 241c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 2429566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatSolve - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n",(double)s1norm,(double)s2norm,bs)); 243c4762a1bSJed Brown } 244c4762a1bSJed Brown } 245c4762a1bSJed Brown 246c4762a1bSJed Brown /* Test MatSolveTranspose() */ 247c4762a1bSJed Brown if (bs < 8) { 248c4762a1bSJed Brown for (i=0; i<10; i++) { 2499566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx,rdm)); 2509566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(B,xx,s2)); 2519566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(A,xx,s1)); 2529566063dSJacob Faibussowitsch PetscCall(VecNorm(s1,NORM_2,&s1norm)); 2539566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_2,&s2norm)); 254c4762a1bSJed Brown rnorm = s2norm-s1norm; 255c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 2569566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatSolveTranspose - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n",(double)s1norm,(double)s2norm,bs)); 257c4762a1bSJed Brown } 258c4762a1bSJed Brown } 259c4762a1bSJed Brown } 260c4762a1bSJed Brown 2619566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2629566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 2639566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 2649566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 2659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xx)); 2669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 2679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 2689566063dSJacob Faibussowitsch PetscCall(VecDestroy(&yy)); 2699566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 2709566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is2)); 2719566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 2729566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 273b122ec5aSJacob Faibussowitsch return 0; 274c4762a1bSJed Brown } 275c4762a1bSJed Brown 276c4762a1bSJed Brown /*TEST 277c4762a1bSJed Brown 278c4762a1bSJed Brown test: 279c4762a1bSJed Brown args: -mat_block_size {{1 2 3 4 5 6 7 8}} 280c4762a1bSJed Brown 281c4762a1bSJed Brown TEST*/ 282