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; 19c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL);CHKERRQ(ierr); 20c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);CHKERRQ(ierr); 2143359b5eSHong Zhang ierr = PetscOptionsGetInt(NULL,NULL,"-prob",&prob,NULL);CHKERRQ(ierr); 22c4762a1bSJed Brown 23ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 24ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 25c4762a1bSJed Brown 2643359b5eSHong Zhang /* Create a BAIJ matrix A */ 27c4762a1bSJed Brown n = mbs*bs; 2843359b5eSHong Zhang ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 2943359b5eSHong Zhang ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);CHKERRQ(ierr); 3043359b5eSHong Zhang ierr = MatSetType(A,MATBAIJ);CHKERRQ(ierr); 3143359b5eSHong Zhang ierr = MatSetFromOptions(A);CHKERRQ(ierr); 3243359b5eSHong Zhang ierr = MatMPIBAIJSetPreallocation(A,bs,d_nz,NULL,o_nz,NULL);CHKERRQ(ierr); 3343359b5eSHong Zhang ierr = MatSeqBAIJSetPreallocation(A,bs,d_nz,NULL);CHKERRQ(ierr); 34c4762a1bSJed Brown ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 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; 41c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 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; 45c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 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; 49c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 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; 55c4762a1bSJed Brown if (i>0) {J = Ii - n1; ierr = MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr);} 56c4762a1bSJed Brown if (i<n1-1) {J = Ii + n1; ierr = MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr);} 57c4762a1bSJed Brown if (j>0) {J = Ii - 1; ierr = MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr);} 58c4762a1bSJed Brown if (j<n1-1) {J = Ii + 1; ierr = MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr);} 59c4762a1bSJed Brown ierr = MatSetValues(A,1,&Ii,1,&Ii,&four,INSERT_VALUES);CHKERRQ(ierr); 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; 70c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 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; 74c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); 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; 78c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); 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; 84c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,1,col,value,INSERT_VALUES);CHKERRQ(ierr); 85c4762a1bSJed Brown col[0]=i; row=i+bs; 86c4762a1bSJed Brown ierr = MatSetValues(A,1,&row,1,col,value,INSERT_VALUES);CHKERRQ(ierr); 87c4762a1bSJed Brown } 88c4762a1bSJed Brown } 89c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 90c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9143359b5eSHong Zhang ierr = MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr); 9243359b5eSHong Zhang 9343359b5eSHong Zhang /* Get SBAIJ matrix sA from A */ 9443359b5eSHong Zhang ierr = MatConvert(A,MATSBAIJ,MAT_INITIAL_MATRIX,&sA);CHKERRQ(ierr); 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* Test MatGetSize(), MatGetLocalSize() */ 97c4762a1bSJed Brown ierr = MatGetSize(sA, &i,&j);CHKERRQ(ierr); 98c4762a1bSJed Brown ierr = MatGetSize(A, &i2,&j2);CHKERRQ(ierr); 99c4762a1bSJed Brown i -= i2; j -= j2; 100c4762a1bSJed Brown if (i || j) { 101c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatGetSize()\n",rank);CHKERRQ(ierr); 102c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr); 103c4762a1bSJed Brown } 104c4762a1bSJed Brown 105c4762a1bSJed Brown ierr = MatGetLocalSize(sA, &i,&j);CHKERRQ(ierr); 106c4762a1bSJed Brown ierr = MatGetLocalSize(A, &i2,&j2);CHKERRQ(ierr); 107c4762a1bSJed Brown i2 -= i; j2 -= j; 108c4762a1bSJed Brown if (i2 || j2) { 109c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatGetLocalSize()\n",rank);CHKERRQ(ierr); 110c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr); 111c4762a1bSJed Brown } 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* vectors */ 114c4762a1bSJed Brown /*--------------------*/ 115c4762a1bSJed Brown /* i is obtained from MatGetLocalSize() */ 116c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr); 117c4762a1bSJed Brown ierr = VecSetSizes(x,i,PETSC_DECIDE);CHKERRQ(ierr); 118c4762a1bSJed Brown ierr = VecSetFromOptions(x);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = VecDuplicate(x,&y);CHKERRQ(ierr); 120c4762a1bSJed Brown ierr = VecDuplicate(x,&u);CHKERRQ(ierr); 121c4762a1bSJed Brown ierr = VecDuplicate(x,&s1);CHKERRQ(ierr); 122c4762a1bSJed Brown ierr = VecDuplicate(x,&s2);CHKERRQ(ierr); 123c4762a1bSJed Brown 124c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rctx);CHKERRQ(ierr); 125c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 126c4762a1bSJed Brown ierr = VecSetRandom(x,rctx);CHKERRQ(ierr); 127c4762a1bSJed Brown ierr = VecSet(u,one);CHKERRQ(ierr); 128c4762a1bSJed Brown 129c4762a1bSJed Brown /* Test MatNorm() */ 130c4762a1bSJed Brown ierr = MatNorm(A,NORM_FROBENIUS,&r1);CHKERRQ(ierr); 131c4762a1bSJed Brown ierr = MatNorm(sA,NORM_FROBENIUS,&r2);CHKERRQ(ierr); 132c4762a1bSJed Brown rnorm = PetscAbsReal(r1-r2)/r2; 133dd400576SPatrick Sanan if (rnorm > tol && rank == 0) { 1348fffc762SJacob Faibussowitsch ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_FROBENIUS(), Anorm=%16.14e, sAnorm=%16.14e bs=%" PetscInt_FMT "\n",(double)r1,(double)r2,bs);CHKERRQ(ierr); 135c4762a1bSJed Brown } 136c4762a1bSJed Brown ierr = MatNorm(A,NORM_INFINITY,&r1);CHKERRQ(ierr); 137c4762a1bSJed Brown ierr = MatNorm(sA,NORM_INFINITY,&r2);CHKERRQ(ierr); 138c4762a1bSJed Brown rnorm = PetscAbsReal(r1-r2)/r2; 139dd400576SPatrick Sanan if (rnorm > tol && rank == 0) { 1408fffc762SJacob Faibussowitsch ierr = PetscPrintf(PETSC_COMM_WORLD,"Error: MatNorm_INFINITY(), Anorm=%16.14e, sAnorm=%16.14e bs=%" PetscInt_FMT "\n",(double)r1,(double)r2,bs);CHKERRQ(ierr); 141c4762a1bSJed Brown } 142c4762a1bSJed Brown ierr = MatNorm(A,NORM_1,&r1);CHKERRQ(ierr); 143c4762a1bSJed Brown ierr = MatNorm(sA,NORM_1,&r2);CHKERRQ(ierr); 144c4762a1bSJed Brown rnorm = PetscAbsReal(r1-r2)/r2; 145dd400576SPatrick Sanan if (rnorm > tol && rank == 0) { 1468fffc762SJacob Faibussowitsch ierr = PetscPrintf(PETSC_COMM_WORLD,"Error: MatNorm_1(), Anorm=%16.14e, sAnorm=%16.14e bs=%" PetscInt_FMT "\n",(double)r1,(double)r2,bs);CHKERRQ(ierr); 147c4762a1bSJed Brown } 148c4762a1bSJed Brown 149c4762a1bSJed Brown /* Test MatGetOwnershipRange() */ 150c4762a1bSJed Brown ierr = MatGetOwnershipRange(sA,&rstart,&rend);CHKERRQ(ierr); 151c4762a1bSJed Brown ierr = MatGetOwnershipRange(A,&i2,&j2);CHKERRQ(ierr); 152c4762a1bSJed Brown i2 -= rstart; j2 -= rend; 153c4762a1bSJed Brown if (i2 || j2) { 154c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MaGetOwnershipRange()\n",rank);CHKERRQ(ierr); 155c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr); 156c4762a1bSJed Brown } 157c4762a1bSJed Brown 158c4762a1bSJed Brown /* Test MatDiagonalScale() */ 159c4762a1bSJed Brown ierr = MatDiagonalScale(A,x,x);CHKERRQ(ierr); 160c4762a1bSJed Brown ierr = MatDiagonalScale(sA,x,x);CHKERRQ(ierr); 161c4762a1bSJed Brown ierr = MatMultEqual(A,sA,10,&flg);CHKERRQ(ierr); 162*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatDiagonalScale"); 163c4762a1bSJed Brown 164c4762a1bSJed Brown /* Test MatGetDiagonal(), MatScale() */ 165c4762a1bSJed Brown ierr = MatGetDiagonal(A,s1);CHKERRQ(ierr); 166c4762a1bSJed Brown ierr = MatGetDiagonal(sA,s2);CHKERRQ(ierr); 167c4762a1bSJed Brown ierr = VecNorm(s1,NORM_1,&r1);CHKERRQ(ierr); 168c4762a1bSJed Brown ierr = VecNorm(s2,NORM_1,&r2);CHKERRQ(ierr); 169c4762a1bSJed Brown r1 -= r2; 170c4762a1bSJed Brown if (r1<-tol || r1>tol) { 171c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatDiagonalScale() or MatGetDiagonal(), r1=%g \n",rank,(double)r1);CHKERRQ(ierr); 172c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr); 173c4762a1bSJed Brown } 174c4762a1bSJed Brown 175c4762a1bSJed Brown ierr = MatScale(A,alpha);CHKERRQ(ierr); 176c4762a1bSJed Brown ierr = MatScale(sA,alpha);CHKERRQ(ierr); 177c4762a1bSJed Brown 178c4762a1bSJed Brown /* Test MatGetRowMaxAbs() */ 179c4762a1bSJed Brown ierr = MatGetRowMaxAbs(A,s1,NULL);CHKERRQ(ierr); 180c4762a1bSJed Brown ierr = MatGetRowMaxAbs(sA,s2,NULL);CHKERRQ(ierr); 181c4762a1bSJed Brown 182c4762a1bSJed Brown ierr = VecNorm(s1,NORM_1,&r1);CHKERRQ(ierr); 183c4762a1bSJed Brown ierr = VecNorm(s2,NORM_1,&r2);CHKERRQ(ierr); 184c4762a1bSJed Brown r1 -= r2; 185c4762a1bSJed Brown if (r1<-tol || r1>tol) { 186c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatGetRowMaxAbs() \n");CHKERRQ(ierr); 187c4762a1bSJed Brown } 188c4762a1bSJed Brown 189c4762a1bSJed Brown /* Test MatMult(), MatMultAdd() */ 190c4762a1bSJed Brown ierr = MatMultEqual(A,sA,10,&flg);CHKERRQ(ierr); 191c4762a1bSJed Brown if (!flg) { 192c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatMult() or MatScale()\n",rank);CHKERRQ(ierr); 193c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr); 194c4762a1bSJed Brown } 195c4762a1bSJed Brown 196c4762a1bSJed Brown ierr = MatMultAddEqual(A,sA,10,&flg);CHKERRQ(ierr); 197c4762a1bSJed Brown if (!flg) { 198c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatMultAdd()\n",rank);CHKERRQ(ierr); 199c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr); 200c4762a1bSJed Brown } 201c4762a1bSJed Brown 202c4762a1bSJed Brown /* Test MatMultTranspose(), MatMultTransposeAdd() */ 203c4762a1bSJed Brown for (i=0; i<10; i++) { 204c4762a1bSJed Brown ierr = VecSetRandom(x,rctx);CHKERRQ(ierr); 205c4762a1bSJed Brown ierr = MatMultTranspose(A,x,s1);CHKERRQ(ierr); 206c4762a1bSJed Brown ierr = MatMultTranspose(sA,x,s2);CHKERRQ(ierr); 207c4762a1bSJed Brown ierr = VecNorm(s1,NORM_1,&r1);CHKERRQ(ierr); 208c4762a1bSJed Brown ierr = VecNorm(s2,NORM_1,&r2);CHKERRQ(ierr); 209c4762a1bSJed Brown r1 -= r2; 210c4762a1bSJed Brown if (r1<-tol || r1>tol) { 211c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatMult() or MatScale(), err=%g\n",rank,(double)r1);CHKERRQ(ierr); 212c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr); 213c4762a1bSJed Brown } 214c4762a1bSJed Brown } 215c4762a1bSJed Brown for (i=0; i<10; i++) { 216c4762a1bSJed Brown ierr = VecSetRandom(x,rctx);CHKERRQ(ierr); 217c4762a1bSJed Brown ierr = VecSetRandom(y,rctx);CHKERRQ(ierr); 218c4762a1bSJed Brown ierr = MatMultTransposeAdd(A,x,y,s1);CHKERRQ(ierr); 219c4762a1bSJed Brown ierr = MatMultTransposeAdd(sA,x,y,s2);CHKERRQ(ierr); 220c4762a1bSJed Brown ierr = VecNorm(s1,NORM_1,&r1);CHKERRQ(ierr); 221c4762a1bSJed Brown ierr = VecNorm(s2,NORM_1,&r2);CHKERRQ(ierr); 222c4762a1bSJed Brown r1 -= r2; 223c4762a1bSJed Brown if (r1<-tol || r1>tol) { 224c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatMultAdd(), err=%g \n",rank,(double)r1);CHKERRQ(ierr); 225c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr); 226c4762a1bSJed Brown } 227c4762a1bSJed Brown } 228c4762a1bSJed Brown 229c4762a1bSJed Brown /* Test MatDuplicate() */ 230c4762a1bSJed Brown ierr = MatDuplicate(sA,MAT_COPY_VALUES,&sB);CHKERRQ(ierr); 231c4762a1bSJed Brown ierr = MatEqual(sA,sB,&flg);CHKERRQ(ierr); 232c4762a1bSJed Brown if (!flg) { 233c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," Error in MatDuplicate(), sA != sB \n");CHKERRQ(ierr); 234c4762a1bSJed Brown } 235c4762a1bSJed Brown ierr = MatMultEqual(sA,sB,5,&flg);CHKERRQ(ierr); 236c4762a1bSJed Brown if (!flg) { 237c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatDuplicate() or MatMult()\n",rank);CHKERRQ(ierr); 238c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr); 239c4762a1bSJed Brown } 240c4762a1bSJed Brown ierr = MatMultAddEqual(sA,sB,5,&flg);CHKERRQ(ierr); 241c4762a1bSJed Brown if (!flg) { 242c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatDuplicate() or MatMultAdd(()\n",rank);CHKERRQ(ierr); 243c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr); 244c4762a1bSJed Brown } 245c4762a1bSJed Brown ierr = MatDestroy(&sB);CHKERRQ(ierr); 246c4762a1bSJed Brown ierr = VecDestroy(&u);CHKERRQ(ierr); 247c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 248c4762a1bSJed Brown ierr = VecDestroy(&y);CHKERRQ(ierr); 249c4762a1bSJed Brown ierr = VecDestroy(&s1);CHKERRQ(ierr); 250c4762a1bSJed Brown ierr = VecDestroy(&s2);CHKERRQ(ierr); 251c4762a1bSJed Brown ierr = MatDestroy(&sA);CHKERRQ(ierr); 252c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 253c4762a1bSJed Brown ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 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