Actual source code: petscpc.h

  1: /*
  2:       Preconditioner module. 
  3: */
 6:  #include petscmat.h

  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: */
 16: #define PCType char*

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

 21:    Level: beginner

 23:   Concepts: preconditioners

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

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

 34:    Level: beginner

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

 38: .seealso: PCSetType(), PC, PCCreate()
 39: E*/
 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 PCKSP             "ksp"
 52: #define PCCOMPOSITE       "composite"
 53: #define PCREDUNDANT       "redundant"
 54: #define PCSPAI            "spai"
 55: #define PCNN              "nn"
 56: #define PCCHOLESKY        "cholesky"
 57: #define PCSAMG            "samg"
 58: #define PCPBJACOBI        "pbjacobi"
 59: #define PCMAT             "mat"
 60: #define PCHYPRE           "hypre"
 61: #define PCFIELDSPLIT      "fieldsplit"
 62: #define PCTFS             "tfs"

 64: /* Logging support */

 69: /*E
 70:     PCSide - If the preconditioner is to be applied to the left, right
 71:      or symmetrically around the operator.

 73:    Level: beginner

 75: .seealso: 
 76: E*/
 77: typedef enum { PC_LEFT,PC_RIGHT,PC_SYMMETRIC } PCSide;

 79: EXTERN PetscErrorCode PCCreate(MPI_Comm,PC*);
 80: EXTERN PetscErrorCode PCSetType(PC,const PCType);
 81: EXTERN PetscErrorCode PCSetUp(PC);
 82: EXTERN PetscErrorCode PCSetUpOnBlocks(PC);
 83: EXTERN PetscErrorCode PCApply(PC,Vec,Vec);
 84: EXTERN PetscErrorCode PCApplySymmetricLeft(PC,Vec,Vec);
 85: EXTERN PetscErrorCode PCApplySymmetricRight(PC,Vec,Vec);
 86: EXTERN PetscErrorCode PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
 87: EXTERN PetscErrorCode PCApplyTranspose(PC,Vec,Vec);
 88: EXTERN PetscErrorCode PCHasApplyTranspose(PC,PetscTruth*);
 89: EXTERN PetscErrorCode PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
 90: EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt);
 91: EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscTruth*);

 93: EXTERN PetscErrorCode        PCRegisterDestroy(void);
 94: EXTERN PetscErrorCode        PCRegisterAll(const char[]);

 97: EXTERN PetscErrorCode PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC));

 99: /*MC
100:    PCRegisterDynamic - Adds a method to the preconditioner package.

102:    Synopsis:
103:    PetscErrorCode PCRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(PC))

105:    Not collective

107:    Input Parameters:
108: +  name_solver - name of a new user-defined solver
109: .  path - path (either absolute or relative) the library containing this solver
110: .  name_create - name of routine to create method context
111: -  routine_create - routine to create method context

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

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

119:    Sample usage:
120: .vb
121:    PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib",
122:               "MySolverCreate",MySolverCreate);
123: .ve

125:    Then, your solver can be chosen with the procedural interface via
126: $     PCSetType(pc,"my_solver")
127:    or at runtime via the option
128: $     -pc_type my_solver

130:    Level: advanced

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

136: .keywords: PC, register

138: .seealso: PCRegisterAll(), PCRegisterDestroy()
139: M*/
140: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
141: #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0)
142: #else
143: #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d)
144: #endif

146: EXTERN PetscErrorCode PCDestroy(PC);
147: EXTERN PetscErrorCode PCSetFromOptions(PC);
148: EXTERN PetscErrorCode PCGetType(PC,PCType*);

150: EXTERN PetscErrorCode PCGetFactoredMatrix(PC,Mat*);
151: EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
152: EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);

154: EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat,MatStructure);
155: EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*,MatStructure*);

157: EXTERN PetscErrorCode PCView(PC,PetscViewer);

159: EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]);
160: EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]);
161: EXTERN PetscErrorCode PCGetOptionsPrefix(PC,char*[]);

163: EXTERN PetscErrorCode PCComputeExplicitOperator(PC,Mat*);

165: /*
166:       These are used to provide extra scaling of preconditioned 
167:    operator for time-stepping schemes like in PVODE 
168: */
169: EXTERN PetscErrorCode PCDiagonalScale(PC,PetscTruth*);
170: EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec);
171: EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec);
172: EXTERN PetscErrorCode PCDiagonalScaleSet(PC,Vec);

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

176: EXTERN PetscErrorCode PCJacobiSetUseRowMax(PC);
177: EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType);
178: EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal);
179: EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt);

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

184: #define USE_PRECONDITIONER_MATRIX 0
185: #define USE_TRUE_MATRIX           1
186: EXTERN PetscErrorCode PCBJacobiSetUseTrueLocal(PC);
187: EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
188: EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);

190: EXTERN PetscErrorCode PCKSPSetUseTrue(PC);

192: EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(void*,Vec,Vec),void*);
193: EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(void*,Vec,Vec));
194: EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(void*));
195: EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(void*,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt),void*);
196: EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(void*,PetscViewer));
197: EXTERN PetscErrorCode PCShellSetName(PC,const char[]);
198: EXTERN PetscErrorCode PCShellGetName(PC,char*[]);

200: EXTERN PetscErrorCode PCLUSetMatOrdering(PC,MatOrderingType);
201: EXTERN PetscErrorCode PCLUSetReuseOrdering(PC,PetscTruth);
202: EXTERN PetscErrorCode PCLUSetReuseFill(PC,PetscTruth);
203: EXTERN PetscErrorCode PCLUSetUseInPlace(PC);
204: EXTERN PetscErrorCode PCLUSetFill(PC,PetscReal);
205: EXTERN PetscErrorCode PCLUSetDamping(PC,PetscReal);
206: EXTERN PetscErrorCode PCLUSetShift(PC,PetscTruth);
207: EXTERN PetscErrorCode PCLUSetPivoting(PC,PetscReal);
208: EXTERN PetscErrorCode PCLUSetPivotInBlocks(PC,PetscTruth);
209: EXTERN PetscErrorCode PCLUSetZeroPivot(PC,PetscReal);

211: EXTERN PetscErrorCode PCCholeskySetMatOrdering(PC,MatOrderingType);
212: EXTERN PetscErrorCode PCCholeskySetReuseOrdering(PC,PetscTruth);
213: EXTERN PetscErrorCode PCCholeskySetReuseFill(PC,PetscTruth);
214: EXTERN PetscErrorCode PCCholeskySetUseInPlace(PC);
215: EXTERN PetscErrorCode PCCholeskySetFill(PC,PetscReal);
216: EXTERN PetscErrorCode PCCholeskySetDamping(PC,PetscReal);
217: EXTERN PetscErrorCode PCCholeskySetShift(PC,PetscTruth);
218: EXTERN PetscErrorCode PCCholeskySetPivotInBlocks(PC,PetscTruth);

220: EXTERN PetscErrorCode PCILUSetMatOrdering(PC,MatOrderingType);
221: EXTERN PetscErrorCode PCILUSetUseInPlace(PC);
222: EXTERN PetscErrorCode PCILUSetFill(PC,PetscReal);
223: EXTERN PetscErrorCode PCILUSetLevels(PC,PetscInt);
224: EXTERN PetscErrorCode PCILUSetReuseOrdering(PC,PetscTruth);
225: EXTERN PetscErrorCode PCILUSetUseDropTolerance(PC,PetscReal,PetscReal,PetscInt);
226: EXTERN PetscErrorCode PCILUDTSetReuseFill(PC,PetscTruth);
227: EXTERN PetscErrorCode PCILUSetAllowDiagonalFill(PC);
228: EXTERN PetscErrorCode PCILUSetDamping(PC,PetscReal);
229: EXTERN PetscErrorCode PCILUSetShift(PC,PetscTruth);
230: EXTERN PetscErrorCode PCILUSetPivotInBlocks(PC,PetscTruth);
231: EXTERN PetscErrorCode PCILUSetZeroPivot(PC,PetscReal);

233: EXTERN PetscErrorCode PCICCSetMatOrdering(PC,MatOrderingType);
234: EXTERN PetscErrorCode PCICCSetFill(PC,PetscReal);
235: EXTERN PetscErrorCode PCICCSetLevels(PC,PetscInt);
236: EXTERN PetscErrorCode PCICCSetDamping(PC,PetscReal);
237: EXTERN PetscErrorCode PCICCSetShift(PC,PetscTruth);
238: EXTERN PetscErrorCode PCICCSetPivotInBlocks(PC,PetscTruth);
239: EXTERN PetscErrorCode PCICCSetZeroPivot(PC,PetscReal);

241: EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[]);
242: EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[]);
243: EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt);
244: /*E
245:     PCASMType - Type of additive Schwarz method to use

247: $  PC_ASM_BASIC - symmetric version where residuals from the ghost points are used
248: $                 and computed values in ghost regions are added together. Classical
249: $                 standard additive Schwarz
250: $  PC_ASM_RESTRICT - residuals from ghost points are used but computed values in ghost
251: $                    region are discarded. Default
252: $  PC_ASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost
253: $                       region are added back in
254: $  PC_ASM_NONE - ghost point residuals are not used, computed ghost values are discarded
255: $                not very good.                

257:    Level: beginner

259: .seealso: PCASMSetType()
260: E*/
261: typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
262: EXTERN PetscErrorCode PCASMSetType(PC,PCASMType);
263: EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt *,IS **);
264: EXTERN PetscErrorCode PCASMSetUseInPlace(PC);
265: EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[]);
266: EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);

268: /*E
269:     PCCompositeType - Determines how two or more preconditioner are composed

271: $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
272: $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
273: $                                computed after the previous preconditioner application
274: $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
275: $                         where first preconditioner is built from alpha I + S and second from
276: $                         alpha I + R

278:    Level: beginner

280: .seealso: PCCompositeSetType()
281: E*/
282: typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL} PCCompositeType;
283: EXTERN PetscErrorCode PCCompositeSetUseTrue(PC);
284: EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType);
285: EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType);
286: EXTERN PetscErrorCode PCCompositeGetPC(PC pc,PetscInt n,PC *);
287: EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar);

289: EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter);
290: EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*);
291: EXTERN PetscErrorCode PCRedundantGetPC(PC,PC*);
292: EXTERN PetscErrorCode MatGetOrderingList(PetscFList *list);

294: EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double);
295: EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt);
296: EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt);
297: EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt);
298: EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt);
299: EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt);
300: EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt);
301: EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt);

303: EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]);
304: EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
305: EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);

307: EXTERN PetscErrorCode PCFieldSplitSetFields(PC,PetscInt,PetscInt*);
308: EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType);

311: #endif /* __PETSCPC_H */