Actual source code: hpddm.cxx

  1: #include <petsc/private/dmimpl.h>
  2: #include <petsc/private/matimpl.h>
  3: #include <petsc/private/petschpddm.h>
  4: #include <petsc/private/pcimpl.h>
  5:                                   /* otherwise, it is assumed that one is compiling libhpddm_petsc => circular dependency */
  6: #if defined(PETSC_HAVE_FORTRAN)
  7: #include <petsc/private/fortranimpl.h>
  8: #endif

 10: static PetscErrorCode (*loadedSym)(HPDDM::Schwarz<PetscScalar>* const, IS, Mat, Mat, Mat, std::vector<Vec>, PC_HPDDM_Level** const) = NULL;

 12: static PetscBool PCHPDDMPackageInitialized = PETSC_FALSE;

 14: PetscLogEvent PC_HPDDM_Strc;
 15: PetscLogEvent PC_HPDDM_PtAP;
 16: PetscLogEvent PC_HPDDM_PtBP;
 17: PetscLogEvent PC_HPDDM_Next;
 18: PetscLogEvent PC_HPDDM_SetUp[PETSC_PCHPDDM_MAXLEVELS];
 19: PetscLogEvent PC_HPDDM_Solve[PETSC_PCHPDDM_MAXLEVELS];

 21: const char *const PCHPDDMCoarseCorrectionTypes[] = { "DEFLATED", "ADDITIVE", "BALANCED", "PCHPDDMCoarseCorrectionType", "PC_HPDDM_COARSE_CORRECTION_", NULL };

 23: static PetscErrorCode PCReset_HPDDM(PC pc)
 24: {
 25:   PC_HPDDM *data = (PC_HPDDM*)pc->data;
 26:   PetscInt i;

 28:   if (data->levels) {
 29:     for (i = 0; i < PETSC_PCHPDDM_MAXLEVELS && data->levels[i]; ++i) {
 30:       KSPDestroy(&data->levels[i]->ksp);
 31:       PCDestroy(&data->levels[i]->pc);
 32:       PetscFree(data->levels[i]);
 33:     }
 34:     PetscFree(data->levels);
 35:   }

 37:   ISDestroy(&data->is);
 38:   MatDestroy(&data->aux);
 39:   MatDestroy(&data->B);
 40:   VecDestroy(&data->normal);
 41:   data->correction = PC_HPDDM_COARSE_CORRECTION_DEFLATED;
 42:   data->Neumann    = PETSC_FALSE;
 43:   data->setup      = NULL;
 44:   data->setup_ctx  = NULL;
 45:   return 0;
 46: }

 48: static PetscErrorCode PCDestroy_HPDDM(PC pc)
 49: {
 50:   PC_HPDDM *data = (PC_HPDDM*)pc->data;

 52:   PCReset_HPDDM(pc);
 53:   PetscFree(data);
 54:   PetscObjectChangeTypeName((PetscObject)pc, NULL);
 55:   PetscObjectComposeFunction((PetscObject)pc, "PCHPDDMSetAuxiliaryMat_C", NULL);
 56:   PetscObjectComposeFunction((PetscObject)pc, "PCHPDDMHasNeumannMat_C", NULL);
 57:   PetscObjectComposeFunction((PetscObject)pc, "PCHPDDMSetRHSMat_C", NULL);
 58:   PetscObjectComposeFunction((PetscObject)pc, "PCHPDDMSetCoarseCorrectionType_C", NULL);
 59:   PetscObjectComposeFunction((PetscObject)pc, "PCHPDDMGetCoarseCorrectionType_C", NULL);
 60:   PetscObjectComposeFunction((PetscObject)pc, "PCHPDDMGetSTShareSubKSP_C", NULL);
 61:   return 0;
 62: }

 64: static PetscErrorCode PCHPDDMSetAuxiliaryMat_HPDDM(PC pc, IS is, Mat A, PetscErrorCode (*setup)(Mat, PetscReal, Vec, Vec, PetscReal, IS, void*), void* setup_ctx)
 65: {
 66:   PC_HPDDM *data = (PC_HPDDM*)pc->data;

 68:   if (is) {
 69:     PetscObjectReference((PetscObject)is);
 70:     if (data->is) { /* new overlap definition resets the PC */
 71:       PCReset_HPDDM(pc);
 72:       pc->setfromoptionscalled = 0;
 73:     }
 74:     ISDestroy(&data->is);
 75:     data->is = is;
 76:   }
 77:   if (A) {
 78:     PetscObjectReference((PetscObject)A);
 79:     MatDestroy(&data->aux);
 80:     data->aux = A;
 81:   }
 82:   if (setup) {
 83:     data->setup = setup;
 84:     data->setup_ctx = setup_ctx;
 85:   }
 86:   return 0;
 87: }

 89: /*@
 90:      PCHPDDMSetAuxiliaryMat - Sets the auxiliary matrix used by PCHPDDM for the concurrent GenEO problems at the finest level. As an example, in a finite element context with nonoverlapping subdomains plus (overlapping) ghost elements, this could be the unassembled (Neumann) local overlapping operator. As opposed to the assembled (Dirichlet) local overlapping operator obtained by summing neighborhood contributions at the interface of ghost elements.

 92:    Input Parameters:
 93: +     pc - preconditioner context
 94: .     is - index set of the local auxiliary, e.g., Neumann, matrix
 95: .     A - auxiliary sequential matrix
 96: .     setup - function for generating the auxiliary matrix
 97: -     setup_ctx - context for setup

 99:    Level: intermediate

101: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCHPDDMSetRHSMat(), MATIS
102: @*/
103: PetscErrorCode PCHPDDMSetAuxiliaryMat(PC pc, IS is, Mat A, PetscErrorCode (*setup)(Mat, PetscReal, Vec, Vec, PetscReal, IS, void*), void* setup_ctx)
104: {
108: #if defined(PETSC_HAVE_FORTRAN)
109:   if (reinterpret_cast<void*>(setup) == reinterpret_cast<void*>(PETSC_NULL_FUNCTION_Fortran)) setup = NULL;
110:   if (setup_ctx == PETSC_NULL_INTEGER_Fortran) setup_ctx = NULL;
111: #endif
112:   PetscTryMethod(pc, "PCHPDDMSetAuxiliaryMat_C", (PC, IS, Mat, PetscErrorCode (*)(Mat, PetscReal, Vec, Vec, PetscReal, IS, void*), void*), (pc, is, A, setup, setup_ctx));
113:   return 0;
114: }

116: static PetscErrorCode PCHPDDMHasNeumannMat_HPDDM(PC pc, PetscBool has)
117: {
118:   PC_HPDDM *data = (PC_HPDDM*)pc->data;

120:   data->Neumann = has;
121:   return 0;
122: }

124: /*@
125:      PCHPDDMHasNeumannMat - Informs PCHPDDM that the Mat passed to PCHPDDMSetAuxiliaryMat() is the local Neumann matrix. This may be used to bypass a call to MatCreateSubMatrices() and to MatConvert() for MATMPISBAIJ matrices. If a DMCreateNeumannOverlap() implementation is available in the DM attached to the Pmat, or the Amat, or the PC, the flag is internally set to PETSC_TRUE. Its default value is otherwise PETSC_FALSE.

127:    Input Parameters:
128: +     pc - preconditioner context
129: -     has - Boolean value

131:    Level: intermediate

133: .seealso:  PCHPDDM, PCHPDDMSetAuxiliaryMat()
134: @*/
135: PetscErrorCode PCHPDDMHasNeumannMat(PC pc, PetscBool has)
136: {
138:   PetscTryMethod(pc, "PCHPDDMHasNeumannMat_C", (PC, PetscBool), (pc, has));
139:   return 0;
140: }

142: static PetscErrorCode PCHPDDMSetRHSMat_HPDDM(PC pc, Mat B)
143: {
144:   PC_HPDDM *data = (PC_HPDDM*)pc->data;

146:   PetscObjectReference((PetscObject)B);
147:   MatDestroy(&data->B);
148:   data->B = B;
149:   return 0;
150: }

152: /*@
153:      PCHPDDMSetRHSMat - Sets the right-hand side matrix used by PCHPDDM for the concurrent GenEO problems at the finest level. Must be used in conjunction with PCHPDDMSetAuxiliaryMat(N), so that Nv = lambda Bv is solved using EPSSetOperators(N, B). It is assumed that N and B are provided using the same numbering. This provides a means to try more advanced methods such as GenEO-II or H-GenEO.

155:    Input Parameters:
156: +     pc - preconditioner context
157: -     B - right-hand side sequential matrix

159:    Level: advanced

161: .seealso:  PCHPDDMSetAuxiliaryMat(), PCHPDDM
162: @*/
163: PetscErrorCode PCHPDDMSetRHSMat(PC pc, Mat B)
164: {
166:   if (B) {
168:     PetscTryMethod(pc, "PCHPDDMSetRHSMat_C", (PC, Mat), (pc, B));
169:   }
170:   return 0;
171: }

173: static PetscErrorCode PCSetFromOptions_HPDDM(PetscOptionItems *PetscOptionsObject, PC pc)
174: {
175:   PC_HPDDM                    *data = (PC_HPDDM*)pc->data;
176:   PC_HPDDM_Level              **levels = data->levels;
177:   char                        prefix[256];
178:   int                         i = 1;
179:   PetscMPIInt                 size, previous;
180:   PetscInt                    n;
181:   PCHPDDMCoarseCorrectionType type;
182:   PetscBool                   flg = PETSC_TRUE;

184:   if (!data->levels) {
185:     PetscCalloc1(PETSC_PCHPDDM_MAXLEVELS, &levels);
186:     data->levels = levels;
187:   }
188:   PetscOptionsHead(PetscOptionsObject, "PCHPDDM options");
189:   MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size);
190:   previous = size;
191:   while (i < PETSC_PCHPDDM_MAXLEVELS) {
192:     PetscInt p = 1;

194:     if (!data->levels[i - 1]) PetscNewLog(pc, data->levels + i - 1);
195:     data->levels[i - 1]->parent = data;
196:     /* if the previous level has a single process, it is not possible to coarsen further */
197:     if (previous == 1 || !flg) break;
198:     data->levels[i - 1]->nu = 0;
199:     data->levels[i - 1]->threshold = -1.0;
200:     PetscSNPrintf(prefix, sizeof(prefix), "-pc_hpddm_levels_%d_eps_nev", i);
201:     PetscOptionsInt(prefix, "Local number of deflation vectors computed by SLEPc", "EPSSetDimensions", data->levels[i - 1]->nu, &data->levels[i - 1]->nu, NULL);
202:     PetscSNPrintf(prefix, sizeof(prefix), "-pc_hpddm_levels_%d_eps_threshold", i);
203:     PetscOptionsReal(prefix, "Local threshold for selecting deflation vectors returned by SLEPc", "PCHPDDM", data->levels[i - 1]->threshold, &data->levels[i - 1]->threshold, NULL);
204:     if (i == 1) {
205:       PetscSNPrintf(prefix, sizeof(prefix), "-pc_hpddm_levels_1_st_share_sub_ksp");
206:       PetscOptionsBool(prefix, "Shared KSP between SLEPc ST and the fine-level subdomain solver", "PCHPDDMGetSTShareSubKSP", PETSC_FALSE, &data->share, NULL);
207:     }
208:     /* if there is no prescribed coarsening, just break out of the loop */
209:     if (data->levels[i - 1]->threshold <= 0.0 && data->levels[i - 1]->nu <= 0) break;
210:     else {
211:       ++i;
212:       PetscSNPrintf(prefix, sizeof(prefix), "-pc_hpddm_levels_%d_eps_nev", i);
213:       PetscOptionsHasName(PetscOptionsObject->options, PetscOptionsObject->prefix, prefix, &flg);
214:       if (!flg) {
215:         PetscSNPrintf(prefix, sizeof(prefix), "-pc_hpddm_levels_%d_eps_threshold", i);
216:         PetscOptionsHasName(PetscOptionsObject->options, PetscOptionsObject->prefix, prefix, &flg);
217:       }
218:       if (flg) {
219:         /* if there are coarsening options for the next level, then register it  */
220:         /* otherwise, don't to avoid having both options levels_N_p and coarse_p */
221:         PetscSNPrintf(prefix, sizeof(prefix), "-pc_hpddm_levels_%d_p", i);
222:         PetscOptionsRangeInt(prefix, "Number of processes used to assemble the coarse operator at this level", "PCHPDDM", p, &p, &flg, 1, PetscMax(1, previous / 2));
223:         previous = p;
224:       }
225:     }
226:   }
227:   data->N = i;
228:   n = 1;
229:   if (i > 1) {
230:     PetscSNPrintf(prefix, sizeof(prefix), "-pc_hpddm_coarse_p");
231:     PetscOptionsRangeInt(prefix, "Number of processes used to assemble the coarsest operator", "PCHPDDM", n, &n, NULL, 1, PetscMax(1, previous / 2));
232: #if defined(PETSC_HAVE_MUMPS)
233:     PetscSNPrintf(prefix, sizeof(prefix), "pc_hpddm_coarse_");
234:     PetscOptionsHasName(NULL, prefix, "-mat_mumps_use_omp_threads", &flg);
235:     if (flg) {
236:       char type[64]; /* same size as in src/ksp/pc/impls/factor/factimpl.c */
237:       if (n == 1) PetscStrcpy(type, MATSOLVERPETSC); /* default solver for a sequential Mat */
238:       PetscOptionsGetString(NULL, prefix, "-pc_factor_mat_solver_type", type, sizeof(type), &flg);
239:       if (flg) PetscStrcmp(type, MATSOLVERMUMPS, &flg);
241:       size = n;
242:       n = -1;
243:       PetscOptionsGetInt(NULL, prefix, "-mat_mumps_use_omp_threads", &n, NULL);
246:     }
247: #endif
248:     PetscOptionsEnum("-pc_hpddm_coarse_correction", "Type of coarse correction applied each iteration", "PCHPDDMSetCoarseCorrectionType", PCHPDDMCoarseCorrectionTypes, (PetscEnum)data->correction, (PetscEnum*)&type, &flg);
249:     if (flg) PCHPDDMSetCoarseCorrectionType(pc, type);
250:     PetscSNPrintf(prefix, sizeof(prefix), "-pc_hpddm_has_neumann");
251:     PetscOptionsBool(prefix, "Is the auxiliary Mat the local Neumann matrix?", "PCHPDDMHasNeumannMat", data->Neumann, &data->Neumann, NULL);
252:     data->log_separate = PETSC_FALSE;
253:     if (PetscDefined(USE_LOG)) {
254:       PetscSNPrintf(prefix, sizeof(prefix), "-pc_hpddm_log_separate");
255:       PetscOptionsBool(prefix, "Log events level by level instead of inside PCSetUp()/KSPSolve()", NULL, data->log_separate, &data->log_separate, NULL);
256:     }
257:   }
258:   PetscOptionsTail();
259:   while (i < PETSC_PCHPDDM_MAXLEVELS && data->levels[i]) PetscFree(data->levels[i++]);
260:   return 0;
261: }

263: static PetscErrorCode PCApply_HPDDM(PC pc, Vec x, Vec y)
264: {
265:   PC_HPDDM *data = (PC_HPDDM*)pc->data;

267:   PetscCitationsRegister(HPDDMCitation, &HPDDMCite);
268:   if (data->levels[0]->ksp) {
269:     if (data->log_separate) PetscLogEventBegin(PC_HPDDM_Solve[0], data->levels[0]->ksp, 0, 0, 0); /* coarser-level events are directly triggered in HPDDM */
270:     KSPSolve(data->levels[0]->ksp, x, y);
271:     if (data->log_separate) PetscLogEventEnd(PC_HPDDM_Solve[0], data->levels[0]->ksp, 0, 0, 0);
272:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "No KSP attached to PCHPDDM");
273:   return 0;
274: }

276: static PetscErrorCode PCMatApply_HPDDM(PC pc, Mat X, Mat Y)
277: {
278:   PC_HPDDM *data = (PC_HPDDM*)pc->data;

280:   PetscCitationsRegister(HPDDMCitation, &HPDDMCite);
281:   if (data->levels[0]->ksp) KSPMatSolve(data->levels[0]->ksp, X, Y);
282:   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "No KSP attached to PCHPDDM");
283:   return 0;
284: }

286: /*@C
287:      PCHPDDMGetComplexities - Computes the grid and operator complexities.

289:    Input Parameter:
290: .     pc - preconditioner context

292:    Output Parameters:
293: +     gc - grid complexity = sum_i(m_i) / m_1
294: -     oc - operator complexity = sum_i(nnz_i) / nnz_1

296:    Notes:
297:      PCGAMG does not follow the usual convention and names the grid complexity what is usually referred to as the operator complexity. PCHPDDM follows what is found in the literature, and in particular, what you get with PCHYPRE and -pc_hypre_boomeramg_print_statistics.

299:    Level: advanced

301: .seealso:  PCMGGetGridComplexity(), PCHPDDM
302: @*/
303: static PetscErrorCode PCHPDDMGetComplexities(PC pc, PetscReal *gc, PetscReal *oc)
304: {
305:   PC_HPDDM       *data = (PC_HPDDM*)pc->data;
306:   MatInfo        info;
307:   PetscInt       n, m;
308:   PetscLogDouble accumulate[2] { }, nnz1 = 1.0, m1 = 1.0;

310:   for (n = 0, *gc = 0, *oc = 0; n < data->N; ++n) {
311:     if (data->levels[n]->ksp) {
312:       Mat P, A;
313:       KSPGetOperators(data->levels[n]->ksp, NULL, &P);
314:       MatGetSize(P, &m, NULL);
315:       accumulate[0] += m;
316:       if (n == 0) {
317:         PetscBool flg;
318:         PetscObjectTypeCompare((PetscObject)P, MATNORMAL, &flg);
319:         if (flg) {
320:           MatConvert(P, MATAIJ, MAT_INITIAL_MATRIX, &A);
321:           P = A;
322:         } else PetscObjectReference((PetscObject)P);
323:       }
324:       if (P->ops->getinfo) {
325:         MatGetInfo(P, MAT_GLOBAL_SUM, &info);
326:         accumulate[1] += info.nz_used;
327:       }
328:       if (n == 0) {
329:         m1 = m;
330:         if (P->ops->getinfo) nnz1 = info.nz_used;
331:         MatDestroy(&P);
332:       }
333:     }
334:   }
335:   *gc = static_cast<PetscReal>(accumulate[0]/m1);
336:   *oc = static_cast<PetscReal>(accumulate[1]/nnz1);
337:   return 0;
338: }

340: static PetscErrorCode PCView_HPDDM(PC pc, PetscViewer viewer)
341: {
342:   PC_HPDDM     *data = (PC_HPDDM*)pc->data;
343:   PetscViewer  subviewer;
344:   PetscSubcomm subcomm;
345:   PetscReal    oc, gc;
346:   PetscInt     i, tabs;
347:   PetscMPIInt  size, color, rank;
348:   PetscBool    ascii;

350:   PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &ascii);
351:   if (ascii) {
352:     PetscViewerASCIIPrintf(viewer, "level%s: %" PetscInt_FMT "\n", data->N > 1 ? "s" : "", data->N);
353:     PCHPDDMGetComplexities(pc, &gc, &oc);
354:     if (data->N > 1) {
355:       PetscViewerASCIIPrintf(viewer, "Neumann matrix attached? %s\n", PetscBools[data->Neumann]);
356:       PetscViewerASCIIPrintf(viewer, "shared subdomain KSP between SLEPc and PETSc? %s\n", PetscBools[data->share]);
357:       PetscViewerASCIIPrintf(viewer, "coarse correction: %s\n", PCHPDDMCoarseCorrectionTypes[data->correction]);
358:       PetscViewerASCIIPrintf(viewer, "on process #0, value%s (+ threshold%s if available) for selecting deflation vectors:", data->N > 2 ? "s" : "", data->N > 2 ? "s" : "");
359:       PetscViewerASCIIGetTab(viewer, &tabs);
360:       PetscViewerASCIISetTab(viewer, 0);
361:       for (i = 1; i < data->N; ++i) {
362:         PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, data->levels[i - 1]->nu);
363:         if (data->levels[i - 1]->threshold > -0.1) PetscViewerASCIIPrintf(viewer, " (%g)", (double)data->levels[i - 1]->threshold);
364:       }
365:       PetscViewerASCIIPrintf(viewer, "\n");
366:       PetscViewerASCIISetTab(viewer, tabs);
367:     }
368:     PetscViewerASCIIPrintf(viewer, "grid and operator complexities: %g %g\n", (double)gc, (double)oc);
369:     if (data->levels[0]->ksp) {
370:       KSPView(data->levels[0]->ksp, viewer);
371:       if (data->levels[0]->pc) PCView(data->levels[0]->pc, viewer);
372:       for (i = 1; i < data->N; ++i) {
373:         if (data->levels[i]->ksp) color = 1;
374:         else color = 0;
375:         MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size);
376:         MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank);
377:         PetscSubcommCreate(PetscObjectComm((PetscObject)pc), &subcomm);
378:         PetscSubcommSetNumber(subcomm, PetscMin(size, 2));
379:         PetscSubcommSetTypeGeneral(subcomm, color, rank);
380:         PetscViewerASCIIPushTab(viewer);
381:         PetscViewerGetSubViewer(viewer, PetscSubcommChild(subcomm), &subviewer);
382:         if (color == 1) {
383:           KSPView(data->levels[i]->ksp, subviewer);
384:           if (data->levels[i]->pc) PCView(data->levels[i]->pc, subviewer);
385:           PetscViewerFlush(subviewer);
386:         }
387:         PetscViewerRestoreSubViewer(viewer, PetscSubcommChild(subcomm), &subviewer);
388:         PetscViewerASCIIPopTab(viewer);
389:         PetscSubcommDestroy(&subcomm);
390:         PetscViewerFlush(viewer);
391:       }
392:     }
393:   }
394:   return 0;
395: }

397: static PetscErrorCode PCPreSolve_HPDDM(PC pc, KSP ksp, Vec, Vec)
398: {
399:   PC_HPDDM  *data = (PC_HPDDM*)pc->data;
400:   PetscBool flg;
401:   Mat       A;

403:   if (ksp) {
404:     PetscObjectTypeCompare((PetscObject)ksp, KSPLSQR, &flg);
405:     if (flg && !data->normal) {
406:       KSPGetOperators(ksp, &A, NULL);
407:       MatCreateVecs(A, NULL, &data->normal); /* temporary Vec used in PCHPDDMShellApply() for coarse grid corrections */
408:     }
409:   }
410:   return 0;
411: }

413: static PetscErrorCode PCHPDDMShellSetUp(PC pc)
414: {
415:   PC_HPDDM_Level *ctx;
416:   Mat            A, P;
417:   Vec            x;
418:   const char     *pcpre;

420:   PCShellGetContext(pc, &ctx);
421:   KSPGetOptionsPrefix(ctx->ksp, &pcpre);
422:   KSPGetOperators(ctx->ksp, &A, &P);
423:   /* smoother */
424:   PCSetOptionsPrefix(ctx->pc, pcpre);
425:   PCSetOperators(ctx->pc, A, P);
426:   if (!ctx->v[0]) {
427:     VecDuplicateVecs(ctx->D, 1, &ctx->v[0]);
428:     if (!std::is_same<PetscScalar, PetscReal>::value) VecDestroy(&ctx->D);
429:     MatCreateVecs(A, &x, NULL);
430:     VecDuplicateVecs(x, 2, &ctx->v[1]);
431:     VecDestroy(&x);
432:   }
433:   std::fill_n(ctx->V, 3, nullptr);
434:   return 0;
435: }

437: template<class Type, typename std::enable_if<std::is_same<Type, Vec>::value>::type* = nullptr>
438: static inline PetscErrorCode PCHPDDMDeflate_Private(PC pc, Type x, Type y)
439: {
440:   PC_HPDDM_Level *ctx;
441:   PetscScalar    *out;

443:   PCShellGetContext(pc, &ctx);
444:   /* going from PETSc to HPDDM numbering */
445:   VecScatterBegin(ctx->scatter, x, ctx->v[0][0], INSERT_VALUES, SCATTER_FORWARD);
446:   VecScatterEnd(ctx->scatter, x, ctx->v[0][0], INSERT_VALUES, SCATTER_FORWARD);
447:   VecGetArrayWrite(ctx->v[0][0], &out);
448:   ctx->P->deflation<false>(NULL, out, 1); /* y = Q x */
449:   VecRestoreArrayWrite(ctx->v[0][0], &out);
450:   /* going from HPDDM to PETSc numbering */
451:   VecScatterBegin(ctx->scatter, ctx->v[0][0], y, INSERT_VALUES, SCATTER_REVERSE);
452:   VecScatterEnd(ctx->scatter, ctx->v[0][0], y, INSERT_VALUES, SCATTER_REVERSE);
453:   return 0;
454: }

456: template<class Type, typename std::enable_if<std::is_same<Type, Mat>::value>::type* = nullptr>
457: static inline PetscErrorCode PCHPDDMDeflate_Private(PC pc, Type X, Type Y)
458: {
459:   PC_HPDDM_Level *ctx;
460:   Vec            vX, vY, vC;
461:   PetscScalar    *out;
462:   PetscInt       i, N;

464:   PCShellGetContext(pc, &ctx);
465:   MatGetSize(X, NULL, &N);
466:   /* going from PETSc to HPDDM numbering */
467:   for (i = 0; i < N; ++i) {
468:     MatDenseGetColumnVecRead(X, i, &vX);
469:     MatDenseGetColumnVecWrite(ctx->V[0], i, &vC);
470:     VecScatterBegin(ctx->scatter, vX, vC, INSERT_VALUES, SCATTER_FORWARD);
471:     VecScatterEnd(ctx->scatter, vX, vC, INSERT_VALUES, SCATTER_FORWARD);
472:     MatDenseRestoreColumnVecWrite(ctx->V[0], i, &vC);
473:     MatDenseRestoreColumnVecRead(X, i, &vX);
474:   }
475:   MatDenseGetArrayWrite(ctx->V[0], &out);
476:   ctx->P->deflation<false>(NULL, out, N); /* Y = Q X */
477:   MatDenseRestoreArrayWrite(ctx->V[0], &out);
478:   /* going from HPDDM to PETSc numbering */
479:   for (i = 0; i < N; ++i) {
480:     MatDenseGetColumnVecRead(ctx->V[0], i, &vC);
481:     MatDenseGetColumnVecWrite(Y, i, &vY);
482:     VecScatterBegin(ctx->scatter, vC, vY, INSERT_VALUES, SCATTER_REVERSE);
483:     VecScatterEnd(ctx->scatter, vC, vY, INSERT_VALUES, SCATTER_REVERSE);
484:     MatDenseRestoreColumnVecWrite(Y, i, &vY);
485:     MatDenseRestoreColumnVecRead(ctx->V[0], i, &vC);
486:   }
487:   return 0;
488: }

490: /*@C
491:      PCHPDDMShellApply - Applies a (2) deflated, (1) additive, or (3) balanced coarse correction. In what follows, E = Z Pmat Z^T and Q = Z^T E^-1 Z.

493: .vb
494:    (1) y =                Pmat^-1              x + Q x,
495:    (2) y =                Pmat^-1 (I - Amat Q) x + Q x (default),
496:    (3) y = (I - Q Amat^T) Pmat^-1 (I - Amat Q) x + Q x.
497: .ve

499:    Input Parameters:
500: +     pc - preconditioner context
501: -     x - input vector

503:    Output Parameter:
504: .     y - output vector

506:    Application Interface Routine: PCApply()

508:    Notes:
509:      The options of Pmat^1 = pc(Pmat) are prefixed by -pc_hpddm_levels_1_pc_. Z is a tall-and-skiny matrix assembled by HPDDM. The number of processes on which (Z Pmat Z^T) is aggregated is set via -pc_hpddm_coarse_p.
510:      The options of (Z Pmat Z^T)^-1 = ksp(Z Pmat Z^T) are prefixed by -pc_hpddm_coarse_ (KSPPREONLY and PCCHOLESKY by default), unless a multilevel correction is turned on, in which case, this function is called recursively at each level except the coarsest one.
511:      (1) and (2) visit the "next" level (in terms of coarsening) once per application, while (3) visits it twice, so it is asymptotically twice costlier. (2) is not symmetric even if both Amat and Pmat are symmetric.

513:    Level: advanced

515: .seealso:  PCHPDDM, PCHPDDMCoarseCorrectionType
516: @*/
517: static PetscErrorCode PCHPDDMShellApply(PC pc, Vec x, Vec y)
518: {
519:   PC_HPDDM_Level *ctx;
520:   Mat            A;

522:   PCShellGetContext(pc, &ctx);
523:   if (ctx->P) {
524:     KSPGetOperators(ctx->ksp, &A, NULL);
525:     PCHPDDMDeflate_Private(pc, x, y);                    /* y = Q x                          */
526:     if (ctx->parent->correction == PC_HPDDM_COARSE_CORRECTION_DEFLATED || ctx->parent->correction == PC_HPDDM_COARSE_CORRECTION_BALANCED) {
527:       if (!ctx->parent->normal || ctx != ctx->parent->levels[0]) MatMult(A, y, ctx->v[1][0]); /* y = A Q x */
528:       else { /* KSPLSQR and finest level */
529:         MatMult(A, y, ctx->parent->normal);              /* y = A Q x                        */
530:         MatMultTranspose(A, ctx->parent->normal, ctx->v[1][0]); /* y = A^T A Q x             */
531:       }
532:       VecWAXPY(ctx->v[1][1], -1.0, ctx->v[1][0], x);     /* y = (I - A Q) x                  */
533:       PCApply(ctx->pc, ctx->v[1][1], ctx->v[1][0]);      /* y = M^-1 (I - A Q) x             */
534:       if (ctx->parent->correction == PC_HPDDM_COARSE_CORRECTION_BALANCED) {
535:         if (!ctx->parent->normal || ctx != ctx->parent->levels[0]) MatMultTranspose(A, ctx->v[1][0], ctx->v[1][1]); /* z = A^T y */
536:         else {
537:           MatMult(A, ctx->v[1][0], ctx->parent->normal);
538:           MatMultTranspose(A, ctx->parent->normal, ctx->v[1][1]); /* z = A^T A y             */
539:         }
540:         PCHPDDMDeflate_Private(pc, ctx->v[1][1], ctx->v[1][1]);
541:         VecAXPBYPCZ(y, -1.0, 1.0, 1.0, ctx->v[1][1], ctx->v[1][0]); /* y = (I - Q A^T) y + Q x */
542:       } else VecAXPY(y, 1.0, ctx->v[1][0]);              /* y = Q M^-1 (I - A Q) x + Q x     */
543:     } else if (ctx->parent->correction == PC_HPDDM_COARSE_CORRECTION_ADDITIVE) {
544:       PCApply(ctx->pc, x, ctx->v[1][0]);
545:       VecAXPY(y, 1.0, ctx->v[1][0]);                     /* y = M^-1 x + Q x                 */
546:     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCSHELL from PCHPDDM called with an unknown PCHPDDMCoarseCorrectionType %d", ctx->parent->correction);
547:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCSHELL from PCHPDDM called with no HPDDM object");
548:   return 0;
549: }

551: /*@C
552:      PCHPDDMShellMatApply - Variant of PCHPDDMShellApply() for blocks of vectors.

554:    Input Parameters:
555: +     pc - preconditioner context
556: -     X - block of input vectors

558:    Output Parameter:
559: .     Y - block of output vectors

561:    Application Interface Routine: PCApply()

563:    Level: advanced

565: .seealso:  PCHPDDM, PCHPDDMShellMatApply(), PCHPDDMCoarseCorrectionType
566: @*/
567: static PetscErrorCode PCHPDDMShellMatApply(PC pc, Mat X, Mat Y)
568: {
569:   PC_HPDDM_Level *ctx;
570:   Mat            A, *ptr;
571:   PetscContainer container = NULL;
572:   PetscScalar    *array;
573:   PetscInt       m, M, N, prev = 0;
574:   PetscBool      reset = PETSC_FALSE;

576:   PCShellGetContext(pc, &ctx);
577:   if (ctx->P) {
578:     MatGetSize(X, NULL, &N);
579:     KSPGetOperators(ctx->ksp, &A, NULL);
580:     PetscObjectQuery((PetscObject)A, "_HPDDM_MatProduct", (PetscObject*)&container);
581:     if (container) { /* MatProduct container already attached */
582:       PetscContainerGetPointer(container, (void**)&ptr);
583:       if (ptr[1] != ctx->V[2]) /* Mat has changed or may have been set first in KSPHPDDM */
584:         for (m = 0; m < 2; ++m) {
585:           MatDestroy(ctx->V + m + 1);
586:           ctx->V[m + 1] = ptr[m];
587:           PetscObjectReference((PetscObject)ctx->V[m + 1]);
588:         }
589:     }
590:     if (ctx->V[1]) MatGetSize(ctx->V[1], NULL, &prev);
591:     if (N != prev || !ctx->V[0]) {
592:       MatDestroy(ctx->V);
593:       VecGetLocalSize(ctx->v[0][0], &m);
594:       MatCreateDense(PetscObjectComm((PetscObject)pc), m, PETSC_DECIDE, PETSC_DECIDE, N, NULL, ctx->V);
595:       if (N != prev) {
596:         MatDestroy(ctx->V + 1);
597:         MatDestroy(ctx->V + 2);
598:         MatGetLocalSize(X, &m, NULL);
599:         MatGetSize(X, &M, NULL);
600:         if (ctx->parent->correction != PC_HPDDM_COARSE_CORRECTION_BALANCED) MatDenseGetArrayWrite(ctx->V[0], &array);
601:         else array = NULL;
602:         MatCreateDense(PetscObjectComm((PetscObject)pc), m, PETSC_DECIDE, M, N, array, ctx->V + 1);
603:         if (ctx->parent->correction != PC_HPDDM_COARSE_CORRECTION_BALANCED) MatDenseRestoreArrayWrite(ctx->V[0], &array);
604:         else {
605:           MatAssemblyBegin(ctx->V[1], MAT_FINAL_ASSEMBLY);
606:           MatAssemblyEnd(ctx->V[1], MAT_FINAL_ASSEMBLY);
607:         }
608:         MatCreateDense(PetscObjectComm((PetscObject)pc), m, PETSC_DECIDE, M, N, NULL, ctx->V + 2);
609:         MatAssemblyBegin(ctx->V[2], MAT_FINAL_ASSEMBLY);
610:         MatAssemblyEnd(ctx->V[2], MAT_FINAL_ASSEMBLY);
611:         MatProductCreateWithMat(A, Y, NULL, ctx->V[1]);
612:         MatProductSetType(ctx->V[1], MATPRODUCT_AB);
613:         MatProductSetFromOptions(ctx->V[1]);
614:         MatProductSymbolic(ctx->V[1]);
615:         if (!container) { /* no MatProduct container attached, create one to be queried in KSPHPDDM or at the next call to PCMatApply() */
616:           PetscContainerCreate(PetscObjectComm((PetscObject)A), &container);
617:           PetscObjectCompose((PetscObject)A, "_HPDDM_MatProduct", (PetscObject)container);
618:         }
619:         PetscContainerSetPointer(container, ctx->V + 1); /* need to compose B and D from MatProductCreateWithMath(A, B, NULL, D), which are stored in the contiguous array ctx->V */
620:       }
621:       if (ctx->parent->correction == PC_HPDDM_COARSE_CORRECTION_BALANCED) {
622:         MatProductCreateWithMat(A, ctx->V[1], NULL, ctx->V[2]);
623:         MatProductSetType(ctx->V[2], MATPRODUCT_AtB);
624:         MatProductSetFromOptions(ctx->V[2]);
625:         MatProductSymbolic(ctx->V[2]);
626:       }
627:       ctx->P->start(N);
628:     }
629:     if (N == prev || container) { /* when MatProduct container is attached, always need to MatProductReplaceMats() since KSPHPDDM may have replaced the Mat as well */
630:       MatProductReplaceMats(NULL, Y, NULL, ctx->V[1]);
631:       if (container && ctx->parent->correction != PC_HPDDM_COARSE_CORRECTION_BALANCED) {
632:         MatDenseGetArrayWrite(ctx->V[0], &array);
633:         MatDensePlaceArray(ctx->V[1], array);
634:         MatDenseRestoreArrayWrite(ctx->V[0], &array);
635:         reset = PETSC_TRUE;
636:       }
637:     }
638:     PCHPDDMDeflate_Private(pc, X, Y);
639:     if (ctx->parent->correction == PC_HPDDM_COARSE_CORRECTION_DEFLATED || ctx->parent->correction == PC_HPDDM_COARSE_CORRECTION_BALANCED) {
640:       MatProductNumeric(ctx->V[1]);
641:       MatCopy(ctx->V[1], ctx->V[2], SAME_NONZERO_PATTERN);
642:       MatAXPY(ctx->V[2], -1.0, X, SAME_NONZERO_PATTERN);
643:       PCMatApply(ctx->pc, ctx->V[2], ctx->V[1]);
644:       if (ctx->parent->correction == PC_HPDDM_COARSE_CORRECTION_BALANCED) {
645:         MatProductNumeric(ctx->V[2]);
646:         PCHPDDMDeflate_Private(pc, ctx->V[2], ctx->V[2]);
647:         MatAXPY(ctx->V[1], -1.0, ctx->V[2], SAME_NONZERO_PATTERN);
648:       }
649:       MatAXPY(Y, -1.0, ctx->V[1], SAME_NONZERO_PATTERN);
650:     } else if (ctx->parent->correction == PC_HPDDM_COARSE_CORRECTION_ADDITIVE) {
651:       PCMatApply(ctx->pc, X, ctx->V[1]);
652:       MatAXPY(Y, 1.0, ctx->V[1], SAME_NONZERO_PATTERN);
653:     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCSHELL from PCHPDDM called with an unknown PCHPDDMCoarseCorrectionType %d", ctx->parent->correction);
654:     if (reset) MatDenseResetArray(ctx->V[1]);
655:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCSHELL from PCHPDDM called with no HPDDM object");
656:   return 0;
657: }

659: static PetscErrorCode PCHPDDMShellDestroy(PC pc)
660: {
661:   PC_HPDDM_Level *ctx;
662:   PetscContainer container;

664:   PCShellGetContext(pc, &ctx);
665:   HPDDM::Schwarz<PetscScalar>::destroy(ctx, PETSC_TRUE);
666:   VecDestroyVecs(1, &ctx->v[0]);
667:   VecDestroyVecs(2, &ctx->v[1]);
668:   PetscObjectQuery((PetscObject)(ctx->pc)->mat, "_HPDDM_MatProduct", (PetscObject*)&container);
669:   PetscContainerDestroy(&container);
670:   PetscObjectCompose((PetscObject)(ctx->pc)->mat, "_HPDDM_MatProduct", NULL);
671:   MatDestroy(ctx->V);
672:   MatDestroy(ctx->V + 1);
673:   MatDestroy(ctx->V + 2);
674:   VecDestroy(&ctx->D);
675:   VecScatterDestroy(&ctx->scatter);
676:   PCDestroy(&ctx->pc);
677:   return 0;
678: }

680: static PetscErrorCode PCHPDDMSolve_Private(const PC_HPDDM_Level *ctx, PetscScalar *rhs, const unsigned short& mu)
681: {
682:   Mat      B, X;
683:   PetscInt n, N, j = 0;

685:   KSPGetOperators(ctx->ksp, &B, NULL);
686:   MatGetLocalSize(B, &n, NULL);
687:   MatGetSize(B, &N, NULL);
688:   if (ctx->parent->log_separate) {
689:     j = std::distance(ctx->parent->levels, std::find(ctx->parent->levels, ctx->parent->levels + ctx->parent->N, ctx));
690:     PetscLogEventBegin(PC_HPDDM_Solve[j], ctx->ksp, 0, 0, 0);
691:   }
692:   if (mu == 1) {
693:     if (!ctx->ksp->vec_rhs) {
694:       VecCreateMPIWithArray(PetscObjectComm((PetscObject)ctx->ksp), 1, n, N, NULL, &ctx->ksp->vec_rhs);
695:       VecCreateMPI(PetscObjectComm((PetscObject)ctx->ksp), n, N, &ctx->ksp->vec_sol);
696:     }
697:     VecPlaceArray(ctx->ksp->vec_rhs, rhs);
698:     KSPSolve(ctx->ksp, NULL, NULL);
699:     VecCopy(ctx->ksp->vec_sol, ctx->ksp->vec_rhs);
700:     VecResetArray(ctx->ksp->vec_rhs);
701:   } else {
702:     MatCreateDense(PetscObjectComm((PetscObject)ctx->ksp), n, PETSC_DECIDE, N, mu, rhs, &B);
703:     MatCreateDense(PetscObjectComm((PetscObject)ctx->ksp), n, PETSC_DECIDE, N, mu, NULL, &X);
704:     KSPMatSolve(ctx->ksp, B, X);
705:     MatCopy(X, B, SAME_NONZERO_PATTERN);
706:     MatDestroy(&X);
707:     MatDestroy(&B);
708:   }
709:   if (ctx->parent->log_separate) PetscLogEventEnd(PC_HPDDM_Solve[j], ctx->ksp, 0, 0, 0);
710:   return 0;
711: }

713: static PetscErrorCode PCHPDDMSetUpNeumannOverlap_Private(PC pc)
714: {
715:   PC_HPDDM *data = (PC_HPDDM*)pc->data;

717:   if (data->setup) {
718:     Mat       P;
719:     Vec       x, xt = NULL;
720:     PetscReal t = 0.0, s = 0.0;

722:     PCGetOperators(pc, NULL, &P);
723:     PetscObjectQuery((PetscObject)P, "__SNES_latest_X", (PetscObject*)&x);
724:     PetscStackPush("PCHPDDM Neumann callback");
725:     (*data->setup)(data->aux, t, x, xt, s, data->is, data->setup_ctx);
726:     PetscStackPop;
727:   }
728:   return 0;
729: }

731: static PetscErrorCode PCHPDDMCreateSubMatrices_Private(Mat mat, PetscInt n, const IS*, const IS*, MatReuse scall, Mat *submat[])
732: {
733:   Mat A;

736:   /* previously composed Mat */
737:   PetscObjectQuery((PetscObject)mat, "_PCHPDDM_SubMatrices", (PetscObject*)&A);
739:   if (scall == MAT_INITIAL_MATRIX) {
740:     PetscCalloc1(1, submat);
741:     MatDuplicate(A, MAT_COPY_VALUES, *submat);
742:   } else MatCopy(A, (*submat)[0], SAME_NONZERO_PATTERN);
743:   return 0;
744: }

746: static PetscErrorCode PCHPDDMCommunicationAvoidingPCASM_Private(PC pc, Mat C, PetscBool sorted)
747: {
748:   void (*op)(void);

750:   /* previously-composed Mat */
751:   PetscObjectCompose((PetscObject)pc->pmat, "_PCHPDDM_SubMatrices", (PetscObject)C);
752:   MatGetOperation(pc->pmat, MATOP_CREATE_SUBMATRICES, &op);
753:   /* trick suggested by Barry https://lists.mcs.anl.gov/pipermail/petsc-dev/2020-January/025491.html */
754:   MatSetOperation(pc->pmat, MATOP_CREATE_SUBMATRICES, (void(*)(void))PCHPDDMCreateSubMatrices_Private);
755:   if (sorted) PCASMSetSortIndices(pc, PETSC_FALSE); /* everything is already sorted */
756:   PCSetFromOptions(pc); /* otherwise -pc_hpddm_levels_1_pc_asm_sub_mat_type is not used */
757:   PCSetUp(pc);
758:   /* reset MatCreateSubMatrices() */
759:   MatSetOperation(pc->pmat, MATOP_CREATE_SUBMATRICES, op);
760:   PetscObjectCompose((PetscObject)pc->pmat, "_PCHPDDM_SubMatrices", NULL);
761:   return 0;
762: }

764: static PetscErrorCode PCHPDDMPermute_Private(IS is, IS in_is, IS *out_is, Mat in_C, Mat *out_C, IS *p)
765: {
766:   IS                           perm;
767:   const PetscInt               *ptr;
768:   PetscInt                     *concatenate, size, n;
769:   std::map<PetscInt, PetscInt> order;
770:   PetscBool                    sorted;

772:   ISSorted(is, &sorted);
773:   if (!sorted) {
774:     ISGetLocalSize(is, &size);
775:     ISGetIndices(is, &ptr);
776:     /* MatCreateSubMatrices(), called by PCASM, follows the global numbering of Pmat */
777:     for (n = 0; n < size; ++n) order.insert(std::make_pair(ptr[n], n));
778:     ISRestoreIndices(is, &ptr);
779:     if (out_C) {
780:       PetscMalloc1(size, &concatenate);
781:       for (const std::pair<const PetscInt, PetscInt>& i : order) *concatenate++ = i.second;
782:       concatenate -= size;
783:       ISCreateGeneral(PetscObjectComm((PetscObject)in_C), size, concatenate, PETSC_OWN_POINTER, &perm);
784:       ISSetPermutation(perm);
785:       /* permute user-provided Mat so that it matches with MatCreateSubMatrices() numbering */
786:       MatPermute(in_C, perm, perm, out_C);
787:       if (p) *p = perm;
788:       else ISDestroy(&perm); /* no need to save the permutation */
789:     }
790:     if (out_is) {
791:       PetscMalloc1(size, &concatenate);
792:       for (const std::pair<const PetscInt, PetscInt>& i : order) *concatenate++ = i.first;
793:       concatenate -= size;
794:       /* permute user-provided IS so that it matches with MatCreateSubMatrices() numbering */
795:       ISCreateGeneral(PetscObjectComm((PetscObject)in_is), size, concatenate, PETSC_OWN_POINTER, out_is);
796:     }
797:   } else { /* input IS is sorted, nothing to permute, simply duplicate inputs when needed */
798:     if (out_C) MatDuplicate(in_C, MAT_COPY_VALUES, out_C);
799:     if (out_is) ISDuplicate(in_is, out_is);
800:   }
801:   return 0;
802: }

804: static PetscErrorCode PCHPDDMDestroySubMatrices_Private(PetscBool flg, PetscBool algebraic, Mat *sub)
805: {
806:   IS is;

808:   if (!flg) {
809:     if (algebraic) {
810:       PetscObjectQuery((PetscObject)sub[0], "_PCHPDDM_Embed", (PetscObject*)&is);
811:       ISDestroy(&is);
812:       PetscObjectCompose((PetscObject)sub[0], "_PCHPDDM_Embed", NULL);
813:       PetscObjectCompose((PetscObject)sub[0], "_PCHPDDM_Compact", NULL);
814:     }
815:     MatDestroySubMatrices(algebraic ? 2 : 1, &sub);
816:   }
817:   return 0;
818: }

820: static PetscErrorCode PCHPDDMAlgebraicAuxiliaryMat_Private(Mat P, IS *is, Mat *sub[], PetscBool block)
821: {
822:   IS        icol[3], irow[2];
823:   Mat       *M, Q;
824:   PetscReal *ptr;
825:   PetscInt  *idx, p = 0, n;
826:   PetscBool flg;

828:   ISCreateStride(PETSC_COMM_SELF, P->cmap->N, 0, 1, icol + 2);
829:   ISSetBlockSize(icol[2], P->cmap->bs);
830:   ISSetIdentity(icol[2]);
831:   PetscObjectTypeCompare((PetscObject)P, MATMPISBAIJ, &flg);
832:   if (flg) {
833:     /* MatCreateSubMatrices() does not handle MATMPISBAIJ properly when iscol != isrow, so convert first to MATMPIBAIJ */
834:     MatConvert(P, MATMPIBAIJ, MAT_INITIAL_MATRIX, &Q);
835:     std::swap(P, Q);
836:   }
837:   MatCreateSubMatrices(P, 1, is, icol + 2, MAT_INITIAL_MATRIX, &M);
838:   if (flg) {
839:     std::swap(P, Q);
840:     MatDestroy(&Q);
841:   }
842:   ISDestroy(icol + 2);
843:   ISCreateStride(PETSC_COMM_SELF, M[0]->rmap->N, 0, 1, irow);
844:   ISSetBlockSize(irow[0], P->cmap->bs);
845:   ISSetIdentity(irow[0]);
846:   if (!block) {
847:     PetscMalloc2(P->cmap->N, &ptr, P->cmap->N, &idx);
848:     MatGetColumnNorms(M[0], NORM_INFINITY, ptr);
849:     /* check for nonzero columns so that M[0] may be expressed in compact form */
850:     for (n = 0; n < P->cmap->N; n += P->cmap->bs)
851:       if (std::find_if(ptr + n, ptr + n + P->cmap->bs, [](PetscReal v) { return v > PETSC_MACHINE_EPSILON; }) != ptr + n + P->cmap->bs) {
852:         std::iota(idx + p, idx + p + P->cmap->bs, n);
853:         p += P->cmap->bs;
854:       }
855:     ISCreateGeneral(PETSC_COMM_SELF, p, idx, PETSC_USE_POINTER, icol + 1);
856:     ISSetBlockSize(icol[1], P->cmap->bs);
857:     ISSetInfo(icol[1], IS_SORTED, IS_GLOBAL, PETSC_TRUE, PETSC_TRUE);
858:     ISEmbed(*is, icol[1], PETSC_FALSE, icol + 2);
859:     irow[1] = irow[0];
860:     /* first Mat will be used in PCASM (if it is used as a PC on this level) and as the left-hand side of GenEO */
861:     icol[0] = is[0];
862:     MatCreateSubMatrices(M[0], 2, irow, icol, MAT_INITIAL_MATRIX, sub);
863:     ISDestroy(icol + 1);
864:     PetscFree2(ptr, idx);
865:     /* IS used to go back and forth between the augmented and the original local linear system, see eq. (3.4) of [2022b] */
866:     PetscObjectCompose((PetscObject)(*sub)[0], "_PCHPDDM_Embed", (PetscObject)icol[2]);
867:     /* Mat used in eq. (3.1) of [2022b] */
868:     PetscObjectCompose((PetscObject)(*sub)[0], "_PCHPDDM_Compact", (PetscObject)(*sub)[1]);
869:   } else {
870:     Mat aux;
871:     MatSetOption(M[0], MAT_SUBMAT_SINGLEIS, PETSC_TRUE);
872:     /* diagonal block of the overlapping rows */
873:     MatCreateSubMatrices(M[0], 1, irow, is, MAT_INITIAL_MATRIX, sub);
874:     MatDuplicate((*sub)[0], MAT_COPY_VALUES, &aux);
875:     MatSetOption(aux, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);
876:     if (P->cmap->bs == 1) { /* scalar case */
877:       Vec sum[2];
878:       MatCreateVecs(aux, sum, sum + 1);
879:       MatGetRowSum(M[0], sum[0]);
880:       MatGetRowSum(aux, sum[1]);
881:       /* off-diagonal block row sum (full rows - diagonal block rows) */
882:       VecAXPY(sum[0], -1.0, sum[1]);
883:       /* subdomain matrix plus off-diagonal block row sum */
884:       MatDiagonalSet(aux, sum[0], ADD_VALUES);
885:       VecDestroy(sum);
886:       VecDestroy(sum + 1);
887:     } else { /* vectorial case */
888:       /* TODO: missing MatGetValuesBlocked(), so the code below is     */
889:       /* an extension of the scalar case for when bs > 1, but it could */
890:       /* be more efficient by avoiding all these MatMatMult()          */
891:       Mat         sum[2], ones;
892:       PetscScalar *ptr;
893:       PetscCalloc1(M[0]->cmap->n * P->cmap->bs, &ptr);
894:       MatCreateDense(PETSC_COMM_SELF, M[0]->cmap->n, P->cmap->bs, M[0]->cmap->n, P->cmap->bs, ptr, &ones);
895:       for (n = 0; n < M[0]->cmap->n; n += P->cmap->bs) {
896:         for (p = 0; p < P->cmap->bs; ++p) ptr[n + p * (M[0]->cmap->n + 1)] = 1.0;
897:       }
898:       MatMatMult(M[0], ones, MAT_INITIAL_MATRIX, PETSC_DEFAULT, sum);
899:       MatDestroy(&ones);
900:       MatCreateDense(PETSC_COMM_SELF, aux->cmap->n, P->cmap->bs, aux->cmap->n, P->cmap->bs, ptr, &ones);
901:       MatDenseSetLDA(ones, M[0]->cmap->n);
902:       MatMatMult(aux, ones, MAT_INITIAL_MATRIX, PETSC_DEFAULT, sum + 1);
903:       MatDestroy(&ones);
904:       PetscFree(ptr);
905:       /* off-diagonal block row sum (full rows - diagonal block rows) */
906:       MatAXPY(sum[0], -1.0, sum[1], SAME_NONZERO_PATTERN);
907:       MatDestroy(sum + 1);
908:       /* re-order values to be consistent with MatSetValuesBlocked()           */
909:       /* equivalent to MatTranspose() which does not truly handle              */
910:       /* MAT_INPLACE_MATRIX in the rectangular case, as it calls PetscMalloc() */
911:       MatDenseGetArrayWrite(sum[0], &ptr);
912:       HPDDM::Wrapper<PetscScalar>::imatcopy<'T'>(P->cmap->bs, sum[0]->rmap->n, ptr, sum[0]->rmap->n, P->cmap->bs);
913:       /* subdomain matrix plus off-diagonal block row sum */
914:       for (n = 0; n < aux->cmap->n / P->cmap->bs; ++n) MatSetValuesBlocked(aux, 1, &n, 1, &n, ptr + n * P->cmap->bs * P->cmap->bs, ADD_VALUES);
915:       MatAssemblyBegin(aux, MAT_FINAL_ASSEMBLY);
916:       MatAssemblyEnd(aux, MAT_FINAL_ASSEMBLY);
917:       MatDenseRestoreArrayWrite(sum[0], &ptr);
918:       MatDestroy(sum);
919:     }
920:     MatSetOption(aux, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE);
921:     /* left-hand side of GenEO, with the same sparsity pattern as PCASM subdomain solvers  */
922:     PetscObjectCompose((PetscObject)(*sub)[0], "_PCHPDDM_Neumann_Mat", (PetscObject)aux);
923:   }
924:   ISDestroy(irow);
925:   MatDestroySubMatrices(1, &M);
926:   return 0;
927: }

929: static PetscErrorCode PCSetUp_HPDDM(PC pc)
930: {
931:   PC_HPDDM                 *data = (PC_HPDDM*)pc->data;
932:   PC                       inner;
933:   KSP                      *ksp;
934:   Mat                      *sub, A, P, N, C = NULL, uaux = NULL, weighted, subA[2];
935:   Vec                      xin, v;
936:   std::vector<Vec>         initial;
937:   IS                       is[1], loc, uis = data->is;
938:   ISLocalToGlobalMapping   l2g;
939:   char                     prefix[256];
940:   const char               *pcpre;
941:   const PetscScalar *const *ev;
942:   PetscInt                 n, requested = data->N, reused = 0;
943:   MatStructure             structure = UNKNOWN_NONZERO_PATTERN;
944:   PetscBool                subdomains = PETSC_FALSE, flg = PETSC_FALSE, ismatis, swap = PETSC_FALSE, algebraic = PETSC_FALSE, block = PETSC_FALSE;
945:   DM                       dm;

948:   PCGetOptionsPrefix(pc, &pcpre);
949:   PCGetOperators(pc, &A, &P);
950:   if (!data->levels[0]->ksp) {
951:     KSPCreate(PetscObjectComm((PetscObject)pc), &data->levels[0]->ksp);
952:     PetscSNPrintf(prefix, sizeof(prefix), "%spc_hpddm_%s_", pcpre ? pcpre : "", data->N > 1 ? "levels_1" : "coarse");
953:     KSPSetOptionsPrefix(data->levels[0]->ksp, prefix);
954:     KSPSetType(data->levels[0]->ksp, KSPPREONLY);
955:   } else if (data->levels[0]->ksp->pc && data->levels[0]->ksp->pc->setupcalled == 1 && data->levels[0]->ksp->pc->reusepreconditioner) {
956:     /* if the fine-level PCSHELL exists, its setup has succeeded, and one wants to reuse it, */
957:     /* then just propagate the appropriate flag to the coarser levels                        */
958:     for (n = 0; n < PETSC_PCHPDDM_MAXLEVELS && data->levels[n]; ++n) {
959:       /* the following KSP and PC may be NULL for some processes, hence the check            */
960:       if (data->levels[n]->ksp) KSPSetReusePreconditioner(data->levels[n]->ksp, PETSC_TRUE);
961:       if (data->levels[n]->pc) PCSetReusePreconditioner(data->levels[n]->pc, PETSC_TRUE);
962:     }
963:     /* early bail out because there is nothing to do */
964:     return 0;
965:   } else {
966:     /* reset coarser levels */
967:     for (n = 1; n < PETSC_PCHPDDM_MAXLEVELS && data->levels[n]; ++n) {
968:       if (data->levels[n]->ksp && data->levels[n]->ksp->pc && data->levels[n]->ksp->pc->setupcalled == 1 && data->levels[n]->ksp->pc->reusepreconditioner && n < data->N) {
969:         reused = data->N - n;
970:         break;
971:       }
972:       KSPDestroy(&data->levels[n]->ksp);
973:       PCDestroy(&data->levels[n]->pc);
974:     }
975:     /* check if some coarser levels are being reused */
976:     MPIU_Allreduce(MPI_IN_PLACE, &reused, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)pc));
977:     const int *addr = data->levels[0]->P ? data->levels[0]->P->getAddrLocal() : &HPDDM::i__0;

979:     if (addr != &HPDDM::i__0 && reused != data->N - 1) {
980:       /* reuse previously computed eigenvectors */
981:       ev = data->levels[0]->P->getVectors();
982:       if (ev) {
983:         initial.reserve(*addr);
984:         VecCreateSeqWithArray(PETSC_COMM_SELF, 1, data->levels[0]->P->getDof(), ev[0], &xin);
985:         for (n = 0; n < *addr; ++n) {
986:           VecDuplicate(xin, &v);
987:           VecPlaceArray(xin, ev[n]);
988:           VecCopy(xin, v);
989:           initial.emplace_back(v);
990:           VecResetArray(xin);
991:         }
992:         VecDestroy(&xin);
993:       }
994:     }
995:   }
996:   data->N -= reused;
997:   KSPSetOperators(data->levels[0]->ksp, A, P);

999:   PetscObjectTypeCompare((PetscObject)P, MATIS, &ismatis);
1000:   if (!data->is && !ismatis) {
1001:     PetscErrorCode (*create)(DM, IS*, Mat*, PetscErrorCode (**)(Mat, PetscReal, Vec, Vec, PetscReal, IS, void*), void**) = NULL;
1002:     PetscErrorCode (*usetup)(Mat, PetscReal, Vec, Vec, PetscReal, IS, void*) = NULL;
1003:     void           *uctx = NULL;

1005:     /* first see if we can get the data from the DM */
1006:     MatGetDM(P, &dm);
1007:     if (!dm) MatGetDM(A, &dm);
1008:     if (!dm) PCGetDM(pc, &dm);
1009:     if (dm) { /* this is the hook for DMPLEX and DMDA for which the auxiliary Mat is the local Neumann matrix */
1010:       PetscObjectQueryFunction((PetscObject)dm, "DMCreateNeumannOverlap_C", &create);
1011:       if (create) {
1012:         (*create)(dm, &uis, &uaux, &usetup, &uctx);
1013:         data->Neumann = PETSC_TRUE;
1014:       }
1015:     }
1016:     if (!create) {
1017:       if (!uis) {
1018:         PetscObjectQuery((PetscObject)pc, "_PCHPDDM_Neumann_IS", (PetscObject*)&uis);
1019:         PetscObjectReference((PetscObject)uis);
1020:       }
1021:       if (!uaux) {
1022:         PetscObjectQuery((PetscObject)pc, "_PCHPDDM_Neumann_Mat", (PetscObject*)&uaux);
1023:         PetscObjectReference((PetscObject)uaux);
1024:       }
1025:       /* look inside the Pmat instead of the PC, needed for MatSchurComplementComputeExplicitOperator() */
1026:       if (!uis) {
1027:         PetscObjectQuery((PetscObject)P, "_PCHPDDM_Neumann_IS", (PetscObject*)&uis);
1028:         PetscObjectReference((PetscObject)uis);
1029:       }
1030:       if (!uaux) {
1031:         PetscObjectQuery((PetscObject)P, "_PCHPDDM_Neumann_Mat", (PetscObject*)&uaux);
1032:         PetscObjectReference((PetscObject)uaux);
1033:       }
1034:     }
1035:     PCHPDDMSetAuxiliaryMat(pc, uis, uaux, usetup, uctx);
1036:     MatDestroy(&uaux);
1037:     ISDestroy(&uis);
1038:   }

1040:   if (!ismatis) {
1041:     PCHPDDMSetUpNeumannOverlap_Private(pc);
1042:     if (!data->is && data->N > 1) {
1043:       char type[256] = { }; /* same size as in src/ksp/pc/interface/pcset.c */
1044:       PetscOptionsGetString(NULL, pcpre, "-pc_hpddm_levels_1_st_pc_type", type, sizeof(type), NULL);
1045:       PetscStrcmp(type, PCMAT, &algebraic);
1046:       PetscOptionsGetBool(NULL, pcpre, "-pc_hpddm_block_splitting", &block, NULL);
1048:       if (block) algebraic = PETSC_TRUE;
1049:       if (algebraic) {
1050:         ISCreateStride(PETSC_COMM_SELF, P->rmap->n, P->rmap->rstart, 1, &data->is);
1051:         MatIncreaseOverlap(P, 1, &data->is, 1);
1052:         ISSort(data->is);
1053:       } else PetscInfo(pc, "Cannot assemble a fully-algebraic coarse operator with an assembled Pmat and -%spc_hpddm_levels_1_st_pc_type != mat and -%spc_hpddm_block_splitting != true\n", pcpre ? pcpre : "", pcpre ? pcpre : "");
1054:     }
1055:   }

1057:   if (data->is || (ismatis && data->N > 1)) {
1058:     if (ismatis) {
1059:       std::initializer_list<std::string> list = { MATSEQBAIJ, MATSEQSBAIJ };
1060:       MatISGetLocalMat(P, &N);
1061:       std::initializer_list<std::string>::const_iterator it = std::find(list.begin(), list.end(), ((PetscObject)N)->type_name);
1062:       MatISRestoreLocalMat(P, &N);
1063:       switch (std::distance(list.begin(), it)) {
1064:       case 0:
1065:         MatConvert(P, MATMPIBAIJ, MAT_INITIAL_MATRIX, &C);
1066:         break;
1067:       case 1:
1068:         /* MatCreateSubMatrices() does not work with MATSBAIJ and unsorted ISes, so convert to MPIBAIJ */
1069:         MatConvert(P, MATMPIBAIJ, MAT_INITIAL_MATRIX, &C);
1070:         MatSetOption(C, MAT_SYMMETRIC, PETSC_TRUE);
1071:         break;
1072:       default:
1073:         MatConvert(P, MATMPIAIJ, MAT_INITIAL_MATRIX, &C);
1074:       }
1075:       MatISGetLocalToGlobalMapping(P, &l2g, NULL);
1076:       PetscObjectReference((PetscObject)P);
1077:       KSPSetOperators(data->levels[0]->ksp, A, C);
1078:       std::swap(C, P);
1079:       ISLocalToGlobalMappingGetSize(l2g, &n);
1080:       ISCreateStride(PETSC_COMM_SELF, n, 0, 1, &loc);
1081:       ISLocalToGlobalMappingApplyIS(l2g, loc, &is[0]);
1082:       ISDestroy(&loc);
1083:       /* the auxiliary Mat is _not_ the local Neumann matrix                                */
1084:       /* it is the local Neumann matrix augmented (with zeros) through MatIncreaseOverlap() */
1085:       data->Neumann = PETSC_FALSE;
1086:       structure = SAME_NONZERO_PATTERN;
1087:       if (data->share) {
1088:         data->share = PETSC_FALSE;
1089:         PetscInfo(pc, "Cannot share subdomain KSP between SLEPc and PETSc with a Pmat of type MATIS\n");
1090:       }
1091:     } else {
1092:       is[0] = data->is;
1093:       if (algebraic) subdomains = PETSC_TRUE;
1094:       PetscOptionsGetBool(NULL, pcpre, "-pc_hpddm_define_subdomains", &subdomains, NULL);
1095:       if (!subdomains && data->share) {
1096:         data->share = PETSC_FALSE;
1097:         PetscInfo(pc, "Cannot share subdomain KSP between SLEPc and PETSc since -%spc_hpddm_define_subdomains is not true\n", pcpre ? pcpre : "");
1098:       }
1099:       if (data->Neumann) {
1102:       }
1103:       if (data->Neumann || block) structure = SAME_NONZERO_PATTERN;
1104:       ISCreateStride(PetscObjectComm((PetscObject)data->is), P->rmap->n, P->rmap->rstart, 1, &loc);
1105:     }
1106:     PetscSNPrintf(prefix, sizeof(prefix), "%spc_hpddm_levels_1_", pcpre ? pcpre : "");
1107:     PetscOptionsGetEnum(NULL, prefix, "-st_matstructure", MatStructures, (PetscEnum*)&structure, &flg); /* if not user-provided, force its value when possible */
1108:     if (!flg && structure == SAME_NONZERO_PATTERN) { /* cannot call STSetMatStructure() yet, insert the appropriate option in the database, parsed by STSetFromOptions() */
1109:       PetscSNPrintf(prefix, sizeof(prefix), "-%spc_hpddm_levels_1_st_matstructure", pcpre ? pcpre : "");
1110:       PetscOptionsSetValue(NULL, prefix, MatStructures[structure]);
1111:     }
1112:     if (data->N > 1 && (data->aux || ismatis || algebraic)) {
1114:       MatSetOption(P, MAT_SUBMAT_SINGLEIS, PETSC_TRUE);
1115:       if (ismatis) {
1116:         /* needed by HPDDM (currently) so that the partition of unity is 0 on subdomain interfaces */
1117:         MatIncreaseOverlap(P, 1, is, 1);
1118:         ISDestroy(&data->is);
1119:         data->is = is[0];
1120:       } else {
1121:         if (PetscDefined(USE_DEBUG)) {
1122:           PetscBool equal;
1123:           IS        intersect;

1125:           ISIntersect(data->is, loc, &intersect);
1126:           ISEqualUnsorted(loc, intersect, &equal);
1127:           ISDestroy(&intersect);
1129:         }
1130:         PetscObjectComposeFunction((PetscObject)pc->pmat, "PCHPDDMAlgebraicAuxiliaryMat_Private_C", PCHPDDMAlgebraicAuxiliaryMat_Private);
1131:         if (!data->Neumann && !algebraic) {
1132:           PetscObjectTypeCompare((PetscObject)P, MATMPISBAIJ, &flg);
1133:           if (flg) {
1134:             /* maybe better to ISSort(is[0]), MatCreateSubMatrices(), and then MatPermute() */
1135:             /* but there is no MatPermute_SeqSBAIJ(), so as before, just use MATMPIBAIJ     */
1136:             MatConvert(P, MATMPIBAIJ, MAT_INITIAL_MATRIX, &uaux);
1137:             flg = PETSC_FALSE;
1138:           }
1139:         }
1140:       }
1141:       if (algebraic) {
1142:         PetscUseMethod(pc->pmat, "PCHPDDMAlgebraicAuxiliaryMat_Private_C", (Mat, IS*, Mat*[], PetscBool), (P, is, &sub, block));
1143:         if (block) {
1144:           PetscObjectQuery((PetscObject)sub[0], "_PCHPDDM_Neumann_Mat", (PetscObject*)&data->aux);
1145:           PetscObjectCompose((PetscObject)sub[0], "_PCHPDDM_Neumann_Mat", NULL);
1146:         }
1147:       } else if (!uaux) {
1148:         if (data->Neumann) sub = &data->aux;
1149:         else MatCreateSubMatrices(P, 1, is, is, MAT_INITIAL_MATRIX, &sub);
1150:       } else {
1151:         MatCreateSubMatrices(uaux, 1, is, is, MAT_INITIAL_MATRIX, &sub);
1152:         MatDestroy(&uaux);
1153:         MatConvert(sub[0], MATSEQSBAIJ, MAT_INPLACE_MATRIX, sub);
1154:       }
1155:       /* Vec holding the partition of unity */
1156:       if (!data->levels[0]->D) {
1157:         ISGetLocalSize(data->is, &n);
1158:         VecCreateMPI(PETSC_COMM_SELF, n, PETSC_DETERMINE, &data->levels[0]->D);
1159:       }
1160:       if (data->share && structure == SAME_NONZERO_PATTERN) { /* share the KSP only when the MatStructure is SAME_NONZERO_PATTERN */
1161:         Mat      D;
1162:         IS       perm = NULL;
1163:         PetscInt size = -1;
1164:         PCHPDDMPermute_Private(*is, data->is, &uis, data->Neumann || block ? sub[0] : data->aux, &C, &perm);
1165:         if (!data->Neumann && !block) {
1166:           MatPermute(sub[0], perm, perm, &D); /* permute since PCASM will call ISSort() */
1167:           MatHeaderReplace(sub[0], &D);
1168:         }
1169:         if (data->B) { /* see PCHPDDMSetRHSMat() */
1170:           MatPermute(data->B, perm, perm, &D);
1171:           MatHeaderReplace(data->B, &D);
1172:         }
1173:         ISDestroy(&perm);
1174:         if (!data->levels[0]->pc) {
1175:           PetscSNPrintf(prefix, sizeof(prefix), "%spc_hpddm_levels_1_", pcpre ? pcpre : "");
1176:           PCCreate(PetscObjectComm((PetscObject)pc), &data->levels[0]->pc);
1177:           PCSetOptionsPrefix(data->levels[0]->pc, prefix);
1178:           PCSetOperators(data->levels[0]->pc, A, P);
1179:         }
1180:         PCSetType(data->levels[0]->pc, PCASM);
1181:         if (!data->levels[0]->pc->setupcalled) PCASMSetLocalSubdomains(data->levels[0]->pc, 1, is, &loc);
1182:         PCSetFromOptions(data->levels[0]->pc);
1183:         if (block) PCHPDDMCommunicationAvoidingPCASM_Private(data->levels[0]->pc, C, algebraic);
1184:         else PCSetUp(data->levels[0]->pc);
1185:         PetscTryMethod(data->levels[0]->pc, "PCASMGetSubKSP_C", (PC, PetscInt*, PetscInt*, KSP**), (data->levels[0]->pc, &size, NULL, &ksp));
1186:         if (size != 1) {
1187:           PCDestroy(&data->levels[0]->pc);
1188:           MatDestroy(&C);
1189:           ISDestroy(&uis);
1190:           data->share = PETSC_FALSE;
1191:           if (size == -1) PetscInfo(pc, "Cannot share subdomain KSP between SLEPc and PETSc since PCASMGetSubKSP() not found in fine-level PC\n");
1192:           else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of subdomain solver %" PetscInt_FMT " != 1", size);
1193:         } else {
1194:           const char *matpre;
1195:           PetscBool  cmp[2];
1196:           KSPGetOperators(ksp[0], subA, subA + 1);
1197:           MatDuplicate(subA[1], MAT_SHARE_NONZERO_PATTERN, &D);
1198:           MatGetOptionsPrefix(subA[1], &matpre);
1199:           MatSetOptionsPrefix(D, matpre);
1200:           PetscObjectTypeCompare((PetscObject)D, MATNORMAL, cmp);
1201:           PetscObjectTypeCompare((PetscObject)C, MATNORMAL, cmp + 1);
1203:           if (!cmp[0]) {
1204:             if (!block) MatAXPY(D, 1.0, C, SUBSET_NONZERO_PATTERN);
1205:             else MatAXPY(D, 1.0, data->aux, SAME_NONZERO_PATTERN);
1206:           } else {
1207:             Mat mat[2];
1208:             MatNormalGetMat(D, mat);
1209:             MatNormalGetMat(C, mat + 1);
1210:             MatAXPY(mat[0], 1.0, mat[1], SUBSET_NONZERO_PATTERN);
1211:           }
1212:           MatPropagateSymmetryOptions(C, D);
1213:           MatDestroy(&C);
1214:           C = D;
1215:           /* swap pointers so that variables stay consistent throughout PCSetUp() */
1216:           std::swap(C, data->aux);
1217:           std::swap(uis, data->is);
1218:           swap = PETSC_TRUE;
1219:         }
1220:       } else if (data->share) {
1221:         data->share = PETSC_FALSE;
1222:         PetscInfo(pc, "Cannot share subdomain KSP between SLEPc and PETSc since -%spc_hpddm_levels_1_st_matstructure %s (!= %s)\n", pcpre ? pcpre : "", MatStructures[structure], MatStructures[SAME_NONZERO_PATTERN]);
1223:       }
1224:       if (!data->levels[0]->scatter) {
1225:         MatCreateVecs(P, &xin, NULL);
1226:         if (ismatis) MatDestroy(&P);
1227:         VecScatterCreate(xin, data->is, data->levels[0]->D, NULL, &data->levels[0]->scatter);
1228:         VecDestroy(&xin);
1229:       }
1230:       if (data->levels[0]->P) {
1231:         /* if the pattern is the same and PCSetUp() has previously succeeded, reuse HPDDM buffers and connectivity */
1232:         HPDDM::Schwarz<PetscScalar>::destroy(data->levels[0], pc->setupcalled < 1 || pc->flag == DIFFERENT_NONZERO_PATTERN ? PETSC_TRUE : PETSC_FALSE);
1233:       }
1234:       if (!data->levels[0]->P) data->levels[0]->P = new HPDDM::Schwarz<PetscScalar>();
1235:       if (data->log_separate) PetscLogEventBegin(PC_HPDDM_SetUp[0], data->levels[0]->ksp, 0, 0, 0);
1236:       else PetscLogEventBegin(PC_HPDDM_Strc, data->levels[0]->ksp, 0, 0, 0);
1237:       /* HPDDM internal data structure */
1238:       data->levels[0]->P->structure(loc, data->is, sub[0], ismatis ? C : data->aux, data->levels);
1239:       if (!data->log_separate) PetscLogEventEnd(PC_HPDDM_Strc, data->levels[0]->ksp, 0, 0, 0);
1240:       /* matrix pencil of the generalized eigenvalue problem on the overlap (GenEO) */
1241:       if (!data->B) {
1242:         MatDuplicate(sub[0], MAT_COPY_VALUES, &weighted);
1243:         PetscObjectTypeCompare((PetscObject)weighted, MATNORMAL, &flg);
1244:         if (!flg) MatDiagonalScale(weighted, data->levels[0]->D, data->levels[0]->D);
1245:         else { /* MATNORMAL applies MatDiagonalScale() in a matrix-free fashion, not what is needed since this won't be passed to SLEPc during the eigensolve */
1246:           MatNormalGetMat(weighted, &data->B);
1247:           MatDiagonalScale(data->B, NULL, data->levels[0]->D);
1248:           data->B = NULL;
1249:           flg = PETSC_FALSE;
1250:         }
1251:         /* neither MatDuplicate() nor MatDiagonaleScale() handles the symmetry options, so propagate the options explicitly */
1252:         /* only useful for -mat_type baij -pc_hpddm_levels_1_st_pc_type cholesky (no problem with MATAIJ or MATSBAIJ)       */
1253:         MatPropagateSymmetryOptions(sub[0], weighted);
1254:       } else weighted = data->B;
1255:       /* SLEPc is used inside the loaded symbol */
1256:       (*loadedSym)(data->levels[0]->P, data->is, ismatis ? C : (algebraic && !block ? sub[0] : data->aux), weighted, data->B, initial, data->levels);
1257:       if (data->share) {
1258:         Mat st[2];
1259:         KSPGetOperators(ksp[0], st, st + 1);
1260:         MatCopy(subA[0], st[0], SAME_NONZERO_PATTERN);
1261:         if (subA[1] != subA[0] || st[1] != st[0]) MatCopy(subA[1], st[1], SAME_NONZERO_PATTERN);
1262:       }
1263:       if (data->log_separate) PetscLogEventEnd(PC_HPDDM_SetUp[0], data->levels[0]->ksp, 0, 0, 0);
1264:       if (ismatis) MatISGetLocalMat(C, &N);
1265:       else N = data->aux;
1266:       P = sub[0];
1267:       /* going through the grid hierarchy */
1268:       for (n = 1; n < data->N; ++n) {
1269:         if (data->log_separate) PetscLogEventBegin(PC_HPDDM_SetUp[n], data->levels[n]->ksp, 0, 0, 0);
1270:         /* method composed in the loaded symbol since there, SLEPc is used as well */
1271:         PetscTryMethod(data->levels[0]->ksp, "PCHPDDMSetUp_Private_C", (Mat*, Mat*, PetscInt, PetscInt* const, PC_HPDDM_Level** const), (&P, &N, n, &data->N, data->levels));
1272:         if (data->log_separate) PetscLogEventEnd(PC_HPDDM_SetUp[n], data->levels[n]->ksp, 0, 0, 0);
1273:       }
1274:       /* reset to NULL to avoid any faulty use */
1275:       PetscObjectComposeFunction((PetscObject)data->levels[0]->ksp, "PCHPDDMSetUp_Private_C", NULL);
1276:       if (!ismatis) PetscObjectComposeFunction((PetscObject)pc->pmat, "PCHPDDMAlgebraicAuxiliaryMat_C", NULL);
1277:       else PetscObjectDereference((PetscObject)C); /* matching PetscObjectReference() above */
1278:       for (n = 0; n < data->N - 1; ++n)
1279:         if (data->levels[n]->P) {
1280:           /* HPDDM internal work buffers */
1281:           data->levels[n]->P->setBuffer();
1282:           data->levels[n]->P->super::start();
1283:         }
1284:       if (ismatis || !subdomains) PCHPDDMDestroySubMatrices_Private(data->Neumann, PetscBool(algebraic && !block), sub);
1285:       if (ismatis) data->is = NULL;
1286:       for (n = 0; n < data->N - 1 + (reused > 0); ++n) {
1287:         if (data->levels[n]->P) {
1288:           PC spc;

1290:           /* force the PC to be PCSHELL to do the coarse grid corrections */
1291:           KSPSetSkipPCSetFromOptions(data->levels[n]->ksp, PETSC_TRUE);
1292:           KSPGetPC(data->levels[n]->ksp, &spc);
1293:           PCSetType(spc, PCSHELL);
1294:           PCShellSetContext(spc, data->levels[n]);
1295:           PCShellSetSetUp(spc, PCHPDDMShellSetUp);
1296:           PCShellSetApply(spc, PCHPDDMShellApply);
1297:           PCShellSetMatApply(spc, PCHPDDMShellMatApply);
1298:           PCShellSetDestroy(spc, PCHPDDMShellDestroy);
1299:           if (!data->levels[n]->pc) PCCreate(PetscObjectComm((PetscObject)data->levels[n]->ksp), &data->levels[n]->pc);
1300:           if (n < reused) {
1301:             PCSetReusePreconditioner(spc, PETSC_TRUE);
1302:             PCSetReusePreconditioner(data->levels[n]->pc, PETSC_TRUE);
1303:           }
1304:           PCSetUp(spc);
1305:         }
1306:       }
1307:     } else flg = reused ? PETSC_FALSE : PETSC_TRUE;
1308:     if (!ismatis && subdomains) {
1309:       if (flg) KSPGetPC(data->levels[0]->ksp, &inner);
1310:       else inner = data->levels[0]->pc;
1311:       if (inner) {
1312:         PCSetType(inner, PCASM); /* inner is the fine-level PC for which one must ensure                       */
1313:                                                       /* PCASMSetLocalSubdomains() has been called when -pc_hpddm_define_subdomains */
1314:         if (!inner->setupcalled) { /* evaluates to PETSC_FALSE when -pc_hpddm_block_splitting */
1315:           PCASMSetLocalSubdomains(inner, 1, is, &loc);
1316:           if (!data->Neumann && data->N > 1) { /* subdomain matrices are already created for the eigenproblem, reuse them for the fine-level PC */
1317:             PCHPDDMPermute_Private(*is, NULL, NULL, sub[0], &C, NULL);
1318:             PCHPDDMCommunicationAvoidingPCASM_Private(inner, C, algebraic);
1319:             MatDestroy(&C);
1320:           }
1321:         }
1322:       }
1323:       if (data->N > 1) PCHPDDMDestroySubMatrices_Private(data->Neumann, PetscBool(algebraic && !block), sub);
1324:     }
1325:     ISDestroy(&loc);
1326:   } else data->N = 1 + reused; /* enforce this value to 1 + reused if there is no way to build another level */
1327:   if (requested != data->N + reused) {
1328:     PetscInfo(pc, "%" PetscInt_FMT " levels requested, only %" PetscInt_FMT " built + %" PetscInt_FMT " reused. Options for level(s) > %" PetscInt_FMT ", including -%spc_hpddm_coarse_ will not be taken into account\n", requested, data->N, reused, data->N, pcpre ? pcpre : "");
1329:     PetscInfo(pc, "It is best to tune parameters, e.g., a higher value for -%spc_hpddm_levels_%" PetscInt_FMT "_eps_threshold so that at least one local deflation vector will be selected\n", pcpre ? pcpre : "", data->N);
1330:     /* cannot use PCHPDDMShellDestroy() because PCSHELL not set for unassembled levels */
1331:     for (n = data->N - 1; n < requested - 1; ++n) {
1332:       if (data->levels[n]->P) {
1333:         HPDDM::Schwarz<PetscScalar>::destroy(data->levels[n], PETSC_TRUE);
1334:         VecDestroyVecs(1, &data->levels[n]->v[0]);
1335:         VecDestroyVecs(2, &data->levels[n]->v[1]);
1336:         MatDestroy(data->levels[n]->V);
1337:         MatDestroy(data->levels[n]->V + 1);
1338:         MatDestroy(data->levels[n]->V + 2);
1339:         VecDestroy(&data->levels[n]->D);
1340:         VecScatterDestroy(&data->levels[n]->scatter);
1341:       }
1342:     }
1343:     if (reused) {
1344:       for (n = reused; n < PETSC_PCHPDDM_MAXLEVELS && data->levels[n]; ++n) {
1345:         KSPDestroy(&data->levels[n]->ksp);
1346:         PCDestroy(&data->levels[n]->pc);
1347:       }
1348:     }
1350:   }
1351:   /* these solvers are created after PCSetFromOptions() is called */
1352:   if (pc->setfromoptionscalled) {
1353:     for (n = 0; n < data->N; ++n) {
1354:       if (data->levels[n]->ksp) KSPSetFromOptions(data->levels[n]->ksp);
1355:       if (data->levels[n]->pc) PCSetFromOptions(data->levels[n]->pc);
1356:     }
1357:     pc->setfromoptionscalled = 0;
1358:   }
1359:   data->N += reused;
1360:   if (data->share && swap) {
1361:     /* swap back pointers so that variables follow the user-provided numbering */
1362:     std::swap(C, data->aux);
1363:     std::swap(uis, data->is);
1364:     MatDestroy(&C);
1365:     ISDestroy(&uis);
1366:   }
1367:   return 0;
1368: }

1370: /*@
1371:      PCHPDDMSetCoarseCorrectionType - Sets the coarse correction type.

1373:    Input Parameters:
1374: +     pc - preconditioner context
1375: -     type - PC_HPDDM_COARSE_CORRECTION_DEFLATED, PC_HPDDM_COARSE_CORRECTION_ADDITIVE, or PC_HPDDM_COARSE_CORRECTION_BALANCED

1377:    Options Database Key:
1378: .   -pc_hpddm_coarse_correction <deflated, additive, balanced> - type of coarse correction to apply

1380:    Level: intermediate

1382: .seealso:  PCHPDDMGetCoarseCorrectionType(), PCHPDDM, PCHPDDMCoarseCorrectionType
1383: @*/
1384: PetscErrorCode PCHPDDMSetCoarseCorrectionType(PC pc, PCHPDDMCoarseCorrectionType type)
1385: {
1388:   PetscTryMethod(pc, "PCHPDDMSetCoarseCorrectionType_C", (PC, PCHPDDMCoarseCorrectionType), (pc, type));
1389:   return 0;
1390: }

1392: /*@
1393:      PCHPDDMGetCoarseCorrectionType - Gets the coarse correction type.

1395:    Input Parameter:
1396: .     pc - preconditioner context

1398:    Output Parameter:
1399: .     type - PC_HPDDM_COARSE_CORRECTION_DEFLATED, PC_HPDDM_COARSE_CORRECTION_ADDITIVE, or PC_HPDDM_COARSE_CORRECTION_BALANCED

1401:    Level: intermediate

1403: .seealso:  PCHPDDMSetCoarseCorrectionType(), PCHPDDM, PCHPDDMCoarseCorrectionType
1404: @*/
1405: PetscErrorCode PCHPDDMGetCoarseCorrectionType(PC pc, PCHPDDMCoarseCorrectionType *type)
1406: {
1408:   if (type) {
1410:     PetscUseMethod(pc, "PCHPDDMGetCoarseCorrectionType_C", (PC, PCHPDDMCoarseCorrectionType*), (pc, type));
1411:   }
1412:   return 0;
1413: }

1415: static PetscErrorCode PCHPDDMSetCoarseCorrectionType_HPDDM(PC pc, PCHPDDMCoarseCorrectionType type)
1416: {
1417:   PC_HPDDM *data = (PC_HPDDM*)pc->data;

1420:   data->correction = type;
1421:   return 0;
1422: }

1424: static PetscErrorCode PCHPDDMGetCoarseCorrectionType_HPDDM(PC pc, PCHPDDMCoarseCorrectionType *type)
1425: {
1426:   PC_HPDDM *data = (PC_HPDDM*)pc->data;

1428:   *type = data->correction;
1429:   return 0;
1430: }

1432: /*@
1433:      PCHPDDMGetSTShareSubKSP - Gets whether the KSP in SLEPc ST and the fine-level subdomain solver is shared.

1435:    Input Parameter:
1436: .     pc - preconditioner context

1438:    Output Parameter:
1439: .     share - whether the KSP is shared or not

1441:    Notes:
1442:      This is not the same as PCGetReusePreconditioner(). The return value is unlikely to be true, but when it is, a symbolic factorization can be skipped
1443:      when using a subdomain PCType such as PCLU or PCCHOLESKY.

1445:    Level: advanced

1447: @*/
1448: PetscErrorCode PCHPDDMGetSTShareSubKSP(PC pc, PetscBool *share)
1449: {
1451:   if (share) {
1453:     PetscUseMethod(pc, "PCHPDDMGetSTShareSubKSP_C", (PC, PetscBool*), (pc, share));
1454:   }
1455:   return 0;
1456: }

1458: static PetscErrorCode PCHPDDMGetSTShareSubKSP_HPDDM(PC pc, PetscBool *share)
1459: {
1460:   PC_HPDDM *data = (PC_HPDDM*)pc->data;

1462:   *share = data->share;
1463:   return 0;
1464: }

1466: PetscErrorCode HPDDMLoadDL_Private(PetscBool *found)
1467: {
1468:   char lib[PETSC_MAX_PATH_LEN], dlib[PETSC_MAX_PATH_LEN], dir[PETSC_MAX_PATH_LEN];

1470:   PetscStrcpy(dir, "${PETSC_LIB_DIR}");
1471:   PetscOptionsGetString(NULL, NULL, "-hpddm_dir", dir, sizeof(dir), NULL);
1472:   PetscSNPrintf(lib, sizeof(lib), "%s/libhpddm_petsc", dir);
1473:   PetscDLLibraryRetrieve(PETSC_COMM_SELF, lib, dlib, 1024, found);
1474:   if (*found) PetscDLLibraryAppend(PETSC_COMM_SELF, &PetscDLLibrariesLoaded, dlib);
1475: #if defined(SLEPC_LIB_DIR) /* this variable is passed during SLEPc ./configure since    */
1476:   else {                   /* slepcconf.h is not yet built (and thus can't be included) */
1477:     PetscStrcpy(dir, HPDDM_STR(SLEPC_LIB_DIR));
1478:     PetscSNPrintf(lib, sizeof(lib), "%s/libhpddm_petsc", dir);
1479:     PetscDLLibraryRetrieve(PETSC_COMM_SELF, lib, dlib, 1024, found);
1480:     if (*found) PetscDLLibraryAppend(PETSC_COMM_SELF, &PetscDLLibrariesLoaded, dlib);
1481: #endif
1482:     else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "%s not found", lib);
1483: #if defined(SLEPC_LIB_DIR)
1484:   }
1485: #endif
1486:   return 0;
1487: }

1489: /*MC
1490:      PCHPDDM - Interface with the HPDDM library.

1492:    This PC may be used to build multilevel spectral domain decomposition methods based on the GenEO framework [2011, 2019]. It may be viewed as an alternative to spectral AMGe or PCBDDC with adaptive selection of constraints. A chronological bibliography of relevant publications linked with PC available in HPDDM through PCHPDDM may be found below. The interface is explained in details in [2021].

1494:    The matrix to be preconditioned (Pmat) may be unassembled (MATIS), assembled (MATMPIAIJ, MATMPIBAIJ, or MATMPISBAIJ), hierarchical (MATHTOOL), or MATNORMAL. For multilevel preconditioning, when using an assembled or hierarchical Pmat, one must provide an auxiliary local Mat (unassembled local operator for GenEO) using PCHPDDMSetAuxiliaryMat(). Calling this routine is not needed when using a MATIS Pmat, assembly done internally using MatConvert().

1496:    Options Database Keys:
1497: +   -pc_hpddm_define_subdomains <true, default=false> - on the finest level, calls PCASMSetLocalSubdomains() with the IS supplied in PCHPDDMSetAuxiliaryMat() (not relevant with an unassembled Pmat)
1498: .   -pc_hpddm_has_neumann <true, default=false> - on the finest level, informs the PC that the local Neumann matrix is supplied in PCHPDDMSetAuxiliaryMat()
1499: -   -pc_hpddm_coarse_correction <type, default=deflated> - determines the PCHPDDMCoarseCorrectionType when calling PCApply

1501:    Options for subdomain solvers, subdomain eigensolvers (for computing deflation vectors), and the coarse solver can be set with
1502: .vb
1503:       -pc_hpddm_levels_%d_pc_
1504:       -pc_hpddm_levels_%d_ksp_
1505:       -pc_hpddm_levels_%d_eps_
1506:       -pc_hpddm_levels_%d_p
1507:       -pc_hpddm_levels_%d_mat_type_
1508:       -pc_hpddm_coarse_
1509:       -pc_hpddm_coarse_p
1510:       -pc_hpddm_coarse_mat_type_
1511: .ve
1512:    e.g., -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_levels_1_eps_nev 10 -pc_hpddm_levels_2_p 4 -pc_hpddm_levels_2_sub_pc_type lu -pc_hpddm_levels_2_eps_nev 10 -pc_hpddm_coarse_p 2 -pc_hpddm_coarse_mat_type baij will use 10 deflation vectors per subdomain on the fine "level 1", aggregate the fine subdomains into 4 "level 2" subdomains, then use 10 deflation vectors per subdomain on "level 2", and assemble the coarse matrix (of dimension 4 x 10 = 40) on two processes as a MATMPIBAIJ (default is MATMPISBAIJ).

1514:    In order to activate a "level N+1" coarse correction, it is mandatory to call -pc_hpddm_levels_N_eps_nev <nu> or -pc_hpddm_levels_N_eps_threshold <val>. The default -pc_hpddm_coarse_p value is 1, meaning that the coarse operator is aggregated on a single process.

1516:    This preconditioner requires that you build PETSc with SLEPc (--download-slepc=1). By default, the underlying concurrent eigenproblems are solved using SLEPc shift-and-invert spectral transformation. This is usually what gives the best performance for GenEO, cf. [2011, 2013]. As stated above, SLEPc options are available through -pc_hpddm_levels_%d_, e.g., -pc_hpddm_levels_1_eps_type arpack -pc_hpddm_levels_1_eps_threshold 0.1 -pc_hpddm_levels_1_st_type sinvert.

1518:    References:
1519: +   2011 - A robust two-level domain decomposition preconditioner for systems of PDEs. Spillane, Dolean, Hauret, Nataf, Pechstein, and Scheichl. Comptes Rendus Mathematique.
1520: .   2013 - Scalable domain decomposition preconditioners for heterogeneous elliptic problems. Jolivet, Hecht, Nataf, and Prud'homme. SC13.
1521: .   2015 - An introduction to domain decomposition methods: algorithms, theory, and parallel implementation. Dolean, Jolivet, and Nataf. SIAM.
1522: .   2019 - A multilevel Schwarz preconditioner based on a hierarchy of robust coarse spaces. Al Daas, Grigori, Jolivet, and Tournier. SIAM Journal on Scientific Computing.
1523: .   2021 - KSPHPDDM and PCHPDDM: extending PETSc with advanced Krylov methods and robust multilevel overlapping Schwarz preconditioners. Jolivet, Roman, and Zampini. Computer & Mathematics with Applications.
1524: .   2022a - A robust algebraic domain decomposition preconditioner for sparse normal equations. Al Daas, Jolivet, and Scott. SIAM Journal on Scientific Computing.
1525: -   2022b - A robust algebraic multilevel domain decomposition preconditioner for sparse symmetric positive definite matrices. Al Daas and Jolivet.

1527:    Level: intermediate

1529: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCHPDDMSetAuxiliaryMat(), MATIS, PCBDDC, PCDEFLATION, PCTELESCOPE
1530: M*/
1531: PETSC_EXTERN PetscErrorCode PCCreate_HPDDM(PC pc)
1532: {
1533:   PC_HPDDM  *data;
1534:   PetscBool found;

1536:   if (!loadedSym) {
1537:     HPDDMLoadDL_Private(&found);
1538:     if (found) PetscDLLibrarySym(PETSC_COMM_SELF, &PetscDLLibrariesLoaded, NULL, "PCHPDDM_Internal", (void**)&loadedSym);
1539:   }
1541:   PetscNewLog(pc, &data);
1542:   pc->data                = data;
1543:   pc->ops->reset          = PCReset_HPDDM;
1544:   pc->ops->destroy        = PCDestroy_HPDDM;
1545:   pc->ops->setfromoptions = PCSetFromOptions_HPDDM;
1546:   pc->ops->setup          = PCSetUp_HPDDM;
1547:   pc->ops->apply          = PCApply_HPDDM;
1548:   pc->ops->matapply       = PCMatApply_HPDDM;
1549:   pc->ops->view           = PCView_HPDDM;
1550:   pc->ops->presolve       = PCPreSolve_HPDDM;
1551:   PetscObjectComposeFunction((PetscObject)pc, "PCHPDDMSetAuxiliaryMat_C", PCHPDDMSetAuxiliaryMat_HPDDM);
1552:   PetscObjectComposeFunction((PetscObject)pc, "PCHPDDMHasNeumannMat_C", PCHPDDMHasNeumannMat_HPDDM);
1553:   PetscObjectComposeFunction((PetscObject)pc, "PCHPDDMSetRHSMat_C", PCHPDDMSetRHSMat_HPDDM);
1554:   PetscObjectComposeFunction((PetscObject)pc, "PCHPDDMSetCoarseCorrectionType_C", PCHPDDMSetCoarseCorrectionType_HPDDM);
1555:   PetscObjectComposeFunction((PetscObject)pc, "PCHPDDMGetCoarseCorrectionType_C", PCHPDDMGetCoarseCorrectionType_HPDDM);
1556:   PetscObjectComposeFunction((PetscObject)pc, "PCHPDDMGetSTShareSubKSP_C", PCHPDDMGetSTShareSubKSP_HPDDM);
1557:   return 0;
1558: }

1560: /*@C
1561:      PCHPDDMInitializePackage - This function initializes everything in the PCHPDDM package. It is called from PCInitializePackage().

1563:    Level: intermediate

1565: .seealso:  PetscInitialize()
1566: @*/
1567: PetscErrorCode PCHPDDMInitializePackage(void)
1568: {
1569:   char     ename[32];
1570:   PetscInt i;

1572:   if (PCHPDDMPackageInitialized) return 0;
1573:   PCHPDDMPackageInitialized = PETSC_TRUE;
1574:   PetscRegisterFinalize(PCHPDDMFinalizePackage);
1575:   /* general events registered once during package initialization */
1576:   /* some of these events are not triggered in libpetsc,          */
1577:   /* but rather directly in libhpddm_petsc,                       */
1578:   /* which is in charge of performing the following operations    */

1580:   /* domain decomposition structure from Pmat sparsity pattern    */
1581:   PetscLogEventRegister("PCHPDDMStrc", PC_CLASSID, &PC_HPDDM_Strc);
1582:   /* Galerkin product, redistribution, and setup (not triggered in libpetsc)                */
1583:   PetscLogEventRegister("PCHPDDMPtAP", PC_CLASSID, &PC_HPDDM_PtAP);
1584:   /* Galerkin product with summation, redistribution, and setup (not triggered in libpetsc) */
1585:   PetscLogEventRegister("PCHPDDMPtBP", PC_CLASSID, &PC_HPDDM_PtBP);
1586:   /* next level construction using PtAP and PtBP (not triggered in libpetsc)                */
1587:   PetscLogEventRegister("PCHPDDMNext", PC_CLASSID, &PC_HPDDM_Next);
1588:   static_assert(PETSC_PCHPDDM_MAXLEVELS <= 9, "PETSC_PCHPDDM_MAXLEVELS value is too high");
1589:   for (i = 1; i < PETSC_PCHPDDM_MAXLEVELS; ++i) {
1590:     PetscSNPrintf(ename, sizeof(ename), "PCHPDDMSetUp L%1d", i);
1591:     /* events during a PCSetUp() at level #i _except_ the assembly */
1592:     /* of the Galerkin operator of the coarser level #(i + 1)      */
1593:     PetscLogEventRegister(ename, PC_CLASSID, &PC_HPDDM_SetUp[i - 1]);
1594:     PetscSNPrintf(ename, sizeof(ename), "PCHPDDMSolve L%1d", i);
1595:     /* events during a PCApply() at level #i _except_              */
1596:     /* the KSPSolve() of the coarser level #(i + 1)                */
1597:     PetscLogEventRegister(ename, PC_CLASSID, &PC_HPDDM_Solve[i - 1]);
1598:   }
1599:   return 0;
1600: }

1602: /*@C
1603:      PCHPDDMFinalizePackage - This function frees everything from the PCHPDDM package. It is called from PetscFinalize().

1605:    Level: intermediate

1607: .seealso:  PetscFinalize()
1608: @*/
1609: PetscErrorCode PCHPDDMFinalizePackage(void)
1610: {
1611:   PCHPDDMPackageInitialized = PETSC_FALSE;
1612:   return 0;
1613: }