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

250: /*TEST

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

257:    test:
258:       requires: hpddm slepc datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
259:       suffix: define_subdomains
260:       nsize: 4
261:       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

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

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

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

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

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

353:    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()
354:       requires: hpddm slepc datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
355:       nsize: 4
356:       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
357:       filter: sed -e "s/Linear solve converged due to CONVERGED_RTOL iterations 1/Linear solve converged due to CONVERGED_RTOL iterations 141/g"
358:       test:
359:         suffix: inertia_petsc
360:         output_file: output/ex76_1.out
361:         args: -pc_hpddm_levels_1_sub_pc_factor_mat_solver_type petsc
362:       test:
363:         suffix: inertia_mumps
364:         output_file: output/ex76_1.out
365:         requires: mumps

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

374: TEST*/