Actual source code: ex27.c
2: static char help[] = "Reads a PETSc matrix and vector from a file and solves the normal equations.\n\n";
4: /*
5: Include "petscksp.h" so that we can use KSP solvers. Note that this file
6: automatically includes:
7: petscsys.h - base PETSc routines petscvec.h - vectors
8: petscmat.h - matrices
9: petscis.h - index sets petscksp.h - Krylov subspace methods
10: petscviewer.h - viewers petscpc.h - preconditioners
11: */
12: #include <petscksp.h>
13: #include <petscviewerhdf5.h>
15: static PetscErrorCode VecLoadIfExists_Private(Vec b, PetscViewer fd, PetscBool *has)
16: {
17: PetscBool hdf5 = PETSC_FALSE;
19: PetscFunctionBeginUser;
20: PetscCall(PetscObjectTypeCompare((PetscObject)fd, PETSCVIEWERHDF5, &hdf5));
21: if (hdf5) {
22: #if defined(PETSC_HAVE_HDF5)
23: PetscCall(PetscViewerHDF5HasObject(fd, (PetscObject)b, has));
24: if (*has) PetscCall(VecLoad(b, fd));
25: #else
26: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "PETSc must be configured with HDF5 to use this feature");
27: #endif
28: } else {
29: PetscErrorCode ierrp;
30: PetscCall(PetscPushErrorHandler(PetscReturnErrorHandler, NULL));
31: ierrp = VecLoad(b, fd);
32: PetscCall(PetscPopErrorHandler());
33: *has = ierrp ? PETSC_FALSE : PETSC_TRUE;
34: }
35: PetscFunctionReturn(PETSC_SUCCESS);
36: }
38: int main(int argc, char **args)
39: {
40: KSP ksp; /* linear solver context */
41: Mat A, N; /* matrix */
42: Vec x, b, r, Ab, v[2]; /* approx solution, RHS, residual */
43: PetscViewer fd; /* viewer */
44: char file[PETSC_MAX_PATH_LEN] = ""; /* input file name */
45: char file_x0[PETSC_MAX_PATH_LEN] = ""; /* name of input file with initial guess */
46: char A_name[128] = "A", b_name[128] = "b", x0_name[128] = "x0"; /* name of the matrix, RHS and initial guess */
47: KSPType ksptype;
48: PetscBool has;
49: PetscInt its, n, m;
50: PetscReal norm;
51: PetscBool nonzero_guess = PETSC_TRUE;
52: PetscBool solve_normal = PETSC_FALSE;
53: PetscBool solve_augmented = PETSC_FALSE;
54: PetscBool truncate = PETSC_FALSE;
55: PetscBool explicit_transpose = PETSC_FALSE;
56: PetscBool hdf5 = PETSC_FALSE;
57: PetscBool test_custom_layout = PETSC_FALSE;
58: PetscMPIInt rank, size;
60: PetscFunctionBeginUser;
61: PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
62: PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
63: PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
64: /*
65: Determine files from which we read the linear system
66: (matrix, right-hand-side and initial guess vector).
67: */
68: PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file, sizeof(file), NULL));
69: PetscCall(PetscOptionsGetBool(NULL, NULL, "-truncate", &truncate, NULL));
70: if (!truncate) PetscCall(PetscOptionsGetString(NULL, NULL, "-f_x0", file_x0, sizeof(file_x0), NULL));
71: PetscCall(PetscOptionsGetString(NULL, NULL, "-A_name", A_name, sizeof(A_name), NULL));
72: PetscCall(PetscOptionsGetString(NULL, NULL, "-b_name", b_name, sizeof(b_name), NULL));
73: PetscCall(PetscOptionsGetString(NULL, NULL, "-x0_name", x0_name, sizeof(x0_name), NULL));
74: /*
75: Decide whether to solve the original system (-solve_normal 0)
76: or the normal equation (-solve_normal 1).
77: */
78: PetscCall(PetscOptionsGetBool(NULL, NULL, "-solve_normal", &solve_normal, NULL));
79: if (!solve_normal) PetscCall(PetscOptionsGetBool(NULL, NULL, "-solve_augmented", &solve_augmented, NULL));
80: if (solve_augmented) PetscCall(PetscOptionsGetBool(NULL, NULL, "-explicit_transpose", &explicit_transpose, NULL));
81: /*
82: Decide whether to use the HDF5 reader.
83: */
84: PetscCall(PetscOptionsGetBool(NULL, NULL, "-hdf5", &hdf5, NULL));
85: /*
86: Decide whether custom matrix layout will be tested.
87: */
88: PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_custom_layout", &test_custom_layout, NULL));
90: /* -----------------------------------------------------------
91: Beginning of linear solver loop
92: ----------------------------------------------------------- */
93: /*
94: Loop through the linear solve 2 times.
95: - The intention here is to preload and solve a small system;
96: then load another (larger) system and solve it as well.
97: This process preloads the instructions with the smaller
98: system so that more accurate performance monitoring (via
99: -log_view) can be done with the larger one (that actually
100: is the system of interest).
101: */
102: PetscPreLoadBegin(PETSC_FALSE, "Load system");
104: /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
105: Load system
106: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108: /*
109: Open binary file. Note that we use FILE_MODE_READ to indicate
110: reading from this file.
111: */
112: if (hdf5) {
113: #if defined(PETSC_HAVE_HDF5)
114: PetscCall(PetscViewerHDF5Open(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd));
115: PetscCall(PetscViewerPushFormat(fd, PETSC_VIEWER_HDF5_MAT));
116: #else
117: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "PETSc must be configured with HDF5 to use this feature");
118: #endif
119: } else {
120: PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd));
121: }
123: /*
124: Load the matrix.
125: Matrix type is set automatically but you can override it by MatSetType() prior to MatLoad().
126: Do that only if you really insist on the given type.
127: */
128: PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
129: PetscCall(PetscObjectSetName((PetscObject)A, A_name));
130: PetscCall(MatSetFromOptions(A));
131: PetscCall(MatLoad(A, fd));
132: if (truncate) {
133: Mat P, B;
134: PetscInt M, N;
135: PetscCall(MatGetLocalSize(A, &m, &n));
136: PetscCall(MatGetSize(A, &M, &N));
137: PetscCall(MatCreateAIJ(PETSC_COMM_WORLD, m, PETSC_DECIDE, M, N / 1.5, 1, NULL, 1, NULL, &P));
138: PetscCall(MatGetOwnershipRangeColumn(P, &m, &n));
139: for (; m < n; ++m) PetscCall(MatSetValue(P, m, m, 1.0, INSERT_VALUES));
140: PetscCall(MatAssemblyBegin(P, MAT_FINAL_ASSEMBLY));
141: PetscCall(MatAssemblyEnd(P, MAT_FINAL_ASSEMBLY));
142: PetscCall(MatShift(P, 1.0));
143: PetscCall(MatMatMult(A, P, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &B));
144: PetscCall(MatDestroy(&P));
145: PetscCall(MatDestroy(&A));
146: A = B;
147: }
148: if (test_custom_layout && size > 1) {
149: /* Perturb the local sizes and create the matrix anew */
150: PetscInt m1, n1;
151: PetscCall(MatGetLocalSize(A, &m, &n));
152: m = rank ? m - 1 : m + size - 1;
153: n = (rank == size - 1) ? n + size - 1 : n - 1;
154: PetscCall(MatDestroy(&A));
155: PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
156: PetscCall(PetscObjectSetName((PetscObject)A, A_name));
157: PetscCall(MatSetSizes(A, m, n, PETSC_DECIDE, PETSC_DECIDE));
158: PetscCall(MatSetFromOptions(A));
159: PetscCall(MatLoad(A, fd));
160: PetscCall(MatGetLocalSize(A, &m1, &n1));
161: PetscCheck(m1 == m && n1 == n, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "resulting sizes differ from requested ones: %" PetscInt_FMT " %" PetscInt_FMT " != %" PetscInt_FMT " %" PetscInt_FMT, m1, n1, m, n);
162: }
163: PetscCall(MatGetLocalSize(A, &m, &n));
165: /*
166: Load the RHS vector if it is present in the file, otherwise use a vector of all ones.
167: */
168: PetscCall(MatCreateVecs(A, &x, &b));
169: PetscCall(PetscObjectSetName((PetscObject)b, b_name));
170: PetscCall(VecSetFromOptions(b));
171: PetscCall(VecLoadIfExists_Private(b, fd, &has));
172: if (!has) {
173: PetscScalar one = 1.0;
174: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Failed to load RHS, so use a vector of all ones.\n"));
175: PetscCall(VecSetFromOptions(b));
176: PetscCall(VecSet(b, one));
177: }
179: /*
180: Load the initial guess vector if it is present in the file, otherwise use a vector of all zeros.
181: */
182: PetscCall(PetscObjectSetName((PetscObject)x, x0_name));
183: PetscCall(VecSetFromOptions(x));
184: if (!truncate) {
185: /* load file_x0 if it is specified, otherwise try to reuse file */
186: if (file_x0[0]) {
187: PetscCall(PetscViewerDestroy(&fd));
188: if (hdf5) {
189: #if defined(PETSC_HAVE_HDF5)
190: PetscCall(PetscViewerHDF5Open(PETSC_COMM_WORLD, file_x0, FILE_MODE_READ, &fd));
191: #endif
192: } else {
193: PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file_x0, FILE_MODE_READ, &fd));
194: }
195: }
196: PetscCall(VecLoadIfExists_Private(x, fd, &has));
197: } else has = PETSC_FALSE;
198: if (truncate || !has) {
199: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Failed to load initial guess, so use a vector of all zeros.\n"));
200: PetscCall(VecSet(x, 0.0));
201: nonzero_guess = PETSC_FALSE;
202: }
203: PetscCall(PetscViewerDestroy(&fd));
205: PetscCall(VecDuplicate(x, &Ab));
207: /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
208: Setup solve for system
209: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
211: /*
212: Conclude profiling last stage; begin profiling next stage.
213: */
214: PetscPreLoadStage("KSPSetUp");
216: PetscCall(MatCreateNormalHermitian(A, &N));
217: PetscCall(MatMultHermitianTranspose(A, b, Ab));
219: /*
220: Create linear solver; set operators; set runtime options.
221: */
222: PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
224: if (solve_normal) {
225: PetscCall(KSPSetOperators(ksp, N, N));
226: } else if (solve_augmented) {
227: Mat array[4], C;
228: Vec view;
229: PetscInt M;
231: PetscCall(MatDestroy(&N));
232: PetscCall(MatGetSize(A, &M, NULL));
233: PetscCall(MatCreateConstantDiagonal(PETSC_COMM_WORLD, m, m, M, M, -1.0, array));
234: array[1] = A;
235: if (!explicit_transpose) PetscCall(MatCreateHermitianTranspose(A, array + 2));
236: else PetscCall(MatHermitianTranspose(A, MAT_INITIAL_MATRIX, array + 2));
237: array[3] = NULL;
238: PetscCall(MatCreateNest(PETSC_COMM_WORLD, 2, NULL, 2, NULL, array, &C));
239: PetscCall(MatNestSetVecType(C, VECNEST));
240: PetscCall(MatCreateVecs(C, v + 1, v));
241: PetscCall(VecSet(v[0], 0.0));
242: PetscCall(VecSet(v[1], 0.0));
243: PetscCall(VecNestGetSubVec(v[0], 0, &view));
244: PetscCall(VecCopy(b, view));
245: PetscCall(VecNestGetSubVec(v[1], 1, &view));
246: PetscCall(VecCopy(x, view));
247: PetscCall(KSPSetOperators(ksp, C, C));
248: PetscCall(MatDestroy(&C));
249: PetscCall(MatDestroy(array));
250: PetscCall(MatDestroy(array + 2));
251: } else {
252: PC pc;
253: PetscCall(KSPSetType(ksp, KSPLSQR));
254: PetscCall(KSPGetPC(ksp, &pc));
255: PetscCall(PCSetType(pc, PCNONE));
256: PetscCall(KSPSetOperators(ksp, A, N));
257: }
258: PetscCall(KSPSetInitialGuessNonzero(ksp, nonzero_guess));
259: PetscCall(KSPSetFromOptions(ksp));
261: /*
262: Here we explicitly call KSPSetUp() and KSPSetUpOnBlocks() to
263: enable more precise profiling of setting up the preconditioner.
264: These calls are optional, since both will be called within
265: KSPSolve() if they haven't been called already.
266: */
267: PetscCall(KSPSetUp(ksp));
268: PetscCall(KSPSetUpOnBlocks(ksp));
270: /*
271: Solve system
272: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
274: /*
275: Begin profiling next stage
276: */
277: PetscPreLoadStage("KSPSolve");
279: /*
280: Solve linear system
281: */
282: if (solve_normal) {
283: PetscCall(KSPSolve(ksp, Ab, x));
284: } else if (solve_augmented) {
285: Vec view;
287: PetscCall(KSPSolve(ksp, v[0], v[1]));
288: PetscCall(VecNestGetSubVec(v[1], 1, &view));
289: PetscCall(VecCopy(view, x));
290: } else {
291: PetscCall(KSPSolve(ksp, b, x));
292: }
293: PetscCall(PetscObjectSetName((PetscObject)x, "x"));
295: /*
296: Conclude profiling this stage
297: */
298: PetscPreLoadStage("Cleanup");
300: /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
301: Check error, print output, free data structures.
302: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
304: /*
305: Check error
306: */
307: PetscCall(VecDuplicate(b, &r));
308: PetscCall(MatMult(A, x, r));
309: PetscCall(VecAXPY(r, -1.0, b));
310: PetscCall(VecNorm(r, NORM_2, &norm));
311: PetscCall(KSPGetIterationNumber(ksp, &its));
312: PetscCall(KSPGetType(ksp, &ksptype));
313: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP type: %s\n", ksptype));
314: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of iterations = %3" PetscInt_FMT "\n", its));
315: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Residual norm %g\n", (double)norm));
317: /*
318: Free work space. All PETSc objects should be destroyed when they
319: are no longer needed.
320: */
321: PetscCall(MatDestroy(&A));
322: PetscCall(VecDestroy(&b));
323: PetscCall(MatDestroy(&N));
324: PetscCall(VecDestroy(&Ab));
325: PetscCall(VecDestroy(&r));
326: PetscCall(VecDestroy(&x));
327: if (solve_augmented) {
328: PetscCall(VecDestroy(v));
329: PetscCall(VecDestroy(v + 1));
330: }
331: PetscCall(KSPDestroy(&ksp));
332: PetscPreLoadEnd();
333: /* -----------------------------------------------------------
334: End of linear solver loop
335: ----------------------------------------------------------- */
337: PetscCall(PetscFinalize());
338: return 0;
339: }
341: /*TEST
343: test:
344: suffix: 1
345: requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES)
346: args: -f ${DATAFILESPATH}/matrices/medium -ksp_view -ksp_monitor_short -ksp_max_it 100 -solve_normal
348: test:
349: suffix: 2
350: nsize: 2
351: requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES)
352: args: -f ${DATAFILESPATH}/matrices/shallow_water1 -ksp_view -ksp_monitor_short -ksp_max_it 100 -solve_normal -pc_type none
354: # Test handling failing VecLoad without abort
355: testset:
356: requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
357: args: -ksp_type cg -ksp_view -ksp_converged_reason -ksp_monitor_short -ksp_max_it 10
358: test:
359: suffix: 3
360: nsize: {{1 2}separate output}
361: args: -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/tiny_system
362: args: -f_x0 ${wPETSC_DIR}/share/petsc/datafiles/matrices/tiny_system_x0
363: test:
364: suffix: 3a
365: nsize: {{1 2}separate output}
366: args: -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/tiny_system
367: args: -f_x0 NONEXISTING_FILE
368: test:
369: suffix: 3b
370: nsize: {{1 2}separate output}
371: args: -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/tiny_system_with_x0 # this file includes all A, b and x0
372: test:
373: # Load square matrix, RHS and initial guess from HDF5 (Version 7.3 MAT-File)
374: suffix: 3b_hdf5
375: requires: hdf5 defined(PETSC_HDF5_HAVE_ZLIB)
376: nsize: {{1 2}separate output}
377: args: -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/tiny_system_with_x0.mat -hdf5
379: # Test least-square algorithms
380: testset:
381: requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES)
382: args: -f ${DATAFILESPATH}/matrices/rectangular_ultrasound_4889x841
383: test:
384: suffix: 4
385: nsize: {{1 2 4}}
386: args: -ksp_converged_reason -ksp_monitor_short -ksp_rtol 1e-5 -ksp_max_it 100
387: args: -solve_normal -ksp_type cg
388: test:
389: suffix: 4a
390: nsize: {{1 2 4}}
391: args: -ksp_converged_reason -ksp_monitor_short -ksp_rtol 1e-5 -ksp_max_it 100
392: args: -ksp_type {{cgls lsqr}separate output}
393: test:
394: # Test KSPLSQR-specific options
395: suffix: 4b
396: nsize: 2
397: args: -ksp_converged_reason -ksp_rtol 1e-3 -ksp_max_it 200 -ksp_view
398: args: -ksp_type lsqr -ksp_convergence_test lsqr -ksp_lsqr_monitor -ksp_lsqr_compute_standard_error -ksp_lsqr_exact_mat_norm {{0 1}separate output}
399: test:
400: suffix: 4c
401: nsize: 4
402: requires: hpddm slepc defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
403: filter: grep -v "shared subdomain KSP between SLEPc and PETSc" | grep -v "total: nonzeros="
404: args: -ksp_converged_reason -ksp_rtol 1e-5 -ksp_max_it 100 -ksp_view
405: args: -ksp_type lsqr -pc_type hpddm -pc_hpddm_define_subdomains -pc_hpddm_levels_1_eps_nev 20 -pc_hpddm_levels_1_st_share_sub_ksp {{false true}shared output}
406: args: -pc_hpddm_levels_1_pc_asm_sub_mat_type aij -pc_hpddm_levels_1_pc_asm_type basic -pc_hpddm_levels_1_sub_pc_type cholesky
407: test:
408: suffix: 4d
409: nsize: 4
410: requires: hpddm slepc suitesparse defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
411: filter: grep -v "shared subdomain KSP between SLEPc and PETSc"
412: args: -ksp_converged_reason -ksp_rtol 1e-5 -ksp_max_it 100 -ksp_view
413: args: -ksp_type lsqr -pc_type hpddm -pc_hpddm_define_subdomains -pc_hpddm_levels_1_eps_nev 20 -pc_hpddm_levels_1_st_share_sub_ksp {{false true}shared output} -pc_hpddm_levels_1_st_pc_type qr
414: args: -pc_hpddm_levels_1_pc_asm_sub_mat_type normalh -pc_hpddm_levels_1_pc_asm_type basic -pc_hpddm_levels_1_sub_pc_type qr
415: test:
416: suffix: 4e
417: nsize: 4
418: requires: hpddm slepc defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
419: args: -solve_augmented -ksp_type gmres
420: args: -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_precondition self -fieldsplit_0_pc_type jacobi -fieldsplit_ksp_type preonly
421: args: -prefix_push fieldsplit_1_ -pc_type hpddm -pc_hpddm_define_subdomains -pc_hpddm_levels_1_eps_nev 20 -pc_hpddm_levels_1_st_share_sub_ksp -pc_hpddm_levels_1_sub_pc_type cholesky -prefix_pop -fieldsplit_1_mat_schur_complement_ainv_type {{diag lump}shared output}
422: test:
423: suffix: 4f
424: nsize: 4
425: requires: hpddm slepc suitesparse defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
426: filter: sed -e "s/(1,0) : type=mpiaij/(1,0) : type=transpose/g" -e "s/hermitiantranspose/transpose/g"
427: args: -solve_augmented -ksp_type gmres -ksp_view -explicit_transpose {{false true}shared output}
428: args: -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_precondition self -fieldsplit_0_pc_type jacobi -fieldsplit_ksp_type preonly
429: args: -prefix_push fieldsplit_1_ -pc_type hpddm -pc_hpddm_define_subdomains -pc_hpddm_levels_1_eps_nev 20 -pc_hpddm_levels_1_st_share_sub_ksp -pc_hpddm_levels_1_sub_pc_type qr -prefix_pop
430: test:
431: suffix: 4g
432: nsize: 4
433: requires: hypre
434: args: -ksp_converged_reason -ksp_monitor_short -ksp_rtol 1e-5 -ksp_max_it 100
435: args: -ksp_type lsqr -pc_type hypre
437: test:
438: # Load rectangular matrix from HDF5 (Version 7.3 MAT-File)
439: suffix: 4a_lsqr_hdf5
440: nsize: {{1 2 4 8}}
441: requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB)
442: args: -f ${DATAFILESPATH}/matrices/matlab/rectangular_ultrasound_4889x841.mat -hdf5
443: args: -ksp_converged_reason -ksp_monitor_short -ksp_rtol 1e-5 -ksp_max_it 100
444: args: -ksp_type lsqr
445: args: -test_custom_layout {{0 1}}
447: # Test for correct cgls convergence reason
448: test:
449: suffix: 5
450: nsize: 1
451: requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES)
452: args: -f ${DATAFILESPATH}/matrices/rectangular_ultrasound_4889x841
453: args: -ksp_converged_reason -ksp_rtol 1e-2 -ksp_max_it 100
454: args: -ksp_type cgls
456: # Load a matrix, RHS and solution from HDF5 (Version 7.3 MAT-File). Test immediate convergence.
457: testset:
458: nsize: {{1 2 4 8}}
459: requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB)
460: args: -ksp_converged_reason -ksp_monitor_short -ksp_rtol 1e-5 -ksp_max_it 10
461: args: -ksp_type lsqr
462: args: -test_custom_layout {{0 1}}
463: args: -hdf5 -x0_name x
464: test:
465: suffix: 6_hdf5
466: args: -f ${DATAFILESPATH}/matrices/matlab/small.mat
467: test:
468: suffix: 6_hdf5_rect
469: args: -f ${DATAFILESPATH}/matrices/matlab/small_rect.mat
470: test:
471: suffix: 6_hdf5_dense
472: args: -f ${DATAFILESPATH}/matrices/matlab/small_dense.mat -mat_type dense
473: test:
474: suffix: 6_hdf5_rect_dense
475: args: -f ${DATAFILESPATH}/matrices/matlab/small_rect_dense.mat -mat_type dense
477: # Test correct handling of local dimensions in PCApply
478: testset:
479: requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES)
480: requires: hdf5 defined(PETSC_HDF5_HAVE_ZLIB)
481: nsize: 3
482: suffix: 7
483: args: -f ${DATAFILESPATH}/matrices/matlab/small.mat -hdf5 -test_custom_layout 1 -ksp_type lsqr -pc_type jacobi
485: # Test complex matrices
486: testset:
487: requires: double complex !defined(PETSC_USE_64BIT_INDICES)
488: args: -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/nh-complex-int32-float64
489: output_file: output/ex27_8.out
490: filter: grep -v "KSP type"
491: test:
492: suffix: 8
493: args: -solve_normal 0 -ksp_type {{lsqr cgls}}
494: test:
495: suffix: 8_normal
496: args: -solve_normal 1 -ksp_type {{cg bicg}}
498: testset:
499: requires: double suitesparse !defined(PETSC_USE_64BIT_INDICES)
500: args: -solve_normal {{0 1}shared output} -pc_type qr
501: output_file: output/ex27_9.out
502: filter: grep -v "KSP type"
503: test:
504: suffix: 9_real
505: requires: !complex
506: args: -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/ns-real-int32-float64
507: test:
508: suffix: 9_complex
509: requires: complex
510: args: -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/nh-complex-int32-float64
512: test:
513: suffix: 10
514: requires: !complex double suitesparse !defined(PETSC_USE_64BIT_INDICES)
515: nsize: 2
516: args: -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/ns-real-int32-float64 -pc_type bjacobi -sub_pc_type qr
518: test:
519: suffix: 11
520: nsize: 4
521: requires: datafilespath double complex !defined(PETSC_USE_64BIT_INDICES) hpddm slepc defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
522: args: -f ${DATAFILESPATH}/matrices/farzad_B_rhs -truncate
523: args: -ksp_converged_reason -ksp_rtol 1e-5 -ksp_max_it 100
524: args: -ksp_type lsqr -pc_type hpddm -pc_hpddm_define_subdomains -pc_hpddm_levels_1_eps_nev 20 -pc_hpddm_levels_1_eps_threshold 1e-6
525: args: -pc_hpddm_levels_1_pc_asm_sub_mat_type aij -pc_hpddm_levels_1_pc_asm_type basic -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_coarse_pc_type lu
527: TEST*/