Actual source code: petscpctypes.h
1: #pragma once
3: /* MANSEC = KSP */
4: /* SUBMANSEC = PC */
6: /*S
7: PC - Abstract PETSc object that manages all preconditioners including direct solvers such as `PCLU`
9: Level: beginner
11: .seealso: [](doc_linsolve), [](sec_pc), `PCCreate()`, `PCSetType()`, `PCType`
12: S*/
13: typedef struct _p_PC *PC;
15: /*J
16: PCType - String with the name of a PETSc preconditioner. These are all the preconditioners and direct solvers that PETSc provides.
18: Level: beginner
20: Notes:
21: Use `PCSetType()` or the options database key `-pc_type` to set the preconditioner to use with a given `PC` object
23: `PCRegister()` is used to register preconditioners that are then accessible via `PCSetType()`
25: .seealso: [](doc_linsolve), [](sec_pc), `PCSetType()`, `PC`, `PCCreate()`, `PCRegister()`, `PCSetFromOptions()`, `PCLU`, `PCJACOBI`, `PCBJACOBI`
26: J*/
27: typedef const char *PCType;
28: #define PCNONE "none"
29: #define PCJACOBI "jacobi"
30: #define PCSOR "sor"
31: #define PCLU "lu"
32: #define PCQR "qr"
33: #define PCSHELL "shell"
34: #define PCAMGX "amgx"
35: #define PCBJACOBI "bjacobi"
36: #define PCMG "mg"
37: #define PCEISENSTAT "eisenstat"
38: #define PCILU "ilu"
39: #define PCICC "icc"
40: #define PCASM "asm"
41: #define PCGASM "gasm"
42: #define PCKSP "ksp"
43: #define PCBJKOKKOS "bjkokkos"
44: #define PCCOMPOSITE "composite"
45: #define PCREDUNDANT "redundant"
46: #define PCSPAI "spai"
47: #define PCNN "nn"
48: #define PCCHOLESKY "cholesky"
49: #define PCPBJACOBI "pbjacobi"
50: #define PCVPBJACOBI "vpbjacobi"
51: #define PCMAT "mat"
52: #define PCHYPRE "hypre"
53: #define PCPARMS "parms"
54: #define PCFIELDSPLIT "fieldsplit"
55: #define PCTFS "tfs"
56: #define PCML "ml"
57: #define PCGALERKIN "galerkin"
58: #define PCEXOTIC "exotic"
59: #define PCCP "cp"
60: #define PCLSC "lsc"
61: #define PCPYTHON "python"
62: #define PCPFMG "pfmg"
63: #define PCSMG "smg"
64: #define PCSYSPFMG "syspfmg"
65: #define PCREDISTRIBUTE "redistribute"
66: #define PCSVD "svd"
67: #define PCGAMG "gamg"
68: #define PCCHOWILUVIENNACL "chowiluviennacl"
69: #define PCROWSCALINGVIENNACL "rowscalingviennacl"
70: #define PCSAVIENNACL "saviennacl"
71: #define PCBDDC "bddc"
72: #define PCKACZMARZ "kaczmarz"
73: #define PCTELESCOPE "telescope"
74: #define PCPATCH "patch"
75: #define PCLMVM "lmvm"
76: #define PCHMG "hmg"
77: #define PCDEFLATION "deflation"
78: #define PCHPDDM "hpddm"
79: #define PCH2OPUS "h2opus"
80: #define PCMPI "mpi"
82: /*E
83: PCSide - Determines if the preconditioner is to be applied to the left, right
84: or symmetrically around the operator in `KSPSolve()`.
86: Values:
87: + `PC_LEFT` - applied after the operator is applied
88: . `PC_RIGHT` - applied before the operator is applied
89: - `PC_SYMMETRIC` - a portion of the preconditioner is applied before the operator and the transpose of this portion is applied after the operator is applied.
91: Level: beginner
93: Note:
94: Certain `KSPType` support only a subset of `PCSide` values
96: .seealso: [](sec_pc), `PC`, `KSPSetPCSide()`, `KSP`, `KSPType`, `KSPGetPCSide()`, `KSPSolve()`
97: E*/
98: typedef enum {
99: PC_SIDE_DEFAULT = -1,
100: PC_LEFT = 0,
101: PC_RIGHT = 1,
102: PC_SYMMETRIC = 2
103: } PCSide;
104: #define PC_SIDE_MAX (PC_SYMMETRIC + 1)
106: /*E
107: PCRichardsonConvergedReason - reason a `PCApplyRichardson()` method terminated
109: Level: advanced
111: .seealso: [](sec_pc), `KSPRICHARDSON`, `PC`, `PCApplyRichardson()`
112: E*/
113: typedef enum {
114: PCRICHARDSON_NOT_SET = 0,
115: PCRICHARDSON_CONVERGED_RTOL = 2,
116: PCRICHARDSON_CONVERGED_ATOL = 3,
117: PCRICHARDSON_CONVERGED_ITS = 4,
118: PCRICHARDSON_DIVERGED_DTOL = -4
119: } PCRichardsonConvergedReason;
121: /*E
122: PCJacobiType - Determines what elements of the matrix are used to form the Jacobi preconditioner, that is with the `PCType` of `PCJACOBI`
124: Values:
125: + `PC_JACOBI_DIAGONAL` - use the diagonal entry, if it is zero use one
126: . `PC_JACOBI_ROWL1` - add sum of absolute values in row i, j != i, to diag_ii
127: . `PC_JACOBI_ROWMAX` - use the maximum absolute value in the row
128: - `PC_JACOBI_ROWSUM` - use the sum of the values in the row (not the absolute values)
130: Level: intermediate
132: .seealso: [](sec_pc), `PCJACOBI`, `PC`
133: E*/
134: typedef enum {
135: PC_JACOBI_DIAGONAL,
136: PC_JACOBI_ROWL1,
137: PC_JACOBI_ROWMAX,
138: PC_JACOBI_ROWSUM
139: } PCJacobiType;
141: /*E
142: PCASMType - Determines the type of additive Schwarz method, `PCASM`, to use
144: Values:
145: + `PC_ASM_BASIC` - Symmetric version where residuals from the ghost points are used
146: and computed values in ghost regions are added together.
147: Classical standard additive Schwarz as introduced in {cite}`dryja1987additive`.
148: . `PC_ASM_RESTRICT` - Residuals from ghost points are used but computed values in ghost
149: region are discarded {cite}`cs99`. Default.
150: . `PC_ASM_INTERPOLATE` - Residuals from ghost points are not used, computed values in ghost
151: region are added back in.
152: - `PC_ASM_NONE` - Residuals from ghost points are not used, computed ghost values are
153: discarded. Not very good.
155: Level: beginner
157: .seealso: [](sec_pc), `PC`, `PCASM`, `PCASMSetType()`, `PCGASMType`
158: E*/
159: typedef enum {
160: PC_ASM_BASIC = 3,
161: PC_ASM_RESTRICT = 1,
162: PC_ASM_INTERPOLATE = 2,
163: PC_ASM_NONE = 0
164: } PCASMType;
166: /*E
167: PCGASMType - Determines the type of generalized additive Schwarz method to use (differs from `PCASM` in allowing multiple processors per subdomain) with the `PCType` of `PCGASM`
169: Values:
170: + `PC_GASM_BASIC` - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied
171: over the outer subdomains. As a result, points in the overlap will receive the sum of the corrections
172: from neighboring subdomains. Classical standard additive Schwarz {cite}`dryja1987additive`.
173: . `PC_GASM_RESTRICT` - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only
174: (i.e., zeroed out over the overlap portion of the outer subdomain before being applied). As a result,
175: each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering
176: assumption) {cite}`cs99`. Default.
177: . `PC_GASM_INTERPOLATE` - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is
178: applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections
179: from neighboring subdomains.
180: - `PC_GASM_NONE` - Residuals and corrections are zeroed out outside the local subdomains. Not very good.
182: Level: beginner
184: Note:
185: Each subdomain has nested inner and outer parts. The inner subdomains are assumed to form a non-overlapping covering of the computational
186: domain, while the outer subdomains contain the inner subdomains and overlap with each other. The `PCGASM` preconditioner will compute
187: a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in
188: (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed.
190: Developer Note:
191: Perhaps better to remove this since it matches `PCASMType`
193: .seealso: [](sec_pc), `PCGASM`, `PCASM`, `PC`, `PCGASMSetType()`, `PCASMType`
194: E*/
195: typedef enum {
196: PC_GASM_BASIC = 3,
197: PC_GASM_RESTRICT = 1,
198: PC_GASM_INTERPOLATE = 2,
199: PC_GASM_NONE = 0
200: } PCGASMType;
202: /*E
203: PCCompositeType - Determines how two or more preconditioner are composed with the `PCType` of `PCCOMPOSITE`
205: Values:
206: + `PC_COMPOSITE_ADDITIVE` - results from application of all preconditioners are added together
207: . `PC_COMPOSITE_MULTIPLICATIVE` - preconditioners are applied sequentially to the residual freshly
208: computed after the previous preconditioner application
209: . `PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE` - preconditioners are applied sequentially to the residual freshly
210: computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditioners)
211: . `PC_COMPOSITE_SPECIAL` - This is very special for a matrix of the form $ \alpha I + R + S$
212: where the first preconditioner is built from $\alpha I + S$ and second from $\alpha I + R$
213: . `PC_COMPOSITE_SCHUR` - composes the Schur complement of the matrix from two blocks, see `PCFIELDSPLIT`
214: - `PC_COMPOSITE_GKB` - the generalized Golub-Kahan bidiagonalization preconditioner, see `PCFIELDSPLIT`
216: Level: beginner
218: .seealso: [](sec_pc), `PCCOMPOSITE`, `PCFIELDSPLIT`, `PC`, `PCCompositeSetType()`, `SNESCompositeType`, `PCCompositeSpecialSetAlpha()`
219: E*/
220: typedef enum {
221: PC_COMPOSITE_ADDITIVE,
222: PC_COMPOSITE_MULTIPLICATIVE,
223: PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,
224: PC_COMPOSITE_SPECIAL,
225: PC_COMPOSITE_SCHUR,
226: PC_COMPOSITE_GKB
227: } PCCompositeType;
229: /*E
230: PCFieldSplitSchurPreType - Determines how to precondition a Schur complement arising with the `PCType` of `PCFIELDSPLIT`
232: Values:
233: + `PC_FIELDSPLIT_SCHUR_PRE_SELF` - the preconditioner for the Schur complement is generated from the symbolic representation of the Schur complement matrix.
234: The only preconditioners that currently work with this symbolic representation matrix object are `PCLSC` and `PCHPDDM`
235: . `PC_FIELDSPLIT_SCHUR_PRE_SELFP` - the preconditioning for the Schur complement is generated from an explicitly-assembled approximation $Sp = A11 - A10 diag(A00)^{-1} A01$.
236: This is only a good preconditioner when $diag(A00)$ is a good preconditioner for $A00$. Optionally, $A00$ can be
237: lumped before extracting the diagonal using the additional option `-fieldsplit_1_mat_schur_complement_ainv_type lump`
238: . `PC_FIELDSPLIT_SCHUR_PRE_A11` - the preconditioner for the Schur complement is generated from $A11$, not the Schur complement matrix
239: . `PC_FIELDSPLIT_SCHUR_PRE_USER` - the preconditioner for the Schur complement is generated from the user provided matrix (pre argument
240: to this function).
241: - `PC_FIELDSPLIT_SCHUR_PRE_FULL` - the preconditioner for the Schur complement is generated from the exact Schur complement matrix representation
242: computed internally by `PCFIELDSPLIT` (this is expensive) useful mostly as a test that the Schur complement approach can work for your problem
244: Level: intermediate
246: .seealso: [](sec_pc), `PCFIELDSPLIT`, `PCFieldSplitSetSchurPre()`, `PC`
247: E*/
248: typedef enum {
249: PC_FIELDSPLIT_SCHUR_PRE_SELF,
250: PC_FIELDSPLIT_SCHUR_PRE_SELFP,
251: PC_FIELDSPLIT_SCHUR_PRE_A11,
252: PC_FIELDSPLIT_SCHUR_PRE_USER,
253: PC_FIELDSPLIT_SCHUR_PRE_FULL
254: } PCFieldSplitSchurPreType;
256: /*E
257: PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use with the `PCType` of `PCFIELDSPLIT`
259: Values:
260: + `PC_FIELDSPLIT_SCHUR_FACT_DIAG` - the preconditioner is solving `D`
261: . `PC_FIELDSPLIT_SCHUR_FACT_LOWER` - the preconditioner is solving `L D`
262: . `PC_FIELDSPLIT_SCHUR_FACT_UPPER` - the preconditioner is solving `D U`
263: - `PC_FIELDSPLIT_SCHUR_FACT_FULL` - the preconditioner is solving `L(D U)`
265: where the matrix is factorized as
266: .vb
267: (A B) = (1 0) (A 0) (1 Ainv*B) = L D U
268: (C E) (C*Ainv 1) (0 S) (0 1)
269: .ve
271: Level: intermediate
273: .seealso: [](sec_pc), `PCFIELDSPLIT`, `PCFieldSplitSetSchurFactType()`, `PC`
274: E*/
275: typedef enum {
276: PC_FIELDSPLIT_SCHUR_FACT_DIAG,
277: PC_FIELDSPLIT_SCHUR_FACT_LOWER,
278: PC_FIELDSPLIT_SCHUR_FACT_UPPER,
279: PC_FIELDSPLIT_SCHUR_FACT_FULL
280: } PCFieldSplitSchurFactType;
282: /*E
283: PCPARMSGlobalType - Determines the global preconditioner method in `PCPARMS`
285: Level: intermediate
287: .seealso: [](sec_pc), `PCPARMS`, `PCPARMSSetGlobal()`, `PC`
288: E*/
289: typedef enum {
290: PC_PARMS_GLOBAL_RAS,
291: PC_PARMS_GLOBAL_SCHUR,
292: PC_PARMS_GLOBAL_BJ
293: } PCPARMSGlobalType;
295: /*E
296: PCPARMSLocalType - Determines the local preconditioner method in `PCPARMS`
298: Level: intermediate
300: .seealso: [](sec_pc), `PCPARMS`, `PCPARMSSetLocal()`, `PC`
301: E*/
302: typedef enum {
303: PC_PARMS_LOCAL_ILU0,
304: PC_PARMS_LOCAL_ILUK,
305: PC_PARMS_LOCAL_ILUT,
306: PC_PARMS_LOCAL_ARMS
307: } PCPARMSLocalType;
309: /*J
310: PCGAMGType - type of generalized algebraic multigrid `PCGAMG` method
312: Values:
313: + `PCGAMGAGG` - (the default) smoothed aggregation algorithm, robust, very well tested
314: . `PCGAMGGEO` - geometric coarsening, uses mesh generator to produce coarser meshes, limited to triangles, not supported, reference implementation (2D)
315: - `PCGAMGCLASSICAL` - classical algebraic multigrid preconditioner, incomplete, not supported, reference implementation
317: Level: intermediate
319: .seealso: [](sec_pc), `PCGAMG`, `PCMG`, `PC`, `PCSetType()`, `PCGAMGSetThreshold()`, `PCGAMGSetThreshold()`, `PCGAMGSetReuseInterpolation()`
320: J*/
321: typedef const char *PCGAMGType;
322: #define PCGAMGAGG "agg"
323: #define PCGAMGGEO "geo"
324: #define PCGAMGCLASSICAL "classical"
326: typedef const char *PCGAMGClassicalType;
327: #define PCGAMGCLASSICALDIRECT "direct"
328: #define PCGAMGCLASSICALSTANDARD "standard"
330: /*E
331: PCMGType - Determines the type of multigrid method that is run with the `PCType` of `PCMG`
333: Values:
334: + `PC_MG_MULTIPLICATIVE` (default) - traditional V or W cycle as determined by `PCMGSetCycleType()`
335: . `PC_MG_ADDITIVE` - the additive multigrid preconditioner where all levels are
336: smoothed before updating the residual. This only uses the
337: down smoother, in the preconditioner the upper smoother is ignored
338: . `PC_MG_FULL` - same as multiplicative except one also performs grid sequencing,
339: that is starts on the coarsest grid, performs a cycle, interpolates
340: to the next, performs a cycle etc. This is much like the F-cycle presented in "Multigrid" by Trottenberg, Oosterlee, Schuller page 49, but that
341: algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm
342: calls the V-cycle only on the coarser level and has a post-smoother instead.
343: - `PC_MG_KASKADE` - Cascadic or Kaskadic multigrid, like full multigrid except one never goes back to a coarser level from a finer
345: Level: beginner
347: .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetType()`, `PCMGSetCycleType()`, `PCMGSetCycleTypeOnLevel()`
348: E*/
349: typedef enum {
350: PC_MG_MULTIPLICATIVE,
351: PC_MG_ADDITIVE,
352: PC_MG_FULL,
353: PC_MG_KASKADE
354: } PCMGType;
355: #define PC_MG_CASCADE PC_MG_KASKADE;
357: /*E
358: PCMGCycleType - Determines which of V-cycle or W-cycle to use with the `PCType` of `PCMG` or `PCGAMG`
360: Values:
361: + `PC_MG_V_CYCLE` - use the V cycle
362: - `PC_MG_W_CYCLE` - use the W cycle
364: Level: beginner
366: .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetCycleType()`
367: E*/
368: typedef enum {
369: PC_MG_CYCLE_V = 1,
370: PC_MG_CYCLE_W = 2
371: } PCMGCycleType;
373: /*E
374: PCMGalerkinType - Determines if the coarse grid operators are computed via the Galerkin process with the `PCType` of `PCMG`
376: Values:
377: + `PC_MG_GALERKIN_PMAT` - computes the `pmat` (matrix from which the preconditioner is built) via the Galerkin process from the finest grid
378: . `PC_MG_GALERKIN_MAT` - computes the `mat` (matrix used to apply the operator) via the Galerkin process from the finest grid
379: . `PC_MG_GALERKIN_BOTH` - computes both the `mat` and `pmat` via the Galerkin process (if pmat == mat the construction is only done once
380: - `PC_MG_GALERKIN_NONE` - neither operator is computed via the Galerkin process, the user must provide the operator
382: Level: beginner
384: Note:
385: Users should never set `PC_MG_GALERKIN_EXTERNAL`, it is used by `PCHYPRE` and `PCML`
387: .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetCycleType()`
388: E*/
389: typedef enum {
390: PC_MG_GALERKIN_BOTH,
391: PC_MG_GALERKIN_PMAT,
392: PC_MG_GALERKIN_MAT,
393: PC_MG_GALERKIN_NONE,
394: PC_MG_GALERKIN_EXTERNAL
395: } PCMGGalerkinType;
397: /*E
398: PCExoticType - Determines which of the face-based or wirebasket-based coarse grid space to use with the `PCType` of `PCEXOTIC`
400: Level: beginner
402: .seealso: [](sec_pc), `PCExoticSetType()`, `PCEXOTIC`
403: E*/
404: typedef enum {
405: PC_EXOTIC_FACE,
406: PC_EXOTIC_WIREBASKET
407: } PCExoticType;
409: /*E
410: PCBDDCInterfaceExtType - Defines how interface balancing is extended into the interior of subdomains with the `PCType` of `PCBDDC`
412: Values:
413: + `PC_BDDC_INTERFACE_EXT_DIRICHLET` - solves Dirichlet interior problem; this is the standard BDDC algorithm
414: - `PC_BDDC_INTERFACE_EXT_LUMP` - skips interior solve; sometimes called $M_1$ and associated with "lumped FETI-DP"
416: Level: intermediate
418: .seealso: [](sec_pc), `PCBDDC`, `PC`
419: E*/
420: typedef enum {
421: PC_BDDC_INTERFACE_EXT_DIRICHLET,
422: PC_BDDC_INTERFACE_EXT_LUMP
423: } PCBDDCInterfaceExtType;
425: /*E
426: PCMGCoarseSpaceType - Function space for coarse space for adaptive interpolation
428: Level: beginner
430: .seealso: [](sec_pc), `PCMGSetAdaptCoarseSpaceType()`, `PCMG`, `PC`
431: E*/
432: typedef enum {
433: PCMG_ADAPT_NONE,
434: PCMG_ADAPT_POLYNOMIAL,
435: PCMG_ADAPT_HARMONIC,
436: PCMG_ADAPT_EIGENVECTOR,
437: PCMG_ADAPT_GENERALIZED_EIGENVECTOR,
438: PCMG_ADAPT_GDSW
439: } PCMGCoarseSpaceType;
441: /*E
442: PCPatchConstructType - Determines the algorithm used to construct patches for the `PCPATCH` preconditioner
444: Level: beginner
446: .seealso: [](sec_pc), `PCPatchSetConstructType()`, `PCPATCH`, `PC`
447: E*/
448: typedef enum {
449: PC_PATCH_STAR,
450: PC_PATCH_VANKA,
451: PC_PATCH_PARDECOMP,
452: PC_PATCH_USER,
453: PC_PATCH_PYTHON
454: } PCPatchConstructType;
456: /*E
457: PCDeflationSpaceType - Type of deflation used by `PCType` `PCDEFLATION`
459: Values:
460: + `PC_DEFLATION_SPACE_HAAR` - directly assembled based on Haar (db2) wavelet with overflowed filter cuted-off
461: . `PC_DEFLATION_SPACE_DB2` - `MATCOMPOSITE` of 1-lvl matices based on db2 (2 coefficient Daubechies / Haar wavelet)
462: . `PC_DEFLATION_SPACE_DB4` - same as above, but with db4 (4 coefficient Daubechies)
463: . `PC_DEFLATION_SPACE_DB8` - same as above, but with db8 (8 coefficient Daubechies)
464: . `PC_DEFLATION_SPACE_DB16` - same as above, but with db16 (16 coefficient Daubechies)
465: . `PC_DEFLATION_SPACE_BIORTH22` - same as above, but with biorthogonal 2.2 (6 coefficients)
466: . `PC_DEFLATION_SPACE_MEYER` - same as above, but with Meyer/FIR (62 coefficients)
467: . `PC_DEFLATION_SPACE_AGGREGATION` - aggregates local indices (given by operator matrix distribution) into a subdomain
468: - `PC_DEFLATION_SPACE_USER` - indicates space set by user
470: Level: intermediate
472: Note:
473: Wavelet-based space (except Haar) can be used in multilevel deflation.
475: .seealso: [](sec_pc), `PCDeflationSetSpaceToCompute()`, `PCDEFLATION`, `PC`
476: E*/
477: typedef enum {
478: PC_DEFLATION_SPACE_HAAR,
479: PC_DEFLATION_SPACE_DB2,
480: PC_DEFLATION_SPACE_DB4,
481: PC_DEFLATION_SPACE_DB8,
482: PC_DEFLATION_SPACE_DB16,
483: PC_DEFLATION_SPACE_BIORTH22,
484: PC_DEFLATION_SPACE_MEYER,
485: PC_DEFLATION_SPACE_AGGREGATION,
486: PC_DEFLATION_SPACE_USER
487: } PCDeflationSpaceType;
489: /*E
490: PCHPDDMCoarseCorrectionType - Type of coarse correction used by `PCType` `PCHPDDM`
492: Values:
493: + `PC_HPDDM_COARSE_CORRECTION_DEFLATED` (default) - eq. (1) in `PCHPDDMShellApply()`
494: . `PC_HPDDM_COARSE_CORRECTION_ADDITIVE` - eq. (2)
495: . `PC_HPDDM_COARSE_CORRECTION_BALANCED` - eq. (3)
496: - `PC_HPDDM_COARSE_CORRECTION_NONE` - no coarse correction (mostly useful for debugging)
498: Level: intermediate
500: .seealso: [](sec_pc), `PCHPDDM`, `PC`, `PCSetType()`, `PCHPDDMShellApply()`
501: E*/
502: typedef enum {
503: PC_HPDDM_COARSE_CORRECTION_DEFLATED,
504: PC_HPDDM_COARSE_CORRECTION_ADDITIVE,
505: PC_HPDDM_COARSE_CORRECTION_BALANCED,
506: PC_HPDDM_COARSE_CORRECTION_NONE
507: } PCHPDDMCoarseCorrectionType;
509: /*E
510: PCHPDDMSchurPreType - Type of `PCHPDDM` preconditioner for a `MATSCHURCOMPLEMENT` generated by `PCFIELDSPLIT` with `PCFieldSplitSchurPreType` set to `PC_FIELDSPLIT_SCHUR_PRE_SELF`
512: Values:
513: + `PC_HPDDM_SCHUR_PRE_LEAST_SQUARES` (default) - only with a near-zero A11 block and A10 = A01^T; a preconditioner for solving A01^T A00^-1 A01 x = b
514: is built by approximating the Schur complement with (inv(sqrt(diag(A00))) A01)^T (inv(sqrt(diag(A00))) A01)
515: and by considering the associated linear least squares problem
516: - `PC_HPDDM_SCHUR_PRE_GENEO` - only with A10 = A01^T, `PCHPDDMSetAuxiliaryMat()` called on the `PC` of the A00 block, and if A11 is nonzero,
517: then `PCHPDDMSetAuxiliaryMat()` must be called on the associated `PC` as well (it is built automatically for the
518: user otherwise); the Schur complement `PC` is set internally to `PCKSP`, with the prefix `-fieldsplit_1_pc_hpddm_`;
519: the operator associated to the `PC` is spectrally equivalent to the original Schur complement
521: Level: advanced
523: .seealso: [](sec_pc), `PCHPDDM`, `PC`, `PCFIELDSPLIT`, `PC_FIELDSPLIT_SCHUR_PRE_SELF`, `PCFieldSplitSetSchurPre()`, `PCHPDDMSetAuxiliaryMat()`
524: E*/
525: typedef enum {
526: PC_HPDDM_SCHUR_PRE_LEAST_SQUARES,
527: PC_HPDDM_SCHUR_PRE_GENEO
528: } PCHPDDMSchurPreType;
530: /*E
531: PCFailedReason - indicates the type of `PC` failure. That is why the construction of the preconditioner, `PCSetUp()`, or its use, `PCApply()`, failed
533: Level: beginner
535: .seealso: [](sec_pc), `PC`, `PCGetFailedReason()`, `PCSetUp()`
536: E*/
537: typedef enum {
538: PC_SETUP_ERROR = -1,
539: PC_NOERROR = 0,
540: PC_FACTOR_STRUCT_ZEROPIVOT = 1,
541: PC_FACTOR_NUMERIC_ZEROPIVOT = 2,
542: PC_FACTOR_OUTMEMORY = 3,
543: PC_FACTOR_OTHER = 4,
544: PC_INCONSISTENT_RHS = 5,
545: PC_SUBPC_ERROR = 6
546: } PCFailedReason;
548: /*E
549: PCGAMGLayoutType - Layout for reduced grids for `PCType` `PCGAMG`
551: Level: intermediate
553: .seealso: [](sec_pc), `PCGAMG`, `PC`, `PCGAMGSetCoarseGridLayoutType()`
554: E*/
555: typedef enum {
556: PCGAMG_LAYOUT_COMPACT,
557: PCGAMG_LAYOUT_SPREAD
558: } PCGAMGLayoutType;