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:   PetscViewer     viewer;
 17:   char            dir[PETSC_MAX_PATH_LEN], name[PETSC_MAX_PATH_LEN], type[256];
 18:   PetscBool3      share = PETSC_BOOL3_UNKNOWN;
 19:   PetscBool       flg, set;
 20: #if defined(PETSC_USE_LOG)
 21:   PetscLogEvent event;
 22: #endif
 23:   PetscEventPerfInfo info1, info2;

 25:   PetscFunctionBeginUser;
 26:   PetscCall(PetscInitialize(&argc, &args, NULL, help));
 27:   PetscCall(PetscLogDefaultBegin());
 28:   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
 29:   PetscCheck(size == 4, PETSC_COMM_WORLD, PETSC_ERR_USER, "This example requires 4 processes");
 30:   PetscCall(PetscOptionsGetInt(NULL, NULL, "-rhs", &N, NULL));
 31:   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
 32:   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
 33:   PetscCall(MatCreate(PETSC_COMM_SELF, &aux));
 34:   PetscCall(ISCreate(PETSC_COMM_SELF, &is));
 35:   PetscCall(PetscStrncpy(dir, ".", sizeof(dir)));
 36:   PetscCall(PetscOptionsGetString(NULL, NULL, "-load_dir", dir, sizeof(dir), NULL));
 37:   /* loading matrices */
 38:   PetscCall(PetscSNPrintf(name, sizeof(name), "%s/sizes_%d_%d.dat", dir, rank, size));
 39:   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, name, FILE_MODE_READ, &viewer));
 40:   PetscCall(ISCreate(PETSC_COMM_SELF, &sizes));
 41:   PetscCall(ISLoad(sizes, viewer));
 42:   PetscCall(ISGetIndices(sizes, &idx));
 43:   PetscCall(MatSetSizes(A, idx[0], idx[1], idx[2], idx[3]));
 44:   PetscCall(MatSetUp(A));
 45:   PetscCall(ISRestoreIndices(sizes, &idx));
 46:   PetscCall(ISDestroy(&sizes));
 47:   PetscCall(PetscViewerDestroy(&viewer));
 48:   PetscCall(PetscSNPrintf(name, sizeof(name), "%s/A.dat", dir));
 49:   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, name, FILE_MODE_READ, &viewer));
 50:   PetscCall(MatLoad(A, viewer));
 51:   PetscCall(PetscViewerDestroy(&viewer));
 52:   PetscCall(PetscSNPrintf(name, sizeof(name), "%s/is_%d_%d.dat", dir, rank, size));
 53:   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, name, FILE_MODE_READ, &viewer));
 54:   PetscCall(ISLoad(is, viewer));
 55:   PetscCall(ISSetBlockSize(is, 2));
 56:   PetscCall(PetscViewerDestroy(&viewer));
 57:   PetscCall(PetscSNPrintf(name, sizeof(name), "%s/Neumann_%d_%d.dat", dir, rank, size));
 58:   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, name, FILE_MODE_READ, &viewer));
 59:   PetscCall(MatLoad(aux, viewer));
 60:   PetscCall(PetscViewerDestroy(&viewer));
 61:   flg = PETSC_FALSE;
 62:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-pc_hpddm_levels_1_st_share_sub_ksp", &flg, &set));
 63:   if (flg) { /* PETSc LU/Cholesky is struggling numerically for bs > 1          */
 64:              /* only set the proper bs for the geneo_share_* tests, 1 otherwise */
 65:     PetscCall(MatSetBlockSizesFromMats(aux, A, A));
 66:     share = PETSC_BOOL3_TRUE;
 67:   } else if (set) share = PETSC_BOOL3_FALSE;
 68:   PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE));
 69:   PetscCall(MatSetOption(aux, MAT_SYMMETRIC, PETSC_TRUE));
 70:   /* ready for testing */
 71:   PetscOptionsBegin(PETSC_COMM_WORLD, "", "", "");
 72:   PetscCall(PetscStrncpy(type, MATAIJ, sizeof(type)));
 73:   PetscCall(PetscOptionsFList("-mat_type", "Matrix type", "MatSetType", MatList, type, type, 256, &flg));
 74:   PetscOptionsEnd();
 75:   PetscCall(MatConvert(A, type, MAT_INPLACE_MATRIX, &A));
 76:   PetscCall(MatConvert(aux, type, MAT_INPLACE_MATRIX, &aux));
 77:   PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
 78:   PetscCall(KSPSetOperators(ksp, A, A));
 79:   PetscCall(KSPGetPC(ksp, &pc));
 80:   PetscCall(PCSetType(pc, PCHPDDM));
 81: #if defined(PETSC_HAVE_HPDDM) && defined(PETSC_HAVE_DYNAMIC_LIBRARIES) && defined(PETSC_USE_SHARED_LIBRARIES)
 82:   flg = PETSC_FALSE;
 83:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-reset", &flg, NULL));
 84:   if (flg) {
 85:     PetscCall(PetscOptionsSetValue(NULL, "-pc_hpddm_block_splitting", "true"));
 86:     PetscCall(PCSetFromOptions(pc));
 87:     PetscCall(PCSetUp(pc));
 88:     PetscCall(PetscOptionsClearValue(NULL, "-pc_hpddm_block_splitting"));
 89:   }
 90:   PetscCall(PCHPDDMSetAuxiliaryMat(pc, is, aux, NULL, NULL));
 91:   PetscCall(PCHPDDMHasNeumannMat(pc, PETSC_FALSE)); /* PETSC_TRUE is fine as well, just testing */
 92:   if (share == PETSC_BOOL3_UNKNOWN) PetscCall(PCHPDDMSetSTShareSubKSP(pc, PetscBool3ToBool(share)));
 93:   flg = PETSC_FALSE;
 94:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-set_rhs", &flg, NULL));
 95:   if (flg) {          /* user-provided RHS for concurrent generalized eigenvalue problems                          */
 96:     Mat      a, c, P; /* usually assembled automatically in PCHPDDM, this is solely for testing PCHPDDMSetRHSMat() */
 97:     PetscInt rstart, rend, location;
 98:     PetscCall(MatDuplicate(aux, MAT_DO_NOT_COPY_VALUES, &B)); /* duplicate so that MatStructure is SAME_NONZERO_PATTERN */
 99:     PetscCall(MatGetDiagonalBlock(A, &a));
100:     PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
101:     PetscCall(ISGetLocalSize(is, &m));
102:     PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, rend - rstart, m, 1, NULL, &P));
103:     for (m = rstart; m < rend; ++m) {
104:       PetscCall(ISLocate(is, m, &location));
105:       PetscCheck(location >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "IS of the auxiliary Mat does not include all local rows of A");
106:       PetscCall(MatSetValue(P, m - rstart, location, 1.0, INSERT_VALUES));
107:     }
108:     PetscCall(MatAssemblyBegin(P, MAT_FINAL_ASSEMBLY));
109:     PetscCall(MatAssemblyEnd(P, MAT_FINAL_ASSEMBLY));
110:     PetscCall(PetscObjectTypeCompare((PetscObject)a, MATSEQAIJ, &flg));
111:     if (flg) PetscCall(MatPtAP(a, P, MAT_INITIAL_MATRIX, 1.0, &X)); /* MatPtAP() is used to extend diagonal blocks with zeros on the overlap */
112:     else { /* workaround for MatPtAP() limitations with some types */ PetscCall(MatConvert(a, MATSEQAIJ, MAT_INITIAL_MATRIX, &c));
113:       PetscCall(MatPtAP(c, P, MAT_INITIAL_MATRIX, 1.0, &X));
114:       PetscCall(MatDestroy(&c));
115:     }
116:     PetscCall(MatDestroy(&P));
117:     PetscCall(MatAXPY(B, 1.0, X, SUBSET_NONZERO_PATTERN));
118:     PetscCall(MatDestroy(&X));
119:     PetscCall(MatSetOption(B, MAT_SYMMETRIC, PETSC_TRUE));
120:     PetscCall(PCHPDDMSetRHSMat(pc, B));
121:     PetscCall(MatDestroy(&B));
122:   }
123: #else
124:   (void)share;
125: #endif
126:   PetscCall(MatDestroy(&aux));
127:   PetscCall(KSPSetFromOptions(ksp));
128:   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &flg));
129:   if (flg) {
130:     flg = PETSC_FALSE;
131:     PetscCall(PetscOptionsGetBool(NULL, NULL, "-pc_hpddm_define_subdomains", &flg, NULL));
132:     if (flg) {
133:       IS rows;
134:       PetscCall(MatGetOwnershipIS(A, &rows, NULL));
135:       PetscCall(PCASMSetLocalSubdomains(pc, 1, &is, &rows));
136:       PetscCall(ISDestroy(&rows));
137:     }
138:   }
139:   PetscCall(ISDestroy(&is));
140:   PetscCall(MatCreateVecs(A, NULL, &b));
141:   PetscCall(VecSet(b, 1.0));
142:   PetscCall(KSPSolve(ksp, b, b));
143:   PetscCall(VecGetLocalSize(b, &m));
144:   PetscCall(VecDestroy(&b));
145:   if (N > 1) {
146:     KSPType type;
147:     PetscCall(PetscOptionsClearValue(NULL, "-ksp_converged_reason"));
148:     PetscCall(KSPSetFromOptions(ksp));
149:     PetscCall(MatCreateDense(PETSC_COMM_WORLD, m, PETSC_DECIDE, PETSC_DECIDE, N, NULL, &B));
150:     PetscCall(MatCreateDense(PETSC_COMM_WORLD, m, PETSC_DECIDE, PETSC_DECIDE, N, NULL, &X));
151:     PetscCall(MatSetRandom(B, NULL));
152:     /* this is algorithmically optimal in the sense that blocks of vectors are coarsened or interpolated using matrix--matrix operations */
153:     /* PCHPDDM however heavily relies on MPI[S]BAIJ format for which there is no efficient MatProduct implementation */
154:     PetscCall(KSPMatSolve(ksp, B, X));
155:     PetscCall(KSPGetType(ksp, &type));
156:     PetscCall(PetscStrcmp(type, KSPHPDDM, &flg));
157: #if defined(PETSC_HAVE_HPDDM)
158:     if (flg) {
159:       PetscReal    norm;
160:       KSPHPDDMType type;
161:       PetscCall(KSPHPDDMGetType(ksp, &type));
162:       if (type == KSP_HPDDM_TYPE_PREONLY || type == KSP_HPDDM_TYPE_CG || type == KSP_HPDDM_TYPE_GMRES || type == KSP_HPDDM_TYPE_GCRODR) {
163:         Mat C;
164:         PetscCall(MatDuplicate(X, MAT_DO_NOT_COPY_VALUES, &C));
165:         PetscCall(KSPSetMatSolveBatchSize(ksp, 1));
166:         PetscCall(KSPMatSolve(ksp, B, C));
167:         PetscCall(MatAYPX(C, -1.0, X, SAME_NONZERO_PATTERN));
168:         PetscCall(MatNorm(C, NORM_INFINITY, &norm));
169:         PetscCall(MatDestroy(&C));
170:         PetscCheck(norm <= 100 * PETSC_MACHINE_EPSILON, PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "KSPMatSolve() and KSPSolve() difference has nonzero norm %g with pseudo-block KSPHPDDMType %s", (double)norm, KSPHPDDMTypes[type]);
171:       }
172:     }
173: #endif
174:     PetscCall(MatDestroy(&X));
175:     PetscCall(MatDestroy(&B));
176:   }
177:   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCHPDDM, &flg));
178: #if defined(PETSC_HAVE_HPDDM) && defined(PETSC_HAVE_DYNAMIC_LIBRARIES) && defined(PETSC_USE_SHARED_LIBRARIES)
179:   if (flg) PetscCall(PCHPDDMGetSTShareSubKSP(pc, &flg));
180: #endif
181:   if (flg && PetscDefined(USE_LOG)) {
182:     PetscCall(PetscLogEventRegister("MatLUFactorSym", PC_CLASSID, &event));
183:     PetscCall(PetscLogEventGetPerfInfo(PETSC_DETERMINE, event, &info1));
184:     PetscCall(PetscLogEventRegister("MatLUFactorNum", PC_CLASSID, &event));
185:     PetscCall(PetscLogEventGetPerfInfo(PETSC_DETERMINE, event, &info2));
186:     if (!info1.count && !info2.count) {
187:       PetscCall(PetscLogEventRegister("MatCholFctrSym", PC_CLASSID, &event));
188:       PetscCall(PetscLogEventGetPerfInfo(PETSC_DETERMINE, event, &info1));
189:       PetscCall(PetscLogEventRegister("MatCholFctrNum", PC_CLASSID, &event));
190:       PetscCall(PetscLogEventGetPerfInfo(PETSC_DETERMINE, event, &info2));
191:       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);
192:     } 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);
193:   }
194: #if defined(PETSC_HAVE_HPDDM) && defined(PETSC_HAVE_DYNAMIC_LIBRARIES) && defined(PETSC_USE_SHARED_LIBRARIES)
195:   if (N == 1) {
196:     flg = PETSC_FALSE;
197:     PetscCall(PetscOptionsGetBool(NULL, NULL, "-successive_solves", &flg, NULL));
198:     if (flg) {
199:       KSPConvergedReason reason[2];
200:       PetscInt           iterations[3];
201:       PetscCall(KSPGetConvergedReason(ksp, reason));
202:       PetscCall(KSPGetTotalIterations(ksp, iterations));
203:       PetscCall(PetscOptionsClearValue(NULL, "-ksp_converged_reason"));
204:       flg = PETSC_FALSE;
205:       PetscCall(PetscOptionsGetBool(NULL, NULL, "-pc_hpddm_block_splitting", &flg, NULL));
206:       if (!flg) {
207:         PetscCall(MatCreate(PETSC_COMM_SELF, &aux));
208:         PetscCall(ISCreate(PETSC_COMM_SELF, &is));
209:         PetscCall(PetscSNPrintf(name, sizeof(name), "%s/is_%d_%d.dat", dir, rank, size));
210:         PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, name, FILE_MODE_READ, &viewer));
211:         PetscCall(ISLoad(is, viewer));
212:         PetscCall(PetscViewerDestroy(&viewer));
213:         PetscCall(PetscSNPrintf(name, sizeof(name), "%s/Neumann_%d_%d.dat", dir, rank, size));
214:         PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, name, FILE_MODE_READ, &viewer));
215:         PetscCall(MatLoad(aux, viewer));
216:         PetscCall(PetscViewerDestroy(&viewer));
217:         PetscCall(MatSetBlockSizesFromMats(aux, A, A));
218:         PetscCall(MatSetOption(aux, MAT_SYMMETRIC, PETSC_TRUE));
219:         PetscCall(MatConvert(aux, type, MAT_INPLACE_MATRIX, &aux));
220:       }
221:       PetscCall(MatCreateVecs(A, NULL, &b));
222:       PetscCall(PetscObjectStateIncrease((PetscObject)A));
223:       if (!flg) PetscCall(PCHPDDMSetAuxiliaryMat(pc, NULL, aux, NULL, NULL));
224:       PetscCall(VecSet(b, 1.0));
225:       PetscCall(KSPSolve(ksp, b, b));
226:       PetscCall(KSPGetConvergedReason(ksp, reason + 1));
227:       PetscCall(KSPGetTotalIterations(ksp, iterations + 1));
228:       iterations[1] -= iterations[0];
229:       PetscCheck(reason[0] == reason[1] && PetscAbs(iterations[0] - iterations[1]) <= 3, PetscObjectComm((PetscObject)ksp), PETSC_ERR_PLIB, "Successive calls to KSPSolve() did not converge for the same reason or with the same number of iterations (+/- 3)");
230:       PetscCall(PetscObjectStateIncrease((PetscObject)A));
231:       if (!flg) PetscCall(PCHPDDMSetAuxiliaryMat(pc, is, aux, NULL, NULL));
232:       PetscCall(PCSetFromOptions(pc));
233:       PetscCall(VecSet(b, 1.0));
234:       PetscCall(KSPSolve(ksp, b, b));
235:       PetscCall(KSPGetConvergedReason(ksp, reason + 1));
236:       PetscCall(KSPGetTotalIterations(ksp, iterations + 2));
237:       iterations[2] -= iterations[0] + iterations[1];
238:       PetscCheck(reason[0] == reason[1] && PetscAbs(iterations[0] - iterations[2]) <= 3, PetscObjectComm((PetscObject)ksp), PETSC_ERR_PLIB, "Successive calls to KSPSolve() did not converge for the same reason or with the same number of iterations (+/- 3)");
239:       PetscCall(VecDestroy(&b));
240:       PetscCall(ISDestroy(&is));
241:       PetscCall(MatDestroy(&aux));
242:     }
243:   }
244: #endif
245:   PetscCall(KSPDestroy(&ksp));
246:   PetscCall(MatDestroy(&A));
247:   PetscCall(PetscFinalize());
248:   return 0;
249: }

251: /*TEST

253:    test:
254:       requires: hpddm slepc datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
255:       nsize: 4
256:       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

258:    test:
259:       requires: hpddm slepc datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
260:       suffix: define_subdomains
261:       nsize: 4
262:       args: -ksp_rtol 1e-3 -ksp_converged_reason -pc_type {{asm hpddm}shared output} -pc_hpddm_coarse_sub_pc_type lu -sub_pc_type lu -pc_hpddm_define_subdomains -options_left no -load_dir ${DATAFILESPATH}/matrices/hpddm/GENEO

264:    testset:
265:       requires: hpddm slepc datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
266:       nsize: 4
267:       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
268:       test:
269:         suffix: geneo
270:         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}
271:       test:
272:         suffix: geneo_block_splitting
273:         output_file: output/ex76_geneo_pc_hpddm_levels_1_eps_nev-15.out
274:         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"
275:         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
276:       test:
277:         suffix: geneo_share
278:         output_file: output/ex76_geneo_pc_hpddm_levels_1_eps_nev-5.out
279:         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}

281:    testset:
282:       requires: hpddm slepc datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
283:       nsize: 4
284:       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
285:       test:
286:         suffix: geneo_share_cholesky
287:         output_file: output/ex76_geneo_share.out
288:         # extra -pc_hpddm_levels_1_eps_gen_non_hermitian needed to avoid failures with PETSc Cholesky
289:         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
290:       test:
291:         suffix: geneo_share_cholesky_matstructure
292:         output_file: output/ex76_geneo_share.out
293:         # extra -pc_hpddm_levels_1_eps_gen_non_hermitian needed to avoid failures with PETSc Cholesky
294:         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}
295:       test:
296:         requires: mumps
297:         suffix: geneo_share_lu
298:         output_file: output/ex76_geneo_share.out
299:         # extra -pc_factor_mat_solver_type mumps needed to avoid failures with PETSc LU
300:         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}
301:       test:
302:         requires: mumps
303:         suffix: geneo_share_lu_matstructure
304:         output_file: output/ex76_geneo_share.out
305:         # extra -pc_factor_mat_solver_type mumps needed to avoid failures with PETSc LU
306:         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
307:       test:
308:         suffix: geneo_share_not_asm
309:         output_file: output/ex76_geneo_pc_hpddm_levels_1_eps_nev-5.out
310:         # extra -pc_hpddm_levels_1_eps_gen_non_hermitian needed to avoid failures with PETSc Cholesky
311:         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

313:    test:
314:       requires: hpddm slepc datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
315:       suffix: fgmres_geneo_20_p_2
316:       nsize: 4
317:       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

319:    testset:
320:       requires: hpddm slepc datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
321:       output_file: output/ex76_fgmres_geneo_20_p_2.out
322:       nsize: 4
323:       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
324:       test:
325:         suffix: fgmres_geneo_20_p_2_geneo
326:         args: -mat_type {{aij sbaij}shared output}
327:       test:
328:         suffix: fgmres_geneo_20_p_2_geneo_algebraic
329:         args: -pc_hpddm_levels_2_st_pc_type mat
330:    # PCHPDDM + KSPHPDDM test to exercise multilevel + multiple RHS in one go
331:    test:
332:       requires: hpddm slepc datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
333:       suffix: fgmres_geneo_20_p_2_geneo_rhs
334:       output_file: output/ex76_fgmres_geneo_20_p_2.out
335:       # for -pc_hpddm_coarse_correction additive
336:       filter: sed -e "s/Linear solve converged due to CONVERGED_RTOL iterations 37/Linear solve converged due to CONVERGED_RTOL iterations 25/g"
337:       nsize: 4
338:       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}

340:    testset:
341:       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)
342:       filter: grep -E -e "Linear solve" -e "      executing" | sed -e "s/MPI =      1/MPI =      2/g" -e "s/OMP =      1/OMP =      2/g"
343:       nsize: 4
344:       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}
345:       test:
346:         suffix: geneo_mumps_use_omp_threads_1
347:         output_file: output/ex76_geneo_mumps_use_omp_threads.out
348:         args: -pc_hpddm_coarse_mat_type {{baij sbaij}shared output}
349:       test:
350:         suffix: geneo_mumps_use_omp_threads_2
351:         output_file: output/ex76_geneo_mumps_use_omp_threads.out
352:         args: -pc_hpddm_coarse_mat_type aij -pc_hpddm_levels_1_eps_threshold 0.3 -pc_hpddm_coarse_pc_type cholesky -pc_hpddm_coarse_mat_chop 1e-12

354:    testset: # converge really poorly because of a tiny -pc_hpddm_levels_1_eps_threshold, but needed for proper code coverage where some subdomains don't call EPSSolve()
355:       requires: hpddm slepc datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
356:       nsize: 4
357:       args: -ksp_converged_reason -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type cholesky -pc_hpddm_levels_1_eps_threshold 0.005 -pc_hpddm_levels_1_eps_use_inertia -load_dir ${DATAFILESPATH}/matrices/hpddm/GENEO -pc_hpddm_coarse_pc_type cholesky -pc_hpddm_levels_1_st_share_sub_ksp -pc_hpddm_define_subdomains -pc_hpddm_has_neumann -ksp_rtol 0.9
358:       filter: sed -e "s/Linear solve converged due to CONVERGED_RTOL iterations 1/Linear solve converged due to CONVERGED_RTOL iterations 141/g"
359:       test:
360:         suffix: inertia_petsc
361:         output_file: output/ex76_1.out
362:         args: -pc_hpddm_levels_1_sub_pc_factor_mat_solver_type petsc
363:       test:
364:         suffix: inertia_mumps
365:         output_file: output/ex76_1.out
366:         requires: mumps

368:    test:
369:       requires: hpddm slepc datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
370:       suffix: reuse_symbolic
371:       output_file: output/ex77_preonly.out
372:       nsize: 4
373:       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

375: TEST*/