Actual source code: petscksp.h
petsc-3.4.0 2013-05-13
1: /*
2: Defines the interface functions for the Krylov subspace accelerators.
3: */
4: #ifndef __PETSCKSP_H
6: #include <petscpc.h>
8: PETSC_EXTERN PetscErrorCode KSPInitializePackage(void);
10: /*S
11: KSP - Abstract PETSc object that manages all Krylov methods
13: Level: beginner
15: Concepts: Krylov methods
17: .seealso: KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP
18: S*/
19: typedef struct _p_KSP* KSP;
21: /*J
22: KSPType - String with the name of a PETSc Krylov method or the creation function
23: with an optional dynamic library name, for example
24: http://www.mcs.anl.gov/petsc/lib.a:mykspcreate()
26: Level: beginner
28: .seealso: KSPSetType(), KSP
29: J*/
30: typedef const char* KSPType;
31: #define KSPRICHARDSON "richardson"
32: #define KSPCHEBYSHEV "chebyshev"
33: #define KSPCG "cg"
34: #define KSPGROPPCG "groppcg"
35: #define KSPPIPECG "pipecg"
36: #define KSPCGNE "cgne"
37: #define KSPNASH "nash"
38: #define KSPSTCG "stcg"
39: #define KSPGLTR "gltr"
40: #define KSPGMRES "gmres"
41: #define KSPFGMRES "fgmres"
42: #define KSPLGMRES "lgmres"
43: #define KSPDGMRES "dgmres"
44: #define KSPPGMRES "pgmres"
45: #define KSPTCQMR "tcqmr"
46: #define KSPBCGS "bcgs"
47: #define KSPIBCGS "ibcgs"
48: #define KSPFBCGS "fbcgs"
49: #define KSPFBCGSR "fbcgsr"
50: #define KSPBCGSL "bcgsl"
51: #define KSPCGS "cgs"
52: #define KSPTFQMR "tfqmr"
53: #define KSPCR "cr"
54: #define KSPPIPECR "pipecr"
55: #define KSPLSQR "lsqr"
56: #define KSPPREONLY "preonly"
57: #define KSPQCG "qcg"
58: #define KSPBICG "bicg"
59: #define KSPMINRES "minres"
60: #define KSPSYMMLQ "symmlq"
61: #define KSPLCD "lcd"
62: #define KSPPYTHON "python"
63: #define KSPGCR "gcr"
64: #define KSPSPECEST "specest"
66: /* Logging support */
67: PETSC_EXTERN PetscClassId KSP_CLASSID;
68: PETSC_EXTERN PetscClassId DMKSP_CLASSID;
70: PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP *);
71: PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType);
72: PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
73: PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
74: PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec);
75: PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec);
76: PETSC_EXTERN PetscErrorCode KSPReset(KSP);
77: PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*);
79: PETSC_EXTERN PetscFunctionList KSPList;
80: PETSC_EXTERN PetscBool KSPRegisterAllCalled;
81: PETSC_EXTERN PetscErrorCode KSPRegisterAll(void);
82: PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP));
83: PETSC_EXTERN PetscErrorCode KSPMatRegisterAll(void);
85: PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType *);
86: PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide);
87: PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*);
88: PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
89: PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
90: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool );
91: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool *);
92: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool );
93: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool *);
94: PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool );
95: PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool *);
96: PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool *);
97: PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool );
98: PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool *);
99: PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool );
100: PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec *);
101: PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec *);
102: PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*);
103: PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*);
104: PETSC_EXTERN PetscErrorCode KSPSetNullSpace(KSP,MatNullSpace);
105: PETSC_EXTERN PetscErrorCode KSPGetNullSpace(KSP,MatNullSpace*);
106: PETSC_EXTERN PetscErrorCode KSPGetVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
108: PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
109: PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
111: PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC);
112: PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*);
114: PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal);
115: PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void**));
116: PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
117: PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **);
118: PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *);
119: PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool );
121: PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*);
122: PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec *);
123: PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
124: PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt);
126: PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*);
127: PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
128: PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
129: PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
130: PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
131: PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*);
132: PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*);
133: PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*);
134: PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*);
135: PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP *);
137: PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *);
138: PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *);
140: PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
141: PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool );
142: PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal);
143: PETSC_EXTERN PetscErrorCode KSPChebyshevSetEstimateEigenvalues(KSP,PetscReal,PetscReal,PetscReal,PetscReal);
144: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetRandom(KSP,PetscRandom);
145: PETSC_EXTERN PetscErrorCode KSPChebyshevSetNewMatrix(KSP);
146: PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
147: PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *);
148: PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*);
150: PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
151: PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*);
152: PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);
154: PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
155: PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
156: PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
157: PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
158: PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
160: PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
161: PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
163: PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
164: PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
165: PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
167: /*E
168: KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
170: Level: advanced
172: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
173: KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
175: E*/
176: typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
177: PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
178: /*MC
179: KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process
181: Level: advanced
183: Note: Possible unstable, but the fastest to compute
185: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
186: KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
187: KSPGMRESModifiedGramSchmidtOrthogonalization()
188: M*/
190: /*MC
191: KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
192: iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
193: poor orthogonality.
195: Level: advanced
197: Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
198: estimate the orthogonality but is more stable.
200: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
201: KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
202: KSPGMRESModifiedGramSchmidtOrthogonalization()
203: M*/
205: /*MC
206: KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
208: Level: advanced
210: Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
211: but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
213: You should only use this if you absolutely know that the iterative refinement is needed.
215: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
216: KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
217: KSPGMRESModifiedGramSchmidtOrthogonalization()
218: M*/
220: PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
221: PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);
223: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
224: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
225: PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
227: PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
228: PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
229: PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);
231: PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
232: PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool );
233: PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
234: PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool);
236: PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
237: PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
239: PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *);
240: PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *);
241: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *);
242: PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *);
243: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
244: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
245: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *);
246: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,void *);
247: PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *);
248: PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *);
249: PETSC_EXTERN PetscErrorCode KSPMonitorAMS(KSP,PetscInt,PetscReal,void*);
250: PETSC_EXTERN PetscErrorCode KSPMonitorAMSCreate(KSP,const char*,void**);
251: PETSC_EXTERN PetscErrorCode KSPMonitorAMSDestroy(void**);
252: PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *);
254: PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
255: PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
257: PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat,MatStructure);
258: PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*,MatStructure*);
259: PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *);
260: PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
261: PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
262: PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);
263: PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt);
264: PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*);
266: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool );
267: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *);
268: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool );
269: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *);
271: PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
272: PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);
274: #define KSP_FILE_CLASSID 1211223
276: PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec);
277: PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
279: PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
280: PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
282: /*E
283: KSPNormType - Norm that is passed in the Krylov convergence
284: test routines.
286: Level: advanced
288: Each solver only supports a subset of these and some may support different ones
289: depending on left or right preconditioning, see KSPSetPCSide()
291: Notes: this must match finclude/petscksp.h
293: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
294: KSPSetConvergenceTest(), KSPSetPCSide()
295: E*/
296: typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
297: #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
298: PETSC_EXTERN const char *const*const KSPNormTypes;
300: /*MC
301: KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
302: possibly save some computation but means the convergence test cannot
303: be based on a norm of a residual etc.
305: Level: advanced
307: Note: Some Krylov methods need to compute a residual norm and then this option is ignored
309: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
310: M*/
312: /*MC
313: KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the
314: convergence test routine.
316: Level: advanced
318: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
319: M*/
321: /*MC
322: KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
323: convergence test routine.
325: Level: advanced
327: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
328: M*/
330: /*MC
331: KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
332: convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS
334: Level: advanced
336: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
337: M*/
339: PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
340: PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
341: PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
342: PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
343: PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool);
345: /*E
346: KSPConvergedReason - reason a Krylov method was said to
347: have converged or diverged
349: Level: beginner
351: Notes: See KSPGetConvergedReason() for explanation of each value
353: Developer notes: this must match finclude/petscksp.h
355: The string versions of these are KSPConvergedReasons; if you change
356: any of the values here also change them that array of names.
358: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
359: E*/
360: typedef enum {/* converged */
361: KSP_CONVERGED_RTOL_NORMAL = 1,
362: KSP_CONVERGED_ATOL_NORMAL = 9,
363: KSP_CONVERGED_RTOL = 2,
364: KSP_CONVERGED_ATOL = 3,
365: KSP_CONVERGED_ITS = 4,
366: KSP_CONVERGED_CG_NEG_CURVE = 5,
367: KSP_CONVERGED_CG_CONSTRAINED = 6,
368: KSP_CONVERGED_STEP_LENGTH = 7,
369: KSP_CONVERGED_HAPPY_BREAKDOWN = 8,
370: /* diverged */
371: KSP_DIVERGED_NULL = -2,
372: KSP_DIVERGED_ITS = -3,
373: KSP_DIVERGED_DTOL = -4,
374: KSP_DIVERGED_BREAKDOWN = -5,
375: KSP_DIVERGED_BREAKDOWN_BICG = -6,
376: KSP_DIVERGED_NONSYMMETRIC = -7,
377: KSP_DIVERGED_INDEFINITE_PC = -8,
378: KSP_DIVERGED_NANORINF = -9,
379: KSP_DIVERGED_INDEFINITE_MAT = -10,
381: KSP_CONVERGED_ITERATING = 0} KSPConvergedReason;
382: PETSC_EXTERN const char *const*KSPConvergedReasons;
384: /*MC
385: KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)
387: Level: beginner
389: See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
390: for left preconditioning it is the 2-norm of the preconditioned residual, and the
391: 2-norm of the residual for right preconditioning
393: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
395: M*/
397: /*MC
398: KSP_CONVERGED_ATOL - norm(r) <= atol
400: Level: beginner
402: See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
403: for left preconditioning it is the 2-norm of the preconditioned residual, and the
404: 2-norm of the residual for right preconditioning
406: Level: beginner
408: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
410: M*/
412: /*MC
413: KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
415: Level: beginner
417: See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
418: for left preconditioning it is the 2-norm of the preconditioned residual, and the
419: 2-norm of the residual for right preconditioning
421: Level: beginner
423: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
425: M*/
427: /*MC
428: KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
429: reached
431: Level: beginner
433: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
435: M*/
437: /*MC
438: KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
439: the preconditioner is applied. Also used when the KSPSkipConverged() convergence
440: test routine is set in KSP.
443: Level: beginner
446: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
448: M*/
450: /*MC
451: KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
452: method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or
453: preconditioner.
455: Level: beginner
457: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
459: M*/
461: /*MC
462: KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
463: method could not continue to enlarge the Krylov space.
466: Level: beginner
469: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
471: M*/
473: /*MC
474: KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
475: symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
477: Level: beginner
479: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
481: M*/
483: /*MC
484: KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
485: positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
486: be positive definite
488: Level: beginner
490: Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
491: the PCICC preconditioner to generate a positive definite preconditioner
493: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
495: M*/
497: /*MC
498: KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
499: while the KSPSolve() is still running.
501: Level: beginner
503: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
505: M*/
507: PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*));
508: PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **);
509: PETSC_EXTERN PetscErrorCode KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
510: PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
511: PETSC_EXTERN PetscErrorCode KSPDefaultConvergedDestroy(void *);
512: PETSC_EXTERN PetscErrorCode KSPDefaultConvergedCreate(void **);
513: PETSC_EXTERN PetscErrorCode KSPDefaultConvergedSetUIRNorm(KSP);
514: PETSC_EXTERN PetscErrorCode KSPDefaultConvergedSetUMIRNorm(KSP);
515: PETSC_EXTERN PetscErrorCode KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
516: PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *);
518: PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *);
520: /*E
521: KSPCGType - Determines what type of CG to use
523: Level: beginner
525: .seealso: KSPCGSetType()
526: E*/
527: typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
528: PETSC_EXTERN const char *const KSPCGTypes[];
530: PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
531: PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool );
533: PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal);
534: PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *);
535: PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *);
537: PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal);
538: PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *);
539: PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *);
541: PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal);
542: PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *);
543: PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *);
544: PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *);
545: PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *);
547: PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]);
549: PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
550: PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);
552: #include <petscdrawtypes.h>
553: PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
554: PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*);
555: PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormDestroy(PetscDrawLG*);
556: PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
557: PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
558: PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG*);
559: PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
561: PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
562: PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
564: /* see src/ksp/ksp/interface/iguess.c */
565: typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool monitor;Mat mat; KSP ksp;}* KSPFischerGuess;
567: PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*);
568: PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*);
569: PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess);
570: PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec);
571: PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec);
572: PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess);
574: PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
575: PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess);
576: PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*);
578: PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
579: PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
580: PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
581: PETSC_EXTERN PetscErrorCode MatSchurComplementSet(Mat,Mat,Mat,Mat,Mat,Mat);
582: PETSC_EXTERN PetscErrorCode MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,Mat,MatStructure);
583: PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
584: PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatReuse,Mat *);
586: PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
587: PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool );
588: PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
589: PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
590: PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
591: PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *);
592: PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*);
593: PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
594: PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*);
595: PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,MatStructure*,void*),void*);
596: PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
597: PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
598: PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
599: PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
601: #endif