Actual source code: ex76.c
1: #include <petscksp.h>
2: #include <petsc/private/petscimpl.h>
4: static char help[] = "Solves a linear system using PCHPDDM.\n\n";
6: int main(int argc, char **args)
7: {
8: Vec b; /* computed solution and RHS */
9: Mat A, aux, X, B; /* linear system matrix */
10: KSP ksp; /* linear solver context */
11: PC pc;
12: IS is, sizes;
13: const PetscInt *idx;
14: PetscMPIInt rank, size;
15: PetscInt m, N = 1;
16: PetscLayout map;
17: PetscViewer viewer;
18: char dir[PETSC_MAX_PATH_LEN], name[PETSC_MAX_PATH_LEN], type[256];
19: PetscBool3 share = PETSC_BOOL3_UNKNOWN;
20: PetscBool flg, set, transpose = PETSC_FALSE;
22: PetscFunctionBeginUser;
23: PetscCall(PetscInitialize(&argc, &args, NULL, help));
24: PetscCall(PetscLogDefaultBegin());
25: PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
26: PetscCheck(size == 4, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This example requires 4 processes");
27: PetscCall(PetscOptionsGetInt(NULL, NULL, "-rhs", &N, NULL));
28: PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
29: PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
30: PetscCall(PetscStrncpy(dir, ".", sizeof(dir)));
31: PetscCall(PetscOptionsGetString(NULL, NULL, "-load_dir", dir, sizeof(dir), NULL));
32: /* loading matrices */
33: PetscCall(PetscSNPrintf(name, sizeof(name), "%s/sizes_%d.dat", dir, size));
34: PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, name, FILE_MODE_READ, &viewer));
35: PetscCall(ISCreate(PETSC_COMM_WORLD, &sizes));
36: PetscCall(ISLoad(sizes, viewer));
37: PetscCall(ISGetIndices(sizes, &idx));
38: PetscCall(MatSetSizes(A, idx[0], idx[1], idx[2], idx[3]));
39: PetscCall(MatCreate(PETSC_COMM_WORLD, &X));
40: PetscCall(MatSetSizes(X, idx[4], idx[4], PETSC_DETERMINE, PETSC_DETERMINE));
41: PetscCall(MatSetUp(X));
42: PetscCall(ISRestoreIndices(sizes, &idx));
43: PetscCall(ISDestroy(&sizes));
44: PetscCall(PetscViewerDestroy(&viewer));
45: PetscCall(PetscSNPrintf(name, sizeof(name), "%s/A.dat", dir));
46: PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, name, FILE_MODE_READ, &viewer));
47: PetscCall(MatLoad(A, viewer));
48: PetscCall(PetscViewerDestroy(&viewer));
49: PetscCall(PetscSNPrintf(name, sizeof(name), "%s/is_%d.dat", dir, size));
50: PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, name, FILE_MODE_READ, &viewer));
51: PetscCall(ISCreate(PETSC_COMM_WORLD, &sizes));
52: PetscCall(MatGetLayouts(X, &map, NULL));
53: PetscCall(ISSetLayout(sizes, map));
54: PetscCall(ISLoad(sizes, viewer));
55: PetscCall(ISGetLocalSize(sizes, &m));
56: PetscCall(ISGetIndices(sizes, &idx));
57: PetscCall(ISCreateGeneral(PETSC_COMM_SELF, m, idx, PETSC_COPY_VALUES, &is));
58: PetscCall(ISRestoreIndices(sizes, &idx));
59: PetscCall(ISDestroy(&sizes));
60: PetscCall(MatGetBlockSize(A, &m));
61: PetscCall(ISSetBlockSize(is, m));
62: PetscCall(PetscViewerDestroy(&viewer));
63: PetscCall(PetscSNPrintf(name, sizeof(name), "%s/Neumann_%d.dat", dir, size));
64: PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, name, FILE_MODE_READ, &viewer));
65: PetscCall(MatLoad(X, viewer));
66: PetscCall(PetscViewerDestroy(&viewer));
67: PetscCall(MatGetDiagonalBlock(X, &B));
68: PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &aux));
69: PetscCall(MatDestroy(&X));
70: flg = PETSC_FALSE;
71: PetscCall(PetscOptionsGetBool(NULL, NULL, "-pc_hpddm_levels_1_st_share_sub_ksp", &flg, &set));
72: if (flg) { /* PETSc LU/Cholesky is struggling numerically for bs > 1 */
73: /* only set the proper bs for the geneo_share_* tests, 1 otherwise */
74: PetscCall(MatSetBlockSizesFromMats(aux, A, A));
75: share = PETSC_BOOL3_TRUE;
76: } else if (set) share = PETSC_BOOL3_FALSE;
77: PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE));
78: PetscCall(MatSetOption(aux, MAT_SYMMETRIC, PETSC_TRUE));
79: /* ready for testing */
80: PetscOptionsBegin(PETSC_COMM_WORLD, "", "", "");
81: PetscCall(PetscStrncpy(type, MATAIJ, sizeof(type)));
82: PetscCall(PetscOptionsFList("-mat_type", "Matrix type", "MatSetType", MatList, type, type, 256, &flg));
83: PetscOptionsEnd();
84: PetscCall(MatConvert(A, type, MAT_INPLACE_MATRIX, &A));
85: PetscCall(MatConvert(aux, type, MAT_INPLACE_MATRIX, &aux));
86: PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
87: PetscCall(KSPSetOperators(ksp, A, A));
88: PetscCall(KSPGetPC(ksp, &pc));
89: PetscCall(PCSetType(pc, PCHPDDM));
90: #if defined(PETSC_HAVE_HPDDM) && defined(PETSC_HAVE_DYNAMIC_LIBRARIES) && defined(PETSC_USE_SHARED_LIBRARIES)
91: flg = PETSC_FALSE;
92: PetscCall(PetscOptionsGetBool(NULL, NULL, "-reset", &flg, NULL));
93: if (flg) {
94: PetscCall(PetscOptionsSetValue(NULL, "-pc_hpddm_block_splitting", "true"));
95: PetscCall(PCSetFromOptions(pc));
96: PetscCall(PCSetUp(pc));
97: PetscCall(PetscOptionsClearValue(NULL, "-pc_hpddm_block_splitting"));
98: }
99: PetscCall(PCHPDDMSetAuxiliaryMat(pc, is, aux, NULL, NULL));
100: PetscCall(PCHPDDMHasNeumannMat(pc, PETSC_FALSE)); /* PETSC_TRUE is fine as well, just testing */
101: if (share == PETSC_BOOL3_UNKNOWN) PetscCall(PCHPDDMSetSTShareSubKSP(pc, PetscBool3ToBool(share)));
102: flg = PETSC_FALSE;
103: PetscCall(PetscOptionsGetBool(NULL, NULL, "-set_rhs", &flg, NULL));
104: if (flg) { /* user-provided RHS for concurrent generalized eigenvalue problems */
105: Mat a, c, P; /* usually assembled automatically in PCHPDDM, this is solely for testing PCHPDDMSetRHSMat() */
106: PetscInt rstart, rend, location;
108: PetscCall(MatDuplicate(aux, MAT_DO_NOT_COPY_VALUES, &B)); /* duplicate so that MatStructure is SAME_NONZERO_PATTERN */
109: PetscCall(MatGetDiagonalBlock(A, &a));
110: PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
111: PetscCall(ISGetLocalSize(is, &m));
112: PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, rend - rstart, m, 1, NULL, &P));
113: for (m = rstart; m < rend; ++m) {
114: PetscCall(ISLocate(is, m, &location));
115: PetscCheck(location >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "IS of the auxiliary Mat does not include all local rows of A");
116: PetscCall(MatSetValue(P, m - rstart, location, 1.0, INSERT_VALUES));
117: }
118: PetscCall(MatAssemblyBegin(P, MAT_FINAL_ASSEMBLY));
119: PetscCall(MatAssemblyEnd(P, MAT_FINAL_ASSEMBLY));
120: PetscCall(PetscObjectTypeCompare((PetscObject)a, MATSEQAIJ, &flg));
121: if (flg) PetscCall(MatPtAP(a, P, MAT_INITIAL_MATRIX, 1.0, &X)); // MatPtAP() is used to extend diagonal blocks with zeros on the overlap
122: else { // workaround for MatPtAP() limitations with some types
123: PetscCall(MatConvert(a, MATSEQAIJ, MAT_INITIAL_MATRIX, &c));
124: PetscCall(MatPtAP(c, P, MAT_INITIAL_MATRIX, 1.0, &X));
125: PetscCall(MatDestroy(&c));
126: }
127: PetscCall(MatDestroy(&P));
128: PetscCall(MatAXPY(B, 1.0, X, SUBSET_NONZERO_PATTERN));
129: PetscCall(MatDestroy(&X));
130: PetscCall(MatSetOption(B, MAT_SYMMETRIC, PETSC_TRUE));
131: PetscCall(PCHPDDMSetRHSMat(pc, B));
132: PetscCall(MatDestroy(&B));
133: }
134: #else
135: (void)share;
136: #endif
137: PetscCall(MatDestroy(&aux));
138: PetscCall(KSPSetFromOptions(ksp));
139: PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &flg));
140: if (flg) {
141: flg = PETSC_FALSE;
142: PetscCall(PetscOptionsGetBool(NULL, NULL, "-pc_hpddm_define_subdomains", &flg, NULL));
143: if (flg) {
144: IS rows;
146: PetscCall(MatGetOwnershipIS(A, &rows, NULL));
147: PetscCall(PCASMSetLocalSubdomains(pc, 1, &is, &rows));
148: PetscCall(ISDestroy(&rows));
149: }
150: }
151: PetscCall(ISDestroy(&is));
152: PetscCall(MatCreateVecs(A, NULL, &b));
153: PetscCall(VecSet(b, 1.0));
154: PetscCall(PetscOptionsGetBool(NULL, NULL, "-transpose", &transpose, NULL));
155: if (!transpose) PetscCall(KSPSolve(ksp, b, b));
156: else PetscCall(KSPSolveTranspose(ksp, b, b));
157: PetscCall(VecGetLocalSize(b, &m));
158: PetscCall(VecDestroy(&b));
159: if (N > 1) {
160: KSPType type;
162: PetscCall(PetscOptionsClearValue(NULL, "-ksp_converged_reason"));
163: PetscCall(KSPSetFromOptions(ksp));
164: PetscCall(MatCreateDense(PETSC_COMM_WORLD, m, PETSC_DECIDE, PETSC_DECIDE, N, NULL, &B));
165: PetscCall(MatCreateDense(PETSC_COMM_WORLD, m, PETSC_DECIDE, PETSC_DECIDE, N, NULL, &X));
166: PetscCall(MatSetRandom(B, NULL));
167: /* this is algorithmically optimal in the sense that blocks of vectors are coarsened or interpolated using matrix--matrix operations */
168: /* PCHPDDM however heavily relies on MPI[S]BAIJ format for which there is no efficient MatProduct implementation */
169: if (!transpose) PetscCall(KSPMatSolve(ksp, B, X));
170: else PetscCall(KSPMatSolveTranspose(ksp, B, X));
171: PetscCall(KSPGetType(ksp, &type));
172: PetscCall(PetscStrcmp(type, KSPHPDDM, &flg));
173: #if defined(PETSC_HAVE_HPDDM)
174: if (flg) {
175: PetscReal norm;
176: KSPHPDDMType type;
178: PetscCall(KSPHPDDMGetType(ksp, &type));
179: if (type == KSP_HPDDM_TYPE_PREONLY || type == KSP_HPDDM_TYPE_CG || type == KSP_HPDDM_TYPE_GMRES || type == KSP_HPDDM_TYPE_GCRODR) {
180: Mat C;
182: PetscCall(MatDuplicate(X, MAT_DO_NOT_COPY_VALUES, &C));
183: PetscCall(KSPSetMatSolveBatchSize(ksp, 1));
184: if (!transpose) PetscCall(KSPMatSolve(ksp, B, C));
185: else PetscCall(KSPMatSolveTranspose(ksp, B, C));
186: PetscCall(MatAYPX(C, -1.0, X, SAME_NONZERO_PATTERN));
187: PetscCall(MatNorm(C, NORM_INFINITY, &norm));
188: PetscCall(MatDestroy(&C));
189: PetscCheck(norm <= 100 * PETSC_MACHINE_EPSILON, PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "KSPMatSolve%s() and KSPSolve%s() difference has nonzero norm %g with pseudo-block KSPHPDDMType %s", (transpose ? "Transpose" : ""), (transpose ? "Transpose" : ""), (double)norm, KSPHPDDMTypes[type]);
190: }
191: }
192: #endif
193: PetscCall(MatDestroy(&X));
194: PetscCall(MatDestroy(&B));
195: }
196: PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCHPDDM, &flg));
197: #if defined(PETSC_HAVE_HPDDM) && defined(PETSC_HAVE_DYNAMIC_LIBRARIES) && defined(PETSC_USE_SHARED_LIBRARIES)
198: if (flg) PetscCall(PCHPDDMGetSTShareSubKSP(pc, &flg));
199: #endif
200: if (flg && PetscDefined(USE_LOG)) {
201: PetscCall(PetscOptionsHasName(NULL, NULL, "-pc_hpddm_harmonic_overlap", &flg));
202: if (!flg) {
203: PetscLogEvent event;
204: PetscEventPerfInfo info1, info2;
206: PetscCall(PetscLogEventRegister("MatLUFactorSym", PC_CLASSID, &event));
207: PetscCall(PetscLogEventGetPerfInfo(PETSC_DETERMINE, event, &info1));
208: PetscCall(PetscLogEventRegister("MatLUFactorNum", PC_CLASSID, &event));
209: PetscCall(PetscLogEventGetPerfInfo(PETSC_DETERMINE, event, &info2));
210: if (!info1.count && !info2.count) {
211: PetscCall(PetscLogEventRegister("MatCholFctrSym", PC_CLASSID, &event));
212: PetscCall(PetscLogEventGetPerfInfo(PETSC_DETERMINE, event, &info1));
213: PetscCall(PetscLogEventRegister("MatCholFctrNum", PC_CLASSID, &event));
214: PetscCall(PetscLogEventGetPerfInfo(PETSC_DETERMINE, event, &info2));
215: PetscCheck(info2.count > info1.count, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cholesky numerical factorization (%d) not called more times than Cholesky symbolic factorization (%d), broken -pc_hpddm_levels_1_st_share_sub_ksp", info2.count, info1.count);
216: } else PetscCheck(info2.count > info1.count, PETSC_COMM_SELF, PETSC_ERR_PLIB, "LU numerical factorization (%d) not called more times than LU symbolic factorization (%d), broken -pc_hpddm_levels_1_st_share_sub_ksp", info2.count, info1.count);
217: }
218: }
219: #if defined(PETSC_HAVE_HPDDM) && defined(PETSC_HAVE_DYNAMIC_LIBRARIES) && defined(PETSC_USE_SHARED_LIBRARIES)
220: if (N == 1) {
221: flg = PETSC_FALSE;
222: PetscCall(PetscOptionsGetBool(NULL, NULL, "-successive_solves", &flg, NULL));
223: if (flg) {
224: KSPConvergedReason reason[2];
225: PetscInt iterations[3];
227: PetscCall(KSPGetConvergedReason(ksp, reason));
228: PetscCall(KSPGetTotalIterations(ksp, iterations));
229: PetscCall(PetscOptionsClearValue(NULL, "-ksp_converged_reason"));
230: PetscCall(KSPSetFromOptions(ksp));
231: flg = PETSC_FALSE;
232: PetscCall(PetscOptionsGetBool(NULL, NULL, "-pc_hpddm_block_splitting", &flg, NULL));
233: if (!flg) {
234: PetscCall(PetscSNPrintf(name, sizeof(name), "%s/sizes_%d.dat", dir, size));
235: PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, name, FILE_MODE_READ, &viewer));
236: PetscCall(ISCreate(PETSC_COMM_WORLD, &sizes));
237: PetscCall(ISLoad(sizes, viewer));
238: PetscCall(ISGetIndices(sizes, &idx));
239: PetscCall(MatCreate(PETSC_COMM_WORLD, &X));
240: PetscCall(MatSetSizes(X, idx[4], idx[4], PETSC_DETERMINE, PETSC_DETERMINE));
241: PetscCall(MatSetUp(X));
242: PetscCall(ISRestoreIndices(sizes, &idx));
243: PetscCall(ISDestroy(&sizes));
244: PetscCall(PetscViewerDestroy(&viewer));
245: PetscCall(PetscSNPrintf(name, sizeof(name), "%s/is_%d.dat", dir, size));
246: PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, name, FILE_MODE_READ, &viewer));
247: PetscCall(ISCreate(PETSC_COMM_WORLD, &sizes));
248: PetscCall(MatGetLayouts(X, &map, NULL));
249: PetscCall(ISSetLayout(sizes, map));
250: PetscCall(ISLoad(sizes, viewer));
251: PetscCall(ISGetLocalSize(sizes, &m));
252: PetscCall(ISGetIndices(sizes, &idx));
253: PetscCall(ISCreateGeneral(PETSC_COMM_SELF, m, idx, PETSC_COPY_VALUES, &is));
254: PetscCall(ISRestoreIndices(sizes, &idx));
255: PetscCall(ISDestroy(&sizes));
256: PetscCall(MatGetBlockSize(A, &m));
257: PetscCall(ISSetBlockSize(is, m));
258: PetscCall(PetscViewerDestroy(&viewer));
259: PetscCall(PetscSNPrintf(name, sizeof(name), "%s/Neumann_%d.dat", dir, size));
260: PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, name, FILE_MODE_READ, &viewer));
261: PetscCall(MatLoad(X, viewer));
262: PetscCall(PetscViewerDestroy(&viewer));
263: PetscCall(MatGetDiagonalBlock(X, &B));
264: PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &aux));
265: PetscCall(MatDestroy(&X));
266: PetscCall(MatSetBlockSizesFromMats(aux, A, A));
267: PetscCall(MatSetOption(aux, MAT_SYMMETRIC, PETSC_TRUE));
268: PetscCall(MatConvert(aux, type, MAT_INPLACE_MATRIX, &aux));
269: }
270: PetscCall(MatCreateVecs(A, NULL, &b));
271: PetscCall(PetscObjectStateIncrease((PetscObject)A));
272: if (!flg) PetscCall(PCHPDDMSetAuxiliaryMat(pc, NULL, aux, NULL, NULL));
273: PetscCall(VecSet(b, 1.0));
274: if (!transpose) PetscCall(KSPSolve(ksp, b, b));
275: else PetscCall(KSPSolveTranspose(ksp, b, b));
276: PetscCall(KSPGetConvergedReason(ksp, reason + 1));
277: PetscCall(KSPGetTotalIterations(ksp, iterations + 1));
278: iterations[1] -= iterations[0];
279: PetscCheck(reason[0] == reason[1] && PetscAbs(iterations[0] - iterations[1]) <= 3, PetscObjectComm((PetscObject)ksp), PETSC_ERR_PLIB, "Successive calls to KSPSolve%s() did not converge for the same reason (%s v. %s) or with the same number of iterations (+/- 3, %" PetscInt_FMT " v. %" PetscInt_FMT ")", (transpose ? "Transpose" : ""), KSPConvergedReasons[reason[0]], KSPConvergedReasons[reason[1]], iterations[0], iterations[1]);
280: PetscCall(PetscObjectStateIncrease((PetscObject)A));
281: if (!flg) PetscCall(PCHPDDMSetAuxiliaryMat(pc, is, aux, NULL, NULL));
282: PetscCall(PCSetFromOptions(pc));
283: PetscCall(VecSet(b, 1.0));
284: if (!transpose) PetscCall(KSPSolve(ksp, b, b));
285: else PetscCall(KSPSolveTranspose(ksp, b, b));
286: PetscCall(KSPGetConvergedReason(ksp, reason + 1));
287: PetscCall(KSPGetTotalIterations(ksp, iterations + 2));
288: iterations[2] -= iterations[0] + iterations[1];
289: PetscCheck(reason[0] == reason[1] && PetscAbs(iterations[0] - iterations[2]) <= 3, PetscObjectComm((PetscObject)ksp), PETSC_ERR_PLIB, "Successive calls to KSPSolve%s() did not converge for the same reason (%s v. %s) or with the same number of iterations (+/- 3, %" PetscInt_FMT " v. %" PetscInt_FMT ")", (transpose ? "Transpose" : ""), KSPConvergedReasons[reason[0]], KSPConvergedReasons[reason[1]], iterations[0], iterations[2]);
290: PetscCall(VecDestroy(&b));
291: PetscCall(ISDestroy(&is));
292: PetscCall(MatDestroy(&aux));
293: }
294: }
295: PetscCall(PetscOptionsGetBool(NULL, NULL, "-viewer", &flg, NULL));
296: if (flg) {
297: PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCHPDDM, &flg));
298: if (flg) {
299: PetscCall(PetscStrncpy(dir, "XXXXXX", sizeof(dir)));
300: if (rank == 0) PetscCall(PetscMkdtemp(dir));
301: PetscCallMPI(MPI_Bcast(dir, 6, MPI_CHAR, 0, PETSC_COMM_WORLD));
302: for (PetscInt i = 0; i < 2; ++i) {
303: PetscCall(PetscSNPrintf(name, sizeof(name), "%s/%s", dir, i == 0 ? "A" : "A.dat"));
304: PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, name, &viewer));
305: PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO_DETAIL));
306: PetscCall(PCView(pc, viewer));
307: PetscCall(PetscViewerPopFormat(viewer));
308: PetscCall(PetscViewerDestroy(&viewer));
309: }
310: PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
311: if (rank == 0) PetscCall(PetscRMTree(dir));
312: }
313: }
314: #endif
315: PetscCall(KSPDestroy(&ksp));
316: PetscCall(MatDestroy(&A));
317: PetscCall(PetscFinalize());
318: return 0;
319: }
321: /*TEST
323: test:
324: requires: hpddm slepc datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
325: nsize: 4
326: args: -ksp_rtol 1e-3 -ksp_converged_reason -pc_type {{bjacobi hpddm}shared output} -pc_hpddm_coarse_sub_pc_type lu -sub_pc_type lu -options_left no -load_dir ${DATAFILESPATH}/matrices/hpddm/GENEO
328: testset:
329: requires: hpddm slepc datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
330: suffix: define_subdomains
331: nsize: 4
332: args: -ksp_rtol 1e-3 -ksp_converged_reason -pc_hpddm_define_subdomains -options_left no -load_dir ${DATAFILESPATH}/matrices/hpddm/GENEO
333: test:
334: args: -pc_type {{asm hpddm}shared output} -pc_hpddm_coarse_sub_pc_type lu -sub_pc_type lu -viewer
335: test:
336: args: -pc_type hpddm -pc_hpddm_levels_1_st_pc_type cholesky -pc_hpddm_levels_1_eps_nev 5 -pc_hpddm_coarse_sub_pc_type lu -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_coarse_correction none
338: testset:
339: requires: hpddm slepc datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
340: nsize: 4
341: args: -ksp_converged_reason -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type cholesky -pc_hpddm_coarse_pc_type redundant -load_dir ${DATAFILESPATH}/matrices/hpddm/GENEO
342: test:
343: suffix: geneo
344: args: -pc_hpddm_coarse_p {{1 2}shared output} -pc_hpddm_levels_1_st_pc_type cholesky -pc_hpddm_levels_1_eps_nev {{5 15}separate output} -mat_type {{aij baij sbaij}shared output}
345: test:
346: suffix: geneo_block_splitting
347: output_file: output/ex76_geneo_pc_hpddm_levels_1_eps_nev-15.out
348: filter: sed -e "s/Linear solve converged due to CONVERGED_RTOL iterations 1[6-9]/Linear solve converged due to CONVERGED_RTOL iterations 11/g"
349: args: -pc_hpddm_coarse_p 2 -pc_hpddm_levels_1_eps_nev 15 -pc_hpddm_block_splitting -pc_hpddm_levels_1_st_pc_type lu -pc_hpddm_levels_1_eps_gen_non_hermitian -mat_type {{aij baij}shared output} -successive_solves
350: test:
351: suffix: geneo_share
352: output_file: output/ex76_geneo_pc_hpddm_levels_1_eps_nev-5.out
353: args: -pc_hpddm_levels_1_st_pc_type cholesky -pc_hpddm_levels_1_eps_nev 5 -pc_hpddm_levels_1_st_share_sub_ksp -reset {{false true}shared output}
354: test:
355: suffix: harmonic_overlap_1_define_false
356: output_file: output/ex76_geneo_share.out
357: filter: sed -e "s/Linear solve converged due to CONVERGED_RTOL iterations 1[0-3]/Linear solve converged due to CONVERGED_RTOL iterations 15/g"
358: args: -pc_hpddm_harmonic_overlap 1 -pc_hpddm_levels_1_eps_nev 30 -pc_hpddm_levels_1_eps_threshold_relative 1e+1 -pc_hpddm_levels_1_st_pc_type lu -pc_hpddm_levels_1_eps_pc_type lu -pc_hpddm_define_subdomains false -pc_hpddm_levels_1_pc_type asm -pc_hpddm_levels_1_pc_asm_overlap 2 -mat_type baij
359: test:
360: suffix: harmonic_overlap_1
361: output_file: output/ex76_geneo_share.out
362: filter: sed -e "s/Linear solve converged due to CONVERGED_RTOL iterations 1[0-3]/Linear solve converged due to CONVERGED_RTOL iterations 15/g"
363: args: -pc_hpddm_harmonic_overlap 1 -pc_hpddm_levels_1_eps_nev 30 -pc_hpddm_levels_1_eps_threshold_relative 1e+1 -pc_hpddm_levels_1_st_pc_type lu -pc_hpddm_levels_1_eps_pc_type lu -mat_type baij
364: test:
365: suffix: harmonic_overlap_1_share_petsc
366: output_file: output/ex76_geneo_share.out
367: filter: sed -e "s/Linear solve converged due to CONVERGED_RTOL iterations 1[0-3]/Linear solve converged due to CONVERGED_RTOL iterations 15/g"
368: args: -pc_hpddm_harmonic_overlap 1 -pc_hpddm_levels_1_eps_nev 30 -pc_hpddm_levels_1_eps_threshold_relative 1e+1 -pc_hpddm_levels_1_st_pc_type lu -pc_hpddm_levels_1_st_share_sub_ksp -pc_hpddm_levels_1_sub_pc_factor_mat_solver_type petsc -pc_hpddm_levels_1_eps_pc_type lu -mat_type baij
369: test:
370: requires: mumps
371: suffix: harmonic_overlap_1_share_mumps
372: output_file: output/ex76_geneo_share.out
373: filter: sed -e "s/Linear solve converged due to CONVERGED_RTOL iterations 1[0-3]/Linear solve converged due to CONVERGED_RTOL iterations 15/g"
374: args: -pc_hpddm_harmonic_overlap 1 -pc_hpddm_levels_1_eps_nev 30 -pc_hpddm_levels_1_eps_threshold_relative 1e+1 -pc_hpddm_levels_1_st_share_sub_ksp -pc_hpddm_levels_1_sub_pc_factor_mat_solver_type mumps -pc_hpddm_coarse_pc_type cholesky -pc_hpddm_coarse_pc_factor_mat_solver_type mumps -pc_hpddm_coarse_mat_mumps_icntl_15 1
375: test:
376: requires: mumps
377: suffix: harmonic_overlap_1_share_mumps_not_set_explicitly
378: output_file: output/ex76_geneo_share.out
379: filter: sed -e "s/Linear solve converged due to CONVERGED_RTOL iterations 1[0-3]/Linear solve converged due to CONVERGED_RTOL iterations 15/g"
380: args: -pc_hpddm_harmonic_overlap 1 -pc_hpddm_levels_1_eps_nev 30 -pc_hpddm_levels_1_eps_threshold_relative 1e+1 -pc_hpddm_levels_1_st_share_sub_ksp -pc_hpddm_levels_1_eps_mat_type baij
381: test:
382: requires: mkl_pardiso
383: suffix: harmonic_overlap_1_share_mkl_pardiso
384: output_file: output/ex76_geneo_share.out
385: filter: sed -e "s/Linear solve converged due to CONVERGED_RTOL iterations [12][0-3]/Linear solve converged due to CONVERGED_RTOL iterations 15/g"
386: args: -pc_hpddm_harmonic_overlap 1 -pc_hpddm_levels_1_eps_nev 30 -pc_hpddm_levels_1_eps_threshold_relative 1e+1 -pc_hpddm_levels_1_st_share_sub_ksp -pc_hpddm_levels_1_eps_mat_type shell -pc_hpddm_levels_1_sub_pc_factor_mat_solver_type mkl_pardiso
387: test:
388: requires: mkl_pardiso !mumps
389: suffix: harmonic_overlap_1_share_mkl_pardiso_no_set_explicitly
390: output_file: output/ex76_geneo_share.out
391: filter: sed -e "s/Linear solve converged due to CONVERGED_RTOL iterations [12][0-3]/Linear solve converged due to CONVERGED_RTOL iterations 15/g"
392: args: -pc_hpddm_harmonic_overlap 1 -pc_hpddm_levels_1_eps_nev 30 -pc_hpddm_levels_1_eps_threshold_relative 1e+1 -pc_hpddm_levels_1_st_share_sub_ksp -pc_hpddm_levels_1_eps_mat_type shell
393: test:
394: suffix: harmonic_overlap_2_threshold_relative
395: output_file: output/ex76_geneo_share.out
396: filter: sed -e "s/Linear solve converged due to CONVERGED_RTOL iterations 9/Linear solve converged due to CONVERGED_RTOL iterations 15/g"
397: args: -pc_hpddm_harmonic_overlap 2 -pc_hpddm_levels_1_svd_nsv 15 -pc_hpddm_levels_1_svd_threshold_relative 1e-1 -pc_hpddm_levels_1_st_share_sub_ksp -mat_type sbaij
398: test:
399: suffix: harmonic_overlap_2
400: output_file: output/ex76_geneo_share.out
401: filter: sed -e "s/Linear solve converged due to CONVERGED_RTOL iterations 9/Linear solve converged due to CONVERGED_RTOL iterations 15/g"
402: args: -pc_hpddm_harmonic_overlap 2 -pc_hpddm_levels_1_svd_nsv 12 -pc_hpddm_levels_1_st_share_sub_ksp -mat_type sbaij
404: testset:
405: requires: hpddm slepc datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
406: nsize: 4
407: args: -ksp_converged_reason -ksp_max_it 150 -pc_type hpddm -pc_hpddm_levels_1_eps_nev 5 -pc_hpddm_coarse_p 1 -pc_hpddm_coarse_pc_type redundant -load_dir ${DATAFILESPATH}/matrices/hpddm/GENEO -pc_hpddm_define_subdomains
408: test:
409: suffix: geneo_share_cholesky
410: output_file: output/ex76_geneo_share.out
411: # extra -pc_hpddm_levels_1_eps_gen_non_hermitian needed to avoid failures with PETSc Cholesky
412: args: -pc_hpddm_levels_1_sub_pc_type cholesky -pc_hpddm_levels_1_st_pc_type cholesky -mat_type {{aij sbaij}shared output} -pc_hpddm_levels_1_eps_gen_non_hermitian -pc_hpddm_has_neumann -pc_hpddm_levels_1_st_share_sub_ksp {{false true}shared output} -successive_solves
413: test:
414: suffix: geneo_share_cholesky_matstructure
415: output_file: output/ex76_geneo_share.out
416: # extra -pc_hpddm_levels_1_eps_gen_non_hermitian needed to avoid failures with PETSc Cholesky
417: filter: sed -e "s/Linear solve converged due to CONVERGED_RTOL iterations 14/Linear solve converged due to CONVERGED_RTOL iterations 15/g"
418: args: -pc_hpddm_levels_1_sub_pc_type cholesky -mat_type {{baij sbaij}shared output} -pc_hpddm_levels_1_eps_gen_non_hermitian -pc_hpddm_levels_1_st_share_sub_ksp -pc_hpddm_levels_1_st_matstructure same -set_rhs {{false true} shared output}
419: test:
420: suffix: geneo_transpose
421: output_file: output/ex76_geneo_share.out
422: filter: sed -e "s/Linear solve converged due to CONVERGED_RTOL iterations 1[234]/Linear solve converged due to CONVERGED_RTOL iterations 15/g" -e "s/Linear solve converged due to CONVERGED_RTOL iterations 26/Linear solve converged due to CONVERGED_RTOL iterations 15/g"
423: args: -pc_hpddm_levels_1_sub_pc_type cholesky -pc_hpddm_levels_1_st_pc_type cholesky -pc_hpddm_levels_1_eps_gen_non_hermitian -pc_hpddm_has_neumann -pc_hpddm_levels_1_st_share_sub_ksp -successive_solves -transpose -pc_hpddm_coarse_correction {{additive deflated balanced}shared output}
424: test:
425: requires: mumps
426: suffix: geneo_share_lu
427: output_file: output/ex76_geneo_share.out
428: # extra -pc_factor_mat_solver_type mumps needed to avoid failures with PETSc LU
429: args: -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_levels_1_st_pc_type lu -mat_type baij -pc_hpddm_levels_1_st_pc_factor_mat_solver_type mumps -pc_hpddm_levels_1_sub_pc_factor_mat_solver_type mumps -pc_hpddm_has_neumann -pc_hpddm_levels_1_st_share_sub_ksp {{false true}shared output}
430: test:
431: requires: mumps
432: suffix: geneo_share_lu_matstructure
433: output_file: output/ex76_geneo_share.out
434: # extra -pc_factor_mat_solver_type mumps needed to avoid failures with PETSc LU
435: args: -pc_hpddm_levels_1_sub_pc_type lu -mat_type aij -pc_hpddm_levels_1_sub_pc_factor_mat_solver_type mumps -pc_hpddm_levels_1_st_share_sub_ksp -pc_hpddm_levels_1_st_matstructure {{same different}shared output} -pc_hpddm_levels_1_st_pc_type lu -pc_hpddm_levels_1_st_pc_factor_mat_solver_type mumps -successive_solves -pc_hpddm_levels_1_eps_target 1e-5
436: test:
437: suffix: geneo_share_not_asm
438: output_file: output/ex76_geneo_pc_hpddm_levels_1_eps_nev-5.out
439: # extra -pc_hpddm_levels_1_eps_gen_non_hermitian needed to avoid failures with PETSc Cholesky
440: args: -pc_hpddm_levels_1_sub_pc_type cholesky -pc_hpddm_levels_1_st_pc_type cholesky -pc_hpddm_levels_1_eps_gen_non_hermitian -pc_hpddm_has_neumann -pc_hpddm_levels_1_st_share_sub_ksp true -pc_hpddm_levels_1_pc_type gasm -successive_solves
442: test:
443: requires: hpddm slepc datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
444: suffix: fgmres_geneo_20_p_2
445: nsize: 4
446: args: -ksp_converged_reason -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_levels_1_eps_nev 20 -pc_hpddm_coarse_p 2 -pc_hpddm_coarse_pc_type redundant -ksp_type fgmres -pc_hpddm_coarse_mat_type {{baij sbaij}shared output} -pc_hpddm_log_separate {{false true}shared output} -load_dir ${DATAFILESPATH}/matrices/hpddm/GENEO
448: testset:
449: requires: hpddm slepc datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
450: output_file: output/ex76_fgmres_geneo_20_p_2.out
451: nsize: 4
452: args: -ksp_converged_reason -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type cholesky -pc_hpddm_levels_1_eps_nev 20 -pc_hpddm_levels_2_p 2 -pc_hpddm_levels_2_mat_type {{baij sbaij}shared output} -pc_hpddm_levels_2_eps_nev {{5 20}shared output} -pc_hpddm_levels_2_sub_pc_type cholesky -pc_hpddm_levels_2_ksp_type gmres -ksp_type fgmres -pc_hpddm_coarse_mat_type {{baij sbaij}shared output} -load_dir ${DATAFILESPATH}/matrices/hpddm/GENEO
453: test:
454: suffix: fgmres_geneo_20_p_2_geneo
455: args: -mat_type {{aij sbaij}shared output}
456: test:
457: suffix: fgmres_geneo_20_p_2_geneo_algebraic
458: args: -pc_hpddm_levels_2_st_pc_type mat
459: # PCHPDDM + KSPHPDDM test to exercise multilevel + multiple RHS in one go
460: test:
461: requires: hpddm slepc datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
462: suffix: fgmres_geneo_20_p_2_geneo_rhs
463: output_file: output/ex76_fgmres_geneo_20_p_2.out
464: # for -pc_hpddm_coarse_correction additive
465: filter: sed -e "s/Linear solve converged due to CONVERGED_RTOL iterations 37/Linear solve converged due to CONVERGED_RTOL iterations 25/g"
466: nsize: 4
467: args: -ksp_converged_reason -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type cholesky -pc_hpddm_levels_1_eps_nev 20 -pc_hpddm_levels_2_p 2 -pc_hpddm_levels_2_mat_type baij -pc_hpddm_levels_2_eps_nev 5 -pc_hpddm_levels_2_sub_pc_type cholesky -pc_hpddm_levels_2_ksp_max_it 10 -pc_hpddm_levels_2_ksp_type hpddm -pc_hpddm_levels_2_ksp_hpddm_type gmres -ksp_type hpddm -ksp_hpddm_variant flexible -pc_hpddm_coarse_mat_type baij -mat_type aij -load_dir ${DATAFILESPATH}/matrices/hpddm/GENEO -rhs 4 -pc_hpddm_coarse_correction {{additive deflated balanced}shared output}
469: testset:
470: requires: hpddm slepc datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES) mumps defined(PETSC_HAVE_OPENMP_SUPPORT)
471: filter: grep -E -e "Linear solve" -e " executing" | sed -e "s/MPI = 1/MPI = 2/g" -e "s/OMP = 1/OMP = 2/g"
472: nsize: 4
473: args: -ksp_converged_reason -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type cholesky -pc_hpddm_levels_1_eps_nev 15 -pc_hpddm_levels_1_st_pc_type cholesky -pc_hpddm_coarse_p {{1 2}shared output} -load_dir ${DATAFILESPATH}/matrices/hpddm/GENEO -pc_hpddm_coarse_pc_factor_mat_solver_type mumps -pc_hpddm_coarse_mat_mumps_icntl_4 2 -pc_hpddm_coarse_mat_mumps_use_omp_threads {{1 2}shared output}
474: test:
475: suffix: geneo_mumps_use_omp_threads_1
476: output_file: output/ex76_geneo_mumps_use_omp_threads.out
477: args: -pc_hpddm_coarse_mat_type {{baij sbaij}shared output}
478: test:
479: suffix: geneo_mumps_use_omp_threads_2
480: output_file: output/ex76_geneo_mumps_use_omp_threads.out
481: args: -pc_hpddm_coarse_mat_type aij -pc_hpddm_levels_1_eps_threshold_absolute 0.4 -pc_hpddm_coarse_pc_type cholesky -pc_hpddm_coarse_mat_filter 1e-12
483: testset: # converge really poorly because of a tiny -pc_hpddm_levels_1_eps_threshold_absolute, but needed for proper code coverage where some subdomains don't call EPSSolve()
484: requires: hpddm slepc datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
485: nsize: 4
486: args: -ksp_converged_reason -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type cholesky -pc_hpddm_levels_1_eps_threshold_absolute 0.005 -pc_hpddm_levels_1_eps_use_inertia -load_dir ${DATAFILESPATH}/matrices/hpddm/GENEO -pc_hpddm_levels_1_st_share_sub_ksp -pc_hpddm_define_subdomains -pc_hpddm_has_neumann -ksp_rtol 0.9
487: filter: sed -e "s/Linear solve converged due to CONVERGED_RTOL iterations 1/Linear solve converged due to CONVERGED_RTOL iterations 141/g"
488: test:
489: suffix: inertia_petsc
490: output_file: output/ex76_1.out
491: args: -pc_hpddm_levels_1_sub_pc_factor_mat_solver_type petsc
492: test:
493: suffix: inertia_mumps
494: output_file: output/ex76_1.out
495: requires: mumps
497: test:
498: requires: hpddm slepc datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
499: suffix: reuse_symbolic
500: output_file: output/empty.out
501: nsize: 4
502: args: -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type cholesky -pc_hpddm_levels_1_eps_nev 20 -rhs 4 -pc_hpddm_coarse_correction {{additive deflated balanced}shared output} -ksp_pc_side {{left right}shared output} -ksp_max_it 20 -ksp_type hpddm -load_dir ${DATAFILESPATH}/matrices/hpddm/GENEO -pc_hpddm_define_subdomains -ksp_error_if_not_converged -transpose {{true false} shared output}
504: TEST*/