1c4762a1bSJed Brown static char help[] = "Tests MatGetRowMax(), MatGetRowMin(), MatGetRowMaxAbs()\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown #define M 5 6c4762a1bSJed Brown #define N 6 7c4762a1bSJed Brown 8d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 9d71ae5a4SJacob Faibussowitsch { 10c4762a1bSJed Brown Mat A; 11fa213d2fSHong Zhang Vec min, max, maxabs, e; 124e879edeSHong Zhang PetscInt m, n, j, imin[M], imax[M], imaxabs[M], indices[N], row, testcase = 0; 13c4762a1bSJed Brown PetscScalar values[N]; 141a254869SHong Zhang PetscMPIInt size, rank; 15fa213d2fSHong Zhang PetscReal enorm; 16c4762a1bSJed Brown 17327415f7SBarry Smith PetscFunctionBeginUser; 18*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help)); 199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 219566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-testcase", &testcase, NULL)); 22c4762a1bSJed Brown 239566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 241a254869SHong Zhang if (testcase == 1) { /* proc[0] holds entire A and other processes have no entry */ 25dd400576SPatrick Sanan if (rank == 0) { 269566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, M, N, PETSC_DECIDE, PETSC_DECIDE)); 271a254869SHong Zhang } else { 289566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, 0, 0, PETSC_DECIDE, PETSC_DECIDE)); 291a254869SHong Zhang } 301a254869SHong Zhang testcase = 0; 311a254869SHong Zhang } else { 329566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, M, N)); 331a254869SHong Zhang } 349566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 359566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 361a254869SHong Zhang 37dd400576SPatrick Sanan if (rank == 0) { /* proc[0] sets matrix A */ 381a254869SHong Zhang for (j = 0; j < N; j++) indices[j] = j; 391a254869SHong Zhang switch (testcase) { 40d71ae5a4SJacob Faibussowitsch case 1: /* see testcast 0 */ 41d71ae5a4SJacob Faibussowitsch break; 421a254869SHong Zhang case 2: 43c4762a1bSJed Brown row = 0; 449371c9d4SSatish Balay values[0] = -2.0; 459371c9d4SSatish Balay values[1] = -2.0; 469371c9d4SSatish Balay values[2] = -2.0; 479371c9d4SSatish Balay values[3] = -4.0; 489371c9d4SSatish Balay values[4] = 1.0; 499371c9d4SSatish Balay values[5] = 1.0; 509566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, N, indices, values, INSERT_VALUES)); 511a254869SHong Zhang row = 2; 529371c9d4SSatish Balay indices[0] = 0; 539371c9d4SSatish Balay indices[1] = 3; 549371c9d4SSatish Balay indices[2] = 5; 559371c9d4SSatish Balay values[0] = -2.0; 569371c9d4SSatish Balay values[1] = -2.0; 579371c9d4SSatish Balay values[2] = -2.0; 589566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 3, indices, values, INSERT_VALUES)); 591a254869SHong Zhang row = 3; 609371c9d4SSatish Balay indices[0] = 0; 619371c9d4SSatish Balay indices[1] = 1; 629371c9d4SSatish Balay indices[2] = 4; 639371c9d4SSatish Balay values[0] = -2.0; 649371c9d4SSatish Balay values[1] = -2.0; 659371c9d4SSatish Balay values[2] = -2.0; 669566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 3, indices, values, INSERT_VALUES)); 67c4762a1bSJed Brown row = 4; 689371c9d4SSatish Balay indices[0] = 0; 699371c9d4SSatish Balay indices[1] = 1; 709371c9d4SSatish Balay indices[2] = 2; 719371c9d4SSatish Balay values[0] = -2.0; 729371c9d4SSatish Balay values[1] = -2.0; 739371c9d4SSatish Balay values[2] = -2.0; 749566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 3, indices, values, INSERT_VALUES)); 751a254869SHong Zhang break; 761a254869SHong Zhang case 3: 771a254869SHong Zhang row = 0; 789371c9d4SSatish Balay values[0] = -2.0; 799371c9d4SSatish Balay values[1] = -2.0; 809371c9d4SSatish Balay values[2] = -2.0; 819566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 3, indices + 1, values, INSERT_VALUES)); 821a254869SHong Zhang row = 1; 839371c9d4SSatish Balay values[0] = -2.0; 849371c9d4SSatish Balay values[1] = -2.0; 859371c9d4SSatish Balay values[2] = -2.0; 869566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 3, indices, values, INSERT_VALUES)); 871a254869SHong Zhang row = 2; 889371c9d4SSatish Balay values[0] = -2.0; 899371c9d4SSatish Balay values[1] = -2.0; 909371c9d4SSatish Balay values[2] = -2.0; 919566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 3, indices, values, INSERT_VALUES)); 921a254869SHong Zhang row = 3; 939371c9d4SSatish Balay values[0] = -2.0; 949371c9d4SSatish Balay values[1] = -2.0; 959371c9d4SSatish Balay values[2] = -2.0; 969371c9d4SSatish Balay values[3] = -1.0; 979566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 4, indices, values, INSERT_VALUES)); 981a254869SHong Zhang row = 4; 999371c9d4SSatish Balay values[0] = -2.0; 1009371c9d4SSatish Balay values[1] = -2.0; 1019371c9d4SSatish Balay values[2] = -2.0; 1029371c9d4SSatish Balay values[3] = -1.0; 1039566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 4, indices, values, INSERT_VALUES)); 1041a254869SHong Zhang break; 1051a254869SHong Zhang 1061a254869SHong Zhang default: 1071a254869SHong Zhang row = 0; 1089371c9d4SSatish Balay values[0] = -1.0; 1099371c9d4SSatish Balay values[1] = 0.0; 1109371c9d4SSatish Balay values[2] = 1.0; 1119371c9d4SSatish Balay values[3] = 3.0; 1129371c9d4SSatish Balay values[4] = 4.0; 1139371c9d4SSatish Balay values[5] = -5.0; 1149566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, N, indices, values, INSERT_VALUES)); 1151a254869SHong Zhang row = 1; 1169566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 3, indices, values, INSERT_VALUES)); 1171a254869SHong Zhang row = 3; 1189566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, indices + 4, values + 4, INSERT_VALUES)); 119c4762a1bSJed Brown row = 4; 1209566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 2, indices + 4, values + 4, INSERT_VALUES)); 1211a254869SHong Zhang } 1221a254869SHong Zhang } 1239566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1249566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 1259566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 126c4762a1bSJed Brown 1279566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, &n)); 1289566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &min)); 1299566063dSJacob Faibussowitsch PetscCall(VecSetSizes(min, m, PETSC_DECIDE)); 1309566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(min)); 1319566063dSJacob Faibussowitsch PetscCall(VecDuplicate(min, &max)); 1329566063dSJacob Faibussowitsch PetscCall(VecDuplicate(min, &maxabs)); 1339566063dSJacob Faibussowitsch PetscCall(VecDuplicate(min, &e)); 134c4762a1bSJed Brown 135f07e67edSHong Zhang /* MatGetRowMax() */ 1369566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n MatGetRowMax\n")); 1379566063dSJacob Faibussowitsch PetscCall(MatGetRowMax(A, max, NULL)); 1389566063dSJacob Faibussowitsch PetscCall(MatGetRowMax(A, max, imax)); 1399566063dSJacob Faibussowitsch PetscCall(VecView(max, PETSC_VIEWER_STDOUT_WORLD)); 1409566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(max, &n)); 1419566063dSJacob Faibussowitsch PetscCall(PetscIntView(n, imax, PETSC_VIEWER_STDOUT_WORLD)); 142fa213d2fSHong Zhang 143f07e67edSHong Zhang /* MatGetRowMin() */ 1449566063dSJacob Faibussowitsch PetscCall(MatScale(A, -1.0)); 1459566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 1469566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n MatGetRowMin\n")); 1479566063dSJacob Faibussowitsch PetscCall(MatGetRowMin(A, min, NULL)); 1489566063dSJacob Faibussowitsch PetscCall(MatGetRowMin(A, min, imin)); 149fa213d2fSHong Zhang 1509566063dSJacob Faibussowitsch PetscCall(VecWAXPY(e, 1.0, max, min)); /* e = max + min */ 1519566063dSJacob Faibussowitsch PetscCall(VecNorm(e, NORM_INFINITY, &enorm)); 152e00437b9SBarry Smith PetscCheck(enorm <= PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "max+min > PETSC_MACHINE_EPSILON "); 153ad540459SPierre Jolivet for (j = 0; j < n; j++) PetscCheck(imin[j] == imax[j], PETSC_COMM_SELF, PETSC_ERR_PLIB, "imin[%" PetscInt_FMT "] %" PetscInt_FMT " != imax %" PetscInt_FMT, j, imin[j], imax[j]); 154fa213d2fSHong Zhang 155f07e67edSHong Zhang /* MatGetRowMaxAbs() */ 1569566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n MatGetRowMaxAbs\n")); 1579566063dSJacob Faibussowitsch PetscCall(MatGetRowMaxAbs(A, maxabs, NULL)); 1589566063dSJacob Faibussowitsch PetscCall(MatGetRowMaxAbs(A, maxabs, imaxabs)); 1599566063dSJacob Faibussowitsch PetscCall(VecView(maxabs, PETSC_VIEWER_STDOUT_WORLD)); 1609566063dSJacob Faibussowitsch PetscCall(PetscIntView(n, imaxabs, PETSC_VIEWER_STDOUT_WORLD)); 1611a254869SHong Zhang 162f07e67edSHong Zhang /* MatGetRowMinAbs() */ 1639566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n MatGetRowMinAbs\n")); 1649566063dSJacob Faibussowitsch PetscCall(MatGetRowMinAbs(A, min, NULL)); 1659566063dSJacob Faibussowitsch PetscCall(MatGetRowMinAbs(A, min, imin)); 1669566063dSJacob Faibussowitsch PetscCall(VecView(min, PETSC_VIEWER_STDOUT_WORLD)); 1679566063dSJacob Faibussowitsch PetscCall(PetscIntView(n, imin, PETSC_VIEWER_STDOUT_WORLD)); 168f07e67edSHong Zhang 169f07e67edSHong Zhang if (size == 1) { 170f07e67edSHong Zhang /* Test MatGetRowMax, MatGetRowMin and MatGetRowMaxAbs for SeqDense and MPIBAIJ matrix */ 1714e879edeSHong Zhang Mat Adense; 1724e879edeSHong Zhang Vec max_d, maxabs_d; 1739566063dSJacob Faibussowitsch PetscCall(VecDuplicate(min, &max_d)); 1749566063dSJacob Faibussowitsch PetscCall(VecDuplicate(min, &maxabs_d)); 175f07e67edSHong Zhang 1769566063dSJacob Faibussowitsch PetscCall(MatScale(A, -1.0)); 1779566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATDENSE, MAT_INITIAL_MATRIX, &Adense)); 1789566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatGetRowMax for seqdense matrix\n")); 1799566063dSJacob Faibussowitsch PetscCall(MatGetRowMax(Adense, max_d, imax)); 1804e879edeSHong Zhang 1819566063dSJacob Faibussowitsch PetscCall(VecWAXPY(e, -1.0, max, max_d)); /* e = -max + max_d */ 1829566063dSJacob Faibussowitsch PetscCall(VecNorm(e, NORM_INFINITY, &enorm)); 183e00437b9SBarry Smith PetscCheck(enorm <= PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "norm(-max + max_d) %g > PETSC_MACHINE_EPSILON", (double)enorm); 1844e879edeSHong Zhang 1859566063dSJacob Faibussowitsch PetscCall(MatScale(Adense, -1.0)); 1869566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatGetRowMin for seqdense matrix\n")); 1879566063dSJacob Faibussowitsch PetscCall(MatGetRowMin(Adense, min, imin)); 188f07e67edSHong Zhang 1899566063dSJacob Faibussowitsch PetscCall(VecWAXPY(e, 1.0, max, min)); /* e = max + min */ 1909566063dSJacob Faibussowitsch PetscCall(VecNorm(e, NORM_INFINITY, &enorm)); 191e00437b9SBarry Smith PetscCheck(enorm <= PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "max+min > PETSC_MACHINE_EPSILON "); 192ad540459SPierre Jolivet for (j = 0; j < n; j++) PetscCheck(imin[j] == imax[j], PETSC_COMM_SELF, PETSC_ERR_PLIB, "imin[%" PetscInt_FMT "] %" PetscInt_FMT " != imax %" PetscInt_FMT " for seqdense matrix", j, imin[j], imax[j]); 193f07e67edSHong Zhang 1949566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatGetRowMaxAbs for seqdense matrix\n")); 1959566063dSJacob Faibussowitsch PetscCall(MatGetRowMaxAbs(Adense, maxabs_d, imaxabs)); 1969566063dSJacob Faibussowitsch PetscCall(VecWAXPY(e, -1.0, maxabs, maxabs_d)); /* e = -maxabs + maxabs_d */ 1979566063dSJacob Faibussowitsch PetscCall(VecNorm(e, NORM_INFINITY, &enorm)); 198e00437b9SBarry Smith PetscCheck(enorm <= PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "norm(-maxabs + maxabs_d) %g > PETSC_MACHINE_EPSILON", (double)enorm); 199c4762a1bSJed Brown 2009566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Adense)); 2019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&max_d)); 2029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&maxabs_d)); 20343359b5eSHong Zhang } 20443359b5eSHong Zhang 20543359b5eSHong Zhang { /* BAIJ matrix */ 2064e879edeSHong Zhang Mat B; 20743359b5eSHong Zhang Vec maxabsB, maxabsB2; 20843359b5eSHong Zhang PetscInt bs = 2, *imaxabsB, *imaxabsB2, rstart, rend, cstart, cend, ncols, col, Brows[2], Bcols[2]; 20943359b5eSHong Zhang const PetscInt *cols; 21043359b5eSHong Zhang const PetscScalar *vals, *vals2; 21143359b5eSHong Zhang PetscScalar Bvals[4]; 21243359b5eSHong Zhang 2139566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(M, &imaxabsB, bs * M, &imaxabsB2)); 21443359b5eSHong Zhang 21543359b5eSHong Zhang /* bs = 1 */ 2169566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATMPIBAIJ, MAT_INITIAL_MATRIX, &B)); 2179566063dSJacob Faibussowitsch PetscCall(VecDuplicate(min, &maxabsB)); 2189566063dSJacob Faibussowitsch PetscCall(MatGetRowMaxAbs(B, maxabsB, NULL)); 2199566063dSJacob Faibussowitsch PetscCall(MatGetRowMaxAbs(B, maxabsB, imaxabsB)); 2209566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n MatGetRowMaxAbs for BAIJ matrix\n")); 2219566063dSJacob Faibussowitsch PetscCall(VecWAXPY(e, -1.0, maxabs, maxabsB)); /* e = -maxabs + maxabsB */ 2229566063dSJacob Faibussowitsch PetscCall(VecNorm(e, NORM_INFINITY, &enorm)); 223e00437b9SBarry Smith PetscCheck(enorm <= PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "norm(-maxabs + maxabs_d) %g > PETSC_MACHINE_EPSILON", (double)enorm); 2244e879edeSHong Zhang 225ad540459SPierre Jolivet for (j = 0; j < n; j++) PetscCheck(imaxabs[j] == imaxabsB[j], PETSC_COMM_SELF, PETSC_ERR_PLIB, "imaxabs[%" PetscInt_FMT "] %" PetscInt_FMT " != imaxabsB %" PetscInt_FMT, j, imin[j], imax[j]); 2269566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 22743359b5eSHong Zhang 22843359b5eSHong Zhang /* Test bs = 2: Create B with bs*bs block structure of A */ 2299566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &maxabsB2)); 2309566063dSJacob Faibussowitsch PetscCall(VecSetSizes(maxabsB2, bs * m, PETSC_DECIDE)); 2319566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(maxabsB2)); 23243359b5eSHong Zhang 2339566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 2349566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(A, &cstart, &cend)); 2359566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &B)); 2369566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, bs * (rend - rstart), bs * (cend - cstart), PETSC_DECIDE, PETSC_DECIDE)); 2379566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 2389566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 23943359b5eSHong Zhang 24043359b5eSHong Zhang for (row = rstart; row < rend; row++) { 2419566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, row, &ncols, &cols, &vals)); 24243359b5eSHong Zhang for (col = 0; col < ncols; col++) { 24343359b5eSHong Zhang for (j = 0; j < bs; j++) { 24443359b5eSHong Zhang Brows[j] = bs * row + j; 24543359b5eSHong Zhang Bcols[j] = bs * cols[col] + j; 24643359b5eSHong Zhang } 24743359b5eSHong Zhang for (j = 0; j < bs * bs; j++) Bvals[j] = vals[col]; 2489566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, bs, Brows, bs, Bcols, Bvals, INSERT_VALUES)); 24943359b5eSHong Zhang } 2509566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, row, &ncols, &cols, &vals)); 25143359b5eSHong Zhang } 2529566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 2539566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 25443359b5eSHong Zhang 2559566063dSJacob Faibussowitsch PetscCall(MatGetRowMaxAbs(B, maxabsB2, imaxabsB2)); 25643359b5eSHong Zhang 25743359b5eSHong Zhang /* Check maxabsB2 and imaxabsB2 */ 2589566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(maxabsB, &vals)); 2599566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(maxabsB2, &vals2)); 260ad540459SPierre Jolivet for (row = 0; row < m; row++) PetscCheck(PetscAbsScalar(vals[row] - vals2[bs * row]) <= PETSC_MACHINE_EPSILON, PETSC_COMM_SELF, PETSC_ERR_PLIB, "row %" PetscInt_FMT " maxabsB != maxabsB2", row); 2619566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(maxabsB, &vals)); 2629566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(maxabsB2, &vals2)); 26343359b5eSHong Zhang 264ad540459SPierre Jolivet for (col = 0; col < n; col++) PetscCheck(imaxabsB[col] == imaxabsB2[bs * col] / bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "col %" PetscInt_FMT " imaxabsB != imaxabsB2", col); 2659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&maxabsB)); 2669566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 2679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&maxabsB2)); 2689566063dSJacob Faibussowitsch PetscCall(PetscFree2(imaxabsB, imaxabsB2)); 269c4762a1bSJed Brown } 270c4762a1bSJed Brown 2719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&min)); 2729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&max)); 2739566063dSJacob Faibussowitsch PetscCall(VecDestroy(&maxabs)); 2749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&e)); 2759566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2769566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 277b122ec5aSJacob Faibussowitsch return 0; 278c4762a1bSJed Brown } 279c4762a1bSJed Brown 280c4762a1bSJed Brown /*TEST 281c4762a1bSJed Brown 282c4762a1bSJed Brown test: 283c4762a1bSJed Brown output_file: output/ex114.out 284c4762a1bSJed Brown 285c4762a1bSJed Brown test: 286c4762a1bSJed Brown suffix: 2 2874e879edeSHong Zhang args: -testcase 1 2884e879edeSHong Zhang output_file: output/ex114.out 289c4762a1bSJed Brown 290c4762a1bSJed Brown test: 291c4762a1bSJed Brown suffix: 3 2924e879edeSHong Zhang args: -testcase 2 2934e879edeSHong Zhang output_file: output/ex114_3.out 2944e879edeSHong Zhang 2954e879edeSHong Zhang test: 2964e879edeSHong Zhang suffix: 4 2974e879edeSHong Zhang args: -testcase 3 2984e879edeSHong Zhang output_file: output/ex114_4.out 2994e879edeSHong Zhang 3004e879edeSHong Zhang test: 3014e879edeSHong Zhang suffix: 5 3024e879edeSHong Zhang nsize: 3 3034e879edeSHong Zhang args: -testcase 0 3044e879edeSHong Zhang output_file: output/ex114_5.out 3054e879edeSHong Zhang 3064e879edeSHong Zhang test: 3074e879edeSHong Zhang suffix: 6 3084e879edeSHong Zhang nsize: 3 3094e879edeSHong Zhang args: -testcase 1 3104e879edeSHong Zhang output_file: output/ex114_6.out 3114e879edeSHong Zhang 3124e879edeSHong Zhang test: 3134e879edeSHong Zhang suffix: 7 3144e879edeSHong Zhang nsize: 3 3154e879edeSHong Zhang args: -testcase 2 3164e879edeSHong Zhang output_file: output/ex114_7.out 3174e879edeSHong Zhang 3184e879edeSHong Zhang test: 3194e879edeSHong Zhang suffix: 8 3204e879edeSHong Zhang nsize: 3 3214e879edeSHong Zhang args: -testcase 3 3224e879edeSHong Zhang output_file: output/ex114_8.out 323c4762a1bSJed Brown 324c4762a1bSJed Brown TEST*/ 325