Actual source code: petscksp.h
petsc-3.6.4 2016-04-12
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. This is the object that manages the
12: linear solves in PETSc (even those such as direct solvers that do no use Krylov accelerators).
14: Level: beginner
16: Concepts: Krylov methods
18: Notes: When a direct solver is used but no Krylov solver is used the KSP object is still used by with a
19: KSPType of KSPPREONLY (meaning application of the preconditioner is only used as the linear solver).
21: .seealso: KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP, KSPDestroy()
22: S*/
23: typedef struct _p_KSP* KSP;
25: /*J
26: KSPType - String with the name of a PETSc Krylov method.
28: Level: beginner
30: .seealso: KSPSetType(), KSP, KSPRegister(), KSPCreate(), KSPSetFromOptions()
31: J*/
32: typedef const char* KSPType;
33: #define KSPRICHARDSON "richardson"
34: #define KSPCHEBYSHEV "chebyshev"
35: #define KSPCG "cg"
36: #define KSPGROPPCG "groppcg"
37: #define KSPPIPECG "pipecg"
38: #define KSPCGNE "cgne"
39: #define KSPNASH "nash"
40: #define KSPSTCG "stcg"
41: #define KSPGLTR "gltr"
42: #define KSPFCG "fcg"
43: #define KSPGMRES "gmres"
44: #define KSPFGMRES "fgmres"
45: #define KSPLGMRES "lgmres"
46: #define KSPDGMRES "dgmres"
47: #define KSPPGMRES "pgmres"
48: #define KSPTCQMR "tcqmr"
49: #define KSPBCGS "bcgs"
50: #define KSPIBCGS "ibcgs"
51: #define KSPFBCGS "fbcgs"
52: #define KSPFBCGSR "fbcgsr"
53: #define KSPBCGSL "bcgsl"
54: #define KSPCGS "cgs"
55: #define KSPTFQMR "tfqmr"
56: #define KSPCR "cr"
57: #define KSPPIPECR "pipecr"
58: #define KSPLSQR "lsqr"
59: #define KSPPREONLY "preonly"
60: #define KSPQCG "qcg"
61: #define KSPBICG "bicg"
62: #define KSPMINRES "minres"
63: #define KSPSYMMLQ "symmlq"
64: #define KSPLCD "lcd"
65: #define KSPPYTHON "python"
66: #define KSPGCR "gcr"
68: /* Logging support */
69: PETSC_EXTERN PetscClassId KSP_CLASSID;
70: PETSC_EXTERN PetscClassId DMKSP_CLASSID;
72: PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP *);
73: PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType);
74: PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType *);
75: PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
76: PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
77: PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec);
78: PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec);
79: PETSC_EXTERN PetscErrorCode KSPReset(KSP);
80: PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*);
81: PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool);
82: PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool);
84: PETSC_EXTERN PetscFunctionList KSPList;
85: PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP));
87: PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide);
88: PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*);
89: PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
90: PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
91: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool );
92: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool *);
93: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool );
94: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool *);
95: PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool );
96: PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool *);
97: PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool );
98: PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool *);
99: PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool );
100: PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool *);
101: PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec *);
102: PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec *);
103: PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*);
104: PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*);
105: PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*);
106: PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
107: PETSC_DEPRECATED("Use KSPCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode KSPGetVecs(KSP ksp,PetscInt n,Vec **x,PetscInt m,Vec **y) {return KSPCreateVecs(ksp,n,x,m,y);}
109: PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
110: PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
112: PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC);
113: PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*);
115: PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal);
116: PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void**));
117: PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
118: PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **);
119: PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *);
120: PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool );
122: PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*);
123: PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec *);
124: PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
125: PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt);
127: PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*);
128: PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
129: PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
130: PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
131: PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
132: PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*);
133: PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*);
134: PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*);
135: PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*);
136: PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP *);
138: PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *);
139: PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *);
141: PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
142: PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool );
143: PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal);
144: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal);
145: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetRandom(KSP,PetscRandom);
146: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*);
147: PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
148: PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt *);
149: PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]);
151: /*E
153: KSPFCGTruncationType - Define how stored directions are used to orthogonalize in FCG
155: KSP_FCG_TRUNC_TYPE_STANDARD uses all (up to mmax) stored directions
156: KSP_FCG_TRUNC_TYPE_NOTAY uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..
158: Level: intermediate
160: .seealso : KSPFCG,KSPFCGSetTruncationType(),KSPFCGGetTruncationType()
162: E*/
163: typedef enum {KSP_FCG_TRUNC_TYPE_STANDARD,KSP_FCG_TRUNC_TYPE_NOTAY} KSPFCGTruncationType;
164: PETSC_EXTERN const char *const KSPFCGTruncationTypes[];
166: PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt);
167: PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*);
168: PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt);
169: PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*);
170: PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCGTruncationType);
171: PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCGTruncationType*);
173: PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
174: PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*);
175: PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);
177: PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
178: PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
179: PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
180: PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
181: PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
183: PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
184: PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
186: PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
187: PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
188: PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
190: /*E
191: KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
193: Level: advanced
195: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
196: KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
198: E*/
199: typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
200: PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
201: /*MC
202: KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
204: Level: advanced
206: Note: Possible unstable, but the fastest to compute
208: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
209: KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
210: KSPGMRESModifiedGramSchmidtOrthogonalization()
211: M*/
213: /*MC
214: KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
215: iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
216: poor orthogonality.
218: Level: advanced
220: Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
221: estimate the orthogonality but is more stable.
223: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
224: KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
225: KSPGMRESModifiedGramSchmidtOrthogonalization()
226: M*/
228: /*MC
229: KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
231: Level: advanced
233: Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
234: but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
236: You should only use this if you absolutely know that the iterative refinement is needed.
238: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
239: KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
240: KSPGMRESModifiedGramSchmidtOrthogonalization()
241: M*/
243: PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
244: PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);
246: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
247: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
248: PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
250: PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
251: PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
252: PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);
254: PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
255: PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool );
256: PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
257: PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool);
259: PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
260: PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
262: PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *);
263: PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *);
264: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *);
265: PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *);
266: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
267: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
268: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *);
269: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,void *);
270: PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *);
271: PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *);
272: PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*);
273: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**);
274: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**);
275: PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *);
277: PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
278: PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
280: PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat);
281: PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*);
282: PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *);
283: PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
284: PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
285: PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);
286: PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt);
287: PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*);
289: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool );
290: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *);
291: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool );
292: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *);
294: PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
295: PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);
296: PETSC_STATIC_INLINE PetscErrorCode KSPViewFromOptions(KSP A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
297: PETSC_EXTERN PetscErrorCode KSPReasonView(KSP,PetscViewer);
298: PETSC_EXTERN PetscErrorCode KSPReasonViewFromOptions(KSP);
300: #define KSP_FILE_CLASSID 1211223
302: PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec);
303: PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
305: PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
306: PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
308: /*E
309: KSPNormType - Norm that is passed in the Krylov convergence
310: test routines.
312: Level: advanced
314: Each solver only supports a subset of these and some may support different ones
315: depending on left or right preconditioning, see KSPSetPCSide()
317: Notes: this must match petsc/finclude/petscksp.h
319: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
320: KSPSetConvergenceTest(), KSPSetPCSide()
321: E*/
322: typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
323: #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
324: PETSC_EXTERN const char *const*const KSPNormTypes;
326: /*MC
327: KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
328: possibly save some computation but means the convergence test cannot
329: be based on a norm of a residual etc.
331: Level: advanced
333: Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored
335: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
336: M*/
338: /*MC
339: KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
340: convergence test routine.
342: Level: advanced
344: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
345: M*/
347: /*MC
348: KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
349: convergence test routine.
351: Level: advanced
353: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
354: M*/
356: /*MC
357: KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
358: convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG
360: Level: advanced
362: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
363: M*/
365: PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
366: PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
367: PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
368: PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
369: PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool);
371: /*E
372: KSPConvergedReason - reason a Krylov method was said to have converged or diverged
374: Level: beginner
376: Notes: See KSPGetConvergedReason() for explanation of each value
378: Developer notes: this must match petsc/finclude/petscksp.h
380: The string versions of these are KSPConvergedReasons; if you change
381: any of the values here also change them that array of names.
383: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
384: E*/
385: typedef enum {/* converged */
386: KSP_CONVERGED_RTOL_NORMAL = 1,
387: KSP_CONVERGED_ATOL_NORMAL = 9,
388: KSP_CONVERGED_RTOL = 2,
389: KSP_CONVERGED_ATOL = 3,
390: KSP_CONVERGED_ITS = 4,
391: KSP_CONVERGED_CG_NEG_CURVE = 5,
392: KSP_CONVERGED_CG_CONSTRAINED = 6,
393: KSP_CONVERGED_STEP_LENGTH = 7,
394: KSP_CONVERGED_HAPPY_BREAKDOWN = 8,
395: /* diverged */
396: KSP_DIVERGED_NULL = -2,
397: KSP_DIVERGED_ITS = -3,
398: KSP_DIVERGED_DTOL = -4,
399: KSP_DIVERGED_BREAKDOWN = -5,
400: KSP_DIVERGED_BREAKDOWN_BICG = -6,
401: KSP_DIVERGED_NONSYMMETRIC = -7,
402: KSP_DIVERGED_INDEFINITE_PC = -8,
403: KSP_DIVERGED_NANORINF = -9,
404: KSP_DIVERGED_INDEFINITE_MAT = -10,
405: KSP_DIVERGED_PCSETUP_FAILED = -11,
407: KSP_CONVERGED_ITERATING = 0} KSPConvergedReason;
408: PETSC_EXTERN const char *const*KSPConvergedReasons;
410: /*MC
411: KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)
413: Level: beginner
415: See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
416: for left preconditioning it is the 2-norm of the preconditioned residual, and the
417: 2-norm of the residual for right preconditioning
419: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
421: M*/
423: /*MC
424: KSP_CONVERGED_ATOL - norm(r) <= atol
426: Level: beginner
428: See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
429: for left preconditioning it is the 2-norm of the preconditioned residual, and the
430: 2-norm of the residual for right preconditioning
432: Level: beginner
434: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
436: M*/
438: /*MC
439: KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
441: Level: beginner
443: See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
444: for left preconditioning it is the 2-norm of the preconditioned residual, and the
445: 2-norm of the residual for right preconditioning
447: Level: beginner
449: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
451: M*/
453: /*MC
454: KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
455: reached
457: Level: beginner
459: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
461: M*/
463: /*MC
464: KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
465: the preconditioner is applied. Also used when the KSPConvergedSkip() convergence
466: test routine is set in KSP.
468: Level: beginner
470: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
472: M*/
474: /*MC
475: KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
476: method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or
477: preconditioner.
479: Level: beginner
481: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
483: M*/
485: /*MC
486: KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
487: method could not continue to enlarge the Krylov space.
489: Level: beginner
491: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
493: M*/
495: /*MC
496: KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
497: symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
499: Level: beginner
501: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
503: M*/
505: /*MC
506: KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
507: positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
508: be positive definite
510: Level: beginner
512: Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
513: the PCICC preconditioner to generate a positive definite preconditioner
515: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
517: M*/
519: /*MC
520: KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
521: while the KSPSolve() is still running.
523: Level: beginner
525: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
527: M*/
529: PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*));
530: PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **);
531: PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
532: PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
533: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *);
534: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **);
535: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
536: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
537: PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
538: PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *);
540: PETSC_DEPRECATED("Use KSPConvergedDefault()") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ }
541: #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
542: PETSC_DEPRECATED("Use KSPConvergedDefaultDestroy()") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ }
543: #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
544: PETSC_DEPRECATED("Use KSPConvergedDefaultCreate()") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ }
545: #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
546: PETSC_DEPRECATED("Use KSPConvergedDefaultSetUIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ }
547: #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
548: PETSC_DEPRECATED("Use KSPConvergedDefaultSetUMIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ }
549: #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
550: PETSC_DEPRECATED("Use KSPConvergedSkip()") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ }
551: #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
553: PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *);
555: /*E
556: KSPCGType - Determines what type of CG to use
558: Level: beginner
560: .seealso: KSPCGSetType()
561: E*/
562: typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
563: PETSC_EXTERN const char *const KSPCGTypes[];
565: PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
566: PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool );
568: PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal);
569: PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *);
570: PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *);
572: PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal);
573: PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *);
574: PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *);
576: PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal);
577: PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *);
578: PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *);
579: PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *);
580: PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *);
582: PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]);
584: PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
585: PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);
587: #include <petscdrawtypes.h>
588: PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(const char[],const char[],int,int,int,int,PetscObject**);
589: PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,PetscObject*);
590: PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormDestroy(PetscObject**);
591: PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(const char[],const char[],int,int,int,int,PetscObject**);
592: PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,PetscObject*);
593: PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormDestroy(PetscObject**);
594: PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
596: PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
597: PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
599: /* see src/ksp/ksp/interface/iguess.c */
600: typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool monitor;Mat mat; KSP ksp;}* KSPFischerGuess;
602: PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*);
603: PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*);
604: PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess);
605: PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec);
606: PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec);
607: PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess);
609: PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
610: PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess);
611: PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*);
613: /*E
614: MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines
616: Level: intermediate
618: .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat()
619: E*/
620: typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP} MatSchurComplementAinvType;
621: PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
623: PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
624: PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
625: PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
626: PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
627: PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
628: PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
629: PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType);
630: PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*);
631: PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*);
632: PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*);
633: PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatSchurComplementAinvType,MatReuse,Mat *);
634: PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*);
636: PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
637: PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool );
638: PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
639: PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
640: PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
641: PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *);
642: PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
643: PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
644: PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
645: PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*);
646: PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
647: PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
648: PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
649: PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
651: PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec);
652: PETSC_EXTERN PetscErrorCode DMPlexProjectField(DM, Vec, void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal [], PetscScalar []), InsertMode, Vec);
653: #endif