1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests various routines in MatMPISBAIJ format.\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,u,s1,s2; 9c4762a1bSJed Brown Mat A,sA,sB; 10c4762a1bSJed Brown PetscRandom rctx; 11c4762a1bSJed Brown PetscReal r1,r2,rnorm,tol = PETSC_SQRT_MACHINE_EPSILON; 12c4762a1bSJed Brown PetscScalar one=1.0, neg_one=-1.0, value[3], four=4.0,alpha=0.1; 1343359b5eSHong Zhang PetscInt n,col[3],n1,block,row,i,j,i2,j2,Ii,J,rstart,rend,bs=1,mbs=16,d_nz=3,o_nz=3,prob=1; 14c4762a1bSJed Brown PetscErrorCode ierr; 15c4762a1bSJed Brown PetscMPIInt size,rank; 16c4762a1bSJed Brown PetscBool flg; 17c4762a1bSJed Brown 18c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL)); 205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL)); 215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-prob",&prob,NULL)); 22c4762a1bSJed Brown 235f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 245f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 25c4762a1bSJed Brown 2643359b5eSHong Zhang /* Create a BAIJ matrix A */ 27c4762a1bSJed Brown n = mbs*bs; 285f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATBAIJ)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIBAIJSetPreallocation(A,bs,d_nz,NULL,o_nz,NULL)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqBAIJSetPreallocation(A,bs,d_nz,NULL)); 345f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 35c4762a1bSJed Brown 36c4762a1bSJed Brown if (bs == 1) { 37c4762a1bSJed Brown if (prob == 1) { /* tridiagonal matrix */ 3843359b5eSHong Zhang value[0] = -1.0; value[1] = 0.0; value[2] = -1.0; 39c4762a1bSJed Brown for (i=1; i<n-1; i++) { 40c4762a1bSJed Brown col[0] = i-1; col[1] = i; col[2] = i+1; 415f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 42c4762a1bSJed Brown } 43c4762a1bSJed Brown i = n - 1; col[0]=0; col[1] = n - 2; col[2] = n - 1; 4443359b5eSHong Zhang value[0]= 0.1; value[1]=-1.0; value[2]=0.0; 455f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 46c4762a1bSJed Brown 47c4762a1bSJed Brown i = 0; col[0] = 0; col[1] = 1; col[2]=n-1; 4843359b5eSHong Zhang value[0] = 0.0; value[1] = -1.0; value[2]=0.1; 495f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 50c4762a1bSJed Brown } else if (prob ==2) { /* matrix for the five point stencil */ 51c4762a1bSJed Brown n1 = (int) PetscSqrtReal((PetscReal)n); 52c4762a1bSJed Brown for (i=0; i<n1; i++) { 53c4762a1bSJed Brown for (j=0; j<n1; j++) { 54c4762a1bSJed Brown Ii = j + n1*i; 555f80ce2aSJacob Faibussowitsch if (i>0) {J = Ii - n1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES));} 565f80ce2aSJacob Faibussowitsch if (i<n1-1) {J = Ii + n1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES));} 575f80ce2aSJacob Faibussowitsch if (j>0) {J = Ii - 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES));} 585f80ce2aSJacob Faibussowitsch if (j<n1-1) {J = Ii + 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES));} 595f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&Ii,&four,INSERT_VALUES)); 60c4762a1bSJed Brown } 61c4762a1bSJed Brown } 62c4762a1bSJed Brown } 63c4762a1bSJed Brown /* end of if (bs == 1) */ 64c4762a1bSJed Brown } else { /* bs > 1 */ 65c4762a1bSJed Brown for (block=0; block<n/bs; block++) { 66c4762a1bSJed Brown /* diagonal blocks */ 67c4762a1bSJed Brown value[0] = -1.0; value[1] = 4.0; value[2] = -1.0; 68c4762a1bSJed Brown for (i=1+block*bs; i<bs-1+block*bs; i++) { 69c4762a1bSJed Brown col[0] = i-1; col[1] = i; col[2] = i+1; 705f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 71c4762a1bSJed Brown } 72c4762a1bSJed Brown i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs; 73c4762a1bSJed Brown value[0]=-1.0; value[1]=4.0; 745f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES)); 75c4762a1bSJed Brown 76c4762a1bSJed Brown i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs; 77c4762a1bSJed Brown value[0]=4.0; value[1] = -1.0; 785f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES)); 79c4762a1bSJed Brown } 80c4762a1bSJed Brown /* off-diagonal blocks */ 81c4762a1bSJed Brown value[0]=-1.0; 82c4762a1bSJed Brown for (i=0; i<(n/bs-1)*bs; i++) { 83c4762a1bSJed Brown col[0]=i+bs; 845f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,1,col,value,INSERT_VALUES)); 85c4762a1bSJed Brown col[0]=i; row=i+bs; 865f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&row,1,col,value,INSERT_VALUES)); 87c4762a1bSJed Brown } 88c4762a1bSJed Brown } 895f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 905f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE)); 9243359b5eSHong Zhang 9343359b5eSHong Zhang /* Get SBAIJ matrix sA from A */ 945f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATSBAIJ,MAT_INITIAL_MATRIX,&sA)); 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* Test MatGetSize(), MatGetLocalSize() */ 975f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(sA, &i,&j)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A, &i2,&j2)); 99c4762a1bSJed Brown i -= i2; j -= j2; 100c4762a1bSJed Brown if (i || j) { 1015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatGetSize()\n",rank)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 103c4762a1bSJed Brown } 104c4762a1bSJed Brown 1055f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(sA, &i,&j)); 1065f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(A, &i2,&j2)); 107c4762a1bSJed Brown i2 -= i; j2 -= j; 108c4762a1bSJed Brown if (i2 || j2) { 1095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatGetLocalSize()\n",rank)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 111c4762a1bSJed Brown } 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* vectors */ 114c4762a1bSJed Brown /*--------------------*/ 115c4762a1bSJed Brown /* i is obtained from MatGetLocalSize() */ 1165f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(x,i,PETSC_DECIDE)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(x)); 1195f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&y)); 1205f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&u)); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&s1)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&s2)); 123c4762a1bSJed Brown 1245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rctx)); 1255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rctx)); 1265f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(x,rctx)); 1275f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(u,one)); 128c4762a1bSJed Brown 129c4762a1bSJed Brown /* Test MatNorm() */ 1305f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(A,NORM_FROBENIUS,&r1)); 1315f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(sA,NORM_FROBENIUS,&r2)); 132c4762a1bSJed Brown rnorm = PetscAbsReal(r1-r2)/r2; 133dd400576SPatrick Sanan if (rnorm > tol && rank == 0) { 1345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_FROBENIUS(), Anorm=%16.14e, sAnorm=%16.14e bs=%" PetscInt_FMT "\n",(double)r1,(double)r2,bs)); 135c4762a1bSJed Brown } 1365f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(A,NORM_INFINITY,&r1)); 1375f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(sA,NORM_INFINITY,&r2)); 138c4762a1bSJed Brown rnorm = PetscAbsReal(r1-r2)/r2; 139dd400576SPatrick Sanan if (rnorm > tol && rank == 0) { 1405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error: MatNorm_INFINITY(), Anorm=%16.14e, sAnorm=%16.14e bs=%" PetscInt_FMT "\n",(double)r1,(double)r2,bs)); 141c4762a1bSJed Brown } 1425f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(A,NORM_1,&r1)); 1435f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(sA,NORM_1,&r2)); 144c4762a1bSJed Brown rnorm = PetscAbsReal(r1-r2)/r2; 145dd400576SPatrick Sanan if (rnorm > tol && rank == 0) { 1465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error: MatNorm_1(), Anorm=%16.14e, sAnorm=%16.14e bs=%" PetscInt_FMT "\n",(double)r1,(double)r2,bs)); 147c4762a1bSJed Brown } 148c4762a1bSJed Brown 149c4762a1bSJed Brown /* Test MatGetOwnershipRange() */ 1505f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(sA,&rstart,&rend)); 1515f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A,&i2,&j2)); 152c4762a1bSJed Brown i2 -= rstart; j2 -= rend; 153c4762a1bSJed Brown if (i2 || j2) { 1545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MaGetOwnershipRange()\n",rank)); 1555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 156c4762a1bSJed Brown } 157c4762a1bSJed Brown 158c4762a1bSJed Brown /* Test MatDiagonalScale() */ 1595f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(A,x,x)); 1605f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(sA,x,x)); 1615f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(A,sA,10,&flg)); 162*28b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatDiagonalScale"); 163c4762a1bSJed Brown 164c4762a1bSJed Brown /* Test MatGetDiagonal(), MatScale() */ 1655f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetDiagonal(A,s1)); 1665f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetDiagonal(sA,s2)); 1675f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s1,NORM_1,&r1)); 1685f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s2,NORM_1,&r2)); 169c4762a1bSJed Brown r1 -= r2; 170c4762a1bSJed Brown if (r1<-tol || r1>tol) { 1715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatDiagonalScale() or MatGetDiagonal(), r1=%g \n",rank,(double)r1)); 1725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 173c4762a1bSJed Brown } 174c4762a1bSJed Brown 1755f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(A,alpha)); 1765f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(sA,alpha)); 177c4762a1bSJed Brown 178c4762a1bSJed Brown /* Test MatGetRowMaxAbs() */ 1795f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowMaxAbs(A,s1,NULL)); 1805f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowMaxAbs(sA,s2,NULL)); 181c4762a1bSJed Brown 1825f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s1,NORM_1,&r1)); 1835f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s2,NORM_1,&r2)); 184c4762a1bSJed Brown r1 -= r2; 185c4762a1bSJed Brown if (r1<-tol || r1>tol) { 1865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatGetRowMaxAbs() \n")); 187c4762a1bSJed Brown } 188c4762a1bSJed Brown 189c4762a1bSJed Brown /* Test MatMult(), MatMultAdd() */ 1905f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(A,sA,10,&flg)); 191c4762a1bSJed Brown if (!flg) { 1925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatMult() or MatScale()\n",rank)); 1935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 194c4762a1bSJed Brown } 195c4762a1bSJed Brown 1965f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAddEqual(A,sA,10,&flg)); 197c4762a1bSJed Brown if (!flg) { 1985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatMultAdd()\n",rank)); 1995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 200c4762a1bSJed Brown } 201c4762a1bSJed Brown 202c4762a1bSJed Brown /* Test MatMultTranspose(), MatMultTransposeAdd() */ 203c4762a1bSJed Brown for (i=0; i<10; i++) { 2045f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(x,rctx)); 2055f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(A,x,s1)); 2065f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(sA,x,s2)); 2075f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s1,NORM_1,&r1)); 2085f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s2,NORM_1,&r2)); 209c4762a1bSJed Brown r1 -= r2; 210c4762a1bSJed Brown if (r1<-tol || r1>tol) { 2115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatMult() or MatScale(), err=%g\n",rank,(double)r1)); 2125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 213c4762a1bSJed Brown } 214c4762a1bSJed Brown } 215c4762a1bSJed Brown for (i=0; i<10; i++) { 2165f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(x,rctx)); 2175f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(y,rctx)); 2185f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeAdd(A,x,y,s1)); 2195f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeAdd(sA,x,y,s2)); 2205f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s1,NORM_1,&r1)); 2215f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s2,NORM_1,&r2)); 222c4762a1bSJed Brown r1 -= r2; 223c4762a1bSJed Brown if (r1<-tol || r1>tol) { 2245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatMultAdd(), err=%g \n",rank,(double)r1)); 2255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 226c4762a1bSJed Brown } 227c4762a1bSJed Brown } 228c4762a1bSJed Brown 229c4762a1bSJed Brown /* Test MatDuplicate() */ 2305f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(sA,MAT_COPY_VALUES,&sB)); 2315f80ce2aSJacob Faibussowitsch CHKERRQ(MatEqual(sA,sB,&flg)); 232c4762a1bSJed Brown if (!flg) { 2335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," Error in MatDuplicate(), sA != sB \n")); 234c4762a1bSJed Brown } 2355f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(sA,sB,5,&flg)); 236c4762a1bSJed Brown if (!flg) { 2375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatDuplicate() or MatMult()\n",rank)); 2385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 239c4762a1bSJed Brown } 2405f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAddEqual(sA,sB,5,&flg)); 241c4762a1bSJed Brown if (!flg) { 2425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatDuplicate() or MatMultAdd(()\n",rank)); 2435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 244c4762a1bSJed Brown } 2455f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&sB)); 2465f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u)); 2475f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 2485f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); 2495f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&s1)); 2505f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&s2)); 2515f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&sA)); 2525f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 2535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rctx)); 254c4762a1bSJed Brown ierr = PetscFinalize(); 255c4762a1bSJed Brown return ierr; 256c4762a1bSJed Brown } 257c4762a1bSJed Brown 258c4762a1bSJed Brown /*TEST 259c4762a1bSJed Brown 260c4762a1bSJed Brown test: 261c4762a1bSJed Brown nsize: {{1 3}} 26243359b5eSHong Zhang args: -bs {{1 2 3 5 7 8}} -mat_ignore_lower_triangular -prob {{1 2}} 263c4762a1bSJed Brown 264c4762a1bSJed Brown TEST*/ 265