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 PCBFBT "bfbt"
61: #define PCLSC "lsc"
62: #define PCPYTHON "python"
63: #define PCPFMG "pfmg"
64: #define PCSMG "smg"
65: #define PCSYSPFMG "syspfmg"
66: #define PCREDISTRIBUTE "redistribute"
67: #define PCSVD "svd"
68: #define PCGAMG "gamg"
69: #define PCCHOWILUVIENNACL "chowiluviennacl"
70: #define PCROWSCALINGVIENNACL "rowscalingviennacl"
71: #define PCSAVIENNACL "saviennacl"
72: #define PCBDDC "bddc"
73: #define PCKACZMARZ "kaczmarz"
74: #define PCTELESCOPE "telescope"
75: #define PCPATCH "patch"
76: #define PCLMVM "lmvm"
77: #define PCHMG "hmg"
78: #define PCDEFLATION "deflation"
79: #define PCHPDDM "hpddm"
80: #define PCH2OPUS "h2opus"
81: #define PCMPI "mpi"
83: /*E
84: PCSide - Determines if the preconditioner is to be applied to the left, right
85: or symmetrically around the operator in `KSPSolve()`.
87: Values:
88: + `PC_LEFT` - applied after the operator is applied
89: . `PC_RIGHT` - applied before the operator is applied
90: - `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.
92: Level: beginner
94: Note:
95: Certain `KSPType` support only a subset of `PCSide` values
97: .seealso: [](sec_pc), `PC`, `KSPSetPCSide()`, `KSP`, `KSPType`, `KSPGetPCSide()`, `KSPSolve()`
98: E*/
99: typedef enum {
100: PC_SIDE_DEFAULT = -1,
101: PC_LEFT = 0,
102: PC_RIGHT = 1,
103: PC_SYMMETRIC = 2
104: } PCSide;
105: #define PC_SIDE_MAX (PC_SYMMETRIC + 1)
107: /*E
108: PCRichardsonConvergedReason - reason a `PCApplyRichardson()` method terminated
110: Level: advanced
112: .seealso: [](sec_pc), `KSPRICHARDSON`, `PC`, `PCApplyRichardson()`
113: E*/
114: typedef enum {
115: PCRICHARDSON_NOT_SET = 0,
116: PCRICHARDSON_CONVERGED_RTOL = 2,
117: PCRICHARDSON_CONVERGED_ATOL = 3,
118: PCRICHARDSON_CONVERGED_ITS = 4,
119: PCRICHARDSON_DIVERGED_DTOL = -4
120: } PCRichardsonConvergedReason;
122: /*E
123: PCJacobiType - Determines what elements of the matrix are used to form the Jacobi preconditioner, that is with the `PCType` of `PCJACOBI`
125: Values:
126: + `PC_JACOBI_DIAGONAL` - use the diagonal entry, if it is zero use one
127: . `PC_JACOBI_ROWL1` - add sum of absolute values in row i, j != i, to diag_ii
128: . `PC_JACOBI_ROWMAX` - use the maximum absolute value in the row
129: - `PC_JACOBI_ROWSUM` - use the sum of the values in the row (not the absolute values)
131: Level: intermediate
133: .seealso: [](sec_pc), `PCJACOBI`, `PC`
134: E*/
135: typedef enum {
136: PC_JACOBI_DIAGONAL,
137: PC_JACOBI_ROWL1,
138: PC_JACOBI_ROWMAX,
139: PC_JACOBI_ROWSUM
140: } PCJacobiType;
142: /*E
143: PCASMType - Determines the type of additive Schwarz method, `PCASM`, to use
145: Values:
146: + `PC_ASM_BASIC` - Symmetric version where residuals from the ghost points are used
147: and computed values in ghost regions are added together.
148: Classical standard additive Schwarz as introduced in {cite}`dryja1987additive`.
149: . `PC_ASM_RESTRICT` - Residuals from ghost points are used but computed values in ghost
150: region are discarded {cite}`cs99`. Default.
151: . `PC_ASM_INTERPOLATE` - Residuals from ghost points are not used, computed values in ghost
152: region are added back in.
153: - `PC_ASM_NONE` - Residuals from ghost points are not used, computed ghost values are
154: discarded. Not very good.
156: Level: beginner
158: .seealso: [](sec_pc), `PC`, `PCASM`, `PCASMSetType()`, `PCGASMType`
159: E*/
160: typedef enum {
161: PC_ASM_BASIC = 3,
162: PC_ASM_RESTRICT = 1,
163: PC_ASM_INTERPOLATE = 2,
164: PC_ASM_NONE = 0
165: } PCASMType;
167: /*E
168: 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`
170: Values:
171: + `PC_GASM_BASIC` - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied
172: over the outer subdomains. As a result, points in the overlap will receive the sum of the corrections
173: from neighboring subdomains. Classical standard additive Schwarz {cite}`dryja1987additive`.
174: . `PC_GASM_RESTRICT` - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only
175: (i.e., zeroed out over the overlap portion of the outer subdomain before being applied). As a result,
176: each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering
177: assumption) {cite}`cs99`. Default.
178: . `PC_GASM_INTERPOLATE` - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is
179: applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections
180: from neighboring subdomains.
181: - `PC_GASM_NONE` - Residuals and corrections are zeroed out outside the local subdomains. Not very good.
183: Level: beginner
185: Note:
186: Each subdomain has nested inner and outer parts. The inner subdomains are assumed to form a non-overlapping covering of the computational
187: domain, while the outer subdomains contain the inner subdomains and overlap with each other. The `PCGASM` preconditioner will compute
188: a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in
189: (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed.
191: Developer Note:
192: Perhaps better to remove this since it matches `PCASMType`
194: .seealso: [](sec_pc), `PCGASM`, `PCASM`, `PC`, `PCGASMSetType()`, `PCASMType`
195: E*/
196: typedef enum {
197: PC_GASM_BASIC = 3,
198: PC_GASM_RESTRICT = 1,
199: PC_GASM_INTERPOLATE = 2,
200: PC_GASM_NONE = 0
201: } PCGASMType;
203: /*E
204: PCCompositeType - Determines how two or more preconditioner are composed with the `PCType` of `PCCOMPOSITE`
206: Values:
207: + `PC_COMPOSITE_ADDITIVE` - results from application of all preconditioners are added together
208: . `PC_COMPOSITE_MULTIPLICATIVE` - preconditioners are applied sequentially to the residual freshly
209: computed after the previous preconditioner application
210: . `PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE` - preconditioners are applied sequentially to the residual freshly
211: computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditioners)
212: . `PC_COMPOSITE_SPECIAL` - This is very special for a matrix of the form $ \alpha I + R + S$
213: where the first preconditioner is built from $\alpha I + S$ and second from $\alpha I + R$
214: . `PC_COMPOSITE_SCHUR` - composes the Schur complement of the matrix from two blocks, see `PCFIELDSPLIT`
215: - `PC_COMPOSITE_GKB` - the generalized Golub-Kahan bidiagonalization preconditioner, see `PCFIELDSPLIT`
217: Level: beginner
219: .seealso: [](sec_pc), `PCCOMPOSITE`, `PCFIELDSPLIT`, `PC`, `PCCompositeSetType()`, `SNESCompositeType`, `PCCompositeSpecialSetAlpha()`
220: E*/
221: typedef enum {
222: PC_COMPOSITE_ADDITIVE,
223: PC_COMPOSITE_MULTIPLICATIVE,
224: PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,
225: PC_COMPOSITE_SPECIAL,
226: PC_COMPOSITE_SCHUR,
227: PC_COMPOSITE_GKB
228: } PCCompositeType;
230: /*E
231: PCFieldSplitSchurPreType - Determines how to precondition a Schur complement arising with the `PCType` of `PCFIELDSPLIT`
233: Values:
234: + `PC_FIELDSPLIT_SCHUR_PRE_SELF` - the preconditioner for the Schur complement is generated from the symbolic representation of the Schur complement matrix.
235: The only preconditioners that currently work with this symbolic representation matrix object are `PCLSC` and `PCHPDDM`
236: . `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$.
237: This is only a good preconditioner when $diag(A00)$ is a good preconditioner for $A00$. Optionally, $A00$ can be
238: lumped before extracting the diagonal using the additional option `-fieldsplit_1_mat_schur_complement_ainv_type lump`
239: . `PC_FIELDSPLIT_SCHUR_PRE_A11` - the preconditioner for the Schur complement is generated from $A11$, not the Schur complement matrix
240: . `PC_FIELDSPLIT_SCHUR_PRE_USER` - the preconditioner for the Schur complement is generated from the user provided matrix (pre argument
241: to this function).
242: - `PC_FIELDSPLIT_SCHUR_PRE_FULL` - the preconditioner for the Schur complement is generated from the exact Schur complement matrix representation
243: computed internally by `PCFIELDSPLIT` (this is expensive) useful mostly as a test that the Schur complement approach can work for your problem
245: Level: intermediate
247: .seealso: [](sec_pc), `PCFIELDSPLIT`, `PCFieldSplitSetSchurPre()`, `PC`
248: E*/
249: typedef enum {
250: PC_FIELDSPLIT_SCHUR_PRE_SELF,
251: PC_FIELDSPLIT_SCHUR_PRE_SELFP,
252: PC_FIELDSPLIT_SCHUR_PRE_A11,
253: PC_FIELDSPLIT_SCHUR_PRE_USER,
254: PC_FIELDSPLIT_SCHUR_PRE_FULL
255: } PCFieldSplitSchurPreType;
257: /*E
258: PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use with the `PCType` of `PCFIELDSPLIT`
260: Values:
261: + `PC_FIELDSPLIT_SCHUR_FACT_DIAG` - the preconditioner is solving `D`
262: . `PC_FIELDSPLIT_SCHUR_FACT_LOWER` - the preconditioner is solving `L D`
263: . `PC_FIELDSPLIT_SCHUR_FACT_UPPER` - the preconditioner is solving `D U`
264: - `PC_FIELDSPLIT_SCHUR_FACT_FULL` - the preconditioner is solving `L(D U)`
266: where the matrix is factorized as
267: .vb
268: (A B) = (1 0) (A 0) (1 Ainv*B) = L D U
269: (C E) (C*Ainv 1) (0 S) (0 1)
270: .ve
272: Level: intermediate
274: .seealso: [](sec_pc), `PCFIELDSPLIT`, `PCFieldSplitSetSchurFactType()`, `PC`
275: E*/
276: typedef enum {
277: PC_FIELDSPLIT_SCHUR_FACT_DIAG,
278: PC_FIELDSPLIT_SCHUR_FACT_LOWER,
279: PC_FIELDSPLIT_SCHUR_FACT_UPPER,
280: PC_FIELDSPLIT_SCHUR_FACT_FULL
281: } PCFieldSplitSchurFactType;
283: /*E
284: PCPARMSGlobalType - Determines the global preconditioner method in `PCPARMS`
286: Level: intermediate
288: .seealso: [](sec_pc), `PCPARMS`, `PCPARMSSetGlobal()`, `PC`
289: E*/
290: typedef enum {
291: PC_PARMS_GLOBAL_RAS,
292: PC_PARMS_GLOBAL_SCHUR,
293: PC_PARMS_GLOBAL_BJ
294: } PCPARMSGlobalType;
296: /*E
297: PCPARMSLocalType - Determines the local preconditioner method in `PCPARMS`
299: Level: intermediate
301: .seealso: [](sec_pc), `PCPARMS`, `PCPARMSSetLocal()`, `PC`
302: E*/
303: typedef enum {
304: PC_PARMS_LOCAL_ILU0,
305: PC_PARMS_LOCAL_ILUK,
306: PC_PARMS_LOCAL_ILUT,
307: PC_PARMS_LOCAL_ARMS
308: } PCPARMSLocalType;
310: /*J
311: PCGAMGType - type of generalized algebraic multigrid `PCGAMG` method
313: Values:
314: + `PCGAMGAGG` - (the default) smoothed aggregation algorithm, robust, very well tested
315: . `PCGAMGGEO` - geometric coarsening, uses mesh generator to produce coarser meshes, limited to triangles, not supported, reference implementation (2D)
316: - `PCGAMGCLASSICAL` - classical algebraic multigrid preconditioner, incomplete, not supported, reference implementation
318: Level: intermediate
320: .seealso: [](sec_pc), `PCGAMG`, `PCMG`, `PC`, `PCSetType()`, `PCGAMGSetThreshold()`, `PCGAMGSetThreshold()`, `PCGAMGSetReuseInterpolation()`
321: J*/
322: typedef const char *PCGAMGType;
323: #define PCGAMGAGG "agg"
324: #define PCGAMGGEO "geo"
325: #define PCGAMGCLASSICAL "classical"
327: typedef const char *PCGAMGClassicalType;
328: #define PCGAMGCLASSICALDIRECT "direct"
329: #define PCGAMGCLASSICALSTANDARD "standard"
331: /*E
332: PCMGType - Determines the type of multigrid method that is run with the `PCType` of `PCMG`
334: Values:
335: + `PC_MG_MULTIPLICATIVE` (default) - traditional V or W cycle as determined by `PCMGSetCycleType()`
336: . `PC_MG_ADDITIVE` - the additive multigrid preconditioner where all levels are
337: smoothed before updating the residual. This only uses the
338: down smoother, in the preconditioner the upper smoother is ignored
339: . `PC_MG_FULL` - same as multiplicative except one also performs grid sequencing,
340: that is starts on the coarsest grid, performs a cycle, interpolates
341: 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
342: algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm
343: calls the V-cycle only on the coarser level and has a post-smoother instead.
344: - `PC_MG_KASKADE` - Cascadic or Kaskadic multigrid, like full multigrid except one never goes back to a coarser level from a finer
346: Level: beginner
348: .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetType()`, `PCMGSetCycleType()`, `PCMGSetCycleTypeOnLevel()`
349: E*/
350: typedef enum {
351: PC_MG_MULTIPLICATIVE,
352: PC_MG_ADDITIVE,
353: PC_MG_FULL,
354: PC_MG_KASKADE
355: } PCMGType;
356: #define PC_MG_CASCADE PC_MG_KASKADE;
358: /*E
359: PCMGCycleType - Determines which of V-cycle or W-cycle to use with the `PCType` of `PCMG` or `PCGAMG`
361: Values:
362: + `PC_MG_V_CYCLE` - use the V cycle
363: - `PC_MG_W_CYCLE` - use the W cycle
365: Level: beginner
367: .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetCycleType()`
368: E*/
369: typedef enum {
370: PC_MG_CYCLE_V = 1,
371: PC_MG_CYCLE_W = 2
372: } PCMGCycleType;
374: /*E
375: PCMGalerkinType - Determines if the coarse grid operators are computed via the Galerkin process with the `PCType` of `PCMG`
377: Values:
378: + `PC_MG_GALERKIN_PMAT` - computes the `pmat` (matrix from which the preconditioner is built) via the Galerkin process from the finest grid
379: . `PC_MG_GALERKIN_MAT` - computes the `mat` (matrix used to apply the operator) via the Galerkin process from the finest grid
380: . `PC_MG_GALERKIN_BOTH` - computes both the `mat` and `pmat` via the Galerkin process (if pmat == mat the construction is only done once
381: - `PC_MG_GALERKIN_NONE` - neither operator is computed via the Galerkin process, the user must provide the operator
383: Level: beginner
385: Note:
386: Users should never set `PC_MG_GALERKIN_EXTERNAL`, it is used by `PCHYPRE` and `PCML`
388: .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetCycleType()`
389: E*/
390: typedef enum {
391: PC_MG_GALERKIN_BOTH,
392: PC_MG_GALERKIN_PMAT,
393: PC_MG_GALERKIN_MAT,
394: PC_MG_GALERKIN_NONE,
395: PC_MG_GALERKIN_EXTERNAL
396: } PCMGGalerkinType;
398: /*E
399: PCExoticType - Determines which of the face-based or wirebasket-based coarse grid space to use with the `PCType` of `PCEXOTIC`
401: Level: beginner
403: .seealso: [](sec_pc), `PCExoticSetType()`, `PCEXOTIC`
404: E*/
405: typedef enum {
406: PC_EXOTIC_FACE,
407: PC_EXOTIC_WIREBASKET
408: } PCExoticType;
410: /*E
411: PCBDDCInterfaceExtType - Defines how interface balancing is extended into the interior of subdomains with the `PCType` of `PCBDDC`
413: Values:
414: + `PC_BDDC_INTERFACE_EXT_DIRICHLET` - solves Dirichlet interior problem; this is the standard BDDC algorithm
415: - `PC_BDDC_INTERFACE_EXT_LUMP` - skips interior solve; sometimes called $M_1$ and associated with "lumped FETI-DP"
417: Level: intermediate
419: .seealso: [](sec_pc), `PCBDDC`, `PC`
420: E*/
421: typedef enum {
422: PC_BDDC_INTERFACE_EXT_DIRICHLET,
423: PC_BDDC_INTERFACE_EXT_LUMP
424: } PCBDDCInterfaceExtType;
426: /*E
427: PCMGCoarseSpaceType - Function space for coarse space for adaptive interpolation
429: Level: beginner
431: .seealso: [](sec_pc), `PCMGSetAdaptCoarseSpaceType()`, `PCMG`, `PC`
432: E*/
433: typedef enum {
434: PCMG_ADAPT_NONE,
435: PCMG_ADAPT_POLYNOMIAL,
436: PCMG_ADAPT_HARMONIC,
437: PCMG_ADAPT_EIGENVECTOR,
438: PCMG_ADAPT_GENERALIZED_EIGENVECTOR,
439: PCMG_ADAPT_GDSW
440: } PCMGCoarseSpaceType;
442: /*E
443: PCPatchConstructType - Determines the algorithm used to construct patches for the `PCPATCH` preconditioner
445: Level: beginner
447: .seealso: [](sec_pc), `PCPatchSetConstructType()`, `PCPATCH`, `PC`
448: E*/
449: typedef enum {
450: PC_PATCH_STAR,
451: PC_PATCH_VANKA,
452: PC_PATCH_PARDECOMP,
453: PC_PATCH_USER,
454: PC_PATCH_PYTHON
455: } PCPatchConstructType;
457: /*E
458: PCDeflationSpaceType - Type of deflation used by `PCType` `PCDEFLATION`
460: Values:
461: + `PC_DEFLATION_SPACE_HAAR` - directly assembled based on Haar (db2) wavelet with overflowed filter cuted-off
462: . `PC_DEFLATION_SPACE_DB2` - `MATCOMPOSITE` of 1-lvl matices based on db2 (2 coefficient Daubechies / Haar wavelet)
463: . `PC_DEFLATION_SPACE_DB4` - same as above, but with db4 (4 coefficient Daubechies)
464: . `PC_DEFLATION_SPACE_DB8` - same as above, but with db8 (8 coefficient Daubechies)
465: . `PC_DEFLATION_SPACE_DB16` - same as above, but with db16 (16 coefficient Daubechies)
466: . `PC_DEFLATION_SPACE_BIORTH22` - same as above, but with biorthogonal 2.2 (6 coefficients)
467: . `PC_DEFLATION_SPACE_MEYER` - same as above, but with Meyer/FIR (62 coefficients)
468: . `PC_DEFLATION_SPACE_AGGREGATION` - aggregates local indices (given by operator matrix distribution) into a subdomain
469: - `PC_DEFLATION_SPACE_USER` - indicates space set by user
471: Level: intermediate
473: Note:
474: Wavelet-based space (except Haar) can be used in multilevel deflation.
476: .seealso: [](sec_pc), `PCDeflationSetSpaceToCompute()`, `PCDEFLATION`, `PC`
477: E*/
478: typedef enum {
479: PC_DEFLATION_SPACE_HAAR,
480: PC_DEFLATION_SPACE_DB2,
481: PC_DEFLATION_SPACE_DB4,
482: PC_DEFLATION_SPACE_DB8,
483: PC_DEFLATION_SPACE_DB16,
484: PC_DEFLATION_SPACE_BIORTH22,
485: PC_DEFLATION_SPACE_MEYER,
486: PC_DEFLATION_SPACE_AGGREGATION,
487: PC_DEFLATION_SPACE_USER
488: } PCDeflationSpaceType;
490: /*E
491: PCHPDDMCoarseCorrectionType - Type of coarse correction used by `PCType` `PCHPDDM`
493: Values:
494: + `PC_HPDDM_COARSE_CORRECTION_DEFLATED` (default) - eq. (1) in `PCHPDDMShellApply()`
495: . `PC_HPDDM_COARSE_CORRECTION_ADDITIVE` - eq. (2)
496: . `PC_HPDDM_COARSE_CORRECTION_BALANCED` - eq. (3)
497: - `PC_HPDDM_COARSE_CORRECTION_NONE` - no coarse correction (mostly useful for debugging)
499: Level: intermediate
501: .seealso: [](sec_pc), `PCHPDDM`, `PC`, `PCSetType()`, `PCHPDDMShellApply()`
502: E*/
503: typedef enum {
504: PC_HPDDM_COARSE_CORRECTION_DEFLATED,
505: PC_HPDDM_COARSE_CORRECTION_ADDITIVE,
506: PC_HPDDM_COARSE_CORRECTION_BALANCED,
507: PC_HPDDM_COARSE_CORRECTION_NONE
508: } PCHPDDMCoarseCorrectionType;
510: /*E
511: PCHPDDMSchurPreType - Type of `PCHPDDM` preconditioner for a `MATSCHURCOMPLEMENT` generated by `PCFIELDSPLIT` with `PCFieldSplitSchurPreType` set to `PC_FIELDSPLIT_SCHUR_PRE_SELF`
513: Values:
514: + `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
515: is built by approximating the Schur complement with (inv(sqrt(diag(A00))) A01)^T (inv(sqrt(diag(A00))) A01)
516: and by considering the associated linear least squares problem
517: - `PC_HPDDM_SCHUR_PRE_GENEO` - only with A10 = A01^T, `PCHPDDMSetAuxiliaryMat()` called on the `PC` of the A00 block, and if A11 is nonzero,
518: then `PCHPDDMSetAuxiliaryMat()` must be called on the associated `PC` as well (it is built automatically for the
519: user otherwise); the Schur complement `PC` is set internally to `PCKSP`, with the prefix `-fieldsplit_1_pc_hpddm_`;
520: the operator associated to the `PC` is spectrally equivalent to the original Schur complement
522: Level: advanced
524: .seealso: [](sec_pc), `PCHPDDM`, `PC`, `PCFIELDSPLIT`, `PC_FIELDSPLIT_SCHUR_PRE_SELF`, `PCFieldSplitSetSchurPre()`, `PCHPDDMSetAuxiliaryMat()`
525: E*/
526: typedef enum {
527: PC_HPDDM_SCHUR_PRE_LEAST_SQUARES,
528: PC_HPDDM_SCHUR_PRE_GENEO
529: } PCHPDDMSchurPreType;
531: /*E
532: PCFailedReason - indicates the type of `PC` failure. That is why the construction of the preconditioner, `PCSetUp()`, or its use, `PCApply()`, failed
534: Level: beginner
536: .seealso: [](sec_pc), `PC`, `PCGetFailedReason()`, `PCSetUp()`
537: E*/
538: typedef enum {
539: PC_SETUP_ERROR = -1,
540: PC_NOERROR = 0,
541: PC_FACTOR_STRUCT_ZEROPIVOT = 1,
542: PC_FACTOR_NUMERIC_ZEROPIVOT = 2,
543: PC_FACTOR_OUTMEMORY = 3,
544: PC_FACTOR_OTHER = 4,
545: PC_INCONSISTENT_RHS = 5,
546: PC_SUBPC_ERROR = 6
547: } PCFailedReason;
549: /*E
550: PCGAMGLayoutType - Layout for reduced grids for `PCType` `PCGAMG`
552: Level: intermediate
554: .seealso: [](sec_pc), `PCGAMG`, `PC`, `PCGAMGSetCoarseGridLayoutType()`
555: E*/
556: typedef enum {
557: PCGAMG_LAYOUT_COMPACT,
558: PCGAMG_LAYOUT_SPREAD
559: } PCGAMGLayoutType;