Actual source code: hpddm.cxx
petsc-3.14.6 2021-03-30
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 const, IS, const Mat, Mat, std::vector<Vec>, PetscInt* const, PC_HPDDM_Level** const) = NULL;
12: static PetscBool PCHPDDMPackageInitialized = PETSC_FALSE;
13: static PetscBool citePC = PETSC_FALSE;
14: static const char hpddmCitationPC[] = "@inproceedings{jolivet2013scalable,\n\tTitle = {{Scalable Domain Decomposition Preconditioners For Heterogeneous Elliptic Problems}},\n\tAuthor = {Jolivet, Pierre and Hecht, Fr\'ed\'eric and Nataf, Fr\'ed\'eric and Prud'homme, Christophe},\n\tOrganization = {ACM},\n\tYear = {2013},\n\tSeries = {SC13},\n\tBooktitle = {Proceedings of the 2013 International Conference for High Performance Computing, Networking, Storage and Analysis}\n}\n";
16: PetscLogEvent PC_HPDDM_Strc;
17: PetscLogEvent PC_HPDDM_PtAP;
18: PetscLogEvent PC_HPDDM_PtBP;
19: PetscLogEvent PC_HPDDM_Next;
21: static const char *PCHPDDMCoarseCorrectionTypes[] = { "deflated", "additive", "balanced" };
23: static PetscErrorCode PCReset_HPDDM(PC pc)
24: {
25: PC_HPDDM *data = (PC_HPDDM*)pc->data;
26: PetscInt i;
30: if (data->levels) {
31: for (i = 0; i < PETSC_HPDDM_MAXLEVELS && data->levels[i]; ++i) {
32: KSPDestroy(&data->levels[i]->ksp);
33: PCDestroy(&data->levels[i]->pc);
34: PetscFree(data->levels[i]);
35: }
36: PetscFree(data->levels);
37: }
39: ISDestroy(&data->is);
40: MatDestroy(&data->aux);
41: MatDestroy(&data->B);
42: VecDestroy(&data->normal);
43: data->correction = PC_HPDDM_COARSE_CORRECTION_DEFLATED;
44: data->Neumann = PETSC_FALSE;
45: data->setup = NULL;
46: data->setup_ctx = NULL;
47: return(0);
48: }
50: static PetscErrorCode PCDestroy_HPDDM(PC pc)
51: {
52: PC_HPDDM *data = (PC_HPDDM*)pc->data;
56: PCReset_HPDDM(pc);
57: PetscFree(data);
58: PetscObjectChangeTypeName((PetscObject)pc, 0);
59: PetscObjectComposeFunction((PetscObject)pc, "PCHPDDMSetAuxiliaryMat_C", NULL);
60: PetscObjectComposeFunction((PetscObject)pc, "PCHPDDMHasNeumannMat_C", NULL);
61: PetscObjectComposeFunction((PetscObject)pc, "PCHPDDMSetRHSMat_C", NULL);
62: PetscObjectComposeFunction((PetscObject)pc, "PCHPDDMSetCoarseCorrectionType_C", NULL);
63: PetscObjectComposeFunction((PetscObject)pc, "PCHPDDMGetCoarseCorrectionType_C", NULL);
64: return(0);
65: }
67: static PetscErrorCode PCHPDDMSetAuxiliaryMat_HPDDM(PC pc, IS is, Mat A, PetscErrorCode (*setup)(Mat, PetscReal, Vec, Vec, PetscReal, IS, void*), void* setup_ctx)
68: {
69: PC_HPDDM *data = (PC_HPDDM*)pc->data;
73: if (is) {
74: PetscObjectReference((PetscObject)is);
75: if (data->is) { /* new overlap definition resets the PC */
76: PCReset_HPDDM(pc);
77: pc->setfromoptionscalled = 0;
78: }
79: ISDestroy(&data->is);
80: data->is = is;
81: }
82: if (A) {
83: PetscObjectReference((PetscObject)A);
84: MatDestroy(&data->aux);
85: data->aux = A;
86: }
87: if (setup) {
88: data->setup = setup;
89: data->setup_ctx = setup_ctx;
90: }
91: return(0);
92: }
94: /*@
95: PCHPDDMSetAuxiliaryMat - Sets the auxiliary matrix used by PCHPDDM for the concurrent GenEO eigenproblems 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.
97: Input Parameters:
98: + pc - preconditioner context
99: . is - index set of the local auxiliary, e.g., Neumann, matrix
100: . A - auxiliary sequential matrix
101: . setup - function for generating the auxiliary matrix
102: - setup_ctx - context for setup
104: Level: intermediate
106: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCHPDDMSetRHSMat(), MATIS
107: @*/
108: PetscErrorCode PCHPDDMSetAuxiliaryMat(PC pc, IS is, Mat A, PetscErrorCode (*setup)(Mat, PetscReal, Vec, Vec, PetscReal, IS, void*), void* setup_ctx)
109: {
116: #if defined(PETSC_HAVE_FORTRAN)
117: if (reinterpret_cast<void*>(setup) == reinterpret_cast<void*>(PETSC_NULL_FUNCTION_Fortran)) setup = NULL;
118: if (setup_ctx == PETSC_NULL_INTEGER_Fortran) setup_ctx = NULL;
119: #endif
120: PetscTryMethod(pc, "PCHPDDMSetAuxiliaryMat_C", (PC, IS, Mat, PetscErrorCode (*)(Mat, PetscReal, Vec, Vec, PetscReal, IS, void*), void*), (pc, is, A, setup, setup_ctx));
121: return(0);
122: }
124: static PetscErrorCode PCHPDDMHasNeumannMat_HPDDM(PC pc, PetscBool has)
125: {
126: PC_HPDDM *data = (PC_HPDDM*)pc->data;
129: data->Neumann = has;
130: return(0);
131: }
133: /*@
134: 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.
136: Input Parameters:
137: + pc - preconditioner context
138: - has - Boolean value
140: Level: intermediate
142: .seealso: PCHPDDM, PCHPDDMSetAuxiliaryMat()
143: @*/
144: PetscErrorCode PCHPDDMHasNeumannMat(PC pc, PetscBool has)
145: {
150: PetscTryMethod(pc, "PCHPDDMHasNeumannMat_C", (PC, PetscBool), (pc, has));
151: return(0);
152: }
154: static PetscErrorCode PCHPDDMSetRHSMat_HPDDM(PC pc, Mat B)
155: {
156: PC_HPDDM *data = (PC_HPDDM*)pc->data;
160: PetscObjectReference((PetscObject)B);
161: MatDestroy(&data->B);
162: data->B = B;
163: return(0);
164: }
166: /*@
167: PCHPDDMSetRHSMat - Sets the right-hand side matrix used by PCHPDDM for the concurrent GenEO eigenproblems at the finest level. Must be used in conjuction 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.
169: Input Parameters:
170: + pc - preconditioner context
171: - B - right-hand side sequential matrix
173: Level: advanced
175: .seealso: PCHPDDMSetAuxiliaryMat(), PCHPDDM
176: @*/
177: PetscErrorCode PCHPDDMSetRHSMat(PC pc, Mat B)
178: {
183: if (B) {
185: PetscTryMethod(pc, "PCHPDDMSetRHSMat_C", (PC, Mat), (pc, B));
186: }
187: return(0);
188: }
190: static PetscErrorCode PCSetFromOptions_HPDDM(PetscOptionItems *PetscOptionsObject, PC pc)
191: {
192: PC_HPDDM *data = (PC_HPDDM*)pc->data;
193: PC_HPDDM_Level **levels = data->levels;
194: char prefix[256];
195: int i = 1;
196: PetscMPIInt size, previous;
197: PetscInt n;
198: PetscBool flg = PETSC_TRUE;
202: if (!data->levels) {
203: PetscCalloc1(PETSC_HPDDM_MAXLEVELS, &levels);
204: data->levels = levels;
205: }
206: PetscOptionsHead(PetscOptionsObject, "PCHPDDM options");
207: MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size);
208: previous = size;
209: while (i < PETSC_HPDDM_MAXLEVELS) {
210: PetscInt p = 1;
212: if (!data->levels[i - 1]) {
213: PetscNewLog(pc, data->levels + i - 1);
214: }
215: data->levels[i - 1]->parent = data;
216: /* if the previous level has a single process, it is not possible to coarsen further */
217: if (previous == 1 || !flg) break;
218: data->levels[i - 1]->nu = 0;
219: data->levels[i - 1]->threshold = -1.0;
220: PetscSNPrintf(prefix, sizeof(prefix), "-pc_hpddm_levels_%d_eps_nev", i);
221: PetscOptionsInt(prefix, "Local number of deflation vectors computed by SLEPc", "none", data->levels[i - 1]->nu, &data->levels[i - 1]->nu, NULL);
222: PetscSNPrintf(prefix, sizeof(prefix), "-pc_hpddm_levels_%d_eps_threshold", i);
223: PetscOptionsReal(prefix, "Local threshold for selecting deflation vectors returned by SLEPc", "none", data->levels[i - 1]->threshold, &data->levels[i - 1]->threshold, NULL);
224: /* if there is no prescribed coarsening, just break out of the loop */
225: if (data->levels[i - 1]->threshold <= 0.0 && data->levels[i - 1]->nu <= 0) break;
226: else {
227: ++i;
228: PetscSNPrintf(prefix, sizeof(prefix), "-pc_hpddm_levels_%d_eps_nev", i);
229: PetscOptionsHasName(PetscOptionsObject->options, PetscOptionsObject->prefix, prefix, &flg);
230: if (!flg) {
231: PetscSNPrintf(prefix, sizeof(prefix), "-pc_hpddm_levels_%d_eps_threshold", i);
232: PetscOptionsHasName(PetscOptionsObject->options, PetscOptionsObject->prefix, prefix, &flg);
233: }
234: if (flg) {
235: /* if there are coarsening options for the next level, then register it */
236: /* otherwise, don't to avoid having both options levels_N_p and coarse_p */
237: PetscSNPrintf(prefix, sizeof(prefix), "-pc_hpddm_levels_%d_p", i);
238: PetscOptionsRangeInt(prefix, "Number of processes used to assemble the coarse operator at this level", "none", p, &p, &flg, 1, PetscMax(1, previous / 2));
239: previous = p;
240: }
241: }
242: }
243: data->N = i;
244: n = 1;
245: if (i > 1) {
246: PetscSNPrintf(prefix, sizeof(prefix), "-pc_hpddm_coarse_p");
247: PetscOptionsRangeInt(prefix, "Number of processes used to assemble the coarsest operator", "none", n, &n, NULL, 1, PetscMax(1, previous / 2));
248: PetscOptionsEList("-pc_hpddm_coarse_correction", "Type of coarse correction applied each iteration", "PCHPDDMSetCoarseCorrectionType", PCHPDDMCoarseCorrectionTypes, 3, PCHPDDMCoarseCorrectionTypes[PC_HPDDM_COARSE_CORRECTION_DEFLATED], &n, &flg);
249: if (flg) data->correction = PCHPDDMCoarseCorrectionType(n);
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: }
253: PetscOptionsTail();
254: while (i < PETSC_HPDDM_MAXLEVELS && data->levels[i]) {
255: PetscFree(data->levels[i++]);
256: }
257: return(0);
258: }
260: static PetscErrorCode PCApply_HPDDM(PC pc, Vec x, Vec y)
261: {
262: PC_HPDDM *data = (PC_HPDDM*)pc->data;
266: PetscCitationsRegister(hpddmCitationPC, &citePC);
267: if (data->levels[0]->ksp) {
268: KSPSolve(data->levels[0]->ksp, x, y);
269: } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "No KSP attached to PCHPDDM");
270: return(0);
271: }
273: static PetscErrorCode PCMatApply_HPDDM(PC pc, Mat X, Mat Y)
274: {
275: PC_HPDDM *data = (PC_HPDDM*)pc->data;
279: PetscCitationsRegister(hpddmCitationPC, &citePC);
280: if (data->levels[0]->ksp) {
281: 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;
312: for (n = 0, *gc = 0, *oc = 0; n < data->N; ++n) {
313: if (data->levels[n]->ksp) {
314: Mat P;
315: KSPGetOperators(data->levels[n]->ksp, NULL, &P);
316: MatGetSize(P, &m, NULL);
317: accumulate[0] += m;
318: MatGetInfo(P, MAT_GLOBAL_SUM, &info);
319: accumulate[1] += info.nz_used;
320: if (n == 0) {
321: m1 = m;
322: nnz1 = info.nz_used;
323: }
324: }
325: }
326: *gc = static_cast<PetscReal>(accumulate[0]/m1);
327: *oc = static_cast<PetscReal>(accumulate[1]/nnz1);
328: return(0);
329: }
331: static PetscErrorCode PCView_HPDDM(PC pc, PetscViewer viewer)
332: {
333: PC_HPDDM *data = (PC_HPDDM*)pc->data;
334: PetscViewer subviewer;
335: PetscSubcomm subcomm;
336: PetscReal oc, gc;
337: PetscInt i, tabs;
338: PetscMPIInt size, color, rank;
339: PetscBool ascii;
343: PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &ascii);
344: if (ascii) {
345: PetscViewerASCIIPrintf(viewer, "level%s: %D\n", data->N > 1 ? "s" : "", data->N);
346: PCHPDDMGetComplexities(pc, &gc, &oc);
347: if (data->N > 1) {
348: PetscViewerASCIIPrintf(viewer, "Neumann matrix attached? %s\n", PetscBools[data->Neumann]);
349: PetscViewerASCIIPrintf(viewer, "coarse correction: %s\n", PCHPDDMCoarseCorrectionTypes[data->correction]);
350: PetscViewerASCIIPrintf(viewer, "on process #0, value%s (+ threshold%s if available) for selecting deflation vectors:", data->N > 2 ? "s" : "", data->N > 2 ? "s" : "");
351: PetscViewerASCIIGetTab(viewer, &tabs);
352: PetscViewerASCIISetTab(viewer, 0);
353: for (i = 1; i < data->N; ++i) {
354: PetscViewerASCIIPrintf(viewer, " %D", data->levels[i - 1]->nu);
355: if (data->levels[i - 1]->threshold > -0.1) {
356: PetscViewerASCIIPrintf(viewer, " (%g)", (double)data->levels[i - 1]->threshold);
357: }
358: }
359: PetscViewerASCIIPrintf(viewer, "\n");
360: PetscViewerASCIISetTab(viewer, tabs);
361: }
362: PetscViewerASCIIPrintf(viewer, "grid and operator complexities: %g %g\n", (double)gc, (double)oc);
363: if (data->levels[0]->ksp) {
364: KSPView(data->levels[0]->ksp, viewer);
365: if (data->levels[0]->pc) {
366: PCView(data->levels[0]->pc, viewer);
367: }
368: for (i = 1; i < data->N; ++i) {
369: if (data->levels[i]->ksp) color = 1;
370: else color = 0;
371: MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size);
372: MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank);
373: PetscSubcommCreate(PetscObjectComm((PetscObject)pc), &subcomm);
374: PetscSubcommSetNumber(subcomm, PetscMin(size, 2));
375: PetscSubcommSetTypeGeneral(subcomm, color, rank);
376: PetscViewerASCIIPushTab(viewer);
377: PetscViewerGetSubViewer(viewer, PetscSubcommChild(subcomm), &subviewer);
378: if (color == 1) {
379: KSPView(data->levels[i]->ksp, subviewer);
380: if (data->levels[i]->pc) {
381: PCView(data->levels[i]->pc, subviewer);
382: }
383: PetscViewerFlush(subviewer);
384: }
385: PetscViewerRestoreSubViewer(viewer, PetscSubcommChild(subcomm), &subviewer);
386: PetscViewerASCIIPopTab(viewer);
387: PetscSubcommDestroy(&subcomm);
388: PetscViewerFlush(viewer);
389: }
390: }
391: }
392: return(0);
393: }
395: static PetscErrorCode PCPreSolve_HPDDM(PC pc, KSP ksp, Vec, Vec)
396: {
397: PC_HPDDM *data = (PC_HPDDM*)pc->data;
398: PetscBool flg;
399: 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;
422: PCShellGetContext(pc, (void**)&ctx);
423: KSPGetOptionsPrefix(ctx->ksp, &pcpre);
424: KSPGetOperators(ctx->ksp, &A, &P);
425: /* smoother */
426: PCSetOptionsPrefix(ctx->pc, pcpre);
427: PCSetOperators(ctx->pc, A, P);
428: if (!ctx->v[0]) {
429: VecDuplicateVecs(ctx->D, 1, &ctx->v[0]);
430: if (!std::is_same<PetscScalar, PetscReal>::value) {
431: VecDestroy(&ctx->D);
432: }
433: MatCreateVecs(A, &x, NULL);
434: VecDuplicateVecs(x, 2, &ctx->v[1]);
435: VecDestroy(&x);
436: }
437: return(0);
438: }
440: template<class Type, typename std::enable_if<std::is_same<Type, Vec>::value>::type* = nullptr>
441: PETSC_STATIC_INLINE PetscErrorCode PCHPDDMDeflate_Private(PC pc, Type x, Type y)
442: {
443: PC_HPDDM_Level *ctx;
444: PetscScalar *out;
448: PCShellGetContext(pc, (void**)&ctx);
449: /* going from PETSc to HPDDM numbering */
450: VecScatterBegin(ctx->scatter, x, ctx->v[0][0], INSERT_VALUES, SCATTER_FORWARD);
451: VecScatterEnd(ctx->scatter, x, ctx->v[0][0], INSERT_VALUES, SCATTER_FORWARD);
452: VecGetArrayWrite(ctx->v[0][0], &out);
453: ctx->P->deflation<false>(NULL, out, 1); /* y = Q x */
454: VecRestoreArrayWrite(ctx->v[0][0], &out);
455: /* going from HPDDM to PETSc numbering */
456: VecScatterBegin(ctx->scatter, ctx->v[0][0], y, INSERT_VALUES, SCATTER_REVERSE);
457: VecScatterEnd(ctx->scatter, ctx->v[0][0], y, INSERT_VALUES, SCATTER_REVERSE);
458: return(0);
459: }
461: template<class Type, typename std::enable_if<std::is_same<Type, Mat>::value>::type* = nullptr>
462: PETSC_STATIC_INLINE PetscErrorCode PCHPDDMDeflate_Private(PC pc, Type X, Type Y)
463: {
464: PC_HPDDM_Level *ctx;
465: Vec vX, vY, vC;
466: PetscScalar *out;
467: PetscInt i, m, N, prev = 0;
471: PCShellGetContext(pc, (void**)&ctx);
472: VecGetLocalSize(ctx->v[0][0], &m);
473: MatGetSize(X, NULL, &N);
474: if (ctx->V) {
475: MatGetSize(ctx->V, NULL, &prev);
476: }
477: if (N != prev) {
478: MatDestroy(&ctx->V);
479: MatCreateDense(PetscObjectComm((PetscObject)pc), m, PETSC_DECIDE, PETSC_DECIDE, N, NULL, &ctx->V);
480: }
481: /* going from PETSc to HPDDM numbering */
482: for (i = 0; i < N; ++i) {
483: MatDenseGetColumnVecRead(X, i, &vX);
484: MatDenseGetColumnVecWrite(ctx->V, i, &vC);
485: VecScatterBegin(ctx->scatter, vX, vC, INSERT_VALUES, SCATTER_FORWARD);
486: VecScatterEnd(ctx->scatter, vX, vC, INSERT_VALUES, SCATTER_FORWARD);
487: MatDenseRestoreColumnVecWrite(ctx->V, i, &vC);
488: MatDenseRestoreColumnVecRead(X, i, &vX);
489: }
490: MatDenseGetArrayWrite(ctx->V, &out);
491: if (N != prev) {
492: ctx->P->start(N);
493: prev = N;
494: }
495: ctx->P->deflation<false>(NULL, out, N); /* Y = Q X */
496: MatDenseRestoreArrayWrite(ctx->V, &out);
497: /* going from HPDDM to PETSc numbering */
498: for (i = 0; i < N; ++i) {
499: MatDenseGetColumnVecRead(ctx->V, i, &vC);
500: MatDenseGetColumnVecWrite(Y, i, &vY);
501: VecScatterBegin(ctx->scatter, vC, vY, INSERT_VALUES, SCATTER_REVERSE);
502: VecScatterEnd(ctx->scatter, vC, vY, INSERT_VALUES, SCATTER_REVERSE);
503: MatDenseRestoreColumnVecWrite(Y, i, &vY);
504: MatDenseRestoreColumnVecRead(ctx->V, i, &vC);
505: }
506: return(0);
507: }
509: /*@C
510: 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.
512: .vb
513: (1) y = Pmat^-1 x + Q x,
514: (2) y = Pmat^-1 (I - Amat Q) x + Q x (default),
515: (3) y = (I - Q Amat^T) Pmat^-1 (I - Amat Q) x + Q x.
516: .ve
518: Input Parameters:
519: + pc - preconditioner context
520: - x - input vector
522: Output Parameter:
523: . y - output vector
525: Application Interface Routine: PCApply()
527: Notes:
528: 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.
529: 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.
530: (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.
532: Level: advanced
534: .seealso: PCHPDDM, PCHPDDMCoarseCorrectionType
535: @*/
536: static PetscErrorCode PCHPDDMShellApply(PC pc, Vec x, Vec y)
537: {
538: PC_HPDDM_Level *ctx;
539: Mat A;
543: PCShellGetContext(pc, (void**)&ctx);
544: if (ctx->P) {
545: KSPGetOperators(ctx->ksp, &A, NULL);
546: PCHPDDMDeflate_Private(pc, x, y); /* y = Q x */
547: if (ctx->parent->correction == PC_HPDDM_COARSE_CORRECTION_DEFLATED || ctx->parent->correction == PC_HPDDM_COARSE_CORRECTION_BALANCED) {
548: if (!ctx->parent->normal || ctx != ctx->parent->levels[0]) {
549: MatMult(A, y, ctx->v[1][0]); /* y = A Q x */
550: } else { /* KSPLSQR and finest level */
551: MatMult(A, y, ctx->parent->normal); /* y = A Q x */
552: MatMultTranspose(A, ctx->parent->normal, ctx->v[1][0]); /* y = A^T A Q x */
553: }
554: VecWAXPY(ctx->v[1][1], -1.0, ctx->v[1][0], x); /* y = (I - A Q) x */
555: PCApply(ctx->pc, ctx->v[1][1], ctx->v[1][0]); /* y = M^-1 (I - A Q) x */
556: if (ctx->parent->correction == PC_HPDDM_COARSE_CORRECTION_BALANCED) {
557: if (!ctx->parent->normal || ctx != ctx->parent->levels[0]) {
558: MatMultTranspose(A, ctx->v[1][0], ctx->v[1][1]); /* z = A^T M^-1 (I - A Q) x */
559: } else {
560: MatMult(A, ctx->v[1][0], ctx->parent->normal);
561: MatMultTranspose(A, ctx->parent->normal, ctx->v[1][1]); /* z = A^T A M^-1 (I - A^T A Q) x */
562: }
563: PCHPDDMDeflate_Private(pc, ctx->v[1][1], ctx->v[1][1]);
564: VecAXPY(ctx->v[1][0], -1.0, ctx->v[1][1]); /* y = (I - Q A^T) M^-1 (I - A Q) x */
565: }
566: VecAXPY(y, 1.0, ctx->v[1][0]); /* y = y + Q x */
567: } else if (ctx->parent->correction == PC_HPDDM_COARSE_CORRECTION_ADDITIVE) {
568: PCApply(ctx->pc, x, ctx->v[1][0]);
569: VecAXPY(y, 1.0, ctx->v[1][0]); /* y = M^-1 x + Q x */
570: } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCSHELL from PCHPDDM called with an unknown PCHPDDMCoarseCorrectionType %d", ctx->parent->correction);
571: } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCSHELL from PCHPDDM called with no HPDDM object");
572: return(0);
573: }
575: /*@C
576: PCHPDDMShellMatApply - Variant of PCHPDDMShellApply() for blocks of vectors
578: Input Parameters:
579: + pc - preconditioner context
580: - X - block of input vectors
582: Output Parameter:
583: . Y - block of output vectors
585: Application Interface Routine: PCApply()
587: Level: advanced
589: .seealso: PCHPDDM, PCHPDDMShellMatApply(), PCHPDDMCoarseCorrectionType
590: @*/
591: static PetscErrorCode PCHPDDMShellMatApply(PC pc, Mat X, Mat Y)
592: {
593: PC_HPDDM_Level *ctx;
594: Mat A, C, D;
598: PCShellGetContext(pc, (void**)&ctx);
599: if (ctx->P) {
600: KSPGetOperators(ctx->ksp, &A, NULL);
601: PCHPDDMDeflate_Private(pc, X, Y);
602: if (ctx->parent->correction == PC_HPDDM_COARSE_CORRECTION_DEFLATED || ctx->parent->correction == PC_HPDDM_COARSE_CORRECTION_BALANCED) {
603: MatMatMult(A, Y, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &C);
604: MatDuplicate(C, MAT_COPY_VALUES, &D);
605: MatAYPX(D, -1.0, X, SAME_NONZERO_PATTERN);
606: PCMatApply(ctx->pc, D, C);
607: if (ctx->parent->correction == PC_HPDDM_COARSE_CORRECTION_BALANCED) {
608: #if 0 // TODO FIXME: there is a bug in MatTransposeMatMult(): results are inconsitent with a column-by-column product
609: MatTransposeMatMult(A, C, MAT_REUSE_MATRIX, PETSC_DEFAULT, &D);
610: #else
611: PetscInt N;
612: MatGetSize(D, NULL, &N);
613: for (PetscInt i = 0; i < N; ++i) {
614: Vec cD, cC;
615: MatDenseGetColumnVecRead(C, i, &cC);
616: MatDenseGetColumnVecWrite(D, i, &cD);
617: MatMultTranspose(A, cC, cD);
618: MatDenseRestoreColumnVecWrite(D, i, &cD);
619: MatDenseRestoreColumnVecRead(C, i, &cC);
620: }
621: #endif
622: }
623: MatAXPY(Y, 1.0, C, SAME_NONZERO_PATTERN);
624: MatDestroy(&D);
625: } else if (ctx->parent->correction == PC_HPDDM_COARSE_CORRECTION_ADDITIVE) {
626: MatDuplicate(X, MAT_DO_NOT_COPY_VALUES, &C);
627: PCMatApply(ctx->pc, X, C);
628: MatAXPY(Y, 1.0, C, SAME_NONZERO_PATTERN);
629: } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCSHELL from PCHPDDM called with an unknown PCHPDDMCoarseCorrectionType %d", ctx->parent->correction);
630: MatDestroy(&C);
631: } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCSHELL from PCHPDDM called with no HPDDM object");
632: return(0);
633: }
635: static PetscErrorCode PCHPDDMShellDestroy(PC pc)
636: {
637: PC_HPDDM_Level *ctx;
641: PCShellGetContext(pc, (void**)&ctx);
642: HPDDM::Schwarz<PetscScalar>::destroy(ctx, PETSC_TRUE);
643: VecDestroyVecs(1, &ctx->v[0]);
644: VecDestroyVecs(2, &ctx->v[1]);
645: MatDestroy(&ctx->V);
646: VecDestroy(&ctx->D);
647: VecScatterDestroy(&ctx->scatter);
648: PCDestroy(&ctx->pc);
649: return(0);
650: }
652: static PetscErrorCode PCHPDDMSetUpNeumannOverlap_Private(PC pc)
653: {
654: PC_HPDDM *data = (PC_HPDDM*)pc->data;
658: if (data->setup) {
659: Mat P;
660: Vec x, xt = NULL;
661: PetscReal t = 0.0, s = 0.0;
663: PCGetOperators(pc, NULL, &P);
664: PetscObjectQuery((PetscObject)P, "__SNES_latest_X", (PetscObject*)&x);
665: PetscStackPush("PCHPDDM Neumann callback");
666: (*data->setup)(data->aux, t, x, xt, s, data->is, data->setup_ctx);
667: PetscStackPop;
668: }
669: return(0);
670: }
672: static PetscErrorCode PCSetUp_HPDDM(PC pc)
673: {
674: PC_HPDDM *data = (PC_HPDDM*)pc->data;
675: PC inner;
676: Mat *sub, A, P, N, C, uaux = NULL;
677: Vec xin, v;
678: std::vector<Vec> initial;
679: IS is[1], loc, uis = data->is;
680: ISLocalToGlobalMapping l2g;
681: char prefix[256];
682: const char *pcpre;
683: const PetscScalar *const *ev;
684: PetscInt n, requested = data->N, reused = 0;
685: PetscBool subdomains = PETSC_FALSE, flag = PETSC_FALSE, ismatis;
686: DM dm;
687: PetscErrorCode ierr;
690: if (!data->levels[0]) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Not a single level allocated");
691: PCGetOptionsPrefix(pc, &pcpre);
692: PCGetOperators(pc, &A, &P);
693: if (!data->levels[0]->ksp) {
694: KSPCreate(PetscObjectComm((PetscObject)pc), &data->levels[0]->ksp);
695: PetscSNPrintf(prefix, sizeof(prefix), "%spc_hpddm_%s_", pcpre ? pcpre : "", data->N > 1 ? "levels_1" : "coarse");
696: KSPSetOptionsPrefix(data->levels[0]->ksp, prefix);
697: KSPSetType(data->levels[0]->ksp, KSPPREONLY);
698: } else if (data->levels[0]->ksp->pc && data->levels[0]->ksp->pc->setupcalled == 1 && data->levels[0]->ksp->pc->reusepreconditioner) {
699: /* if the fine level PCSHELL exists, its setup has succeeded, and one wants to reuse it, */
700: /* then just propagate the appropriate flag to the coarser levels */
701: for (n = 0; n < PETSC_HPDDM_MAXLEVELS && data->levels[n]; ++n) {
702: /* the following KSP and PC may be NULL for some processes, hence the check */
703: if (data->levels[n]->ksp) {
704: KSPSetReusePreconditioner(data->levels[n]->ksp, PETSC_TRUE);
705: }
706: if (data->levels[n]->pc) {
707: PCSetReusePreconditioner(data->levels[n]->pc, PETSC_TRUE);
708: }
709: }
710: /* early bail out because there is nothing to do */
711: return(0);
712: } else {
713: /* reset coarser levels */
714: for (n = 1; n < PETSC_HPDDM_MAXLEVELS && data->levels[n]; ++n) {
715: 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) {
716: reused = data->N - n;
717: break;
718: }
719: KSPDestroy(&data->levels[n]->ksp);
720: PCDestroy(&data->levels[n]->pc);
721: }
722: /* check if some coarser levels are being reused */
723: MPI_Allreduce(MPI_IN_PLACE, &reused, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)pc));
724: const int *addr = data->levels[0]->P ? data->levels[0]->P->getAddrLocal() : &HPDDM::i__0;
726: if (addr != &HPDDM::i__0 && reused != data->N - 1) {
727: /* reuse previously computed eigenvectors */
728: ev = data->levels[0]->P->getVectors();
729: if (ev) {
730: initial.reserve(*addr);
731: VecCreateSeqWithArray(PETSC_COMM_SELF, 1, data->levels[0]->P->getDof(), ev[0], &xin);
732: for (n = 0; n < *addr; ++n) {
733: VecDuplicate(xin, &v);
734: VecPlaceArray(xin, ev[n]);
735: VecCopy(xin, v);
736: initial.emplace_back(v);
737: VecResetArray(xin);
738: }
739: VecDestroy(&xin);
740: }
741: }
742: }
743: data->N -= reused;
744: KSPSetOperators(data->levels[0]->ksp, A, P);
746: PetscObjectTypeCompare((PetscObject)P, MATIS, &ismatis);
747: if (!data->is && !ismatis) {
748: PetscErrorCode (*create)(DM, IS*, Mat*, PetscErrorCode (**)(Mat, PetscReal, Vec, Vec, PetscReal, IS, void*), void**) = NULL;
749: PetscErrorCode (*usetup)(Mat, PetscReal, Vec, Vec, PetscReal, IS, void*) = NULL;
750: void *uctx = NULL;
752: /* first see if we can get the data from the DM */
753: MatGetDM(P, &dm);
754: if (!dm) {
755: MatGetDM(A, &dm);
756: }
757: if (!dm) {
758: PCGetDM(pc, &dm);
759: }
760: if (dm) { /* this is the hook for DMPLEX and DMDA for which the auxiliary Mat is the local Neumann matrix */
761: PetscObjectQueryFunction((PetscObject)dm, "DMCreateNeumannOverlap_C", &create);
762: if (create) {
763: (*create)(dm, &uis, &uaux, &usetup, &uctx);
764: data->Neumann = PETSC_TRUE;
765: }
766: }
767: if (!create) {
768: if (!uis) {
769: PetscObjectQuery((PetscObject)pc, "_PCHPDDM_Neumann_IS", (PetscObject*)&uis);
770: PetscObjectReference((PetscObject)uis);
771: }
772: if (!uaux) {
773: PetscObjectQuery((PetscObject)pc, "_PCHPDDM_Neumann_Mat", (PetscObject*)&uaux);
774: PetscObjectReference((PetscObject)uaux);
775: }
776: }
777: PCHPDDMSetAuxiliaryMat(pc, uis, uaux, usetup, uctx);
778: MatDestroy(&uaux);
779: ISDestroy(&uis);
780: }
782: if (!ismatis) {
783: PCHPDDMSetUpNeumannOverlap_Private(pc);
784: }
786: if (data->is || (ismatis && data->N > 1)) {
787: if (ismatis) {
788: std::initializer_list<std::string> list = { MATSEQBAIJ, MATSEQSBAIJ };
789: MatISGetLocalMat(P, &N);
790: std::initializer_list<std::string>::const_iterator it = std::find(list.begin(), list.end(), ((PetscObject)N)->type_name);
791: MatISRestoreLocalMat(P, &N);
792: switch (std::distance(list.begin(), it)) {
793: case 0:
794: MatConvert(P, MATMPIBAIJ, MAT_INITIAL_MATRIX, &C);
795: break;
796: case 1:
797: /* MatCreateSubMatrices does not work with MATSBAIJ and unsorted ISes, so convert to MPIBAIJ */
798: MatConvert(P, MATMPIBAIJ, MAT_INITIAL_MATRIX, &C);
799: MatSetOption(C, MAT_SYMMETRIC, PETSC_TRUE);
800: break;
801: default:
802: MatConvert(P, MATMPIAIJ, MAT_INITIAL_MATRIX, &C);
803: }
804: MatGetLocalToGlobalMapping(P, &l2g, NULL);
805: PetscObjectReference((PetscObject)P);
806: KSPSetOperators(data->levels[0]->ksp, A, C);
807: std::swap(C, P);
808: ISLocalToGlobalMappingGetSize(l2g, &n);
809: ISCreateStride(PETSC_COMM_SELF, n, 0, 1, &loc);
810: ISLocalToGlobalMappingApplyIS(l2g, loc, &is[0]);
811: ISDestroy(&loc);
812: /* the auxiliary Mat is _not_ the local Neumann matrix */
813: /* it is the local Neumann matrix augmented (with zeros) through MatIncreaseOverlap() */
814: data->Neumann = PETSC_FALSE;
815: } else {
816: is[0] = data->is;
817: PetscOptionsGetBool(NULL, pcpre, "-pc_hpddm_define_subdomains", &subdomains, NULL);
818: PetscOptionsGetBool(NULL, pcpre, "-pc_hpddm_has_neumann", &data->Neumann, NULL);
819: ISCreateStride(PetscObjectComm((PetscObject)data->is), P->rmap->n, P->rmap->rstart, 1, &loc);
820: }
821: if (data->N > 1 && (data->aux || ismatis)) {
822: if (!loadedSym) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "HPDDM library not loaded, cannot use more than one level");
823: MatSetOption(P, MAT_SUBMAT_SINGLEIS, PETSC_TRUE);
824: if (ismatis) {
825: /* needed by HPDDM (currently) so that the partition of unity is 0 on subdomain interfaces */
826: MatIncreaseOverlap(P, 1, is, 1);
827: ISDestroy(&data->is);
828: data->is = is[0];
829: } else {
830: if (PetscDefined(USE_DEBUG)) {
831: PetscBool equal;
832: IS intersect;
834: ISIntersect(data->is, loc, &intersect);
835: ISEqualUnsorted(loc, intersect, &equal);
836: ISDestroy(&intersect);
837: if (!equal) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "IS of the auxiliary Mat does not include all local rows of A");
838: }
839: if (!data->Neumann) {
840: PetscObjectTypeCompare((PetscObject)P, MATMPISBAIJ, &flag);
841: if (flag) {
842: /* maybe better to ISSort(is[0]), MatCreateSubMatrices, and then MatPermute */
843: /* but there is no MatPermute_SeqSBAIJ, so as before, just use MATMPIBAIJ */
844: MatConvert(P, MATMPIBAIJ, MAT_INITIAL_MATRIX, &uaux);
845: flag = PETSC_FALSE;
846: }
847: }
848: }
849: if (!uaux) {
850: if (data->Neumann) sub = &data->aux;
851: else {
852: MatCreateSubMatrices(P, 1, is, is, MAT_INITIAL_MATRIX, &sub);
853: }
854: } else {
855: MatCreateSubMatrices(uaux, 1, is, is, MAT_INITIAL_MATRIX, &sub);
856: MatDestroy(&uaux);
857: MatConvert(sub[0], MATSEQSBAIJ, MAT_INPLACE_MATRIX, sub);
858: }
859: /* Vec holding the partition of unity */
860: if (!data->levels[0]->D) {
861: ISGetLocalSize(data->is, &n);
862: VecCreateMPI(PETSC_COMM_SELF, n, PETSC_DETERMINE, &data->levels[0]->D);
863: }
864: if (!data->levels[0]->scatter) {
865: MatCreateVecs(P, &xin, NULL);
866: if (ismatis) {
867: MatDestroy(&P);
868: }
869: VecScatterCreate(xin, data->is, data->levels[0]->D, NULL, &data->levels[0]->scatter);
870: VecDestroy(&xin);
871: }
872: if (data->levels[0]->P) {
873: /* if the pattern is the same and PCSetUp has previously succeeded, reuse HPDDM buffers and connectivity */
874: HPDDM::Schwarz<PetscScalar>::destroy(data->levels[0], pc->setupcalled < 1 || pc->flag == DIFFERENT_NONZERO_PATTERN ? PETSC_TRUE : PETSC_FALSE);
875: }
876: if (!data->levels[0]->P) data->levels[0]->P = new HPDDM::Schwarz<PetscScalar>();
877: (*loadedSym)(data->levels[0]->P, loc, data->is, sub[0], ismatis ? C : data->aux, initial, &data->N, data->levels);
878: if (!data->Neumann)
879: MatDestroySubMatrices(1, &sub);
880: if (ismatis) data->is = NULL;
881: for (n = 0; n < data->N - 1 + (reused > 0); ++n) {
882: if (data->levels[n]->P) {
883: PC spc;
885: /* force the PC to be PCSHELL to do the coarse grid corrections */
886: KSPSetSkipPCSetFromOptions(data->levels[n]->ksp, PETSC_TRUE);
887: KSPGetPC(data->levels[n]->ksp, &spc);
888: PCSetType(spc, PCSHELL);
889: PCShellSetContext(spc, data->levels[n]);
890: PCShellSetSetUp(spc, PCHPDDMShellSetUp);
891: PCShellSetApply(spc, PCHPDDMShellApply);
892: PCShellSetMatApply(spc, PCHPDDMShellMatApply);
893: PCShellSetDestroy(spc, PCHPDDMShellDestroy);
894: if (!data->levels[n]->pc) {
895: PCCreate(PetscObjectComm((PetscObject)data->levels[n]->ksp), &data->levels[n]->pc);
896: }
897: if (n < reused) {
898: PCSetReusePreconditioner(spc, PETSC_TRUE);
899: PCSetReusePreconditioner(data->levels[n]->pc, PETSC_TRUE);
900: }
901: PCSetUp(spc);
902: }
903: }
904: } else flag = reused ? PETSC_FALSE : PETSC_TRUE;
905: if (!ismatis && subdomains) {
906: if (flag) {
907: KSPGetPC(data->levels[0]->ksp, &inner);
908: } else inner = data->levels[0]->pc;
909: if (inner) {
910: PCSetType(inner, PCASM);
911: if (!inner->setupcalled) {
912: PCASMSetLocalSubdomains(inner, 1, is, &loc);
913: }
914: }
915: }
916: ISDestroy(&loc);
917: } else data->N = 1 + reused; /* enforce this value to 1 + reused if there is no way to build another level */
918: if (requested != data->N + reused) {
919: PetscInfo5(pc, "%D levels requested, only %D built + %D reused. Options for level(s) > %D, including -%spc_hpddm_coarse_ will not be taken into account\n", requested, data->N, reused, data->N, pcpre ? pcpre : "");
920: PetscInfo2(pc, "It is best to tune parameters, e.g., a higher value for -%spc_hpddm_levels_%D_eps_threshold so that at least one local deflation vector will be selected\n", pcpre ? pcpre : "", data->N);
921: /* cannot use PCHPDDMShellDestroy because PCSHELL not set for unassembled levels */
922: for (n = data->N - 1; n < requested - 1; ++n) {
923: if (data->levels[n]->P) {
924: HPDDM::Schwarz<PetscScalar>::destroy(data->levels[n], PETSC_TRUE);
925: VecDestroyVecs(1, &data->levels[n]->v[0]);
926: VecDestroyVecs(2, &data->levels[n]->v[1]);
927: MatDestroy(&data->levels[n]->V);
928: VecDestroy(&data->levels[n]->D);
929: VecScatterDestroy(&data->levels[n]->scatter);
930: }
931: }
932: if (reused) {
933: for (n = reused; n < PETSC_HPDDM_MAXLEVELS && data->levels[n]; ++n) {
934: KSPDestroy(&data->levels[n]->ksp);
935: PCDestroy(&data->levels[n]->pc);
936: }
937: }
938: if (PetscDefined(USE_DEBUG)) SETERRQ7(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "%D levels requested, only %D built + %D reused. Options for level(s) > %D, including -%spc_hpddm_coarse_ will not be taken into account. It is best to tune parameters, e.g., a higher value for -%spc_hpddm_levels_%D_eps_threshold so that at least one local deflation vector will be selected. If you don't want this to error out, compile --with-debugging=0", requested, data->N, reused, data->N, pcpre ? pcpre : "", pcpre ? pcpre : "", data->N);
939: }
941: /* these solvers are created after PCSetFromOptions is called */
942: if (pc->setfromoptionscalled) {
943: for (n = 0; n < data->N; ++n) {
944: if (data->levels[n]->ksp) {
945: KSPSetFromOptions(data->levels[n]->ksp);
946: }
947: if (data->levels[n]->pc) {
948: PCSetFromOptions(data->levels[n]->pc);
949: }
950: }
951: pc->setfromoptionscalled = 0;
952: }
953: data->N += reused;
954: return(0);
955: }
957: /*@
958: PCHPDDMSetCoarseCorrectionType - Sets the coarse correction type.
960: Input Parameters:
961: + pc - preconditioner context
962: - type - PC_HPDDM_COARSE_CORRECTION_DEFLATED, PC_HPDDM_COARSE_CORRECTION_ADDITIVE, or PC_HPDDM_COARSE_CORRECTION_BALANCED
964: Options Database Key:
965: . -pc_hpddm_coarse_correction <deflated, additive, balanced> - type of coarse correction to apply
967: Level: intermediate
969: .seealso: PCHPDDMGetCoarseCorrectionType(), PCHPDDM, PCHPDDMCoarseCorrectionType
970: @*/
971: PetscErrorCode PCHPDDMSetCoarseCorrectionType(PC pc, PCHPDDMCoarseCorrectionType type)
972: {
977: PetscTryMethod(pc, "PCHPDDMSetCoarseCorrectionType_C", (PC, PCHPDDMCoarseCorrectionType), (pc, type));
978: return(0);
979: }
981: /*@
982: PCHPDDMGetCoarseCorrectionType - Gets the coarse correction type.
984: Input Parameter:
985: . pc - preconditioner context
987: Output Parameter:
988: . type - PC_HPDDM_COARSE_CORRECTION_DEFLATED, PC_HPDDM_COARSE_CORRECTION_ADDITIVE, or PC_HPDDM_COARSE_CORRECTION_BALANCED
990: Level: intermediate
992: .seealso: PCHPDDMSetCoarseCorrectionType(), PCHPDDM, PCHPDDMCoarseCorrectionType
993: @*/
994: PetscErrorCode PCHPDDMGetCoarseCorrectionType(PC pc, PCHPDDMCoarseCorrectionType *type)
995: {
1000: PetscUseMethod(pc, "PCHPDDMGetCoarseCorrectionType_C", (PC, PCHPDDMCoarseCorrectionType*), (pc, type));
1001: return(0);
1002: }
1004: static PetscErrorCode PCHPDDMSetCoarseCorrectionType_HPDDM(PC pc, PCHPDDMCoarseCorrectionType type)
1005: {
1006: PC_HPDDM *data = (PC_HPDDM*)pc->data;
1009: if (type < 0 || type > 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PCHPDDMCoarseCorrectionType %d", type);
1010: data->correction = type;
1011: return(0);
1012: }
1014: static PetscErrorCode PCHPDDMGetCoarseCorrectionType_HPDDM(PC pc, PCHPDDMCoarseCorrectionType *type)
1015: {
1016: PC_HPDDM *data = (PC_HPDDM*)pc->data;
1019: *type = data->correction;
1020: return(0);
1021: }
1023: PetscErrorCode HPDDMLoadDL_Private(PetscBool *found) {
1024: char lib[PETSC_MAX_PATH_LEN], dlib[PETSC_MAX_PATH_LEN], dir[PETSC_MAX_PATH_LEN];
1028: PetscStrcpy(dir, "${PETSC_LIB_DIR}");
1029: PetscOptionsGetString(NULL, NULL, "-hpddm_dir", dir, sizeof(dir), NULL);
1030: PetscSNPrintf(lib, sizeof(lib), "%s/libhpddm_petsc", dir);
1031: PetscDLLibraryRetrieve(PETSC_COMM_SELF, lib, dlib, 1024, found);
1032: if (*found) {
1033: PetscDLLibraryAppend(PETSC_COMM_SELF, &PetscDLLibrariesLoaded, dlib);
1034: #if defined(SLEPC_LIB_DIR) /* this variable is passed during SLEPc ./configure since */
1035: } else { /* slepcconf.h is not yet build (and thus can't be included) */
1036: PetscStrcpy(dir, HPDDM_STR(SLEPC_LIB_DIR));
1037: PetscSNPrintf(lib, sizeof(lib), "%s/libhpddm_petsc", dir);
1038: PetscDLLibraryRetrieve(PETSC_COMM_SELF, lib, dlib, 1024, found);
1039: if (*found) {
1040: PetscDLLibraryAppend(PETSC_COMM_SELF, &PetscDLLibrariesLoaded, dlib);
1041: #endif
1042: } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "%s not found", lib);
1043: #if defined(SLEPC_LIB_DIR)
1044: }
1045: #endif
1046: return(0);
1047: }
1049: /*MC
1050: PCHPDDM - Interface with the HPDDM library.
1052: 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.
1054: The matrix to be preconditioned (Pmat) may be unassembled (MATIS) or assembled (MATMPIAIJ, MATMPIBAIJ, or MATMPISBAIJ). For multilevel preconditioning, when using an assembled 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).
1056: Options Database Keys:
1057: + -pc_hpddm_define_subdomains <true, default=false> - on the finest level, calls PCASMSetLocalSubdomains() with the IS supplied in PCHPDDMSetAuxiliaryMat() (only relevant with an assembled Pmat)
1058: . -pc_hpddm_has_neumann <true, default=false> - on the finest level, informs the PC that the local Neumann matrix is supplied in PCHPDDMSetAuxiliaryMat()
1059: - -pc_hpddm_coarse_correction <type, default=deflated> - determines the PCHPDDMCoarseCorrectionType when calling PCApply
1061: Options for subdomain solvers, subdomain eigensolvers (for computing deflation vectors), and the coarse solver can be set with
1062: .vb
1063: -pc_hpddm_levels_%d_pc_
1064: -pc_hpddm_levels_%d_ksp_
1065: -pc_hpddm_levels_%d_eps_
1066: -pc_hpddm_levels_%d_p
1067: -pc_hpddm_levels_%d_mat_type_
1068: -pc_hpddm_coarse_
1069: -pc_hpddm_coarse_p
1070: -pc_hpddm_coarse_mat_type_
1071: .ve
1072: 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).
1074: 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.
1076: 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.
1078: References:
1079: + 2011 - A robust two-level domain decomposition preconditioner for systems of PDEs. Spillane, Dolean, Hauret, Nataf, Pechstein, and Scheichl. Comptes Rendus Mathematique.
1080: . 2013 - Scalable Domain Decomposition Preconditioners For Heterogeneous Elliptic Problems. Jolivet, Hecht, Nataf, and Prud'homme. SC13.
1081: . 2015 - An Introduction to Domain Decomposition Methods: Algorithms, Theory, and Parallel Implementation. Dolean, Jolivet, and Nataf. SIAM.
1082: - 2019 - A Multilevel Schwarz Preconditioner Based on a Hierarchy of Robust Coarse Spaces. Al Daas, Grigori, Jolivet, and Tournier.
1084: Level: intermediate
1086: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCHPDDMSetAuxiliaryMat(), MATIS, PCBDDC, PCDEFLATION, PCTELESCOPE
1087: M*/
1088: PETSC_EXTERN PetscErrorCode PCCreate_HPDDM(PC pc)
1089: {
1090: PC_HPDDM *data;
1091: PetscBool found;
1095: if (!loadedSym) {
1096: HPDDMLoadDL_Private(&found);
1097: if (found) {
1098: PetscDLLibrarySym(PETSC_COMM_SELF, &PetscDLLibrariesLoaded, NULL, "PCHPDDM_Internal", (void**)&loadedSym);
1099: }
1100: }
1101: if (!loadedSym) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCHPDDM_Internal symbol not found in loaded libhpddm_petsc");
1102: PetscNewLog(pc, &data);
1103: pc->data = data;
1104: pc->ops->reset = PCReset_HPDDM;
1105: pc->ops->destroy = PCDestroy_HPDDM;
1106: pc->ops->setfromoptions = PCSetFromOptions_HPDDM;
1107: pc->ops->setup = PCSetUp_HPDDM;
1108: pc->ops->apply = PCApply_HPDDM;
1109: pc->ops->matapply = PCMatApply_HPDDM;
1110: pc->ops->view = PCView_HPDDM;
1111: pc->ops->applytranspose = 0;
1112: pc->ops->applysymmetricleft = 0;
1113: pc->ops->applysymmetricright = 0;
1114: pc->ops->presolve = PCPreSolve_HPDDM;
1115: PetscObjectComposeFunction((PetscObject)pc, "PCHPDDMSetAuxiliaryMat_C", PCHPDDMSetAuxiliaryMat_HPDDM);
1116: PetscObjectComposeFunction((PetscObject)pc, "PCHPDDMHasNeumannMat_C", PCHPDDMHasNeumannMat_HPDDM);
1117: PetscObjectComposeFunction((PetscObject)pc, "PCHPDDMSetRHSMat_C", PCHPDDMSetRHSMat_HPDDM);
1118: PetscObjectComposeFunction((PetscObject)pc, "PCHPDDMSetCoarseCorrectionType_C", PCHPDDMSetCoarseCorrectionType_HPDDM);
1119: PetscObjectComposeFunction((PetscObject)pc, "PCHPDDMGetCoarseCorrectionType_C", PCHPDDMGetCoarseCorrectionType_HPDDM);
1120: return(0);
1121: }
1123: /*@C
1124: PCHPDDMInitializePackage - This function initializes everything in the PCHPDDM package. It is called from PCInitializePackage().
1126: Level: intermediate
1128: .seealso: PetscInitialize()
1129: @*/
1130: PetscErrorCode PCHPDDMInitializePackage(void)
1131: {
1135: if (PCHPDDMPackageInitialized) return(0);
1136: PCHPDDMPackageInitialized = PETSC_TRUE;
1137: PetscRegisterFinalize(PCHPDDMFinalizePackage);
1138: /* general events registered once during package initialization */
1139: /* these events are not triggered in libpetsc, */
1140: /* but rather directly in libhpddm_petsc, */
1141: /* which is in charge of performing the following operations */
1143: /* domain decomposition structure from Pmat sparsity pattern */
1144: PetscLogEventRegister("PCHPDDMStrc", PC_CLASSID, &PC_HPDDM_Strc);
1145: /* Galerkin product, redistribution, and setup */
1146: PetscLogEventRegister("PCHPDDMPtAP", PC_CLASSID, &PC_HPDDM_PtAP);
1147: /* Galerkin product with summation, redistribution, and setup */
1148: PetscLogEventRegister("PCHPDDMPtBP", PC_CLASSID, &PC_HPDDM_PtBP);
1149: /* next level construction using PtAP and PtBP */
1150: PetscLogEventRegister("PCHPDDMNext", PC_CLASSID, &PC_HPDDM_Next);
1151: return(0);
1152: }
1154: /*@C
1155: PCHPDDMFinalizePackage - This function frees everything from the PCHPDDM package. It is called from PetscFinalize().
1157: Level: intermediate
1159: .seealso: PetscFinalize()
1160: @*/
1161: PetscErrorCode PCHPDDMFinalizePackage(void)
1162: {
1164: PCHPDDMPackageInitialized = PETSC_FALSE;
1165: return(0);
1166: }