1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests various routines in MatMPISBAIJ format.\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 69371c9d4SSatish Balay int main(int argc, char **args) { 7c4762a1bSJed Brown Vec x, y, u, s1, s2; 8c4762a1bSJed Brown Mat A, sA, sB; 9c4762a1bSJed Brown PetscRandom rctx; 10c4762a1bSJed Brown PetscReal r1, r2, rnorm, tol = PETSC_SQRT_MACHINE_EPSILON; 11c4762a1bSJed Brown PetscScalar one = 1.0, neg_one = -1.0, value[3], four = 4.0, alpha = 0.1; 1243359b5eSHong 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; 13c4762a1bSJed Brown PetscMPIInt size, rank; 14c4762a1bSJed Brown PetscBool flg; 15c4762a1bSJed Brown 16327415f7SBarry Smith PetscFunctionBeginUser; 179566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 189566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mbs", &mbs, NULL)); 199566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-bs", &bs, NULL)); 209566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-prob", &prob, NULL)); 21c4762a1bSJed Brown 229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 24c4762a1bSJed Brown 2543359b5eSHong Zhang /* Create a BAIJ matrix A */ 26c4762a1bSJed Brown n = mbs * bs; 279566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 289566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, n, n)); 299566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATBAIJ)); 309566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 319566063dSJacob Faibussowitsch PetscCall(MatMPIBAIJSetPreallocation(A, bs, d_nz, NULL, o_nz, NULL)); 329566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(A, bs, d_nz, NULL)); 339566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 34c4762a1bSJed Brown 35c4762a1bSJed Brown if (bs == 1) { 36c4762a1bSJed Brown if (prob == 1) { /* tridiagonal matrix */ 379371c9d4SSatish Balay value[0] = -1.0; 389371c9d4SSatish Balay value[1] = 0.0; 399371c9d4SSatish Balay value[2] = -1.0; 40c4762a1bSJed Brown for (i = 1; i < n - 1; i++) { 419371c9d4SSatish Balay col[0] = i - 1; 429371c9d4SSatish Balay col[1] = i; 439371c9d4SSatish Balay col[2] = i + 1; 449566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES)); 45c4762a1bSJed Brown } 469371c9d4SSatish Balay i = n - 1; 479371c9d4SSatish Balay col[0] = 0; 489371c9d4SSatish Balay col[1] = n - 2; 499371c9d4SSatish Balay col[2] = n - 1; 509371c9d4SSatish Balay value[0] = 0.1; 519371c9d4SSatish Balay value[1] = -1.0; 529371c9d4SSatish Balay value[2] = 0.0; 539566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES)); 54c4762a1bSJed Brown 559371c9d4SSatish Balay i = 0; 569371c9d4SSatish Balay col[0] = 0; 579371c9d4SSatish Balay col[1] = 1; 589371c9d4SSatish Balay col[2] = n - 1; 599371c9d4SSatish Balay value[0] = 0.0; 609371c9d4SSatish Balay value[1] = -1.0; 619371c9d4SSatish Balay value[2] = 0.1; 629566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES)); 63c4762a1bSJed Brown } else if (prob == 2) { /* matrix for the five point stencil */ 64c4762a1bSJed Brown n1 = (int)PetscSqrtReal((PetscReal)n); 65c4762a1bSJed Brown for (i = 0; i < n1; i++) { 66c4762a1bSJed Brown for (j = 0; j < n1; j++) { 67c4762a1bSJed Brown Ii = j + n1 * i; 689371c9d4SSatish Balay if (i > 0) { 699371c9d4SSatish Balay J = Ii - n1; 709371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES)); 719371c9d4SSatish Balay } 729371c9d4SSatish Balay if (i < n1 - 1) { 739371c9d4SSatish Balay J = Ii + n1; 749371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES)); 759371c9d4SSatish Balay } 769371c9d4SSatish Balay if (j > 0) { 779371c9d4SSatish Balay J = Ii - 1; 789371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES)); 799371c9d4SSatish Balay } 809371c9d4SSatish Balay if (j < n1 - 1) { 819371c9d4SSatish Balay J = Ii + 1; 829371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES)); 839371c9d4SSatish Balay } 849566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &Ii, 1, &Ii, &four, INSERT_VALUES)); 85c4762a1bSJed Brown } 86c4762a1bSJed Brown } 87c4762a1bSJed Brown } 88c4762a1bSJed Brown /* end of if (bs == 1) */ 89c4762a1bSJed Brown } else { /* bs > 1 */ 90c4762a1bSJed Brown for (block = 0; block < n / bs; block++) { 91c4762a1bSJed Brown /* diagonal blocks */ 929371c9d4SSatish Balay value[0] = -1.0; 939371c9d4SSatish Balay value[1] = 4.0; 949371c9d4SSatish Balay value[2] = -1.0; 95c4762a1bSJed Brown for (i = 1 + block * bs; i < bs - 1 + block * bs; i++) { 969371c9d4SSatish Balay col[0] = i - 1; 979371c9d4SSatish Balay col[1] = i; 989371c9d4SSatish Balay col[2] = i + 1; 999566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES)); 100c4762a1bSJed Brown } 1019371c9d4SSatish Balay i = bs - 1 + block * bs; 1029371c9d4SSatish Balay col[0] = bs - 2 + block * bs; 1039371c9d4SSatish Balay col[1] = bs - 1 + block * bs; 1049371c9d4SSatish Balay value[0] = -1.0; 1059371c9d4SSatish Balay value[1] = 4.0; 1069566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES)); 107c4762a1bSJed Brown 1089371c9d4SSatish Balay i = 0 + block * bs; 1099371c9d4SSatish Balay col[0] = 0 + block * bs; 1109371c9d4SSatish Balay col[1] = 1 + block * bs; 1119371c9d4SSatish Balay value[0] = 4.0; 1129371c9d4SSatish Balay value[1] = -1.0; 1139566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES)); 114c4762a1bSJed Brown } 115c4762a1bSJed Brown /* off-diagonal blocks */ 116c4762a1bSJed Brown value[0] = -1.0; 117c4762a1bSJed Brown for (i = 0; i < (n / bs - 1) * bs; i++) { 118c4762a1bSJed Brown col[0] = i + bs; 1199566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 1, col, value, INSERT_VALUES)); 1209371c9d4SSatish Balay col[0] = i; 1219371c9d4SSatish Balay row = i + bs; 1229566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, col, value, INSERT_VALUES)); 123c4762a1bSJed Brown } 124c4762a1bSJed Brown } 1259566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1269566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 1279566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE)); 12843359b5eSHong Zhang 12943359b5eSHong Zhang /* Get SBAIJ matrix sA from A */ 1309566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATSBAIJ, MAT_INITIAL_MATRIX, &sA)); 131c4762a1bSJed Brown 132c4762a1bSJed Brown /* Test MatGetSize(), MatGetLocalSize() */ 1339566063dSJacob Faibussowitsch PetscCall(MatGetSize(sA, &i, &j)); 1349566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &i2, &j2)); 1359371c9d4SSatish Balay i -= i2; 1369371c9d4SSatish Balay j -= j2; 137c4762a1bSJed Brown if (i || j) { 1389566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d], Error: MatGetSize()\n", rank)); 1399566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 140c4762a1bSJed Brown } 141c4762a1bSJed Brown 1429566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(sA, &i, &j)); 1439566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &i2, &j2)); 1449371c9d4SSatish Balay i2 -= i; 1459371c9d4SSatish Balay j2 -= j; 146c4762a1bSJed Brown if (i2 || j2) { 1479566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d], Error: MatGetLocalSize()\n", rank)); 1489566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 149c4762a1bSJed Brown } 150c4762a1bSJed Brown 151c4762a1bSJed Brown /* vectors */ 152c4762a1bSJed Brown /*--------------------*/ 153c4762a1bSJed Brown /* i is obtained from MatGetLocalSize() */ 1549566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &x)); 1559566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x, i, PETSC_DECIDE)); 1569566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 1579566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &y)); 1589566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &u)); 1599566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &s1)); 1609566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &s2)); 161c4762a1bSJed Brown 1629566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rctx)); 1639566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rctx)); 1649566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rctx)); 1659566063dSJacob Faibussowitsch PetscCall(VecSet(u, one)); 166c4762a1bSJed Brown 167c4762a1bSJed Brown /* Test MatNorm() */ 1689566063dSJacob Faibussowitsch PetscCall(MatNorm(A, NORM_FROBENIUS, &r1)); 1699566063dSJacob Faibussowitsch PetscCall(MatNorm(sA, NORM_FROBENIUS, &r2)); 170c4762a1bSJed Brown rnorm = PetscAbsReal(r1 - r2) / r2; 171*48a46eb9SPierre 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)); 1729566063dSJacob Faibussowitsch PetscCall(MatNorm(A, NORM_INFINITY, &r1)); 1739566063dSJacob Faibussowitsch PetscCall(MatNorm(sA, NORM_INFINITY, &r2)); 174c4762a1bSJed Brown rnorm = PetscAbsReal(r1 - r2) / r2; 175*48a46eb9SPierre 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)); 1769566063dSJacob Faibussowitsch PetscCall(MatNorm(A, NORM_1, &r1)); 1779566063dSJacob Faibussowitsch PetscCall(MatNorm(sA, NORM_1, &r2)); 178c4762a1bSJed Brown rnorm = PetscAbsReal(r1 - r2) / r2; 179*48a46eb9SPierre 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)); 180c4762a1bSJed Brown 181c4762a1bSJed Brown /* Test MatGetOwnershipRange() */ 1829566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(sA, &rstart, &rend)); 1839566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &i2, &j2)); 1849371c9d4SSatish Balay i2 -= rstart; 1859371c9d4SSatish Balay j2 -= rend; 186c4762a1bSJed Brown if (i2 || j2) { 1879566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d], Error: MaGetOwnershipRange()\n", rank)); 1889566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 189c4762a1bSJed Brown } 190c4762a1bSJed Brown 191c4762a1bSJed Brown /* Test MatDiagonalScale() */ 1929566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A, x, x)); 1939566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(sA, x, x)); 1949566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A, sA, 10, &flg)); 19528b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMETYPE, "Error in MatDiagonalScale"); 196c4762a1bSJed Brown 197c4762a1bSJed Brown /* Test MatGetDiagonal(), MatScale() */ 1989566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(A, s1)); 1999566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(sA, s2)); 2009566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_1, &r1)); 2019566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_1, &r2)); 202c4762a1bSJed Brown r1 -= r2; 203c4762a1bSJed Brown if (r1 < -tol || r1 > tol) { 2049566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d], Error: MatDiagonalScale() or MatGetDiagonal(), r1=%g \n", rank, (double)r1)); 2059566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 206c4762a1bSJed Brown } 207c4762a1bSJed Brown 2089566063dSJacob Faibussowitsch PetscCall(MatScale(A, alpha)); 2099566063dSJacob Faibussowitsch PetscCall(MatScale(sA, alpha)); 210c4762a1bSJed Brown 211c4762a1bSJed Brown /* Test MatGetRowMaxAbs() */ 2129566063dSJacob Faibussowitsch PetscCall(MatGetRowMaxAbs(A, s1, NULL)); 2139566063dSJacob Faibussowitsch PetscCall(MatGetRowMaxAbs(sA, s2, NULL)); 214c4762a1bSJed Brown 2159566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_1, &r1)); 2169566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_1, &r2)); 217c4762a1bSJed Brown r1 -= r2; 218*48a46eb9SPierre Jolivet if (r1 < -tol || r1 > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatGetRowMaxAbs() \n")); 219c4762a1bSJed Brown 220c4762a1bSJed Brown /* Test MatMult(), MatMultAdd() */ 2219566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A, sA, 10, &flg)); 222c4762a1bSJed Brown if (!flg) { 2239566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d], Error: MatMult() or MatScale()\n", rank)); 2249566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 225c4762a1bSJed Brown } 226c4762a1bSJed Brown 2279566063dSJacob Faibussowitsch PetscCall(MatMultAddEqual(A, sA, 10, &flg)); 228c4762a1bSJed Brown if (!flg) { 2299566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d], Error: MatMultAdd()\n", rank)); 2309566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 231c4762a1bSJed Brown } 232c4762a1bSJed Brown 233c4762a1bSJed Brown /* Test MatMultTranspose(), MatMultTransposeAdd() */ 234c4762a1bSJed Brown for (i = 0; i < 10; i++) { 2359566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rctx)); 2369566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, x, s1)); 2379566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(sA, x, s2)); 2389566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_1, &r1)); 2399566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_1, &r2)); 240c4762a1bSJed Brown r1 -= r2; 241c4762a1bSJed Brown if (r1 < -tol || r1 > tol) { 2429566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d], Error: MatMult() or MatScale(), err=%g\n", rank, (double)r1)); 2439566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 244c4762a1bSJed Brown } 245c4762a1bSJed Brown } 246c4762a1bSJed Brown for (i = 0; i < 10; i++) { 2479566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rctx)); 2489566063dSJacob Faibussowitsch PetscCall(VecSetRandom(y, rctx)); 2499566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(A, x, y, s1)); 2509566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(sA, x, y, s2)); 2519566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_1, &r1)); 2529566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_1, &r2)); 253c4762a1bSJed Brown r1 -= r2; 254c4762a1bSJed Brown if (r1 < -tol || r1 > tol) { 2559566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d], Error: MatMultAdd(), err=%g \n", rank, (double)r1)); 2569566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 257c4762a1bSJed Brown } 258c4762a1bSJed Brown } 259c4762a1bSJed Brown 260c4762a1bSJed Brown /* Test MatDuplicate() */ 2619566063dSJacob Faibussowitsch PetscCall(MatDuplicate(sA, MAT_COPY_VALUES, &sB)); 2629566063dSJacob Faibussowitsch PetscCall(MatEqual(sA, sB, &flg)); 263*48a46eb9SPierre Jolivet if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Error in MatDuplicate(), sA != sB \n")); 2649566063dSJacob Faibussowitsch PetscCall(MatMultEqual(sA, sB, 5, &flg)); 265c4762a1bSJed Brown if (!flg) { 2669566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d], Error: MatDuplicate() or MatMult()\n", rank)); 2679566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 268c4762a1bSJed Brown } 2699566063dSJacob Faibussowitsch PetscCall(MatMultAddEqual(sA, sB, 5, &flg)); 270c4762a1bSJed Brown if (!flg) { 2719566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d], Error: MatDuplicate() or MatMultAdd(()\n", rank)); 2729566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 273c4762a1bSJed Brown } 2749566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sB)); 2759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 2769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 2789566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 2799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 2809566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sA)); 2819566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2829566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rctx)); 2839566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 284b122ec5aSJacob Faibussowitsch return 0; 285c4762a1bSJed Brown } 286c4762a1bSJed Brown 287c4762a1bSJed Brown /*TEST 288c4762a1bSJed Brown 289c4762a1bSJed Brown test: 290c4762a1bSJed Brown nsize: {{1 3}} 29143359b5eSHong Zhang args: -bs {{1 2 3 5 7 8}} -mat_ignore_lower_triangular -prob {{1 2}} 292c4762a1bSJed Brown 293c4762a1bSJed Brown TEST*/ 294