Actual source code: petscksp.h
1: /*
2: Defines the interface functions for the Krylov subspace accelerators.
3: */
4: #ifndef PETSCKSP_H
5: #define 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: Notes:
17: When a direct solver is used, but no Krylov solver is used, the KSP object is still used but with a
18: KSPType of KSPPREONLY, meaning that only application of the preconditioner is used as the linear solver.
20: .seealso: KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP, KSPDestroy(), KSPCG, KSPGMRES
21: S*/
22: typedef struct _p_KSP* KSP;
24: /*J
25: KSPType - String with the name of a PETSc Krylov method.
27: Level: beginner
29: .seealso: KSPSetType(), KSP, KSPRegister(), KSPCreate(), KSPSetFromOptions()
30: J*/
31: typedef const char* KSPType;
32: #define KSPRICHARDSON "richardson"
33: #define KSPCHEBYSHEV "chebyshev"
34: #define KSPCG "cg"
35: #define KSPGROPPCG "groppcg"
36: #define KSPPIPECG "pipecg"
37: #define KSPPIPECGRR "pipecgrr"
38: #define KSPPIPELCG "pipelcg"
39: #define KSPPIPEPRCG "pipeprcg"
40: #define KSPPIPECG2 "pipecg2"
41: #define KSPCGNE "cgne"
42: #define KSPNASH "nash"
43: #define KSPSTCG "stcg"
44: #define KSPGLTR "gltr"
45: #define KSPCGNASH PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGNASH macro is deprecated use KSPNASH (since version 3.11)\"") "nash"
46: #define KSPCGSTCG PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGSTCG macro is deprecated use KSPSTCG (since version 3.11)\"") "stcg"
47: #define KSPCGGLTR PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGGLTR macro is deprecated use KSPSGLTR (since version 3.11)\"") "gltr"
48: #define KSPFCG "fcg"
49: #define KSPPIPEFCG "pipefcg"
50: #define KSPGMRES "gmres"
51: #define KSPPIPEFGMRES "pipefgmres"
52: #define KSPFGMRES "fgmres"
53: #define KSPLGMRES "lgmres"
54: #define KSPDGMRES "dgmres"
55: #define KSPPGMRES "pgmres"
56: #define KSPTCQMR "tcqmr"
57: #define KSPBCGS "bcgs"
58: #define KSPIBCGS "ibcgs"
59: #define KSPFBCGS "fbcgs"
60: #define KSPFBCGSR "fbcgsr"
61: #define KSPBCGSL "bcgsl"
62: #define KSPPIPEBCGS "pipebcgs"
63: #define KSPCGS "cgs"
64: #define KSPTFQMR "tfqmr"
65: #define KSPCR "cr"
66: #define KSPPIPECR "pipecr"
67: #define KSPLSQR "lsqr"
68: #define KSPPREONLY "preonly"
69: #define KSPQCG "qcg"
70: #define KSPBICG "bicg"
71: #define KSPMINRES "minres"
72: #define KSPSYMMLQ "symmlq"
73: #define KSPLCD "lcd"
74: #define KSPPYTHON "python"
75: #define KSPGCR "gcr"
76: #define KSPPIPEGCR "pipegcr"
77: #define KSPTSIRM "tsirm"
78: #define KSPCGLS "cgls"
79: #define KSPFETIDP "fetidp"
80: #define KSPHPDDM "hpddm"
82: /* Logging support */
83: PETSC_EXTERN PetscClassId KSP_CLASSID;
84: PETSC_EXTERN PetscClassId KSPGUESS_CLASSID;
85: PETSC_EXTERN PetscClassId DMKSP_CLASSID;
87: PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP*);
88: PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType);
89: PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType*);
90: PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
91: PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
92: PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec);
93: PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec);
94: PETSC_EXTERN PetscErrorCode KSPSetUseExplicitTranspose(KSP,PetscBool);
95: PETSC_EXTERN PetscErrorCode KSPMatSolve(KSP,Mat,Mat);
96: PETSC_EXTERN PetscErrorCode KSPSetMatSolveBatchSize(KSP,PetscInt);
97: PETSC_DEPRECATED_FUNCTION("Use KSPSetMatSolveBatchSize() (since version 3.15)") PETSC_STATIC_INLINE PetscErrorCode KSPSetMatSolveBlockSize(KSP ksp,PetscInt n) {return KSPSetMatSolveBatchSize(ksp,n);}
98: PETSC_EXTERN PetscErrorCode KSPGetMatSolveBatchSize(KSP,PetscInt*);
99: PETSC_DEPRECATED_FUNCTION("Use KSPGetMatSolveBatchSize() (since version 3.15)") PETSC_STATIC_INLINE PetscErrorCode KSPGetMatSolveBlockSize(KSP ksp,PetscInt *n) {return KSPGetMatSolveBatchSize(ksp,n);}
100: PETSC_EXTERN PetscErrorCode KSPReset(KSP);
101: PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP);
102: PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*);
103: PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool);
104: PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP,PetscBool*);
105: PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool);
106: PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP,PC,Vec);
108: PETSC_EXTERN PetscFunctionList KSPList;
109: PETSC_EXTERN PetscFunctionList KSPGuessList;
110: PETSC_EXTERN PetscFunctionList KSPMonitorList;
111: PETSC_EXTERN PetscFunctionList KSPMonitorCreateList;
112: PETSC_EXTERN PetscFunctionList KSPMonitorDestroyList;
113: PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP));
114: PETSC_EXTERN PetscErrorCode KSPMonitorRegister(const char[], PetscViewerType, PetscViewerFormat, PetscErrorCode (*)(KSP, PetscInt, PetscReal, PetscViewerAndFormat *), PetscErrorCode (*)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **), PetscErrorCode (*)(PetscViewerAndFormat **));
116: PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide);
117: PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*);
118: PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
119: PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
120: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool);
121: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool*);
122: PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool);
123: PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool*);
124: PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool);
125: PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool);
126: PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool*);
127: PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool);
128: PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool*);
129: PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec*);
130: PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec*);
131: PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*);
132: PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*);
133: PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*);
134: PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
135: PETSC_DEPRECATED_FUNCTION("Use KSPCreateVecs() (since version 3.6)") PETSC_STATIC_INLINE PetscErrorCode KSPGetVecs(KSP ksp,PetscInt n,Vec **x,PetscInt m,Vec **y) {return KSPCreateVecs(ksp,n,x,m,y);}
137: PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
138: PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
140: PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC);
141: PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*);
143: PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal);
144: PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode (*)(void**));
145: PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
146: PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void**);
147: PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,const PetscReal*[],PetscInt*);
148: PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool);
149: PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP,const PetscReal*[],PetscInt*);
150: PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP,PetscReal[],PetscInt,PetscBool);
152: PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*);
153: PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec*);
154: PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
155: PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt);
157: PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*);
158: PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC,KSP);
159: PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
160: PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
161: PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
162: PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
163: PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC,PetscInt*,KSP*[]);
164: PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*);
165: PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*);
166: PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*);
167: PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*);
168: PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP*);
169: PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC,KSP*);
170: /*
171: PCMGCoarseList contains the list of coarse space constructor currently registered
172: These are added with PCMGRegisterCoarseSpaceConstructor()
173: */
174: PETSC_EXTERN PetscFunctionList PCMGCoarseList;
175: PETSC_EXTERN PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char[], PetscErrorCode (*)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **));
176: PETSC_EXTERN PetscErrorCode PCMGGetCoarseSpaceConstructor(const char[], PetscErrorCode (**)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **));
179: PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec*);
180: PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec*);
182: PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
183: PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool);
184: PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal);
185: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal);
186: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP,PetscBool);
187: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*);
188: PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
189: PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt*);
190: PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]);
191: PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt*,Vec[],PetscReal[],PetscReal[]);
193: /*E
195: KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods
197: KSP_FCD_TRUNC_TYPE_STANDARD uses all (up to mmax) stored directions
198: KSP_FCD_TRUNC_TYPE_NOTAY uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..
200: Level: intermediate
201: .seealso : KSPFCG,KSPPIPEFCG,KSPPIPEGCR,KSPFCGSetTruncationType(),KSPFCGGetTruncationType()
203: E*/
204: typedef enum {KSP_FCD_TRUNC_TYPE_STANDARD,KSP_FCD_TRUNC_TYPE_NOTAY} KSPFCDTruncationType;
205: PETSC_EXTERN const char *const KSPFCDTruncationTypes[];
207: PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt);
208: PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*);
209: PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt);
210: PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*);
211: PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType);
212: PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*);
214: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt);
215: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*);
216: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt);
217: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*);
218: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType);
219: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*);
221: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt);
222: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*);
223: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt);
224: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*);
225: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType);
226: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*);
227: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool);
228: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*);
230: PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
231: PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*);
232: PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);
233: PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP,PetscReal);
235: PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
236: PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
237: PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
238: PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
239: PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
241: PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
242: PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
244: PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP,PetscScalar);
246: PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
247: PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
248: PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
250: PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP,PC*);
251: PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP,PC);
252: PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP,KSP*);
253: PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP,Mat);
255: PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationSpace(KSP,Mat);
256: PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationSpace(KSP,Mat*);
257: PETSC_DEPRECATED_FUNCTION("Use KSPMatSolve() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X) { return KSPMatSolve(ksp, B, X); }
258: /*E
259: KSPHPDDMType - Type of Krylov method used by KSPHPDDM
261: Level: intermediate
263: Values:
264: + KSP_HPDDM_TYPE_GMRES (default)
265: . KSP_HPDDM_TYPE_BGMRES
266: . KSP_HPDDM_TYPE_CG
267: . KSP_HPDDM_TYPE_BCG
268: . KSP_HPDDM_TYPE_GCRODR
269: . KSP_HPDDM_TYPE_BGCRODR
270: . KSP_HPDDM_TYPE_BFBCG
271: - KSP_HPDDM_TYPE_PREONLY
273: .seealso: KSPHPDDM, KSPHPDDMSetType()
274: E*/
275: typedef enum {
276: KSP_HPDDM_TYPE_GMRES = 0,
277: KSP_HPDDM_TYPE_BGMRES = 1,
278: KSP_HPDDM_TYPE_CG = 2,
279: KSP_HPDDM_TYPE_BCG = 3,
280: KSP_HPDDM_TYPE_GCRODR = 4,
281: KSP_HPDDM_TYPE_BGCRODR = 5,
282: KSP_HPDDM_TYPE_BFBCG = 6,
283: KSP_HPDDM_TYPE_PREONLY = 7
284: } KSPHPDDMType;
285: PETSC_EXTERN const char *const KSPHPDDMTypes[];
286: PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP,KSPHPDDMType);
287: PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP,KSPHPDDMType*);
289: /*E
290: KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
292: Level: advanced
294: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
295: KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
297: E*/
298: typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
299: PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
300: /*MC
301: KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
303: Level: advanced
305: Note: Possible unstable, but the fastest to compute
307: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
308: KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
309: KSPGMRESModifiedGramSchmidtOrthogonalization()
310: M*/
312: /*MC
313: KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
314: iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
315: poor orthogonality.
317: Level: advanced
319: Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
320: estimate the orthogonality but is more stable.
322: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
323: KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
324: KSPGMRESModifiedGramSchmidtOrthogonalization()
325: M*/
327: /*MC
328: KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
330: Level: advanced
332: Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
333: but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
335: You should only use this if you absolutely know that the iterative refinement is needed.
337: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
338: KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
339: KSPGMRESModifiedGramSchmidtOrthogonalization()
340: M*/
342: PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
343: PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);
345: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
346: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
347: PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
349: PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
350: PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
351: PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);
353: PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
354: PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool);
355: PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
356: PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool);
358: PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
359: PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);
360: PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
362: PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP,const char[],const char[],void*);
363: PETSC_EXTERN PetscErrorCode KSPMonitorLGCreate(MPI_Comm,const char[],const char[],const char[],PetscInt,const char *[],int,int,int,int,PetscDrawLG *);
364: PETSC_EXTERN PetscErrorCode KSPMonitorResidual(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
365: PETSC_EXTERN PetscErrorCode KSPMonitorResidualDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
366: PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
367: PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **);
368: PETSC_EXTERN PetscErrorCode KSPMonitorResidualShort(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
369: PETSC_EXTERN PetscErrorCode KSPMonitorResidualRange(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
370: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidual(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
371: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
372: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
373: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **);
374: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMax(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
375: PETSC_EXTERN PetscErrorCode KSPMonitorError(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
376: PETSC_EXTERN PetscErrorCode KSPMonitorErrorDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
377: PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
378: PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **);
379: PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
380: PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
381: PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
382: PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **);
383: PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
384: PETSC_EXTERN PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **);
385: PETSC_DEPRECATED_FUNCTION("Use KSPMonitorResidual() (since version 3.15)") PETSC_STATIC_INLINE PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {return KSPMonitorResidual(ksp,n,rnorm,vf);}
386: PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidual() (since version 3.15)") PETSC_STATIC_INLINE PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {return KSPMonitorTrueResidual(ksp,n,rnorm,vf);}
387: PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidualMax() (since version 3.15)") PETSC_STATIC_INLINE PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {return KSPMonitorTrueResidualMax(ksp,n,rnorm,vf);}
389: PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void*);
390: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
391: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
392: PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*);
393: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**);
394: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**);
396: PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
397: PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
399: PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat);
400: PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*);
401: PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool*,PetscBool*);
402: PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
403: PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
404: PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);
406: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool);
407: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool*);
408: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool);
409: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool*);
411: PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
412: PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);
413: PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP,PetscObject,const char[]);
414: PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP,PetscViewer);
415: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP,PetscErrorCode (*)(KSP,void*),void *vctx,PetscErrorCode (*)(void**));
416: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP);
417: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP);
418: PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP,PetscViewer);
420: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonView() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode KSPReasonView(KSP ksp,PetscViewer v) {return KSPConvergedReasonView(ksp,v);}
421: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonViewFromOptions() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode KSPReasonViewFromOptions(KSP ksp) {return KSPConvergedReasonViewFromOptions(ksp);}
423: #define KSP_FILE_CLASSID 1211223
425: PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP,PetscBool);
426: PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP,PetscBool);
427: PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
428: PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP,PetscReal*,PetscReal*);
429: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidual(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
430: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
431: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **);
433: PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
434: PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
435: PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*);
437: /*E
438: KSPNormType - Norm that is passed in the Krylov convergence
439: test routines.
441: Level: advanced
443: Each solver only supports a subset of these and some may support different ones
444: depending on left or right preconditioning, see KSPSetPCSide()
446: Notes:
447: this must match petsc/finclude/petscksp.h
449: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
450: KSPSetConvergenceTest(), KSPSetPCSide()
451: E*/
452: typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
453: #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
454: PETSC_EXTERN const char *const*const KSPNormTypes;
456: /*MC
457: KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
458: possibly save some computation but means the convergence test cannot
459: be based on a norm of a residual etc.
461: Level: advanced
463: Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored
465: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
466: M*/
468: /*MC
469: KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
470: convergence test routine.
472: Level: advanced
474: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
475: M*/
477: /*MC
478: KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
479: convergence test routine.
481: Level: advanced
483: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
484: M*/
486: /*MC
487: KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
488: convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG, KSPPIPEFCG, KSPPIPEGCR
490: Level: advanced
492: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
493: M*/
495: PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
496: PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
497: PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
498: PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
499: PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool);
501: #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM("Use KSP_DIVERGED_PC_FAILED (since version 3.11)")
502: /*E
503: KSPConvergedReason - reason a Krylov method was said to have converged or diverged
505: Level: beginner
507: Notes:
508: See KSPGetConvergedReason() for explanation of each value
510: Developer Notes:
511: this must match petsc/finclude/petscksp.h
513: The string versions of these are KSPConvergedReasons; if you change
514: any of the values here also change them that array of names.
516: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances(), KSPConvergedReasonView()
517: E*/
518: typedef enum {/* converged */
519: KSP_CONVERGED_RTOL_NORMAL = 1,
520: KSP_CONVERGED_ATOL_NORMAL = 9,
521: KSP_CONVERGED_RTOL = 2,
522: KSP_CONVERGED_ATOL = 3,
523: KSP_CONVERGED_ITS = 4,
524: KSP_CONVERGED_CG_NEG_CURVE = 5,
525: KSP_CONVERGED_CG_CONSTRAINED = 6,
526: KSP_CONVERGED_STEP_LENGTH = 7,
527: KSP_CONVERGED_HAPPY_BREAKDOWN = 8,
528: /* diverged */
529: KSP_DIVERGED_NULL = -2,
530: KSP_DIVERGED_ITS = -3,
531: KSP_DIVERGED_DTOL = -4,
532: KSP_DIVERGED_BREAKDOWN = -5,
533: KSP_DIVERGED_BREAKDOWN_BICG = -6,
534: KSP_DIVERGED_NONSYMMETRIC = -7,
535: KSP_DIVERGED_INDEFINITE_PC = -8,
536: KSP_DIVERGED_NANORINF = -9,
537: KSP_DIVERGED_INDEFINITE_MAT = -10,
538: KSP_DIVERGED_PC_FAILED = -11,
539: KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11,
541: KSP_CONVERGED_ITERATING = 0} KSPConvergedReason;
542: PETSC_EXTERN const char *const*KSPConvergedReasons;
544: /*MC
545: KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if KSPConvergedDefaultSetUIRNorm() was called
547: Level: beginner
549: See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
550: for left preconditioning it is the 2-norm of the preconditioned residual, and the
551: 2-norm of the residual for right preconditioning
553: See also KSP_CONVERGED_ATOL which may apply before this tolerance.
555: .seealso: KSP_CONVERGED_ATOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
557: M*/
559: /*MC
560: KSP_CONVERGED_ATOL - norm(r) <= atol
562: Level: beginner
564: See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
565: for left preconditioning it is the 2-norm of the preconditioned residual, and the
566: 2-norm of the residual for right preconditioning
568: See also KSP_CONVERGED_RTOL which may apply before this tolerance.
570: .seealso: KSP_CONVERGED_RTOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
572: M*/
574: /*MC
575: KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
577: Level: beginner
579: See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
580: for left preconditioning it is the 2-norm of the preconditioned residual, and the
581: 2-norm of the residual for right preconditioning
583: Level: beginner
585: .seealso: KSP_CONVERGED_ATOL, KSP_DIVERGED_RTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
587: M*/
589: /*MC
590: KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
591: reached
593: Level: beginner
595: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
597: M*/
599: /*MC
600: KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
601: the preconditioner is applied. Also used when the KSPConvergedSkip() convergence
602: test routine is set in KSP.
604: Level: beginner
606: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
608: M*/
610: /*MC
611: KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
612: method could not continue to enlarge the Krylov space. Could be due to a singular matrix or
613: preconditioner. In KSPHPDDM, this is also returned when some search directions within a block
614: are colinear.
616: Level: beginner
618: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
620: M*/
622: /*MC
623: KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
624: method could not continue to enlarge the Krylov space.
626: Level: beginner
628: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
630: M*/
632: /*MC
633: KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
634: symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
636: Level: beginner
638: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
640: M*/
642: /*MC
643: KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
644: positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
645: be positive definite
647: Level: beginner
649: Notes:
650: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
651: the PCICC preconditioner to generate a positive definite preconditioner
653: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
655: M*/
657: /*MC
658: KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a
659: zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner
660: such as PCFIELDSPLIT.
662: Level: beginner
664: Notes:
665: Run with -ksp_error_if_not_converged to stop the program when the error is detected and print an error message with details.
668: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
670: M*/
672: /*MC
673: KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
674: while the KSPSolve() is still running.
676: Level: beginner
678: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
680: M*/
682: PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void*,PetscErrorCode (*)(void*));
683: PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*));
684: PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*));
685: PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void**);
686: PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
687: PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
688: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void*);
689: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void**);
690: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
691: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
692: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP,PetscBool);
693: PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
694: PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason*);
695: PETSC_EXTERN PetscErrorCode KSPGetConvergedReasonString(KSP,const char**);
696: PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP,PetscReal*,PetscReal*,PetscReal*,PetscReal*);
698: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ }
699: #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
700: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ }
701: #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
702: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ }
703: #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
704: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ }
705: #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
706: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ }
707: #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
708: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ }
709: #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
711: PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP,MatType,Mat*);
712: PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPComputeExplicitOperator(KSP A,Mat* B) { return KSPComputeOperator(A,NULL,B); }
714: /*E
715: KSPCGType - Determines what type of CG to use
717: Level: beginner
719: .seealso: KSPCGSetType()
720: E*/
721: typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
722: PETSC_EXTERN const char *const KSPCGTypes[];
724: PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
725: PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool);
727: PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP,PetscReal);
728: PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP,PetscReal*);
729: PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP,PetscReal*);
731: PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal*);
732: PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal*);
733: PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp,PetscReal *x) {return KSPGLTRGetMinEig(ksp,x);}
734: PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetLambda(KSP ksp,PetscReal *x) {return KSPGLTRGetLambda(ksp,x);}
736: PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]);
738: PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
739: PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);
741: #include <petscdrawtypes.h>
742: PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
744: PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
745: PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
747: /*S
748: KSPGuess - Abstract PETSc object that manages all initial guess methods in Krylov methods.
750: Level: beginner
752: .seealso: KSPCreate(), KSPSetGuessType(), KSPGuessType
753: S*/
754: typedef struct _p_KSPGuess* KSPGuess;
755: /*J
756: KSPGuessType - String with the name of a PETSc initial guess for Krylov methods.
758: Level: beginner
760: .seealso: KSPGuess
761: J*/
762: typedef const char* KSPGuessType;
763: #define KSPGUESSFISCHER "fischer"
764: #define KSPGUESSPOD "pod"
765: PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[],PetscErrorCode (*)(KSPGuess));
766: PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP,KSPGuess);
767: PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP,KSPGuess*);
768: PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess,PetscViewer);
769: PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess*);
770: PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm,KSPGuess*);
771: PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess,KSPGuessType);
772: PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess,KSPGuessType*);
773: PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
774: PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess,Vec,Vec);
775: PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess,Vec,Vec);
776: PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
777: PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess,PetscInt,PetscInt);
778: PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
779: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool);
780: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool*);
782: /*E
783: MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines
785: Level: intermediate
787: .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat()
788: E*/
789: typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG} MatSchurComplementAinvType;
790: PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
792: typedef enum {MAT_LMVM_SYMBROYDEN_SCALE_NONE = 0,
793: MAT_LMVM_SYMBROYDEN_SCALE_SCALAR = 1,
794: MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2,
795: MAT_LMVM_SYMBROYDEN_SCALE_USER = 3} MatLMVMSymBroydenScaleType;
796: PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[];
798: PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
799: PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
800: PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
801: PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
802: PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
803: PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
804: PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType);
805: PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*);
806: PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*);
807: PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*);
808: PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatSchurComplementAinvType,MatReuse,Mat*);
809: PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*);
811: PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm,PetscInt,PetscInt,Mat*);
812: PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm,PetscInt,PetscInt,Mat*);
813: PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm,PetscInt,PetscInt,Mat*);
814: PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
815: PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
816: PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
817: PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
818: PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
820: PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
821: PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool*);
822: PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
823: PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
824: PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
825: PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
826: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
827: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
828: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
829: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
830: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
831: PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
832: PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
833: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat*);
834: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC*);
835: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP*);
836: PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt);
837: PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt*);
838: PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt*);
839: PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar);
840: PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType);
842: PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
843: PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool);
844: PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
845: PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
846: PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
847: PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void*);
848: PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
849: PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
850: PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
851: PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*);
852: PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
853: PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
854: PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
855: PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
857: PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec);
858: PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec,
859: void (**)(PetscInt, PetscInt, PetscInt,
860: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
861: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
862: PetscReal, const PetscReal [], PetscInt, const PetscScalar[], PetscScalar []), InsertMode, Vec);
865: PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, PetscInt, Vec[], Vec[], Mat *, void *);
866: PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, PetscInt, Vec[], Vec[], PetscReal);
867: #endif