Actual source code: petscpc.h

petsc-3.5.4 2015-05-23
Report Typos and Errors
  1: /*
  2:       Preconditioner module.
  3: */
  6: #include <petscmat.h>
  7: #include <petscdmtypes.h>

  9: PETSC_EXTERN PetscErrorCode PCInitializePackage(void);

 11: /*
 12:     PCList contains the list of preconditioners currently registered
 13:    These are added with PCRegister()
 14: */
 15: PETSC_EXTERN PetscFunctionList PCList;

 17: /*S
 18:      PC - Abstract PETSc object that manages all preconditioners including direct solvers such as PCLU

 20:    Level: beginner

 22:   Concepts: preconditioners

 24: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types)
 25: S*/
 26: typedef struct _p_PC* PC;

 28: /*J
 29:     PCType - String with the name of a PETSc preconditioner method.

 31:    Level: beginner

 33:    Notes: Click on the links below to see details on a particular solver

 35:           PCRegister() is used to register preconditioners that are then accessible via PCSetType()

 37: .seealso: PCSetType(), PC, PCCreate(), PCRegister(), PCSetFromOptions()
 38: J*/
 39: typedef const char* PCType;
 40: #define PCNONE            "none"
 41: #define PCJACOBI          "jacobi"
 42: #define PCSOR             "sor"
 43: #define PCLU              "lu"
 44: #define PCSHELL           "shell"
 45: #define PCBJACOBI         "bjacobi"
 46: #define PCMG              "mg"
 47: #define PCEISENSTAT       "eisenstat"
 48: #define PCILU             "ilu"
 49: #define PCICC             "icc"
 50: #define PCASM             "asm"
 51: #define PCGASM            "gasm"
 52: #define PCKSP             "ksp"
 53: #define PCCOMPOSITE       "composite"
 54: #define PCREDUNDANT       "redundant"
 55: #define PCSPAI            "spai"
 56: #define PCNN              "nn"
 57: #define PCCHOLESKY        "cholesky"
 58: #define PCPBJACOBI        "pbjacobi"
 59: #define PCMAT             "mat"
 60: #define PCHYPRE           "hypre"
 61: #define PCPARMS           "parms"
 62: #define PCFIELDSPLIT      "fieldsplit"
 63: #define PCTFS             "tfs"
 64: #define PCML              "ml"
 65: #define PCGALERKIN        "galerkin"
 66: #define PCEXOTIC          "exotic"
 67: #define PCCP              "cp"
 68: #define PCBFBT            "bfbt"
 69: #define PCLSC             "lsc"
 70: #define PCPYTHON          "python"
 71: #define PCPFMG            "pfmg"
 72: #define PCSYSPFMG         "syspfmg"
 73: #define PCREDISTRIBUTE    "redistribute"
 74: #define PCSVD             "svd"
 75: #define PCGAMG            "gamg"
 76: #define PCSACUSP          "sacusp"        /* these four run on NVIDIA GPUs using CUSP */
 77: #define PCSACUSPPOLY      "sacusppoly"
 78: #define PCBICGSTABCUSP    "bicgstabcusp"
 79: #define PCAINVCUSP        "ainvcusp"
 80: #define PCBDDC            "bddc"
 81: #define PCKACZMARZ        "kaczmarz"

 83: /* Logging support */
 84: PETSC_EXTERN PetscClassId PC_CLASSID;

 86: /*E
 87:     PCSide - If the preconditioner is to be applied to the left, right
 88:      or symmetrically around the operator.

 90:    Level: beginner

 92: .seealso:
 93: E*/
 94: typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide;
 95: #define PC_SIDE_MAX (PC_SYMMETRIC + 1)
 96: PETSC_EXTERN const char *const *const PCSides;

 98: PETSC_EXTERN PetscErrorCode PCCreate(MPI_Comm,PC*);
 99: PETSC_EXTERN PetscErrorCode PCSetType(PC,PCType);
100: PETSC_EXTERN PetscErrorCode PCSetUp(PC);
101: PETSC_EXTERN PetscErrorCode PCSetUpOnBlocks(PC);
102: PETSC_EXTERN PetscErrorCode PCApply(PC,Vec,Vec);
103: PETSC_EXTERN PetscErrorCode PCApplySymmetricLeft(PC,Vec,Vec);
104: PETSC_EXTERN PetscErrorCode PCApplySymmetricRight(PC,Vec,Vec);
105: PETSC_EXTERN PetscErrorCode PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
106: PETSC_EXTERN PetscErrorCode PCApplyTranspose(PC,Vec,Vec);
107: PETSC_EXTERN PetscErrorCode PCApplyTransposeExists(PC,PetscBool *);
108: PETSC_EXTERN PetscErrorCode PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
109: PETSC_EXTERN PetscErrorCode PCSetReusePreconditioner(PC,PetscBool);
110: PETSC_EXTERN PetscErrorCode PCGetReusePreconditioner(PC,PetscBool*);

112: #define PC_FILE_CLASSID 1211222

114: /*E
115:     PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates

117:    Level: advanced

119:    Notes: this must match finclude/petscpc.h and the KSPConvergedReason values in petscksp.h

121: .seealso: PCApplyRichardson()
122: E*/
123: typedef enum {
124:               PCRICHARDSON_CONVERGED_RTOL               =  2,
125:               PCRICHARDSON_CONVERGED_ATOL               =  3,
126:               PCRICHARDSON_CONVERGED_ITS                =  4,
127:               PCRICHARDSON_DIVERGED_DTOL                = -4} PCRichardsonConvergedReason;

129: PETSC_EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*);
130: PETSC_EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscBool *);
131: PETSC_EXTERN PetscErrorCode PCSetInitialGuessNonzero(PC,PetscBool );
132: PETSC_EXTERN PetscErrorCode PCSetUseAmat(PC,PetscBool);
133: PETSC_EXTERN PetscErrorCode PCGetUseAmat(PC,PetscBool*);

135: PETSC_EXTERN PetscErrorCode PCRegisterAll(void);
136: PETSC_EXTERN PetscBool PCRegisterAllCalled;

138: PETSC_EXTERN PetscErrorCode PCRegister(const char[],PetscErrorCode(*)(PC));

140: PETSC_EXTERN PetscErrorCode PCReset(PC);
141: PETSC_EXTERN PetscErrorCode PCDestroy(PC*);
142: PETSC_EXTERN PetscErrorCode PCSetFromOptions(PC);
143: PETSC_EXTERN PetscErrorCode PCGetType(PC,PCType*);

145: PETSC_EXTERN PetscErrorCode PCFactorGetMatrix(PC,Mat*);
146: PETSC_EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
147: PETSC_EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);

149: PETSC_EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat);
150: PETSC_EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*);
151: PETSC_EXTERN PetscErrorCode PCGetOperatorsSet(PC,PetscBool *,PetscBool *);

153: PETSC_EXTERN PetscErrorCode PCView(PC,PetscViewer);
154: PETSC_EXTERN PetscErrorCode PCLoad(PC,PetscViewer);
155: PETSC_STATIC_INLINE PetscErrorCode PCViewFromOptions(PC A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);}

157: PETSC_EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]);
158: PETSC_EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]);
159: PETSC_EXTERN PetscErrorCode PCGetOptionsPrefix(PC,const char*[]);

161: PETSC_EXTERN PetscErrorCode PCComputeExplicitOperator(PC,Mat*);

163: /*
164:       These are used to provide extra scaling of preconditioned
165:    operator for time-stepping schemes like in SUNDIALS
166: */
167: PETSC_EXTERN PetscErrorCode PCGetDiagonalScale(PC,PetscBool *);
168: PETSC_EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec);
169: PETSC_EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec);
170: PETSC_EXTERN PetscErrorCode PCSetDiagonalScale(PC,Vec);

172: /* ------------- options specific to particular preconditioners --------- */

174: PETSC_EXTERN PetscErrorCode PCJacobiSetUseRowMax(PC);
175: PETSC_EXTERN PetscErrorCode PCJacobiSetUseRowSum(PC);
176: PETSC_EXTERN PetscErrorCode PCJacobiSetUseAbs(PC);
177: PETSC_EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType);
178: PETSC_EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal);
179: PETSC_EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt);

181: PETSC_EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal);
182: PETSC_EXTERN PetscErrorCode PCEisenstatNoDiagonalScaling(PC);

184: PETSC_EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
185: PETSC_EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);

187: PETSC_EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec));
188: PETSC_EXTERN PetscErrorCode PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec));
189: PETSC_EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec));
190: PETSC_EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(PC));
191: PETSC_EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*));
192: PETSC_EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer));
193: PETSC_EXTERN PetscErrorCode PCShellSetDestroy(PC,PetscErrorCode (*)(PC));
194: PETSC_EXTERN PetscErrorCode PCShellGetContext(PC,void**);
195: PETSC_EXTERN PetscErrorCode PCShellSetContext(PC,void*);
196: PETSC_EXTERN PetscErrorCode PCShellSetName(PC,const char[]);
197: PETSC_EXTERN PetscErrorCode PCShellGetName(PC,const char*[]);

199: PETSC_EXTERN PetscErrorCode PCFactorSetZeroPivot(PC,PetscReal);

201: PETSC_EXTERN PetscErrorCode PCFactorSetShiftType(PC,MatFactorShiftType);
202: PETSC_EXTERN PetscErrorCode PCFactorSetShiftAmount(PC,PetscReal);

204: PETSC_EXTERN PetscErrorCode PCFactorSetMatSolverPackage(PC,const MatSolverPackage);
205: PETSC_EXTERN PetscErrorCode PCFactorGetMatSolverPackage(PC,const MatSolverPackage*);
206: PETSC_EXTERN PetscErrorCode PCFactorSetUpMatSolverPackage(PC);

208: PETSC_EXTERN PetscErrorCode PCFactorSetFill(PC,PetscReal);
209: PETSC_EXTERN PetscErrorCode PCFactorSetColumnPivot(PC,PetscReal);
210: PETSC_EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC,PetscReal);

212: PETSC_EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC,MatOrderingType);
213: PETSC_EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC,PetscBool );
214: PETSC_EXTERN PetscErrorCode PCFactorSetReuseFill(PC,PetscBool );
215: PETSC_EXTERN PetscErrorCode PCFactorSetUseInPlace(PC);
216: PETSC_EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC);
217: PETSC_EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscBool );

219: PETSC_EXTERN PetscErrorCode PCFactorGetLevels(PC,PetscInt*);
220: PETSC_EXTERN PetscErrorCode PCFactorSetLevels(PC,PetscInt);
221: PETSC_EXTERN PetscErrorCode PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt);

223: PETSC_EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
224: PETSC_EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]);
225: PETSC_EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt);
226: PETSC_EXTERN PetscErrorCode PCASMSetDMSubdomains(PC,PetscBool);
227: PETSC_EXTERN PetscErrorCode PCASMGetDMSubdomains(PC,PetscBool*);
228: PETSC_EXTERN PetscErrorCode PCASMSetSortIndices(PC,PetscBool );

230: /*E
231:     PCASMType - Type of additive Schwarz method to use

233: $  PC_ASM_BASIC        - Symmetric version where residuals from the ghost points are used
234: $                        and computed values in ghost regions are added together.
235: $                        Classical standard additive Schwarz.
236: $  PC_ASM_RESTRICT     - Residuals from ghost points are used but computed values in ghost
237: $                        region are discarded.
238: $                        Default.
239: $  PC_ASM_INTERPOLATE  - Residuals from ghost points are not used, computed values in ghost
240: $                        region are added back in.
241: $  PC_ASM_NONE         - Residuals from ghost points are not used, computed ghost values are
242: $                        discarded.
243: $                        Not very good.

245:    Level: beginner

247: .seealso: PCASMSetType()
248: E*/
249: typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
250: PETSC_EXTERN const char *const PCASMTypes[];

252: PETSC_EXTERN PetscErrorCode PCASMSetType(PC,PCASMType);
253: PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains(Mat,PetscInt,IS*[]);
254: PETSC_EXTERN PetscErrorCode PCASMDestroySubdomains(PetscInt,IS[],IS[]);
255: PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
256: PETSC_EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
257: PETSC_EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);

259: /*E
260:     PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain).

262:    Each subdomain has nested inner and outer parts.  The inner subdomains are assumed to form a non-overlapping covering of the computational
263:    domain, while the outer subdomains contain the inner subdomains and overlap with each other.  This preconditioner will compute
264:    a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in
265:    (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed.

267: $  PC_GASM_BASIC       - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied
268: $                        over the outer subdomains.  As a result, points in the overlap will receive the sum of the corrections
269: $                        from neighboring subdomains.
270: $                        Classical standard additive Schwarz.
271: $  PC_GASM_RESTRICT    - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only
272: $                        (i.e., zeroed out over the overlap portion of the outer subdomain before being applied).  As a result,
273: $                        each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering
274: $                        assumption).
275: $                        Default.
276: $  PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is
277: $                        applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections
278: $                        from neighboring subdomains.
279: $
280: $  PC_GASM_NONE        - Residuals and corrections are zeroed out outside the local subdomains.
281: $                        Not very good.

283:    Level: beginner

285: .seealso: PCGASMSetType()
286: E*/
287: typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;
288: PETSC_EXTERN const char *const PCGASMTypes[];

290: PETSC_EXTERN PetscErrorCode PCGASMSetSubdomains(PC,PetscInt,IS[],IS[]);
291: PETSC_EXTERN PetscErrorCode PCGASMSetTotalSubdomains(PC,PetscInt,PetscBool);
292: PETSC_EXTERN PetscErrorCode PCGASMSetOverlap(PC,PetscInt);
293: PETSC_EXTERN PetscErrorCode PCGASMSetDMSubdomains(PC,PetscBool);
294: PETSC_EXTERN PetscErrorCode PCGASMGetDMSubdomains(PC,PetscBool*);
295: PETSC_EXTERN PetscErrorCode PCGASMSetSortIndices(PC,PetscBool );

297: PETSC_EXTERN PetscErrorCode PCGASMSetType(PC,PCGASMType);
298: PETSC_EXTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat,PetscInt,PetscInt,IS*[],IS*[]);
299: PETSC_EXTERN PetscErrorCode PCGASMDestroySubdomains(PetscInt,IS[],IS[]);
300: PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
301: PETSC_EXTERN PetscErrorCode PCGASMGetSubdomains(PC,PetscInt*,IS*[],IS*[]);
302: PETSC_EXTERN PetscErrorCode PCGASMGetSubmatrices(PC,PetscInt*,Mat*[]);

304: /*E
305:     PCCompositeType - Determines how two or more preconditioner are composed

307: $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
308: $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
309: $                                computed after the previous preconditioner application
310: $  PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
311: $                                computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions)
312: $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
313: $                         where first preconditioner is built from alpha I + S and second from
314: $                         alpha I + R

316:    Level: beginner

318: .seealso: PCCompositeSetType()
319: E*/
320: typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType;
321: PETSC_EXTERN const char *const PCCompositeTypes[];

323: PETSC_EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType);
324: PETSC_EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType);
325: PETSC_EXTERN PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *);
326: PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar);

328: PETSC_EXTERN PetscErrorCode PCRedundantSetNumber(PC,PetscInt);
329: PETSC_EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter);
330: PETSC_EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*);

332: PETSC_EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double);
333: PETSC_EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt);
334: PETSC_EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt);
335: PETSC_EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt);
336: PETSC_EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt);
337: PETSC_EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt);
338: PETSC_EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt);
339: PETSC_EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt);

341: PETSC_EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]);
342: PETSC_EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]);
343: PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
344: PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);

346: PETSC_EXTERN PetscErrorCode PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*,const PetscInt*);
347: PETSC_EXTERN PetscErrorCode PCFieldSplitGetType(PC,PCCompositeType*);
348: PETSC_EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType);
349: PETSC_EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt);
350: PETSC_EXTERN PetscErrorCode PCFieldSplitSetIS(PC,const char[],IS);
351: PETSC_EXTERN PetscErrorCode PCFieldSplitGetIS(PC,const char[],IS*);
352: PETSC_EXTERN PetscErrorCode PCFieldSplitSetDMSplits(PC,PetscBool);
353: PETSC_EXTERN PetscErrorCode PCFieldSplitGetDMSplits(PC,PetscBool*);
354: PETSC_EXTERN PetscErrorCode PCFieldSplitSetDiagUseAmat(PC,PetscBool);
355: PETSC_EXTERN PetscErrorCode PCFieldSplitGetDiagUseAmat(PC,PetscBool*);
356: PETSC_EXTERN PetscErrorCode PCFieldSplitSetOffDiagUseAmat(PC,PetscBool);
357: PETSC_EXTERN PetscErrorCode PCFieldSplitGetOffDiagUseAmat(PC,PetscBool*);


360: /*E
361:     PCFieldSplitSchurPreType - Determines how to precondition Schur complement

363:     Level: intermediate

365: .seealso: PCFieldSplitSetSchurPre()
366: E*/
367: typedef enum {PC_FIELDSPLIT_SCHUR_PRE_SELF,PC_FIELDSPLIT_SCHUR_PRE_SELFP,PC_FIELDSPLIT_SCHUR_PRE_A11,PC_FIELDSPLIT_SCHUR_PRE_USER,PC_FIELDSPLIT_SCHUR_PRE_FULL} PCFieldSplitSchurPreType;
368: PETSC_EXTERN const char *const PCFieldSplitSchurPreTypes[];

370: /*E
371:     PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use

373:     Level: intermediate

375: .seealso: PCFieldSplitSetSchurFactType()
376: E*/
377: typedef enum {
378:   PC_FIELDSPLIT_SCHUR_FACT_DIAG,
379:   PC_FIELDSPLIT_SCHUR_FACT_LOWER,
380:   PC_FIELDSPLIT_SCHUR_FACT_UPPER,
381:   PC_FIELDSPLIT_SCHUR_FACT_FULL
382: } PCFieldSplitSchurFactType;
383: PETSC_EXTERN const char *const PCFieldSplitSchurFactTypes[];

385: PETSC_EXTERN PETSC_DEPRECATED("Use PCFieldSplitSetSchurPre") PetscErrorCode PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat);
386: PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurPre(PC,PCFieldSplitSchurPreType,Mat);
387: PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurPre(PC,PCFieldSplitSchurPreType*,Mat*);
388: PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC,PCFieldSplitSchurFactType);
389: PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*);
390: PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetS(PC,Mat *S);
391: PETSC_EXTERN PetscErrorCode PCFieldSplitSchurRestoreS(PC,Mat *S);

393: PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat);
394: PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat);

396: PETSC_EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscInt,PetscReal*);
397: PETSC_EXTERN PetscErrorCode PCSASetVectors(PC,PetscInt,PetscReal *);

399: PETSC_EXTERN PetscErrorCode PCPythonSetType(PC,const char[]);

401: PETSC_EXTERN PetscErrorCode PCSetDM(PC,DM);
402: PETSC_EXTERN PetscErrorCode PCGetDM(PC,DM*);

404: PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC,void*);
405: PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC,void*);

407: PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetTolerance(PC,PetscReal);
408: PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetIterations(PC,PetscInt);
409: PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool);

411: PETSC_EXTERN PetscErrorCode PCAINVCUSPSetDropTolerance(PC,PetscReal);
412: PETSC_EXTERN PetscErrorCode PCAINVCUSPUseScaling(PC,PetscBool);
413: PETSC_EXTERN PetscErrorCode PCAINVCUSPSetNonzeros(PC,PetscInt);
414: PETSC_EXTERN PetscErrorCode PCAINVCUSPSetLinParameter(PC,PetscInt);
415: /*E
416:     PCPARMSGlobalType - Determines the global preconditioner method in PARMS

418:     Level: intermediate

420: .seealso: PCPARMSSetGlobal()
421: E*/
422: typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType;
423: PETSC_EXTERN const char *const PCPARMSGlobalTypes[];
424: /*E
425:     PCPARMSLocalType - Determines the local preconditioner method in PARMS

427:     Level: intermediate

429: .seealso: PCPARMSSetLocal()
430: E*/
431: typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType;
432: PETSC_EXTERN const char *const PCPARMSLocalTypes[];

434: PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC pc,PCPARMSGlobalType type);
435: PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC pc,PCPARMSLocalType type);
436: PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC pc,PetscReal tol,PetscInt maxits);
437: PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC pc,PetscInt restart);
438: PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC pc,PetscBool nonsym);
439: PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2);

441: /*E
442:     PCGAMGType - type of generalized algebraic multigrid (PCGAMG) method

444:     Level: intermediate

446: .seealso: PCMG, PCSetType(), PCGAMGSetThreshold(), PCGAMGSetThreshold(), PCGAMGSetReuseProl()
447: E*/
448: typedef const char *PCGAMGType;
449: #define PCGAMGAGG         "agg"
450: #define PCGAMGGEO         "geo"
451: #define PCGAMGCLASSICAL   "classical"
452: PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC,PetscInt);
453: PETSC_EXTERN PetscErrorCode PCGAMGSetRepartitioning(PC,PetscBool);
454: PETSC_EXTERN PetscErrorCode PCGAMGSetUseASMAggs(PC,PetscBool);
455: PETSC_EXTERN PetscErrorCode PCGAMGSetSolverType(PC,char[],PetscInt);
456: PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC,PetscReal);
457: PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC,PetscInt);
458: PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC,PetscInt);
459: PETSC_EXTERN PetscErrorCode PCGAMGSetType( PC,PCGAMGType );
460: PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n);
461: PETSC_EXTERN PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n);
462: PETSC_EXTERN PetscErrorCode PCGAMGSetSquareGraph(PC,PetscBool);
463: PETSC_EXTERN PetscErrorCode PCGAMGSetReuseProl(PC,PetscBool);
464: PETSC_EXTERN PetscErrorCode PCGAMGFinalizePackage(void);
465: PETSC_EXTERN PetscErrorCode PCGAMGInitializePackage(void);

467: typedef const char *PCGAMGClassicalType;
468: #define PCGAMGCLASSICALDIRECT   "direct"
469: #define PCGAMGCLASSICALSTANDARD "standard"
470: PETSC_EXTERN PetscErrorCode PCGAMGClassicalSetType(PC,PCGAMGClassicalType);

472: #if defined(PETSC_HAVE_PCBDDC)
473: PETSC_EXTERN PetscErrorCode PCBDDCSetChangeOfBasisLocalMat(PC,Mat);
474: PETSC_EXTERN PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC,IS);
475: PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseningRatio(PC,PetscInt);
476: PETSC_EXTERN PetscErrorCode PCBDDCSetLevels(PC,PetscInt);
477: PETSC_EXTERN PetscErrorCode PCBDDCSetNullSpace(PC,MatNullSpace);
478: PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC,IS);
479: PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC,IS);
480: PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC,IS*);
481: PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC,IS*);
482: PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC,IS);
483: PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC,IS);
484: PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC,IS*);
485: PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC,IS*);
486: PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(PC,PetscInt,IS[]);
487: PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplittingLocal(PC,PetscInt,IS[]);
488: PETSC_EXTERN PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode);
489: PETSC_EXTERN PetscErrorCode PCBDDCCreateFETIDPOperators(PC,Mat*,PC*);
490: PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat,Vec,Vec);
491: PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat,Vec,Vec);
492: #endif

494: PETSC_EXTERN PetscErrorCode PCISSetUseStiffnessScaling(PC,PetscBool);
495: PETSC_EXTERN PetscErrorCode PCISSetSubdomainScalingFactor(PC,PetscScalar);
496: PETSC_EXTERN PetscErrorCode PCISSetSubdomainDiagonalScaling(PC,Vec);

498: /*E
499:     PCMGType - Determines the type of multigrid method that is run.

501:    Level: beginner

503:    Values:
504: +  PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycles()
505: .  PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are
506:                 smoothed before updating the residual. This only uses the
507:                 down smoother, in the preconditioner the upper smoother is ignored
508: .  PC_MG_FULL - same as multiplicative except one also performs grid sequencing,
509:             that is starts on the coarsest grid, performs a cycle, interpolates
510:             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
511:             algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm
512:             calls the V-cycle only on the coarser level and has a post-smoother instead.
513: -  PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level
514:                from a finer

516: .seealso: PCMGSetType()

518: E*/
519: typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType;
520: PETSC_EXTERN const char *const PCMGTypes[];
521: #define PC_MG_CASCADE PC_MG_KASKADE;

523: /*E
524:     PCMGCycleType - Use V-cycle or W-cycle

526:    Level: beginner

528:    Values:
529: +  PC_MG_V_CYCLE
530: -  PC_MG_W_CYCLE

532: .seealso: PCMGSetCycleType()

534: E*/
535: typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType;
536: PETSC_EXTERN const char *const PCMGCycleTypes[];

538: PETSC_EXTERN PetscErrorCode PCMGSetType(PC,PCMGType);
539: PETSC_EXTERN PetscErrorCode PCMGSetLevels(PC,PetscInt,MPI_Comm*);
540: PETSC_EXTERN PetscErrorCode PCMGGetLevels(PC,PetscInt*);

542: PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothUp(PC,PetscInt);
543: PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothDown(PC,PetscInt);
544: PETSC_EXTERN PetscErrorCode PCMGSetCycleType(PC,PCMGCycleType);
545: PETSC_EXTERN PetscErrorCode PCMGSetCycleTypeOnLevel(PC,PetscInt,PCMGCycleType);
546: PETSC_EXTERN PetscErrorCode PCMGSetCyclesOnLevel(PC,PetscInt,PetscInt);
547: PETSC_EXTERN PetscErrorCode PCMGMultiplicativeSetCycles(PC,PetscInt);
548: PETSC_EXTERN PetscErrorCode PCMGSetGalerkin(PC,PetscBool);
549: PETSC_EXTERN PetscErrorCode PCMGGetGalerkin(PC,PetscBool*);


552: PETSC_EXTERN PetscErrorCode PCMGSetRhs(PC,PetscInt,Vec);
553: PETSC_EXTERN PetscErrorCode PCMGSetX(PC,PetscInt,Vec);
554: PETSC_EXTERN PetscErrorCode PCMGSetR(PC,PetscInt,Vec);

556: PETSC_EXTERN PetscErrorCode PCMGSetRestriction(PC,PetscInt,Mat);
557: PETSC_EXTERN PetscErrorCode PCMGGetRestriction(PC,PetscInt,Mat*);
558: PETSC_EXTERN PetscErrorCode PCMGSetInterpolation(PC,PetscInt,Mat);
559: PETSC_EXTERN PetscErrorCode PCMGGetInterpolation(PC,PetscInt,Mat*);
560: PETSC_EXTERN PetscErrorCode PCMGSetRScale(PC,PetscInt,Vec);
561: PETSC_EXTERN PetscErrorCode PCMGGetRScale(PC,PetscInt,Vec*);
562: PETSC_EXTERN PetscErrorCode PCMGSetResidual(PC,PetscInt,PetscErrorCode (*)(Mat,Vec,Vec,Vec),Mat);
563: PETSC_EXTERN PetscErrorCode PCMGResidualDefault(Mat,Vec,Vec,Vec);

565: /*E
566:     PCExoticType - Face based or wirebasket based coarse grid space

568:    Level: beginner

570: .seealso: PCExoticSetType(), PCEXOTIC
571: E*/
572: typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType;
573: PETSC_EXTERN const char *const PCExoticTypes[];
574: PETSC_EXTERN PetscErrorCode PCExoticSetType(PC,PCExoticType);

576: #endif /* __PETSCPC_H */