Actual source code: ex52.c
1: static char help[] = "Solves a linear system in parallel with KSP. Modified from ex2.c \n\
2: Illustrate how to use external packages MUMPS, SUPERLU and STRUMPACK \n\
3: Input parameters include:\n\
4: -random_exact_sol : use a random exact solution vector\n\
5: -view_exact_sol : write exact solution vector to stdout\n\
6: -m <mesh_x> : number of mesh points in x-direction\n\
7: -n <mesh_y> : number of mesh points in y-direction\n\n";
9: #include <petscksp.h>
11: #if defined(PETSC_HAVE_MUMPS)
12: /* Subroutine contributed by Varun Hiremath */
13: PetscErrorCode printMumpsMemoryInfo(Mat F)
14: {
15: PetscInt maxMem, sumMem;
17: PetscFunctionBeginUser;
18: PetscCall(MatMumpsGetInfog(F, 16, &maxMem));
19: PetscCall(MatMumpsGetInfog(F, 17, &sumMem));
20: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n MUMPS INFOG(16) :: Max memory in MB = %" PetscInt_FMT, maxMem));
21: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n MUMPS INFOG(17) :: Sum memory in MB = %" PetscInt_FMT "\n", sumMem));
22: PetscFunctionReturn(PETSC_SUCCESS);
23: }
24: #endif
26: int main(int argc, char **args)
27: {
28: Vec x, b, u; /* approx solution, RHS, exact solution */
29: Mat A, F;
30: KSP ksp; /* linear solver context */
31: PC pc;
32: PetscRandom rctx; /* random number generator context */
33: PetscReal norm; /* norm of solution error */
34: PetscInt i, j, Ii, J, Istart, Iend, m = 8, n = 7, its;
35: PetscBool flg = PETSC_FALSE, flg_ilu = PETSC_FALSE, flg_ch = PETSC_FALSE;
36: #if defined(PETSC_HAVE_MUMPS)
37: char tmpdir[PETSC_MAX_PATH_LEN];
38: PetscBool flg_mumps = PETSC_FALSE, flg_mumps_ch = PETSC_FALSE, test_mumps_ooc_api = PETSC_FALSE;
39: #endif
40: #if defined(PETSC_HAVE_SUPERLU) || defined(PETSC_HAVE_SUPERLU_DIST)
41: PetscBool flg_superlu = PETSC_FALSE;
42: #endif
43: #if defined(PETSC_HAVE_STRUMPACK)
44: PetscBool flg_strumpack = PETSC_FALSE;
45: #endif
46: PetscScalar v;
47: PetscMPIInt rank, size;
48: PetscLogStage stage;
50: #if defined(PETSC_HAVE_STRUMPACK) && defined(PETSC_HAVE_SLATE)
51: PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_MULTIPLE;
52: #endif
53: PetscFunctionBeginUser;
54: PetscCall(PetscInitialize(&argc, &args, NULL, help));
55: PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
56: PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
57: PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
58: PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
59: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
60: Compute the matrix and right-hand-side vector that define
61: the linear system, Ax = b.
62: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
63: PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
64: PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n));
65: PetscCall(MatSetFromOptions(A));
66: PetscCall(MatMPIAIJSetPreallocation(A, 5, NULL, 5, NULL));
67: PetscCall(MatSeqAIJSetPreallocation(A, 5, NULL));
68: PetscCall(MatSetUp(A));
70: /*
71: Currently, all PETSc parallel matrix formats are partitioned by
72: contiguous chunks of rows across the processors. Determine which
73: rows of the matrix are locally owned.
74: */
75: PetscCall(MatGetOwnershipRange(A, &Istart, &Iend));
77: /*
78: Set matrix elements for the 2-D, five-point stencil in parallel.
79: - Each processor needs to insert only elements that it owns
80: locally (but any non-local elements will be sent to the
81: appropriate processor during matrix assembly).
82: - Always specify global rows and columns of matrix entries.
84: Note: this uses the less common natural ordering that orders first
85: all the unknowns for x = h then for x = 2h etc; Hence you see J = Ii +- n
86: instead of J = I +- m as you might expect. The more standard ordering
87: would first do all variables for y = h, then y = 2h etc.
89: */
90: PetscCall(PetscLogStageRegister("Assembly", &stage));
91: PetscCall(PetscLogStagePush(stage));
92: for (Ii = Istart; Ii < Iend; Ii++) {
93: v = -1.0;
94: i = Ii / n;
95: j = Ii - i * n;
96: if (i > 0) {
97: J = Ii - n;
98: PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES));
99: }
100: if (i < m - 1) {
101: J = Ii + n;
102: PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES));
103: }
104: if (j > 0) {
105: J = Ii - 1;
106: PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES));
107: }
108: if (j < n - 1) {
109: J = Ii + 1;
110: PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES));
111: }
112: v = 4.0;
113: PetscCall(MatSetValues(A, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
114: }
116: /*
117: Assemble matrix, using the 2-step process:
118: MatAssemblyBegin(), MatAssemblyEnd()
119: Computations can be done while messages are in transition
120: by placing code between these two statements.
121: */
122: PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
123: PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
124: PetscCall(PetscLogStagePop());
126: /* A is symmetric. Set symmetric flag to enable ICC/Cholesky preconditioner */
127: PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE));
129: /* Create parallel vectors */
130: PetscCall(MatCreateVecs(A, &u, &b));
131: PetscCall(VecDuplicate(u, &x));
133: /*
134: Set exact solution; then compute right-hand-side vector.
135: By default we use an exact solution of a vector with all
136: elements of 1.0; Alternatively, using the runtime option
137: -random_sol forms a solution vector with random components.
138: */
139: PetscCall(PetscOptionsGetBool(NULL, NULL, "-random_exact_sol", &flg, NULL));
140: if (flg) {
141: PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rctx));
142: PetscCall(PetscRandomSetFromOptions(rctx));
143: PetscCall(VecSetRandom(u, rctx));
144: PetscCall(PetscRandomDestroy(&rctx));
145: } else {
146: PetscCall(VecSet(u, 1.0));
147: }
148: PetscCall(MatMult(A, u, b));
150: /*
151: View the exact solution vector if desired
152: */
153: flg = PETSC_FALSE;
154: PetscCall(PetscOptionsGetBool(NULL, NULL, "-view_exact_sol", &flg, NULL));
155: if (flg) PetscCall(VecView(u, PETSC_VIEWER_STDOUT_WORLD));
157: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158: Create the linear solver and set various options
159: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
161: /*
162: Create linear solver context
163: */
164: PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
165: PetscCall(KSPSetOperators(ksp, A, A));
167: /*
168: Example of how to use external package MUMPS
169: Note: runtime options
170: '-ksp_type preonly -pc_type lu -pc_factor_mat_solver_type mumps -mat_mumps_icntl_7 3 -mat_mumps_icntl_1 0.0'
171: are equivalent to these procedural calls
172: */
173: #if defined(PETSC_HAVE_MUMPS)
174: flg_mumps = PETSC_FALSE;
175: flg_mumps_ch = PETSC_FALSE;
176: PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_mumps_lu", &flg_mumps, NULL));
177: PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_mumps_ch", &flg_mumps_ch, NULL));
178: if (flg_mumps || flg_mumps_ch) {
179: PetscCall(KSPSetType(ksp, KSPPREONLY));
180: PetscInt ival, icntl;
181: PetscReal val;
182: PetscCall(KSPGetPC(ksp, &pc));
183: if (flg_mumps) PetscCall(PCSetType(pc, PCLU));
184: else if (flg_mumps_ch) {
185: PetscCall(MatSetOption(A, MAT_SPD, PETSC_TRUE)); /* set MUMPS id%SYM=1 */
186: PetscCall(PCSetType(pc, PCCHOLESKY));
187: }
188: PetscCall(PCFactorSetMatSolverType(pc, MATSOLVERMUMPS));
189: PetscCall(PCFactorSetUpMatSolverType(pc)); /* call MatGetFactor() to create F */
190: PetscCall(PCFactorGetMatrix(pc, &F));
191: PetscCall(MatMumpsSetIcntl(F, 24, 1));
192: PetscCall(MatMumpsGetIcntl(F, 24, &ival));
193: PetscCheck(ival == 1, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "ICNTL(24) = %" PetscInt_FMT " (!= 1)", ival);
194: PetscCall(MatMumpsSetCntl(F, 3, 1e-6));
195: PetscCall(MatMumpsGetCntl(F, 3, &val));
196: PetscCheck(PetscEqualReal(val, 1e-6), PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CNTL(3) = %g (!= %g)", (double)val, 1e-6);
197: if (flg_mumps) {
198: /* Zero the first and last rows in the rank, they should then show up in corresponding null pivot rows output via
199: MatMumpsGetNullPivots */
200: flg = PETSC_FALSE;
201: PetscCall(PetscOptionsGetBool(NULL, NULL, "-zero_first_and_last_rows", &flg, NULL));
202: if (flg) {
203: PetscInt rows[2];
204: rows[0] = Istart; /* first row of the rank */
205: rows[1] = Iend - 1; /* last row of the rank */
206: PetscCall(MatZeroRows(A, 2, rows, 0.0, NULL, NULL));
207: }
208: /* Get memory estimates from MUMPS' MatLUFactorSymbolic(), e.g. INFOG(16), INFOG(17).
209: KSPSetUp() below will do nothing inside MatLUFactorSymbolic() */
210: MatFactorInfo info;
211: PetscCall(MatLUFactorSymbolic(F, A, NULL, NULL, &info));
212: flg = PETSC_FALSE;
213: PetscCall(PetscOptionsGetBool(NULL, NULL, "-print_mumps_memory", &flg, NULL));
214: if (flg) PetscCall(printMumpsMemoryInfo(F));
215: }
217: /* sequential ordering */
218: icntl = 7;
219: ival = 2;
220: PetscCall(MatMumpsSetIcntl(F, icntl, ival));
222: /* threshold for row pivot detection */
223: PetscCall(MatMumpsGetIcntl(F, 24, &ival));
224: PetscCheck(ival == 1, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "ICNTL(24) = %" PetscInt_FMT " (!= 1)", ival);
225: icntl = 3;
226: PetscCall(MatMumpsGetCntl(F, icntl, &val));
227: PetscCheck(PetscEqualReal(val, 1e-6), PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CNTL(3) = %g (!= %g)", (double)val, 1e-6);
229: /* compute determinant of A */
230: PetscCall(MatMumpsSetIcntl(F, 33, 1));
231: }
232: #endif
234: /*
235: Example of how to use external package SuperLU
236: Note: runtime options
237: '-ksp_type preonly -pc_type ilu -pc_factor_mat_solver_type superlu -mat_superlu_ilu_droptol 1.e-8'
238: are equivalent to these procedual calls
239: */
240: #if defined(PETSC_HAVE_SUPERLU) || defined(PETSC_HAVE_SUPERLU_DIST)
241: flg_ilu = PETSC_FALSE;
242: flg_superlu = PETSC_FALSE;
243: PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_superlu_lu", &flg_superlu, NULL));
244: PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_superlu_ilu", &flg_ilu, NULL));
245: if (flg_superlu || flg_ilu) {
246: PetscCall(KSPSetType(ksp, KSPPREONLY));
247: PetscCall(KSPGetPC(ksp, &pc));
248: if (flg_superlu) PetscCall(PCSetType(pc, PCLU));
249: else if (flg_ilu) PetscCall(PCSetType(pc, PCILU));
250: if (size == 1) {
251: #if !defined(PETSC_HAVE_SUPERLU)
252: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "This test requires SUPERLU");
253: #else
254: PetscCall(PCFactorSetMatSolverType(pc, MATSOLVERSUPERLU));
255: #endif
256: } else {
257: #if !defined(PETSC_HAVE_SUPERLU_DIST)
258: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "This test requires SUPERLU_DIST");
259: #else
260: PetscCall(PCFactorSetMatSolverType(pc, MATSOLVERSUPERLU_DIST));
261: #endif
262: }
263: PetscCall(PCFactorSetUpMatSolverType(pc)); /* call MatGetFactor() to create F */
264: PetscCall(PCFactorGetMatrix(pc, &F));
265: #if defined(PETSC_HAVE_SUPERLU)
266: if (size == 1) PetscCall(MatSuperluSetILUDropTol(F, 1.e-8));
267: #endif
268: }
269: #endif
271: /*
272: Example of how to use external package STRUMPACK
273: Note: runtime options
274: '-pc_type lu/ilu \
275: -pc_factor_mat_solver_type strumpack \
276: -mat_strumpack_reordering GEOMETRIC \
277: -mat_strumpack_geometric_xyz n,m \
278: -mat_strumpack_colperm 0 \
279: -mat_strumpack_compression_rel_tol 1.e-3 \
280: -mat_strumpack_compression_min_sep_size 15 \
281: -mat_strumpack_leaf_size 4'
282: are equivalent to these procedural calls
284: We refer to the STRUMPACK manual for more info on
285: how to tune the preconditioner, see for instance:
286: https://portal.nersc.gov/project/sparse/strumpack/master/prec.html
287: */
288: #if defined(PETSC_HAVE_STRUMPACK)
289: flg_ilu = PETSC_FALSE;
290: flg_strumpack = PETSC_FALSE;
291: PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_strumpack_lu", &flg_strumpack, NULL));
292: PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_strumpack_ilu", &flg_ilu, NULL));
293: if (flg_strumpack || flg_ilu) {
294: PetscCall(KSPSetType(ksp, KSPPREONLY));
295: PetscCall(KSPGetPC(ksp, &pc));
296: if (flg_strumpack) PetscCall(PCSetType(pc, PCLU));
297: else if (flg_ilu) PetscCall(PCSetType(pc, PCILU));
298: #if !defined(PETSC_HAVE_STRUMPACK)
299: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "This test requires STRUMPACK");
300: #endif
301: PetscCall(PCFactorSetMatSolverType(pc, MATSOLVERSTRUMPACK));
302: PetscCall(PCFactorSetUpMatSolverType(pc)); /* call MatGetFactor() to create F */
303: PetscCall(PCFactorGetMatrix(pc, &F));
305: /* Set the fill-reducing reordering, MAT_STRUMPACK_METIS is */
306: /* always supported, but is sequential. Parallel alternatives are */
307: /* MAT_STRUMPACK_PARMETIS and MAT_STRUMPACK_PTSCOTCH, but */
308: /* strumpack needs to be configured with support for these. */
309: /*PetscCall(MatSTRUMPACKSetReordering(F, MAT_STRUMPACK_METIS)); */
310: /* However, since this is a problem on a regular grid, we can use */
311: /* a simple geometric nested dissection implementation, which */
312: /* requires passing the grid dimensions to strumpack. */
313: PetscCall(MatSTRUMPACKSetReordering(F, MAT_STRUMPACK_GEOMETRIC));
314: PetscCall(MatSTRUMPACKSetGeometricNxyz(F, n, m, PETSC_DECIDE));
315: /* These are optional, defaults are 1 */
316: PetscCall(MatSTRUMPACKSetGeometricComponents(F, 1));
317: PetscCall(MatSTRUMPACKSetGeometricWidth(F, 1));
319: /* Since this is a simple discretization, the diagonal is always */
320: /* nonzero, and there is no need for the extra MC64 permutation. */
321: PetscCall(MatSTRUMPACKSetColPerm(F, PETSC_FALSE));
323: if (flg_ilu) {
324: /* The compression tolerance used when doing low-rank compression */
325: /* in the preconditioner. This is problem specific! */
326: PetscCall(MatSTRUMPACKSetCompRelTol(F, 1.e-3));
328: /* Set a small minimum (dense) matrix size for compression to */
329: /* demonstrate the preconditioner on small problems. */
330: /* For performance the default value should be better. */
331: /* This size corresponds to the size of separators in the graph. */
332: /* For instance on an m x n mesh, the top level separator is of */
333: /* size m (if m <= n) */
334: /*PetscCall(MatSTRUMPACKSetCompMinSepSize(F,15));*/
336: /* Set the size of the diagonal blocks (the leafs) in the HSS */
337: /* approximation. The default value should be better for real */
338: /* problems. This is mostly for illustration on a small problem. */
339: /*PetscCall(MatSTRUMPACKSetCompLeafSize(F,4));*/
340: }
341: }
342: #endif
344: /*
345: Example of how to use procedural calls that are equivalent to
346: '-ksp_type preonly -pc_type lu/ilu -pc_factor_mat_solver_type petsc'
347: */
348: flg = PETSC_FALSE;
349: flg_ilu = PETSC_FALSE;
350: flg_ch = PETSC_FALSE;
351: PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_petsc_lu", &flg, NULL));
352: PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_petsc_ilu", &flg_ilu, NULL));
353: PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_petsc_ch", &flg_ch, NULL));
354: if (flg || flg_ilu || flg_ch) {
355: Vec diag;
357: PetscCall(KSPSetType(ksp, KSPPREONLY));
358: PetscCall(KSPGetPC(ksp, &pc));
359: if (flg) PetscCall(PCSetType(pc, PCLU));
360: else if (flg_ilu) PetscCall(PCSetType(pc, PCILU));
361: else if (flg_ch) PetscCall(PCSetType(pc, PCCHOLESKY));
362: PetscCall(PCFactorSetMatSolverType(pc, MATSOLVERPETSC));
363: PetscCall(PCFactorSetUpMatSolverType(pc)); /* call MatGetFactor() to create F */
364: PetscCall(PCFactorGetMatrix(pc, &F));
366: /* Test MatGetDiagonal() */
367: PetscCall(KSPSetUp(ksp));
368: PetscCall(VecDuplicate(x, &diag));
369: PetscCall(MatGetDiagonal(F, diag));
370: /* PetscCall(VecView(diag,PETSC_VIEWER_STDOUT_WORLD)); */
371: PetscCall(VecDestroy(&diag));
372: }
374: PetscCall(KSPSetFromOptions(ksp));
376: #if defined(PETSC_HAVE_MUMPS)
377: // The OOC options must be set before PCSetUp()/KSPSetUp(), as MUMPS requires them be set after the initialization phase and before the (numeric) factorization phase.
378: if (flg_mumps || flg_mumps_ch) {
379: PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_mumps_ooc_api", &test_mumps_ooc_api, NULL));
380: if (test_mumps_ooc_api) {
381: PetscCall(PetscStrncpy(tmpdir, "OOC_XXXXXX", sizeof(tmpdir)));
382: if (rank == 0) PetscCall(PetscMkdtemp(tmpdir));
383: PetscCallMPI(MPI_Bcast(tmpdir, 11, MPI_CHAR, 0, PETSC_COMM_WORLD));
384: PetscCall(MatMumpsSetIcntl(F, 22, 1)); // enable out-of-core
385: PetscCall(MatMumpsSetOocTmpDir(F, tmpdir));
386: }
387: }
388: #endif
390: /* Get info from matrix factors */
391: PetscCall(KSPSetUp(ksp));
393: #if defined(PETSC_HAVE_MUMPS)
394: if (flg_mumps || flg_mumps_ch) {
395: PetscInt icntl, infog34, num_null_pivots, *null_pivots;
396: PetscReal cntl, rinfo12, rinfo13;
397: icntl = 3;
398: PetscCall(MatMumpsGetCntl(F, icntl, &cntl));
400: /* compute determinant and check for any null pivots*/
401: if (rank == 0) {
402: PetscCall(MatMumpsGetInfog(F, 34, &infog34));
403: PetscCall(MatMumpsGetRinfog(F, 12, &rinfo12));
404: PetscCall(MatMumpsGetRinfog(F, 13, &rinfo13));
405: PetscCall(MatMumpsGetNullPivots(F, &num_null_pivots, &null_pivots));
406: PetscCall(PetscPrintf(PETSC_COMM_SELF, " Mumps row pivot threshold = %g\n", (double)cntl));
407: PetscCall(PetscPrintf(PETSC_COMM_SELF, " Mumps determinant = (%g, %g) * 2^%" PetscInt_FMT " \n", (double)rinfo12, (double)rinfo13, infog34));
408: if (num_null_pivots > 0) {
409: PetscCall(PetscPrintf(PETSC_COMM_SELF, " Mumps num of null pivots detected = %" PetscInt_FMT "\n", num_null_pivots));
410: PetscCall(PetscSortInt(num_null_pivots, null_pivots)); /* just make the printf deterministic */
411: for (j = 0; j < num_null_pivots; j++) PetscCall(PetscPrintf(PETSC_COMM_SELF, " Mumps row with null pivots is = %" PetscInt_FMT "\n", null_pivots[j]));
412: }
413: PetscCall(PetscFree(null_pivots));
414: }
415: }
416: #endif
418: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
419: Solve the linear system
420: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
421: PetscCall(KSPSolve(ksp, b, x));
423: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
424: Check solution and clean up
425: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
426: PetscCall(VecAXPY(x, -1.0, u));
427: PetscCall(VecNorm(x, NORM_2, &norm));
428: PetscCall(KSPGetIterationNumber(ksp, &its));
430: /*
431: Print convergence information. PetscPrintf() produces a single
432: print statement from all processes that share a communicator.
433: An alternative is PetscFPrintf(), which prints to a file.
434: */
435: if (norm < 1.e-12) {
436: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of error < 1.e-12 iterations %" PetscInt_FMT "\n", its));
437: } else {
438: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of error %g iterations %" PetscInt_FMT "\n", (double)norm, its));
439: }
441: #if defined(PETSC_HAVE_MUMPS)
442: // Get the OCC tmpdir via F before it is destroyed in KSPDestroy().
443: if (test_mumps_ooc_api) {
444: const char *dir;
446: PetscCall(MatMumpsGetOocTmpDir(F, &dir));
447: PetscCall(PetscPrintf(PETSC_COMM_WORLD, " OOC_TMPDIR = %s\n", dir));
448: }
449: #endif
451: /*
452: Free work space. All PETSc objects should be destroyed when they
453: are no longer needed.
454: */
455: PetscCall(KSPDestroy(&ksp));
457: #if defined(PETSC_HAVE_MUMPS)
458: // Files created by MUMPS under the OOC tmpdir were automatically deleted by MUMPS with JOB_END (-2)
459: // during KSP/PCDestroy(), but we need to remove the tmpdir created by us.
460: if (test_mumps_ooc_api && rank == 0) PetscCall(PetscRMTree(tmpdir));
461: #endif
463: PetscCall(VecDestroy(&u));
464: PetscCall(VecDestroy(&x));
465: PetscCall(VecDestroy(&b));
466: PetscCall(MatDestroy(&A));
468: /*
469: Always call PetscFinalize() before exiting a program. This routine
470: - finalizes the PETSc libraries as well as MPI
471: - provides summary and diagnostic information if certain runtime
472: options are chosen (e.g., -log_view).
473: */
474: PetscCall(PetscFinalize());
475: return 0;
476: }
478: /*TEST
480: test:
481: args: -use_petsc_lu
482: output_file: output/ex52_2.out
484: testset:
485: suffix: mumps
486: nsize: 3
487: requires: mumps
488: output_file: output/ex52_1.out
489: args: -use_mumps_lu
491: test:
492: suffix: mumps
494: test:
495: suffix: mumps_ooc
496: args: -mat_mumps_icntl_22 1 -mat_mumps_ooc_tmpdir /tmp
498: test:
499: suffix: mumps_ooc_api
500: args: -mat_mumps_icntl_22 1 -test_mumps_ooc_api
501: filter: grep -v "OOC_TMPDIR"
503: test:
504: suffix: mumps_2
505: nsize: 3
506: requires: mumps
507: args: -use_mumps_ch
508: output_file: output/ex52_1.out
510: test:
511: suffix: mumps_3
512: nsize: 3
513: requires: mumps
514: args: -use_mumps_ch -mat_type sbaij
515: output_file: output/ex52_1.out
517: test:
518: suffix: mumps_4
519: nsize: 3
520: requires: mumps !complex !single
521: args: -use_mumps_lu -m 50 -n 50 -print_mumps_memory
522: output_file: output/ex52_4.out
524: test:
525: suffix: mumps_5
526: nsize: 3
527: requires: mumps !complex !single
528: args: -use_mumps_lu -m 50 -n 50 -zero_first_and_last_rows
529: output_file: output/ex52_5.out
531: test:
532: suffix: mumps_omp_2
533: nsize: 4
534: requires: mumps hwloc openmp pthread defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
535: args: -use_mumps_lu -mat_mumps_use_omp_threads 2
536: output_file: output/ex52_1.out
538: test:
539: suffix: mumps_omp_3
540: nsize: 4
541: requires: mumps hwloc openmp pthread defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
542: args: -use_mumps_ch -mat_mumps_use_omp_threads 3 -mat_mumps_icntl_48 0
543: # Ignore the warning since we are intentionally testing the imbalanced case
544: filter: grep -v "Warning: number of OpenMP threads"
545: output_file: output/ex52_1.out
547: test:
548: suffix: mumps_omp_4
549: nsize: 4
550: requires: mumps hwloc openmp pthread defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
551: # let PETSc guess a proper number for threads
552: args: -use_mumps_ch -mat_type sbaij -mat_mumps_use_omp_threads
553: output_file: output/ex52_1.out
555: testset:
556: suffix: strumpack_2
557: nsize: {{1 2}}
558: requires: strumpack
559: args: -use_strumpack_lu
560: output_file: output/ex52_3.out
562: test:
563: suffix: aij
564: args: -mat_type aij
566: test:
567: requires: kokkos_kernels
568: suffix: kok
569: args: -mat_type aijkokkos
571: test:
572: requires: cuda
573: suffix: cuda
574: args: -mat_type aijcusparse
576: test:
577: requires: hip
578: suffix: hip
579: args: -mat_type aijhipsparse
581: test:
582: suffix: strumpack_ilu
583: nsize: {{1 2}}
584: requires: strumpack
585: args: -use_strumpack_ilu
586: output_file: output/ex52_3.out
588: testset:
589: suffix: superlu_dist
590: nsize: {{1 2}}
591: requires: superlu superlu_dist
592: args: -use_superlu_lu
593: output_file: output/ex52_2.out
595: test:
596: suffix: aij
597: args: -mat_type aij
599: test:
600: requires: kokkos_kernels
601: suffix: kok
602: args: -mat_type aijkokkos
604: test:
605: requires: cuda
606: suffix: cuda
607: args: -mat_type aijcusparse
609: test:
610: requires: hip
611: suffix: hip
612: args: -mat_type aijhipsparse
614: test:
615: suffix: superlu_ilu
616: requires: superlu superlu_dist
617: args: -use_superlu_ilu
618: output_file: output/ex52_2.out
620: TEST*/