Actual source code: petscpc.h

petsc-3.3-p7 2013-05-11
  1: /*
  2:       Preconditioner module. 
  3: */
 6:  #include petscdm.h

  8: PETSC_EXTERN PetscErrorCode PCInitializePackage(const char[]);

 10: /*
 11:     PCList contains the list of preconditioners currently registered
 12:    These are added with the PCRegisterDynamic() macro
 13: */
 14: PETSC_EXTERN PetscFList PCList;

 16: /*S
 17:      PC - Abstract PETSc object that manages all preconditioners

 19:    Level: beginner

 21:   Concepts: preconditioners

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

 27: /*J
 28:     PCType - String with the name of a PETSc preconditioner method or the creation function
 29:        with an optional dynamic library name, for example
 30:        http://www.mcs.anl.gov/petsc/lib.a:mypccreate()

 32:    Level: beginner

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

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

 85: /* Logging support */
 86: PETSC_EXTERN PetscClassId PC_CLASSID;

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

 92:    Level: beginner

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

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

112: /*E
113:     PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates

115:    Level: advanced

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

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

127: PETSC_EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*);
128: PETSC_EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscBool *);
129: PETSC_EXTERN PetscErrorCode PCSetInitialGuessNonzero(PC,PetscBool );

131: PETSC_EXTERN PetscErrorCode PCRegisterDestroy(void);
132: PETSC_EXTERN PetscErrorCode PCRegisterAll(const char[]);
133: PETSC_EXTERN PetscBool PCRegisterAllCalled;

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

137: /*MC
138:    PCRegisterDynamic - Adds a method to the preconditioner package.

140:    Synopsis:
141:    PetscErrorCode PCRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(PC))

143:    Not collective

145:    Input Parameters:
146: +  name_solver - name of a new user-defined solver
147: .  path - path (either absolute or relative) the library containing this solver
148: .  name_create - name of routine to create method context
149: -  routine_create - routine to create method context

151:    Notes:
152:    PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners.

154:    If dynamic libraries are used, then the fourth input argument (routine_create)
155:    is ignored.

157:    Sample usage:
158: .vb
159:    PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib",
160:               "MySolverCreate",MySolverCreate);
161: .ve

163:    Then, your solver can be chosen with the procedural interface via
164: $     PCSetType(pc,"my_solver")
165:    or at runtime via the option
166: $     -pc_type my_solver

168:    Level: advanced

170:    Notes: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},  or ${any environmental variable}
171:            occuring in pathname will be replaced with appropriate values.
172:          If your function is not being put into a shared library then use PCRegister() instead

174: .keywords: PC, register

176: .seealso: PCRegisterAll(), PCRegisterDestroy()
177: M*/
178: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
179: #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0)
180: #else
181: #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d)
182: #endif

184: PETSC_EXTERN PetscErrorCode PCReset(PC);
185: PETSC_EXTERN PetscErrorCode PCDestroy(PC*);
186: PETSC_EXTERN PetscErrorCode PCSetFromOptions(PC);
187: PETSC_EXTERN PetscErrorCode PCGetType(PC,const PCType*);

189: PETSC_EXTERN PetscErrorCode PCFactorGetMatrix(PC,Mat*);
190: PETSC_EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
191: PETSC_EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);

193: PETSC_EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat,MatStructure);
194: PETSC_EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*,MatStructure*);
195: PETSC_EXTERN PetscErrorCode PCGetOperatorsSet(PC,PetscBool *,PetscBool *);

197: PETSC_EXTERN PetscErrorCode PCView(PC,PetscViewer);

199: PETSC_EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]);
200: PETSC_EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]);
201: PETSC_EXTERN PetscErrorCode PCGetOptionsPrefix(PC,const char*[]);

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

205: /*
206:       These are used to provide extra scaling of preconditioned 
207:    operator for time-stepping schemes like in SUNDIALS 
208: */
209: PETSC_EXTERN PetscErrorCode PCGetDiagonalScale(PC,PetscBool *);
210: PETSC_EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec);
211: PETSC_EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec);
212: PETSC_EXTERN PetscErrorCode PCSetDiagonalScale(PC,Vec);

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

216: PETSC_EXTERN PetscErrorCode PCJacobiSetUseRowMax(PC);
217: PETSC_EXTERN PetscErrorCode PCJacobiSetUseRowSum(PC);
218: PETSC_EXTERN PetscErrorCode PCJacobiSetUseAbs(PC);
219: PETSC_EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType);
220: PETSC_EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal);
221: PETSC_EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt);

223: PETSC_EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal);
224: PETSC_EXTERN PetscErrorCode PCEisenstatNoDiagonalScaling(PC);

226: #define USE_PRECONDITIONER_MATRIX 0
227: #define USE_TRUE_MATRIX           1
228: PETSC_EXTERN PetscErrorCode PCBJacobiSetUseTrueLocal(PC);
229: PETSC_EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
230: PETSC_EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);

232: PETSC_EXTERN PetscErrorCode PCKSPSetUseTrue(PC);

234: PETSC_EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec));
235: PETSC_EXTERN PetscErrorCode PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec));
236: PETSC_EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec));
237: PETSC_EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(PC));
238: PETSC_EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*));
239: PETSC_EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer));
240: PETSC_EXTERN PetscErrorCode PCShellSetDestroy(PC,PetscErrorCode (*)(PC));
241: PETSC_EXTERN PetscErrorCode PCShellGetContext(PC,void**);
242: PETSC_EXTERN PetscErrorCode PCShellSetContext(PC,void*);
243: PETSC_EXTERN PetscErrorCode PCShellSetName(PC,const char[]);
244: PETSC_EXTERN PetscErrorCode PCShellGetName(PC,const char*[]);

246: PETSC_EXTERN PetscErrorCode PCFactorSetZeroPivot(PC,PetscReal);

248: PETSC_EXTERN PetscErrorCode PCFactorSetShiftType(PC,MatFactorShiftType);
249: PETSC_EXTERN PetscErrorCode PCFactorSetShiftAmount(PC,PetscReal);

251: PETSC_EXTERN PetscErrorCode PCFactorSetMatSolverPackage(PC,const MatSolverPackage);
252: PETSC_EXTERN PetscErrorCode PCFactorGetMatSolverPackage(PC,const MatSolverPackage*);
253: PETSC_EXTERN PetscErrorCode PCFactorSetUpMatSolverPackage(PC);

255: PETSC_EXTERN PetscErrorCode PCFactorSetFill(PC,PetscReal);
256: PETSC_EXTERN PetscErrorCode PCFactorSetColumnPivot(PC,PetscReal);
257: PETSC_EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC,PetscReal);

259: PETSC_EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC,const MatOrderingType);
260: PETSC_EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC,PetscBool );
261: PETSC_EXTERN PetscErrorCode PCFactorSetReuseFill(PC,PetscBool );
262: PETSC_EXTERN PetscErrorCode PCFactorSetUseInPlace(PC);
263: PETSC_EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC);
264: PETSC_EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscBool );

266: PETSC_EXTERN PetscErrorCode PCFactorSetLevels(PC,PetscInt);
267: PETSC_EXTERN PetscErrorCode PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt);

269: PETSC_EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
270: PETSC_EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]);
271: PETSC_EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt);
272: PETSC_EXTERN PetscErrorCode PCASMSetSortIndices(PC,PetscBool );

274: /*E
275:     PCASMType - Type of additive Schwarz method to use

277: $  PC_ASM_BASIC        - Symmetric version where residuals from the ghost points are used
278: $                        and computed values in ghost regions are added together. 
279: $                        Classical standard additive Schwarz.
280: $  PC_ASM_RESTRICT     - Residuals from ghost points are used but computed values in ghost
281: $                        region are discarded. 
282: $                        Default.
283: $  PC_ASM_INTERPOLATE  - Residuals from ghost points are not used, computed values in ghost
284: $                        region are added back in.
285: $  PC_ASM_NONE         - Residuals from ghost points are not used, computed ghost values are 
286: $                        discarded.
287: $                        Not very good.

289:    Level: beginner

291: .seealso: PCASMSetType()
292: E*/
293: typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
294: PETSC_EXTERN const char *PCASMTypes[];

296: PETSC_EXTERN PetscErrorCode PCASMSetType(PC,PCASMType);
297: PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains(Mat,PetscInt,IS*[]);
298: PETSC_EXTERN PetscErrorCode PCASMDestroySubdomains(PetscInt,IS[],IS[]);
299: PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
300: PETSC_EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
301: PETSC_EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);

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

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

311: $  PC_GASM_BASIC       - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied
312: $                        over the outer subdomains.  As a result, points in the overlap will receive the sum of the corrections 
313: $                        from neighboring subdomains.
314: $                        Classical standard additive Schwarz.
315: $  PC_GASM_RESTRICT    - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only
316: $                        (i.e., zeroed out over the overlap portion of the outer subdomain before being applied).  As a result, 
317: $                        each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering 
318: $                        assumption).
319: $                        Default.
320: $  PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is 
321: $                        applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections 
322: $                        from neighboring subdomains.
323: $
324: $  PC_GASM_NONE        - Residuals and corrections are zeroed out outside the local subdomains.
325: $                        Not very good.

327:    Level: beginner

329: .seealso: PCGASMSetType()
330: E*/
331: typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;
332: PETSC_EXTERN const char *PCGASMTypes[];

334: PETSC_EXTERN PetscErrorCode PCGASMSetSubdomains(PC,PetscInt,IS[],IS[]);
335: PETSC_EXTERN PetscErrorCode PCGASMSetTotalSubdomains(PC,PetscInt,PetscBool);
336: PETSC_EXTERN PetscErrorCode PCGASMSetOverlap(PC,PetscInt);
337: PETSC_EXTERN PetscErrorCode PCGASMSetSortIndices(PC,PetscBool );

339: PETSC_EXTERN PetscErrorCode PCGASMSetType(PC,PCGASMType);
340: PETSC_EXTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat,PetscInt,PetscInt,IS*[],IS*[]);
341: PETSC_EXTERN PetscErrorCode PCGASMDestroySubdomains(PetscInt,IS[],IS[]);
342: PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
343: PETSC_EXTERN PetscErrorCode PCGASMGetSubdomains(PC,PetscInt*,IS*[],IS*[]);
344: PETSC_EXTERN PetscErrorCode PCGASMGetSubmatrices(PC,PetscInt*,Mat*[]);

346: /*E
347:     PCCompositeType - Determines how two or more preconditioner are composed

349: $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
350: $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
351: $                                computed after the previous preconditioner application
352: $  PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 
353: $                                computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions)
354: $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
355: $                         where first preconditioner is built from alpha I + S and second from
356: $                         alpha I + R

358:    Level: beginner

360: .seealso: PCCompositeSetType()
361: E*/
362: typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType;
363: PETSC_EXTERN const char *PCCompositeTypes[];

365: PETSC_EXTERN PetscErrorCode PCCompositeSetUseTrue(PC);
366: PETSC_EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType);
367: PETSC_EXTERN PetscErrorCode PCCompositeAddPC(PC,const PCType);
368: PETSC_EXTERN PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *);
369: PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar);

371: PETSC_EXTERN PetscErrorCode PCRedundantSetNumber(PC,PetscInt);
372: PETSC_EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter);
373: PETSC_EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*);

375: PETSC_EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double);
376: PETSC_EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt);
377: PETSC_EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt);
378: PETSC_EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt);
379: PETSC_EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt);
380: PETSC_EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt);
381: PETSC_EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt);
382: PETSC_EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt);

384: PETSC_EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]);
385: PETSC_EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]);
386: PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
387: PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);

389: PETSC_EXTERN PetscErrorCode PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*,const PetscInt*);
390: PETSC_EXTERN PetscErrorCode PCFieldSplitGetType(PC,PCCompositeType*);
391: PETSC_EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType);
392: PETSC_EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt);
393: PETSC_EXTERN PetscErrorCode PCFieldSplitSetIS(PC,const char[],IS);
394: PETSC_EXTERN PetscErrorCode PCFieldSplitGetIS(PC,const char[],IS*);

396: /*E
397:     PCFieldSplitSchurPreType - Determines how to precondition Schur complement

399:     Level: intermediate

401: .seealso: PCFieldSplitSchurPrecondition()
402: E*/
403: typedef enum {PC_FIELDSPLIT_SCHUR_PRE_SELF,PC_FIELDSPLIT_SCHUR_PRE_DIAG,PC_FIELDSPLIT_SCHUR_PRE_USER} PCFieldSplitSchurPreType;
404: PETSC_EXTERN const char *const PCFieldSplitSchurPreTypes[];

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

409:     Level: intermediate

411: .seealso: PCFieldSplitSetSchurFactType()
412: E*/
413: typedef enum {
414:   PC_FIELDSPLIT_SCHUR_FACT_DIAG,
415:   PC_FIELDSPLIT_SCHUR_FACT_LOWER,
416:   PC_FIELDSPLIT_SCHUR_FACT_UPPER,
417:   PC_FIELDSPLIT_SCHUR_FACT_FULL
418: } PCFieldSplitSchurFactType;
419: PETSC_EXTERN const char *const PCFieldSplitSchurFactTypes[];

421: PETSC_EXTERN PetscErrorCode PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat);
422: PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC,PCFieldSplitSchurFactType);
423: PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*);

425: PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat);
426: PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat);

428: PETSC_EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscInt,PetscReal*);
429: PETSC_EXTERN PetscErrorCode PCSASetVectors(PC,PetscInt,PetscReal *);

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

433: PETSC_EXTERN PetscErrorCode PCSetDM(PC,DM);
434: PETSC_EXTERN PetscErrorCode PCGetDM(PC,DM*);

436: PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC,void*);
437: PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC,void*);

439: PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetTolerance(PC,PetscReal);
440: PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetIterations(PC,PetscInt);
441: PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool);

443: PETSC_EXTERN PetscErrorCode PCAINVCUSPSetDropTolerance(PC,PetscReal);
444: PETSC_EXTERN PetscErrorCode PCAINVCUSPUseScaling(PC,PetscBool);
445: PETSC_EXTERN PetscErrorCode PCAINVCUSPSetNonzeros(PC,PetscInt);
446: PETSC_EXTERN PetscErrorCode PCAINVCUSPSetLinParameter(PC,PetscInt);
447: /*E
448:     PCPARMSGlobalType - Determines the global preconditioner method in PARMS

450:     Level: intermediate

452: .seealso: PCPARMSSetGlobal()
453: E*/
454: typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType;
455: PETSC_EXTERN const char *PCPARMSGlobalTypes[];
456: /*E
457:     PCPARMSLocalType - Determines the local preconditioner method in PARMS

459:     Level: intermediate

461: .seealso: PCPARMSSetLocal()
462: E*/
463: typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType;
464: PETSC_EXTERN const char *PCPARMSLocalTypes[];

466: PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC pc,PCPARMSGlobalType type);
467: PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC pc,PCPARMSLocalType type);
468: PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC pc,PetscReal tol,PetscInt maxits);
469: PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC pc,PetscInt restart);
470: PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC pc,PetscBool nonsym);
471: PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2);

473: PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC,PetscInt);
474: PETSC_EXTERN PetscErrorCode PCGAMGSetRepartitioning(PC,PetscBool);
475: PETSC_EXTERN PetscErrorCode PCGAMGSetUseASMAggs(PC,PetscBool);
476: PETSC_EXTERN PetscErrorCode PCGAMGSetSolverType(PC,char[],PetscInt);
477: PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC,PetscReal);
478: PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC,PetscInt);
479: PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC,PetscInt);
480: #define PCGAMGType char*
481: PETSC_EXTERN PetscErrorCode PCGAMGSetType( PC,const PCGAMGType );
482: PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n);
483: PETSC_EXTERN PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n);
484: PETSC_EXTERN PetscErrorCode PCGAMGSetSquareGraph(PC,PetscBool);

486: #if defined(PETSC_HAVE_PCBDDC)
487: /* Enum defining how to treat the coarse problem */
488: typedef enum {SEQUENTIAL_BDDC,REPLICATED_BDDC,PARALLEL_BDDC,MULTILEVEL_BDDC} CoarseProblemType;
489: PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC,IS);
490: PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC,IS*);
491: PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC,IS);
492: PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC,IS*);
493: PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseProblemType(PC,CoarseProblemType);
494: PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(PC,PetscInt,IS[]);
495: #endif

497: PETSC_EXTERN PetscErrorCode PCISSetSubdomainScalingFactor(PC,PetscScalar);


500: #endif /* __PETSCPC_H */