Actual source code: petscpc.h

  1: /*
  2:       Preconditioner module. 
  3: */
 6:  #include petscdm.h
  7: PETSC_EXTERN_CXX_BEGIN

  9: extern PetscErrorCode   PCInitializePackage(const char[]);

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

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

 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 or the creation function
 30:        with an optional dynamic library name, for example
 31:        http://www.mcs.anl.gov/petsc/lib.a:mypccreate()

 33:    Level: beginner

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

 37: .seealso: PCSetType(), PC, PCCreate()
 38: J*/
 39: #define PCType char*
 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 PCPROMETHEUS      "prometheus"
 66: #define PCGALERKIN        "galerkin"
 67: #define PCEXOTIC          "exotic"
 68: #define PCHMPI            "hmpi"
 69: #define PCSUPPORTGRAPH    "supportgraph"
 70: #define PCASA             "asa"
 71: #define PCCP              "cp"
 72: #define PCBFBT            "bfbt"
 73: #define PCLSC             "lsc"
 74: #define PCPYTHON          "python"
 75: #define PCPFMG            "pfmg"
 76: #define PCSYSPFMG         "syspfmg"
 77: #define PCREDISTRIBUTE    "redistribute"
 78: #define PCSVD             "svd"
 79: #define PCGAMG            "gamg"
 80: #define PCSACUSP          "sacusp"        /* these four run on NVIDIA GPUs using CUSP */
 81: #define PCSACUSPPOLY      "sacusppoly"
 82: #define PCBICGSTABCUSP    "bicgstabcusp"
 83: #define PCAINVCUSP        "ainvcusp"
 84: #define PCBDDC            "bddc"

 86: /* Logging support */
 87: extern PetscClassId  PC_CLASSID;

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

 93:    Level: beginner

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

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

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

116:    Level: advanced

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

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

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

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

136: extern PetscErrorCode  PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC));

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

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

144:    Not collective

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

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

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

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

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

169:    Level: advanced

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

175: .keywords: PC, register

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

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

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

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

198: extern PetscErrorCode  PCView(PC,PetscViewer);

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

204: extern PetscErrorCode  PCComputeExplicitOperator(PC,Mat*);

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

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

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

224: extern PetscErrorCode  PCEisenstatSetOmega(PC,PetscReal);
225: extern PetscErrorCode  PCEisenstatNoDiagonalScaling(PC);

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

233: extern PetscErrorCode  PCKSPSetUseTrue(PC);

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

247: extern PetscErrorCode  PCFactorSetZeroPivot(PC,PetscReal);

249: extern PetscErrorCode  PCFactorSetShiftType(PC,MatFactorShiftType);
250: extern PetscErrorCode  PCFactorSetShiftAmount(PC,PetscReal);

252: extern PetscErrorCode  PCFactorSetMatSolverPackage(PC,const MatSolverPackage);
253: extern PetscErrorCode  PCFactorGetMatSolverPackage(PC,const MatSolverPackage*);
254: extern PetscErrorCode  PCFactorSetUpMatSolverPackage(PC);

256: extern PetscErrorCode  PCFactorSetFill(PC,PetscReal);
257: extern PetscErrorCode  PCFactorSetColumnPivot(PC,PetscReal);
258: extern PetscErrorCode  PCFactorReorderForNonzeroDiagonal(PC,PetscReal);

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

267: extern PetscErrorCode  PCFactorSetLevels(PC,PetscInt);
268: extern PetscErrorCode  PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt);

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

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

278: $  PC_ASM_BASIC - symmetric version where residuals from the ghost points are used
279: $                 and computed values in ghost regions are added together. Classical
280: $                 standard additive Schwarz
281: $  PC_ASM_RESTRICT - residuals from ghost points are used but computed values in ghost
282: $                    region are discarded. 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 - ghost point residuals are not used, computed ghost values are discarded
286: $                not very good.                

288:    Level: beginner

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

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

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

305: $  PC_GASM_BASIC - symmetric version where residuals from the ghost points are used
306: $                 and computed values in ghost regions are added together. Classical
307: $                 standard additive Schwarz
308: $  PC_GASM_RESTRICT - residuals from ghost points are used but computed values in ghost
309: $                    region are discarded. Default
310: $  PC_GASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost
311: $                       region are added back in
312: $  PC_GASM_NONE - ghost point residuals are not used, computed ghost values are discarded
313: $                not very good.                

315:    Level: beginner

317: .seealso: PCGASMSetType()
318: E*/
319: typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;
320: extern const char *PCGASMTypes[];

322: extern PetscErrorCode  PCGASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
323: extern PetscErrorCode  PCGASMSetTotalSubdomains(PC,PetscInt);
324: extern PetscErrorCode  PCGASMSetOverlap(PC,PetscInt);
325: extern PetscErrorCode  PCGASMSetSortIndices(PC,PetscBool );

327: extern PetscErrorCode  PCGASMSetType(PC,PCGASMType);
328: extern PetscErrorCode  PCGASMCreateSubdomains(Mat,PetscInt,IS*[]);
329: extern PetscErrorCode  PCGASMDestroySubdomains(PetscInt,IS[],IS[]);
330: extern PetscErrorCode  PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
331: extern PetscErrorCode  PCGASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
332: extern PetscErrorCode  PCGASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);

334: /*E
335:     PCCompositeType - Determines how two or more preconditioner are composed

337: $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
338: $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
339: $                                computed after the previous preconditioner application
340: $  PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 
341: $                                computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions)
342: $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
343: $                         where first preconditioner is built from alpha I + S and second from
344: $                         alpha I + R

346:    Level: beginner

348: .seealso: PCCompositeSetType()
349: E*/
350: typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType;
351: extern const char *PCCompositeTypes[];

353: extern PetscErrorCode  PCCompositeSetUseTrue(PC);
354: extern PetscErrorCode  PCCompositeSetType(PC,PCCompositeType);
355: extern PetscErrorCode  PCCompositeAddPC(PC,PCType);
356: extern PetscErrorCode  PCCompositeGetPC(PC,PetscInt,PC *);
357: extern PetscErrorCode  PCCompositeSpecialSetAlpha(PC,PetscScalar);

359: extern PetscErrorCode  PCRedundantSetNumber(PC,PetscInt);
360: extern PetscErrorCode  PCRedundantSetScatter(PC,VecScatter,VecScatter);
361: extern PetscErrorCode  PCRedundantGetOperators(PC,Mat*,Mat*);

363: extern PetscErrorCode  PCSPAISetEpsilon(PC,double);
364: extern PetscErrorCode  PCSPAISetNBSteps(PC,PetscInt);
365: extern PetscErrorCode  PCSPAISetMax(PC,PetscInt);
366: extern PetscErrorCode  PCSPAISetMaxNew(PC,PetscInt);
367: extern PetscErrorCode  PCSPAISetBlockSize(PC,PetscInt);
368: extern PetscErrorCode  PCSPAISetCacheSize(PC,PetscInt);
369: extern PetscErrorCode  PCSPAISetVerbose(PC,PetscInt);
370: extern PetscErrorCode  PCSPAISetSp(PC,PetscInt);

372: extern PetscErrorCode  PCHYPRESetType(PC,const char[]);
373: extern PetscErrorCode  PCHYPREGetType(PC,const char*[]);
374: extern PetscErrorCode  PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
375: extern PetscErrorCode  PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);

377: extern PetscErrorCode  PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*);
378: extern PetscErrorCode  PCFieldSplitSetType(PC,PCCompositeType);
379: extern PetscErrorCode  PCFieldSplitSetBlockSize(PC,PetscInt);
380: extern PetscErrorCode  PCFieldSplitSetIS(PC,const char[],IS);
381: extern PetscErrorCode  PCFieldSplitGetIS(PC,const char[],IS*);

383: /*E
384:     PCFieldSplitSchurPreType - Determines how to precondition Schur complement

386:     Level: intermediate

388: .seealso: PCFieldSplitSchurPrecondition()
389: E*/
390: typedef enum {PC_FIELDSPLIT_SCHUR_PRE_SELF,PC_FIELDSPLIT_SCHUR_PRE_DIAG,PC_FIELDSPLIT_SCHUR_PRE_USER} PCFieldSplitSchurPreType;
391: extern const char *const PCFieldSplitSchurPreTypes[];

393: extern PetscErrorCode  PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat);
394: extern PetscErrorCode  PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*);

396: extern PetscErrorCode  PCGalerkinSetRestriction(PC,Mat);
397: extern PetscErrorCode  PCGalerkinSetInterpolation(PC,Mat);

399: extern PetscErrorCode  PCSetCoordinates(PC,PetscInt,PetscReal*);
400: extern PetscErrorCode  PCSASetVectors(PC,PetscInt,PetscReal *);

402: extern PetscErrorCode  PCPythonSetType(PC,const char[]);

404: extern PetscErrorCode  PCSetDM(PC,DM);
405: extern PetscErrorCode  PCGetDM(PC,DM*);

407: extern PetscErrorCode  PCSetApplicationContext(PC,void*);
408: extern PetscErrorCode  PCGetApplicationContext(PC,void*);

410: extern PetscErrorCode  PCBiCGStabCUSPSetTolerance(PC,PetscReal);
411: extern PetscErrorCode  PCBiCGStabCUSPSetIterations(PC,PetscInt);
412: extern PetscErrorCode  PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool);

414: extern PetscErrorCode  PCAINVCUSPSetDropTolerance(PC,PetscReal);
415: extern PetscErrorCode  PCAINVCUSPUseScaling(PC,PetscBool);
416: extern PetscErrorCode  PCAINVCUSPSetNonzeros(PC,PetscInt);
417: extern PetscErrorCode  PCAINVCUSPSetLinParameter(PC,PetscInt);
418: /*E
419:     PCPARMSGlobalType - Determines the global preconditioner method in PARMS

421:     Level: intermediate

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

430:     Level: intermediate

432: .seealso: PCPARMSSetLocal()
433: E*/
434: typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType;
435: extern const char *PCPARMSLocalTypes[];

437: extern PetscErrorCode PCPARMSSetGlobal(PC pc,PCPARMSGlobalType type);
438: extern PetscErrorCode PCPARMSSetLocal(PC pc,PCPARMSLocalType type);
439: extern PetscErrorCode PCPARMSSetSolveTolerances(PC pc,PetscReal tol,PetscInt maxits);
440: extern PetscErrorCode PCPARMSSetSolveRestart(PC pc,PetscInt restart);
441: extern PetscErrorCode PCPARMSSetNonsymPerm(PC pc,PetscBool nonsym);
442: extern PetscErrorCode PCPARMSSetFill(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2);

444: extern PetscErrorCode PCGAMGSetProcEqLim(PC,PetscInt);
445: extern PetscErrorCode PCGAMGSetRepartitioning(PC,PetscBool);
446: extern PetscErrorCode PCGAMGSetSolverType(PC,char[],PetscInt);
447: extern PetscErrorCode PCGAMGSetThreshold(PC,PetscReal);
448: extern PetscErrorCode PCGAMGSetCoarseEqLim(PC,PetscInt);
449: extern PetscErrorCode PCGAMGSetNlevels(PC,PetscInt);
450: #define PCGAMGType char*
451: extern PetscErrorCode PCGAMGSetType( PC,const PCGAMGType );
452: extern PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n);
453: extern PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n);

455: #if defined(PETSC_HAVE_PCBDDC)
456: /* Enum defining how to treat the coarse problem */
457: typedef enum {SEQUENTIAL_BDDC,REPLICATED_BDDC,PARALLEL_BDDC,MULTILEVEL_BDDC} CoarseProblemType;
458: extern PetscErrorCode PCBDDCSetNeumannBoundaries(PC,IS);
459: extern PetscErrorCode PCBDDCGetNeumannBoundaries(PC,IS*);
460: extern PetscErrorCode PCBDDCSetCoarseProblemType(PC,CoarseProblemType);
461: extern PetscErrorCode PCBDDCSetDofsSplitting(PC,PetscInt,IS[]);
462: #endif

464: PETSC_EXTERN_CXX_END

466: #endif /* __PETSCPC_H */