1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests various routines in MatMPISBAIJ format.\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 7*d71ae5a4SJacob Faibussowitsch { 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 PetscMPIInt size, rank; 15c4762a1bSJed Brown PetscBool flg; 16c4762a1bSJed Brown 17327415f7SBarry Smith PetscFunctionBeginUser; 189566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 199566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mbs", &mbs, NULL)); 209566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-bs", &bs, NULL)); 219566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-prob", &prob, NULL)); 22c4762a1bSJed Brown 239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 25c4762a1bSJed Brown 2643359b5eSHong Zhang /* Create a BAIJ matrix A */ 27c4762a1bSJed Brown n = mbs * bs; 289566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 299566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, n, n)); 309566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATBAIJ)); 319566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 329566063dSJacob Faibussowitsch PetscCall(MatMPIBAIJSetPreallocation(A, bs, d_nz, NULL, o_nz, NULL)); 339566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(A, bs, d_nz, NULL)); 349566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 35c4762a1bSJed Brown 36c4762a1bSJed Brown if (bs == 1) { 37c4762a1bSJed Brown if (prob == 1) { /* tridiagonal matrix */ 389371c9d4SSatish Balay value[0] = -1.0; 399371c9d4SSatish Balay value[1] = 0.0; 409371c9d4SSatish Balay value[2] = -1.0; 41c4762a1bSJed Brown for (i = 1; i < n - 1; i++) { 429371c9d4SSatish Balay col[0] = i - 1; 439371c9d4SSatish Balay col[1] = i; 449371c9d4SSatish Balay col[2] = i + 1; 459566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES)); 46c4762a1bSJed Brown } 479371c9d4SSatish Balay i = n - 1; 489371c9d4SSatish Balay col[0] = 0; 499371c9d4SSatish Balay col[1] = n - 2; 509371c9d4SSatish Balay col[2] = n - 1; 519371c9d4SSatish Balay value[0] = 0.1; 529371c9d4SSatish Balay value[1] = -1.0; 539371c9d4SSatish Balay value[2] = 0.0; 549566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES)); 55c4762a1bSJed Brown 569371c9d4SSatish Balay i = 0; 579371c9d4SSatish Balay col[0] = 0; 589371c9d4SSatish Balay col[1] = 1; 599371c9d4SSatish Balay col[2] = n - 1; 609371c9d4SSatish Balay value[0] = 0.0; 619371c9d4SSatish Balay value[1] = -1.0; 629371c9d4SSatish Balay value[2] = 0.1; 639566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES)); 64c4762a1bSJed Brown } else if (prob == 2) { /* matrix for the five point stencil */ 65c4762a1bSJed Brown n1 = (int)PetscSqrtReal((PetscReal)n); 66c4762a1bSJed Brown for (i = 0; i < n1; i++) { 67c4762a1bSJed Brown for (j = 0; j < n1; j++) { 68c4762a1bSJed Brown Ii = j + n1 * i; 699371c9d4SSatish Balay if (i > 0) { 709371c9d4SSatish Balay J = Ii - n1; 719371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES)); 729371c9d4SSatish Balay } 739371c9d4SSatish Balay if (i < n1 - 1) { 749371c9d4SSatish Balay J = Ii + n1; 759371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES)); 769371c9d4SSatish Balay } 779371c9d4SSatish Balay if (j > 0) { 789371c9d4SSatish Balay J = Ii - 1; 799371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES)); 809371c9d4SSatish Balay } 819371c9d4SSatish Balay if (j < n1 - 1) { 829371c9d4SSatish Balay J = Ii + 1; 839371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES)); 849371c9d4SSatish Balay } 859566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &Ii, 1, &Ii, &four, INSERT_VALUES)); 86c4762a1bSJed Brown } 87c4762a1bSJed Brown } 88c4762a1bSJed Brown } 89c4762a1bSJed Brown /* end of if (bs == 1) */ 90c4762a1bSJed Brown } else { /* bs > 1 */ 91c4762a1bSJed Brown for (block = 0; block < n / bs; block++) { 92c4762a1bSJed Brown /* diagonal blocks */ 939371c9d4SSatish Balay value[0] = -1.0; 949371c9d4SSatish Balay value[1] = 4.0; 959371c9d4SSatish Balay value[2] = -1.0; 96c4762a1bSJed Brown for (i = 1 + block * bs; i < bs - 1 + block * bs; i++) { 979371c9d4SSatish Balay col[0] = i - 1; 989371c9d4SSatish Balay col[1] = i; 999371c9d4SSatish Balay col[2] = i + 1; 1009566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES)); 101c4762a1bSJed Brown } 1029371c9d4SSatish Balay i = bs - 1 + block * bs; 1039371c9d4SSatish Balay col[0] = bs - 2 + block * bs; 1049371c9d4SSatish Balay col[1] = bs - 1 + block * bs; 1059371c9d4SSatish Balay value[0] = -1.0; 1069371c9d4SSatish Balay value[1] = 4.0; 1079566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES)); 108c4762a1bSJed Brown 1099371c9d4SSatish Balay i = 0 + block * bs; 1109371c9d4SSatish Balay col[0] = 0 + block * bs; 1119371c9d4SSatish Balay col[1] = 1 + block * bs; 1129371c9d4SSatish Balay value[0] = 4.0; 1139371c9d4SSatish Balay value[1] = -1.0; 1149566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES)); 115c4762a1bSJed Brown } 116c4762a1bSJed Brown /* off-diagonal blocks */ 117c4762a1bSJed Brown value[0] = -1.0; 118c4762a1bSJed Brown for (i = 0; i < (n / bs - 1) * bs; i++) { 119c4762a1bSJed Brown col[0] = i + bs; 1209566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 1, col, value, INSERT_VALUES)); 1219371c9d4SSatish Balay col[0] = i; 1229371c9d4SSatish Balay row = i + bs; 1239566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, col, value, INSERT_VALUES)); 124c4762a1bSJed Brown } 125c4762a1bSJed Brown } 1269566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1279566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 1289566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE)); 12943359b5eSHong Zhang 13043359b5eSHong Zhang /* Get SBAIJ matrix sA from A */ 1319566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATSBAIJ, MAT_INITIAL_MATRIX, &sA)); 132c4762a1bSJed Brown 133c4762a1bSJed Brown /* Test MatGetSize(), MatGetLocalSize() */ 1349566063dSJacob Faibussowitsch PetscCall(MatGetSize(sA, &i, &j)); 1359566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &i2, &j2)); 1369371c9d4SSatish Balay i -= i2; 1379371c9d4SSatish Balay j -= j2; 138c4762a1bSJed Brown if (i || j) { 1399566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d], Error: MatGetSize()\n", rank)); 1409566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 141c4762a1bSJed Brown } 142c4762a1bSJed Brown 1439566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(sA, &i, &j)); 1449566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &i2, &j2)); 1459371c9d4SSatish Balay i2 -= i; 1469371c9d4SSatish Balay j2 -= j; 147c4762a1bSJed Brown if (i2 || j2) { 1489566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d], Error: MatGetLocalSize()\n", rank)); 1499566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 150c4762a1bSJed Brown } 151c4762a1bSJed Brown 152c4762a1bSJed Brown /* vectors */ 153c4762a1bSJed Brown /*--------------------*/ 154c4762a1bSJed Brown /* i is obtained from MatGetLocalSize() */ 1559566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &x)); 1569566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x, i, PETSC_DECIDE)); 1579566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 1589566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &y)); 1599566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &u)); 1609566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &s1)); 1619566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &s2)); 162c4762a1bSJed Brown 1639566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rctx)); 1649566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rctx)); 1659566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rctx)); 1669566063dSJacob Faibussowitsch PetscCall(VecSet(u, one)); 167c4762a1bSJed Brown 168c4762a1bSJed Brown /* Test MatNorm() */ 1699566063dSJacob Faibussowitsch PetscCall(MatNorm(A, NORM_FROBENIUS, &r1)); 1709566063dSJacob Faibussowitsch PetscCall(MatNorm(sA, NORM_FROBENIUS, &r2)); 171c4762a1bSJed Brown rnorm = PetscAbsReal(r1 - r2) / r2; 17248a46eb9SPierre Jolivet if (rnorm > tol && rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatNorm_FROBENIUS(), Anorm=%16.14e, sAnorm=%16.14e bs=%" PetscInt_FMT "\n", (double)r1, (double)r2, bs)); 1739566063dSJacob Faibussowitsch PetscCall(MatNorm(A, NORM_INFINITY, &r1)); 1749566063dSJacob Faibussowitsch PetscCall(MatNorm(sA, NORM_INFINITY, &r2)); 175c4762a1bSJed Brown rnorm = PetscAbsReal(r1 - r2) / r2; 17648a46eb9SPierre Jolivet if (rnorm > tol && rank == 0) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error: MatNorm_INFINITY(), Anorm=%16.14e, sAnorm=%16.14e bs=%" PetscInt_FMT "\n", (double)r1, (double)r2, bs)); 1779566063dSJacob Faibussowitsch PetscCall(MatNorm(A, NORM_1, &r1)); 1789566063dSJacob Faibussowitsch PetscCall(MatNorm(sA, NORM_1, &r2)); 179c4762a1bSJed Brown rnorm = PetscAbsReal(r1 - r2) / r2; 18048a46eb9SPierre Jolivet if (rnorm > tol && rank == 0) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error: MatNorm_1(), Anorm=%16.14e, sAnorm=%16.14e bs=%" PetscInt_FMT "\n", (double)r1, (double)r2, bs)); 181c4762a1bSJed Brown 182c4762a1bSJed Brown /* Test MatGetOwnershipRange() */ 1839566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(sA, &rstart, &rend)); 1849566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &i2, &j2)); 1859371c9d4SSatish Balay i2 -= rstart; 1869371c9d4SSatish Balay j2 -= rend; 187c4762a1bSJed Brown if (i2 || j2) { 1889566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d], Error: MaGetOwnershipRange()\n", rank)); 1899566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 190c4762a1bSJed Brown } 191c4762a1bSJed Brown 192c4762a1bSJed Brown /* Test MatDiagonalScale() */ 1939566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A, x, x)); 1949566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(sA, x, x)); 1959566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A, sA, 10, &flg)); 19628b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMETYPE, "Error in MatDiagonalScale"); 197c4762a1bSJed Brown 198c4762a1bSJed Brown /* Test MatGetDiagonal(), MatScale() */ 1999566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(A, s1)); 2009566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(sA, s2)); 2019566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_1, &r1)); 2029566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_1, &r2)); 203c4762a1bSJed Brown r1 -= r2; 204c4762a1bSJed Brown if (r1 < -tol || r1 > tol) { 2059566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d], Error: MatDiagonalScale() or MatGetDiagonal(), r1=%g \n", rank, (double)r1)); 2069566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 207c4762a1bSJed Brown } 208c4762a1bSJed Brown 2099566063dSJacob Faibussowitsch PetscCall(MatScale(A, alpha)); 2109566063dSJacob Faibussowitsch PetscCall(MatScale(sA, alpha)); 211c4762a1bSJed Brown 212c4762a1bSJed Brown /* Test MatGetRowMaxAbs() */ 2139566063dSJacob Faibussowitsch PetscCall(MatGetRowMaxAbs(A, s1, NULL)); 2149566063dSJacob Faibussowitsch PetscCall(MatGetRowMaxAbs(sA, s2, NULL)); 215c4762a1bSJed Brown 2169566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_1, &r1)); 2179566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_1, &r2)); 218c4762a1bSJed Brown r1 -= r2; 21948a46eb9SPierre Jolivet if (r1 < -tol || r1 > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatGetRowMaxAbs() \n")); 220c4762a1bSJed Brown 221c4762a1bSJed Brown /* Test MatMult(), MatMultAdd() */ 2229566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A, sA, 10, &flg)); 223c4762a1bSJed Brown if (!flg) { 2249566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d], Error: MatMult() or MatScale()\n", rank)); 2259566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 226c4762a1bSJed Brown } 227c4762a1bSJed Brown 2289566063dSJacob Faibussowitsch PetscCall(MatMultAddEqual(A, sA, 10, &flg)); 229c4762a1bSJed Brown if (!flg) { 2309566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d], Error: MatMultAdd()\n", rank)); 2319566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 232c4762a1bSJed Brown } 233c4762a1bSJed Brown 234c4762a1bSJed Brown /* Test MatMultTranspose(), MatMultTransposeAdd() */ 235c4762a1bSJed Brown for (i = 0; i < 10; i++) { 2369566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rctx)); 2379566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, x, s1)); 2389566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(sA, x, s2)); 2399566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_1, &r1)); 2409566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_1, &r2)); 241c4762a1bSJed Brown r1 -= r2; 242c4762a1bSJed Brown if (r1 < -tol || r1 > tol) { 2439566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d], Error: MatMult() or MatScale(), err=%g\n", rank, (double)r1)); 2449566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 245c4762a1bSJed Brown } 246c4762a1bSJed Brown } 247c4762a1bSJed Brown for (i = 0; i < 10; i++) { 2489566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rctx)); 2499566063dSJacob Faibussowitsch PetscCall(VecSetRandom(y, rctx)); 2509566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(A, x, y, s1)); 2519566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(sA, x, y, s2)); 2529566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_1, &r1)); 2539566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_1, &r2)); 254c4762a1bSJed Brown r1 -= r2; 255c4762a1bSJed Brown if (r1 < -tol || r1 > tol) { 2569566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d], Error: MatMultAdd(), err=%g \n", rank, (double)r1)); 2579566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 258c4762a1bSJed Brown } 259c4762a1bSJed Brown } 260c4762a1bSJed Brown 261c4762a1bSJed Brown /* Test MatDuplicate() */ 2629566063dSJacob Faibussowitsch PetscCall(MatDuplicate(sA, MAT_COPY_VALUES, &sB)); 2639566063dSJacob Faibussowitsch PetscCall(MatEqual(sA, sB, &flg)); 26448a46eb9SPierre Jolivet if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Error in MatDuplicate(), sA != sB \n")); 2659566063dSJacob Faibussowitsch PetscCall(MatMultEqual(sA, sB, 5, &flg)); 266c4762a1bSJed Brown if (!flg) { 2679566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d], Error: MatDuplicate() or MatMult()\n", rank)); 2689566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 269c4762a1bSJed Brown } 2709566063dSJacob Faibussowitsch PetscCall(MatMultAddEqual(sA, sB, 5, &flg)); 271c4762a1bSJed Brown if (!flg) { 2729566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d], Error: MatDuplicate() or MatMultAdd(()\n", rank)); 2739566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 274c4762a1bSJed Brown } 2759566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sB)); 2769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 2779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2789566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 2799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 2809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 2819566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sA)); 2829566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2839566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rctx)); 2849566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 285b122ec5aSJacob Faibussowitsch return 0; 286c4762a1bSJed Brown } 287c4762a1bSJed Brown 288c4762a1bSJed Brown /*TEST 289c4762a1bSJed Brown 290c4762a1bSJed Brown test: 291c4762a1bSJed Brown nsize: {{1 3}} 29243359b5eSHong Zhang args: -bs {{1 2 3 5 7 8}} -mat_ignore_lower_triangular -prob {{1 2}} 293c4762a1bSJed Brown 294c4762a1bSJed Brown TEST*/ 295