1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatGetRowMax(), MatGetRowMin(), MatGetRowMaxAbs()\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown #define M 5 7c4762a1bSJed Brown #define N 6 8c4762a1bSJed Brown 9*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 10*d71ae5a4SJacob Faibussowitsch { 11c4762a1bSJed Brown Mat A; 12fa213d2fSHong Zhang Vec min, max, maxabs, e; 134e879edeSHong Zhang PetscInt m, n, j, imin[M], imax[M], imaxabs[M], indices[N], row, testcase = 0; 14c4762a1bSJed Brown PetscScalar values[N]; 151a254869SHong Zhang PetscMPIInt size, rank; 16fa213d2fSHong Zhang PetscReal enorm; 17c4762a1bSJed Brown 18327415f7SBarry Smith PetscFunctionBeginUser; 199566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 229566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-testcase", &testcase, NULL)); 23c4762a1bSJed Brown 249566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 251a254869SHong Zhang if (testcase == 1) { /* proc[0] holds entire A and other processes have no entry */ 26dd400576SPatrick Sanan if (rank == 0) { 279566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, M, N, PETSC_DECIDE, PETSC_DECIDE)); 281a254869SHong Zhang } else { 299566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, 0, 0, PETSC_DECIDE, PETSC_DECIDE)); 301a254869SHong Zhang } 311a254869SHong Zhang testcase = 0; 321a254869SHong Zhang } else { 339566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, M, N)); 341a254869SHong Zhang } 359566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 369566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 371a254869SHong Zhang 38dd400576SPatrick Sanan if (rank == 0) { /* proc[0] sets matrix A */ 391a254869SHong Zhang for (j = 0; j < N; j++) indices[j] = j; 401a254869SHong Zhang switch (testcase) { 41*d71ae5a4SJacob Faibussowitsch case 1: /* see testcast 0 */ 42*d71ae5a4SJacob Faibussowitsch break; 431a254869SHong Zhang case 2: 44c4762a1bSJed Brown row = 0; 459371c9d4SSatish Balay values[0] = -2.0; 469371c9d4SSatish Balay values[1] = -2.0; 479371c9d4SSatish Balay values[2] = -2.0; 489371c9d4SSatish Balay values[3] = -4.0; 499371c9d4SSatish Balay values[4] = 1.0; 509371c9d4SSatish Balay values[5] = 1.0; 519566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, N, indices, values, INSERT_VALUES)); 521a254869SHong Zhang row = 2; 539371c9d4SSatish Balay indices[0] = 0; 549371c9d4SSatish Balay indices[1] = 3; 559371c9d4SSatish Balay indices[2] = 5; 569371c9d4SSatish Balay values[0] = -2.0; 579371c9d4SSatish Balay values[1] = -2.0; 589371c9d4SSatish Balay values[2] = -2.0; 599566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 3, indices, values, INSERT_VALUES)); 601a254869SHong Zhang row = 3; 619371c9d4SSatish Balay indices[0] = 0; 629371c9d4SSatish Balay indices[1] = 1; 639371c9d4SSatish Balay indices[2] = 4; 649371c9d4SSatish Balay values[0] = -2.0; 659371c9d4SSatish Balay values[1] = -2.0; 669371c9d4SSatish Balay values[2] = -2.0; 679566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 3, indices, values, INSERT_VALUES)); 68c4762a1bSJed Brown row = 4; 699371c9d4SSatish Balay indices[0] = 0; 709371c9d4SSatish Balay indices[1] = 1; 719371c9d4SSatish Balay indices[2] = 2; 729371c9d4SSatish Balay values[0] = -2.0; 739371c9d4SSatish Balay values[1] = -2.0; 749371c9d4SSatish Balay values[2] = -2.0; 759566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 3, indices, values, INSERT_VALUES)); 761a254869SHong Zhang break; 771a254869SHong Zhang case 3: 781a254869SHong Zhang row = 0; 799371c9d4SSatish Balay values[0] = -2.0; 809371c9d4SSatish Balay values[1] = -2.0; 819371c9d4SSatish Balay values[2] = -2.0; 829566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 3, indices + 1, values, INSERT_VALUES)); 831a254869SHong Zhang row = 1; 849371c9d4SSatish Balay values[0] = -2.0; 859371c9d4SSatish Balay values[1] = -2.0; 869371c9d4SSatish Balay values[2] = -2.0; 879566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 3, indices, values, INSERT_VALUES)); 881a254869SHong Zhang row = 2; 899371c9d4SSatish Balay values[0] = -2.0; 909371c9d4SSatish Balay values[1] = -2.0; 919371c9d4SSatish Balay values[2] = -2.0; 929566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 3, indices, values, INSERT_VALUES)); 931a254869SHong Zhang row = 3; 949371c9d4SSatish Balay values[0] = -2.0; 959371c9d4SSatish Balay values[1] = -2.0; 969371c9d4SSatish Balay values[2] = -2.0; 979371c9d4SSatish Balay values[3] = -1.0; 989566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 4, indices, values, INSERT_VALUES)); 991a254869SHong Zhang row = 4; 1009371c9d4SSatish Balay values[0] = -2.0; 1019371c9d4SSatish Balay values[1] = -2.0; 1029371c9d4SSatish Balay values[2] = -2.0; 1039371c9d4SSatish Balay values[3] = -1.0; 1049566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 4, indices, values, INSERT_VALUES)); 1051a254869SHong Zhang break; 1061a254869SHong Zhang 1071a254869SHong Zhang default: 1081a254869SHong Zhang row = 0; 1099371c9d4SSatish Balay values[0] = -1.0; 1109371c9d4SSatish Balay values[1] = 0.0; 1119371c9d4SSatish Balay values[2] = 1.0; 1129371c9d4SSatish Balay values[3] = 3.0; 1139371c9d4SSatish Balay values[4] = 4.0; 1149371c9d4SSatish Balay values[5] = -5.0; 1159566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, N, indices, values, INSERT_VALUES)); 1161a254869SHong Zhang row = 1; 1179566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 3, indices, values, INSERT_VALUES)); 1181a254869SHong Zhang row = 3; 1199566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, indices + 4, values + 4, INSERT_VALUES)); 120c4762a1bSJed Brown row = 4; 1219566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 2, indices + 4, values + 4, INSERT_VALUES)); 1221a254869SHong Zhang } 1231a254869SHong Zhang } 1249566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1259566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 1269566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 127c4762a1bSJed Brown 1289566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, &n)); 1299566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &min)); 1309566063dSJacob Faibussowitsch PetscCall(VecSetSizes(min, m, PETSC_DECIDE)); 1319566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(min)); 1329566063dSJacob Faibussowitsch PetscCall(VecDuplicate(min, &max)); 1339566063dSJacob Faibussowitsch PetscCall(VecDuplicate(min, &maxabs)); 1349566063dSJacob Faibussowitsch PetscCall(VecDuplicate(min, &e)); 135c4762a1bSJed Brown 136f07e67edSHong Zhang /* MatGetRowMax() */ 1379566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n MatGetRowMax\n")); 1389566063dSJacob Faibussowitsch PetscCall(MatGetRowMax(A, max, NULL)); 1399566063dSJacob Faibussowitsch PetscCall(MatGetRowMax(A, max, imax)); 1409566063dSJacob Faibussowitsch PetscCall(VecView(max, PETSC_VIEWER_STDOUT_WORLD)); 1419566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(max, &n)); 1429566063dSJacob Faibussowitsch PetscCall(PetscIntView(n, imax, PETSC_VIEWER_STDOUT_WORLD)); 143fa213d2fSHong Zhang 144f07e67edSHong Zhang /* MatGetRowMin() */ 1459566063dSJacob Faibussowitsch PetscCall(MatScale(A, -1.0)); 1469566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 1479566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n MatGetRowMin\n")); 1489566063dSJacob Faibussowitsch PetscCall(MatGetRowMin(A, min, NULL)); 1499566063dSJacob Faibussowitsch PetscCall(MatGetRowMin(A, min, imin)); 150fa213d2fSHong Zhang 1519566063dSJacob Faibussowitsch PetscCall(VecWAXPY(e, 1.0, max, min)); /* e = max + min */ 1529566063dSJacob Faibussowitsch PetscCall(VecNorm(e, NORM_INFINITY, &enorm)); 153e00437b9SBarry Smith PetscCheck(enorm <= PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "max+min > PETSC_MACHINE_EPSILON "); 154ad540459SPierre 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]); 155fa213d2fSHong Zhang 156f07e67edSHong Zhang /* MatGetRowMaxAbs() */ 1579566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n MatGetRowMaxAbs\n")); 1589566063dSJacob Faibussowitsch PetscCall(MatGetRowMaxAbs(A, maxabs, NULL)); 1599566063dSJacob Faibussowitsch PetscCall(MatGetRowMaxAbs(A, maxabs, imaxabs)); 1609566063dSJacob Faibussowitsch PetscCall(VecView(maxabs, PETSC_VIEWER_STDOUT_WORLD)); 1619566063dSJacob Faibussowitsch PetscCall(PetscIntView(n, imaxabs, PETSC_VIEWER_STDOUT_WORLD)); 1621a254869SHong Zhang 163f07e67edSHong Zhang /* MatGetRowMinAbs() */ 1649566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n MatGetRowMinAbs\n")); 1659566063dSJacob Faibussowitsch PetscCall(MatGetRowMinAbs(A, min, NULL)); 1669566063dSJacob Faibussowitsch PetscCall(MatGetRowMinAbs(A, min, imin)); 1679566063dSJacob Faibussowitsch PetscCall(VecView(min, PETSC_VIEWER_STDOUT_WORLD)); 1689566063dSJacob Faibussowitsch PetscCall(PetscIntView(n, imin, PETSC_VIEWER_STDOUT_WORLD)); 169f07e67edSHong Zhang 170f07e67edSHong Zhang if (size == 1) { 171f07e67edSHong Zhang /* Test MatGetRowMax, MatGetRowMin and MatGetRowMaxAbs for SeqDense and MPIBAIJ matrix */ 1724e879edeSHong Zhang Mat Adense; 1734e879edeSHong Zhang Vec max_d, maxabs_d; 1749566063dSJacob Faibussowitsch PetscCall(VecDuplicate(min, &max_d)); 1759566063dSJacob Faibussowitsch PetscCall(VecDuplicate(min, &maxabs_d)); 176f07e67edSHong Zhang 1779566063dSJacob Faibussowitsch PetscCall(MatScale(A, -1.0)); 1789566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATDENSE, MAT_INITIAL_MATRIX, &Adense)); 1799566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatGetRowMax for seqdense matrix\n")); 1809566063dSJacob Faibussowitsch PetscCall(MatGetRowMax(Adense, max_d, imax)); 1814e879edeSHong Zhang 1829566063dSJacob Faibussowitsch PetscCall(VecWAXPY(e, -1.0, max, max_d)); /* e = -max + max_d */ 1839566063dSJacob Faibussowitsch PetscCall(VecNorm(e, NORM_INFINITY, &enorm)); 184e00437b9SBarry Smith PetscCheck(enorm <= PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "norm(-max + max_d) %g > PETSC_MACHINE_EPSILON", (double)enorm); 1854e879edeSHong Zhang 1869566063dSJacob Faibussowitsch PetscCall(MatScale(Adense, -1.0)); 1879566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatGetRowMin for seqdense matrix\n")); 1889566063dSJacob Faibussowitsch PetscCall(MatGetRowMin(Adense, min, imin)); 189f07e67edSHong Zhang 1909566063dSJacob Faibussowitsch PetscCall(VecWAXPY(e, 1.0, max, min)); /* e = max + min */ 1919566063dSJacob Faibussowitsch PetscCall(VecNorm(e, NORM_INFINITY, &enorm)); 192e00437b9SBarry Smith PetscCheck(enorm <= PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "max+min > PETSC_MACHINE_EPSILON "); 193ad540459SPierre 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]); 194f07e67edSHong Zhang 1959566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatGetRowMaxAbs for seqdense matrix\n")); 1969566063dSJacob Faibussowitsch PetscCall(MatGetRowMaxAbs(Adense, maxabs_d, imaxabs)); 1979566063dSJacob Faibussowitsch PetscCall(VecWAXPY(e, -1.0, maxabs, maxabs_d)); /* e = -maxabs + maxabs_d */ 1989566063dSJacob Faibussowitsch PetscCall(VecNorm(e, NORM_INFINITY, &enorm)); 199e00437b9SBarry Smith PetscCheck(enorm <= PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "norm(-maxabs + maxabs_d) %g > PETSC_MACHINE_EPSILON", (double)enorm); 200c4762a1bSJed Brown 2019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Adense)); 2029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&max_d)); 2039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&maxabs_d)); 20443359b5eSHong Zhang } 20543359b5eSHong Zhang 20643359b5eSHong Zhang { /* BAIJ matrix */ 2074e879edeSHong Zhang Mat B; 20843359b5eSHong Zhang Vec maxabsB, maxabsB2; 20943359b5eSHong Zhang PetscInt bs = 2, *imaxabsB, *imaxabsB2, rstart, rend, cstart, cend, ncols, col, Brows[2], Bcols[2]; 21043359b5eSHong Zhang const PetscInt *cols; 21143359b5eSHong Zhang const PetscScalar *vals, *vals2; 21243359b5eSHong Zhang PetscScalar Bvals[4]; 21343359b5eSHong Zhang 2149566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(M, &imaxabsB, bs * M, &imaxabsB2)); 21543359b5eSHong Zhang 21643359b5eSHong Zhang /* bs = 1 */ 2179566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATMPIBAIJ, MAT_INITIAL_MATRIX, &B)); 2189566063dSJacob Faibussowitsch PetscCall(VecDuplicate(min, &maxabsB)); 2199566063dSJacob Faibussowitsch PetscCall(MatGetRowMaxAbs(B, maxabsB, NULL)); 2209566063dSJacob Faibussowitsch PetscCall(MatGetRowMaxAbs(B, maxabsB, imaxabsB)); 2219566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n MatGetRowMaxAbs for BAIJ matrix\n")); 2229566063dSJacob Faibussowitsch PetscCall(VecWAXPY(e, -1.0, maxabs, maxabsB)); /* e = -maxabs + maxabsB */ 2239566063dSJacob Faibussowitsch PetscCall(VecNorm(e, NORM_INFINITY, &enorm)); 224e00437b9SBarry Smith PetscCheck(enorm <= PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "norm(-maxabs + maxabs_d) %g > PETSC_MACHINE_EPSILON", (double)enorm); 2254e879edeSHong Zhang 226ad540459SPierre 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]); 2279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 22843359b5eSHong Zhang 22943359b5eSHong Zhang /* Test bs = 2: Create B with bs*bs block structure of A */ 2309566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &maxabsB2)); 2319566063dSJacob Faibussowitsch PetscCall(VecSetSizes(maxabsB2, bs * m, PETSC_DECIDE)); 2329566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(maxabsB2)); 23343359b5eSHong Zhang 2349566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 2359566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(A, &cstart, &cend)); 2369566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &B)); 2379566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, bs * (rend - rstart), bs * (cend - cstart), PETSC_DECIDE, PETSC_DECIDE)); 2389566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 2399566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 24043359b5eSHong Zhang 24143359b5eSHong Zhang for (row = rstart; row < rend; row++) { 2429566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, row, &ncols, &cols, &vals)); 24343359b5eSHong Zhang for (col = 0; col < ncols; col++) { 24443359b5eSHong Zhang for (j = 0; j < bs; j++) { 24543359b5eSHong Zhang Brows[j] = bs * row + j; 24643359b5eSHong Zhang Bcols[j] = bs * cols[col] + j; 24743359b5eSHong Zhang } 24843359b5eSHong Zhang for (j = 0; j < bs * bs; j++) Bvals[j] = vals[col]; 2499566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, bs, Brows, bs, Bcols, Bvals, INSERT_VALUES)); 25043359b5eSHong Zhang } 2519566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, row, &ncols, &cols, &vals)); 25243359b5eSHong Zhang } 2539566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 2549566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 25543359b5eSHong Zhang 2569566063dSJacob Faibussowitsch PetscCall(MatGetRowMaxAbs(B, maxabsB2, imaxabsB2)); 25743359b5eSHong Zhang 25843359b5eSHong Zhang /* Check maxabsB2 and imaxabsB2 */ 2599566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(maxabsB, &vals)); 2609566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(maxabsB2, &vals2)); 261ad540459SPierre 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); 2629566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(maxabsB, &vals)); 2639566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(maxabsB2, &vals2)); 26443359b5eSHong Zhang 265ad540459SPierre 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); 2669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&maxabsB)); 2679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 2689566063dSJacob Faibussowitsch PetscCall(VecDestroy(&maxabsB2)); 2699566063dSJacob Faibussowitsch PetscCall(PetscFree2(imaxabsB, imaxabsB2)); 270c4762a1bSJed Brown } 271c4762a1bSJed Brown 2729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&min)); 2739566063dSJacob Faibussowitsch PetscCall(VecDestroy(&max)); 2749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&maxabs)); 2759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&e)); 2769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2779566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 278b122ec5aSJacob Faibussowitsch return 0; 279c4762a1bSJed Brown } 280c4762a1bSJed Brown 281c4762a1bSJed Brown /*TEST 282c4762a1bSJed Brown 283c4762a1bSJed Brown test: 284c4762a1bSJed Brown output_file: output/ex114.out 285c4762a1bSJed Brown 286c4762a1bSJed Brown test: 287c4762a1bSJed Brown suffix: 2 2884e879edeSHong Zhang args: -testcase 1 2894e879edeSHong Zhang output_file: output/ex114.out 290c4762a1bSJed Brown 291c4762a1bSJed Brown test: 292c4762a1bSJed Brown suffix: 3 2934e879edeSHong Zhang args: -testcase 2 2944e879edeSHong Zhang output_file: output/ex114_3.out 2954e879edeSHong Zhang 2964e879edeSHong Zhang test: 2974e879edeSHong Zhang suffix: 4 2984e879edeSHong Zhang args: -testcase 3 2994e879edeSHong Zhang output_file: output/ex114_4.out 3004e879edeSHong Zhang 3014e879edeSHong Zhang test: 3024e879edeSHong Zhang suffix: 5 3034e879edeSHong Zhang nsize: 3 3044e879edeSHong Zhang args: -testcase 0 3054e879edeSHong Zhang output_file: output/ex114_5.out 3064e879edeSHong Zhang 3074e879edeSHong Zhang test: 3084e879edeSHong Zhang suffix: 6 3094e879edeSHong Zhang nsize: 3 3104e879edeSHong Zhang args: -testcase 1 3114e879edeSHong Zhang output_file: output/ex114_6.out 3124e879edeSHong Zhang 3134e879edeSHong Zhang test: 3144e879edeSHong Zhang suffix: 7 3154e879edeSHong Zhang nsize: 3 3164e879edeSHong Zhang args: -testcase 2 3174e879edeSHong Zhang output_file: output/ex114_7.out 3184e879edeSHong Zhang 3194e879edeSHong Zhang test: 3204e879edeSHong Zhang suffix: 8 3214e879edeSHong Zhang nsize: 3 3224e879edeSHong Zhang args: -testcase 3 3234e879edeSHong Zhang output_file: output/ex114_8.out 324c4762a1bSJed Brown 325c4762a1bSJed Brown TEST*/ 326