Actual source code: petscksp.h
petsc-3.7.3 2016-08-01
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 KSPPIPECGRR "pipecgrr"
39: #define KSPCGNE "cgne"
40: #define KSPNASH "nash"
41: #define KSPSTCG "stcg"
42: #define KSPGLTR "gltr"
43: #define KSPFCG "fcg"
44: #define KSPPIPEFCG "pipefcg"
45: #define KSPGMRES "gmres"
46: #define KSPPIPEFGMRES "pipefgmres"
47: #define KSPFGMRES "fgmres"
48: #define KSPLGMRES "lgmres"
49: #define KSPDGMRES "dgmres"
50: #define KSPPGMRES "pgmres"
51: #define KSPTCQMR "tcqmr"
52: #define KSPBCGS "bcgs"
53: #define KSPIBCGS "ibcgs"
54: #define KSPFBCGS "fbcgs"
55: #define KSPFBCGSR "fbcgsr"
56: #define KSPBCGSL "bcgsl"
57: #define KSPCGS "cgs"
58: #define KSPTFQMR "tfqmr"
59: #define KSPCR "cr"
60: #define KSPPIPECR "pipecr"
61: #define KSPLSQR "lsqr"
62: #define KSPPREONLY "preonly"
63: #define KSPQCG "qcg"
64: #define KSPBICG "bicg"
65: #define KSPMINRES "minres"
66: #define KSPSYMMLQ "symmlq"
67: #define KSPLCD "lcd"
68: #define KSPPYTHON "python"
69: #define KSPGCR "gcr"
70: #define KSPPIPEGCR "pipegcr"
71: #define KSPTSIRM "tsirm"
72: #define KSPCGLS "cgls"
74: /* Logging support */
75: PETSC_EXTERN PetscClassId KSP_CLASSID;
76: PETSC_EXTERN PetscClassId DMKSP_CLASSID;
78: PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP *);
79: PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType);
80: PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType *);
81: PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
82: PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
83: PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec);
84: PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec);
85: PETSC_EXTERN PetscErrorCode KSPReset(KSP);
86: PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*);
87: PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool);
88: PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool);
90: PETSC_EXTERN PetscFunctionList KSPList;
91: PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP));
93: PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide);
94: PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*);
95: PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
96: PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
97: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool );
98: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool *);
99: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool );
100: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool *);
101: PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool );
102: PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool *);
103: PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool );
104: PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool );
105: PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool *);
106: PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool );
107: PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool *);
108: PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec *);
109: PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec *);
110: PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*);
111: PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*);
112: PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*);
113: PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
114: 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);}
116: PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
117: PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
119: PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC);
120: PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*);
122: PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal);
123: PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void**));
124: PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
125: PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **);
126: PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *);
127: PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool );
129: PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*);
130: PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec *);
131: PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
132: PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt);
134: PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*);
135: PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
136: PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
137: PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
138: PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
139: PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*);
140: PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*);
141: PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*);
142: PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*);
143: PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP *);
145: PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *);
146: PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *);
148: PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
149: PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool );
150: PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal);
151: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal);
152: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseRandom(KSP,PetscBool);
153: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetRandom(KSP,PetscRandom);
154: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*);
155: PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
156: PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt *);
157: PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]);
158: PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt *,Vec[],PetscReal[],PetscReal[]);
160: /*E
162: KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods
164: KSP_FCD_TRUNC_TYPE_STANDARD uses all (up to mmax) stored directions
165: KSP_FCD_TRUNC_TYPE_NOTAY uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..
167: Level: intermediate
168: .seealso : KSPFCG,KSPPIPEFCG,KSPPIPEGCR,KSPFCGSetTruncationType(),KSPFCGGetTruncationType()
170: E*/
171: typedef enum {KSP_FCD_TRUNC_TYPE_STANDARD,KSP_FCD_TRUNC_TYPE_NOTAY} KSPFCDTruncationType;
172: PETSC_EXTERN const char *const KSPFCDTruncationTypes[];
174: PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt);
175: PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*);
176: PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt);
177: PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*);
178: PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType);
179: PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*);
181: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt);
182: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*);
183: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt);
184: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*);
185: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType);
186: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*);
188: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt);
189: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*);
190: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt);
191: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*);
192: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType);
193: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*);
194: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool);
195: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*);
197: PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
198: PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*);
199: PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);
201: PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
202: PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
203: PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
204: PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
205: PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
207: PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
208: PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
210: PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP,PetscScalar);
212: PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
213: PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
214: PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
216: /*E
217: KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
219: Level: advanced
221: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
222: KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
224: E*/
225: typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
226: PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
227: /*MC
228: KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
230: Level: advanced
232: Note: Possible unstable, but the fastest to compute
234: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
235: KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
236: KSPGMRESModifiedGramSchmidtOrthogonalization()
237: M*/
239: /*MC
240: KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
241: iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
242: poor orthogonality.
244: Level: advanced
246: Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
247: estimate the orthogonality but is more stable.
249: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
250: KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
251: KSPGMRESModifiedGramSchmidtOrthogonalization()
252: M*/
254: /*MC
255: KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
257: Level: advanced
259: Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
260: but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
262: You should only use this if you absolutely know that the iterative refinement is needed.
264: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
265: KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
266: KSPGMRESModifiedGramSchmidtOrthogonalization()
267: M*/
269: PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
270: PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);
272: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
273: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
274: PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
276: PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
277: PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
278: PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);
280: PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
281: PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool );
282: PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
283: PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool);
285: PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
286: PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
288: PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,PetscViewerAndFormat *);
289: PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat *);
290: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat *);
291: PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,PetscViewerAndFormat *);
292: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
293: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
294: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat *);
295: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat *);
296: PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,PetscViewerAndFormat *);
297: PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,PetscViewerAndFormat *);
298: PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*);
299: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**);
300: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**);
301: PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *);
302: PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP,const char[],const char[],const char [],PetscErrorCode (*)(KSP,PetscInt,PetscReal,PetscViewerAndFormat*));
304: PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
305: PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
307: PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat);
308: PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*);
309: PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *);
310: PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
311: PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
312: PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);
313: PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt);
314: PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*);
316: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool );
317: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *);
318: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool );
319: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *);
321: PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
322: PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);
323: PETSC_STATIC_INLINE PetscErrorCode KSPViewFromOptions(KSP A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
324: PETSC_EXTERN PetscErrorCode KSPReasonView(KSP,PetscViewer);
325: PETSC_EXTERN PetscErrorCode KSPReasonViewFromOptions(KSP);
327: #define KSP_FILE_CLASSID 1211223
329: PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec);
330: PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
332: PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
333: PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
334: PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*);
336: /*E
337: KSPNormType - Norm that is passed in the Krylov convergence
338: test routines.
340: Level: advanced
342: Each solver only supports a subset of these and some may support different ones
343: depending on left or right preconditioning, see KSPSetPCSide()
345: Notes: this must match petsc/finclude/petscksp.h
347: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
348: KSPSetConvergenceTest(), KSPSetPCSide()
349: E*/
350: typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
351: #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
352: PETSC_EXTERN const char *const*const KSPNormTypes;
354: /*MC
355: KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
356: possibly save some computation but means the convergence test cannot
357: be based on a norm of a residual etc.
359: Level: advanced
361: Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored
363: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
364: M*/
366: /*MC
367: KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
368: convergence test routine.
370: Level: advanced
372: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
373: M*/
375: /*MC
376: KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
377: convergence test routine.
379: Level: advanced
381: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
382: M*/
384: /*MC
385: KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
386: convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG, KSPPIPEFCG, KSPPIPEGCR
388: Level: advanced
390: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
391: M*/
393: PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
394: PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
395: PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
396: PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
397: PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool);
399: /*E
400: KSPConvergedReason - reason a Krylov method was said to have converged or diverged
402: Level: beginner
404: Notes: See KSPGetConvergedReason() for explanation of each value
406: Developer notes: this must match petsc/finclude/petscksp.h
408: The string versions of these are KSPConvergedReasons; if you change
409: any of the values here also change them that array of names.
411: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
412: E*/
413: typedef enum {/* converged */
414: KSP_CONVERGED_RTOL_NORMAL = 1,
415: KSP_CONVERGED_ATOL_NORMAL = 9,
416: KSP_CONVERGED_RTOL = 2,
417: KSP_CONVERGED_ATOL = 3,
418: KSP_CONVERGED_ITS = 4,
419: KSP_CONVERGED_CG_NEG_CURVE = 5,
420: KSP_CONVERGED_CG_CONSTRAINED = 6,
421: KSP_CONVERGED_STEP_LENGTH = 7,
422: KSP_CONVERGED_HAPPY_BREAKDOWN = 8,
423: /* diverged */
424: KSP_DIVERGED_NULL = -2,
425: KSP_DIVERGED_ITS = -3,
426: KSP_DIVERGED_DTOL = -4,
427: KSP_DIVERGED_BREAKDOWN = -5,
428: KSP_DIVERGED_BREAKDOWN_BICG = -6,
429: KSP_DIVERGED_NONSYMMETRIC = -7,
430: KSP_DIVERGED_INDEFINITE_PC = -8,
431: KSP_DIVERGED_NANORINF = -9,
432: KSP_DIVERGED_INDEFINITE_MAT = -10,
433: KSP_DIVERGED_PCSETUP_FAILED = -11,
435: KSP_CONVERGED_ITERATING = 0} KSPConvergedReason;
436: PETSC_EXTERN const char *const*KSPConvergedReasons;
438: /*MC
439: KSP_CONVERGED_RTOL - norm(r) <= rtol*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: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
449: M*/
451: /*MC
452: KSP_CONVERGED_ATOL - norm(r) <= atol
454: Level: beginner
456: See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
457: for left preconditioning it is the 2-norm of the preconditioned residual, and the
458: 2-norm of the residual for right preconditioning
460: Level: beginner
462: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
464: M*/
466: /*MC
467: KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
469: Level: beginner
471: See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
472: for left preconditioning it is the 2-norm of the preconditioned residual, and the
473: 2-norm of the residual for right preconditioning
475: Level: beginner
477: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
479: M*/
481: /*MC
482: KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
483: reached
485: Level: beginner
487: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
489: M*/
491: /*MC
492: KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
493: the preconditioner is applied. Also used when the KSPConvergedSkip() convergence
494: test routine is set in KSP.
496: Level: beginner
498: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
500: M*/
502: /*MC
503: KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
504: method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or
505: preconditioner.
507: Level: beginner
509: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
511: M*/
513: /*MC
514: KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
515: method could not continue to enlarge the Krylov space.
517: Level: beginner
519: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
521: M*/
523: /*MC
524: KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
525: symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
527: Level: beginner
529: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
531: M*/
533: /*MC
534: KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
535: positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
536: be positive definite
538: Level: beginner
540: Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
541: the PCICC preconditioner to generate a positive definite preconditioner
543: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
545: M*/
547: /*MC
548: KSP_DIVERGED_PCSETUP_FAILED - It was not possible to build the requested preconditioner. This is usually due to a
549: zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner
550: such as PCFIELDSPLIT.
552: Level: beginner
554: Notes: Run with -ksp_error_if_not_converged to stop the program when the error is detected and print an error message with details.
557: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
559: M*/
561: /*MC
562: KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
563: while the KSPSolve() is still running.
565: Level: beginner
567: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
569: M*/
571: PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*));
572: PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **);
573: PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
574: PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
575: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *);
576: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **);
577: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
578: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
579: PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
580: PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *);
582: PETSC_DEPRECATED("Use KSPConvergedDefault()") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ }
583: #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
584: PETSC_DEPRECATED("Use KSPConvergedDefaultDestroy()") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ }
585: #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
586: PETSC_DEPRECATED("Use KSPConvergedDefaultCreate()") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ }
587: #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
588: PETSC_DEPRECATED("Use KSPConvergedDefaultSetUIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ }
589: #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
590: PETSC_DEPRECATED("Use KSPConvergedDefaultSetUMIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ }
591: #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
592: PETSC_DEPRECATED("Use KSPConvergedSkip()") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ }
593: #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
595: PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *);
597: /*E
598: KSPCGType - Determines what type of CG to use
600: Level: beginner
602: .seealso: KSPCGSetType()
603: E*/
604: typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
605: PETSC_EXTERN const char *const KSPCGTypes[];
607: PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
608: PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool );
610: PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal);
611: PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *);
612: PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *);
614: PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal);
615: PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *);
616: PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *);
618: PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal);
619: PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *);
620: PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *);
621: PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *);
622: PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *);
624: PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]);
626: PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
627: PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);
629: #include <petscdrawtypes.h>
630: PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
631: PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*);
632: PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
633: PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
634: PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
636: PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
637: PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
639: /* see src/ksp/ksp/interface/iguess.c */
640: typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool monitor;Mat mat; KSP ksp;}* KSPFischerGuess;
642: PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*);
643: PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*);
644: PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess);
645: PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec);
646: PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec);
647: PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess);
649: PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
650: PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess);
651: PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*);
653: /*E
654: MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines
656: Level: intermediate
658: .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat()
659: E*/
660: typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP} MatSchurComplementAinvType;
661: PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
663: PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
664: PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
665: PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
666: PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
667: PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
668: PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
669: PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType);
670: PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*);
671: PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*);
672: PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*);
673: PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatSchurComplementAinvType,MatReuse,Mat *);
674: PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*);
676: PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
677: PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool );
678: PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
679: PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
680: PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
681: PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *);
682: PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
683: PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
684: PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
685: PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*);
686: PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
687: PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
688: PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
689: PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
691: PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec);
692: PETSC_EXTERN PetscErrorCode DMProjectField(DM, Vec,
693: void (**)(PetscInt, PetscInt, PetscInt,
694: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
695: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
696: PetscReal, const PetscReal [], PetscScalar []), InsertMode, Vec);
697: #endif