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

  7: #include <petscpc.h>

  9: /* SUBMANSEC = KSP */

 11: PETSC_EXTERN PetscErrorCode KSPInitializePackage(void);

 13: /*S
 14:      KSP - Abstract PETSc object that manages all Krylov methods. This is the object that manages the
 15:          linear solves in PETSc (even those such as direct solvers that do no use Krylov accelerators).

 17:    Level: beginner

 19:    Note:
 20:     When a direct solver is used, but no Krylov solver is used, the KSP object is still used but with a
 21:     `KSPType` of `KSPPREONLY`, meaning that only application of the preconditioner is used as the linear solver.

 23: .seealso: `KSPCreate()`, `KSPSetType()`, `KSPType`, `SNES`, `TS`, `PC`, `KSP`, `KSPDestroy()`, `KSPCG`, `KSPGMRES`
 24: S*/
 25: typedef struct _p_KSP *KSP;

 27: /*J
 28:     KSPType - String with the name of a PETSc Krylov method.

 30:    Level: beginner

 32: .seealso: `KSPSetType()`, `KSP`, `KSPRegister()`, `KSPCreate()`, `KSPSetFromOptions()`
 33: J*/
 34: typedef const char *KSPType;
 35: #define KSPRICHARDSON "richardson"
 36: #define KSPCHEBYSHEV  "chebyshev"
 37: #define KSPCG         "cg"
 38: #define KSPGROPPCG    "groppcg"
 39: #define KSPPIPECG     "pipecg"
 40: #define KSPPIPECGRR   "pipecgrr"
 41: #define KSPPIPELCG    "pipelcg"
 42: #define KSPPIPEPRCG   "pipeprcg"
 43: #define KSPPIPECG2    "pipecg2"
 44: #define KSPCGNE       "cgne"
 45: #define KSPNASH       "nash"
 46: #define KSPSTCG       "stcg"
 47: #define KSPGLTR       "gltr"
 48: #define KSPCGNASH     PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGNASH macro is deprecated use KSPNASH (since version 3.11)\"") "nash"
 49: #define KSPCGSTCG     PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGSTCG macro is deprecated use KSPSTCG (since version 3.11)\"") "stcg"
 50: #define KSPCGGLTR     PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGGLTR macro is deprecated use KSPSGLTR (since version 3.11)\"") "gltr"
 51: #define KSPFCG        "fcg"
 52: #define KSPPIPEFCG    "pipefcg"
 53: #define KSPGMRES      "gmres"
 54: #define KSPPIPEFGMRES "pipefgmres"
 55: #define KSPFGMRES     "fgmres"
 56: #define KSPLGMRES     "lgmres"
 57: #define KSPDGMRES     "dgmres"
 58: #define KSPPGMRES     "pgmres"
 59: #define KSPTCQMR      "tcqmr"
 60: #define KSPBCGS       "bcgs"
 61: #define KSPIBCGS      "ibcgs"
 62: #define KSPQMRCGS     "qmrcgs"
 63: #define KSPFBCGS      "fbcgs"
 64: #define KSPFBCGSR     "fbcgsr"
 65: #define KSPBCGSL      "bcgsl"
 66: #define KSPPIPEBCGS   "pipebcgs"
 67: #define KSPCGS        "cgs"
 68: #define KSPTFQMR      "tfqmr"
 69: #define KSPCR         "cr"
 70: #define KSPPIPECR     "pipecr"
 71: #define KSPLSQR       "lsqr"
 72: #define KSPPREONLY    "preonly"
 73: #define KSPNONE       "none"
 74: #define KSPQCG        "qcg"
 75: #define KSPBICG       "bicg"
 76: #define KSPMINRES     "minres"
 77: #define KSPSYMMLQ     "symmlq"
 78: #define KSPLCD        "lcd"
 79: #define KSPPYTHON     "python"
 80: #define KSPGCR        "gcr"
 81: #define KSPPIPEGCR    "pipegcr"
 82: #define KSPTSIRM      "tsirm"
 83: #define KSPCGLS       "cgls"
 84: #define KSPFETIDP     "fetidp"
 85: #define KSPHPDDM      "hpddm"

 87: /* Logging support */
 88: PETSC_EXTERN PetscClassId KSP_CLASSID;
 89: PETSC_EXTERN PetscClassId KSPGUESS_CLASSID;
 90: PETSC_EXTERN PetscClassId DMKSP_CLASSID;

 92: PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm, KSP *);
 93: PETSC_EXTERN PetscErrorCode KSPSetType(KSP, KSPType);
 94: PETSC_EXTERN PetscErrorCode KSPGetType(KSP, KSPType *);
 95: PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
 96: PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
 97: PETSC_EXTERN PetscErrorCode KSPSolve(KSP, Vec, Vec);
 98: PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP, Vec, Vec);
 99: PETSC_EXTERN PetscErrorCode KSPSetUseExplicitTranspose(KSP, PetscBool);
100: PETSC_EXTERN PetscErrorCode KSPMatSolve(KSP, Mat, Mat);
101: PETSC_EXTERN PetscErrorCode KSPSetMatSolveBatchSize(KSP, PetscInt);
102: PETSC_DEPRECATED_FUNCTION("Use KSPSetMatSolveBatchSize() (since version 3.15)") static inline PetscErrorCode KSPSetMatSolveBlockSize(KSP ksp, PetscInt n) {
103:   return KSPSetMatSolveBatchSize(ksp, n);
104: }
105: PETSC_EXTERN PetscErrorCode KSPGetMatSolveBatchSize(KSP, PetscInt *);
106: PETSC_DEPRECATED_FUNCTION("Use KSPGetMatSolveBatchSize() (since version 3.15)") static inline PetscErrorCode KSPGetMatSolveBlockSize(KSP ksp, PetscInt *n) {
107:   return KSPGetMatSolveBatchSize(ksp, n);
108: }
109: PETSC_EXTERN PetscErrorCode KSPReset(KSP);
110: PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP);
111: PETSC_EXTERN PetscErrorCode KSPDestroy(KSP *);
112: PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP, PetscBool);
113: PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP, PetscBool *);
114: PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP, PetscBool);
115: PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP, PC, Vec);

117: PETSC_EXTERN PetscFunctionList KSPList;
118: PETSC_EXTERN PetscFunctionList KSPGuessList;
119: PETSC_EXTERN PetscFunctionList KSPMonitorList;
120: PETSC_EXTERN PetscFunctionList KSPMonitorCreateList;
121: PETSC_EXTERN PetscFunctionList KSPMonitorDestroyList;
122: PETSC_EXTERN PetscErrorCode    KSPRegister(const char[], PetscErrorCode (*)(KSP));
123: PETSC_EXTERN PetscErrorCode    KSPMonitorRegister(const char[], PetscViewerType, PetscViewerFormat, PetscErrorCode (*)(KSP, PetscInt, PetscReal, PetscViewerAndFormat *), PetscErrorCode (*)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **), PetscErrorCode (*)(PetscViewerAndFormat **));

125: PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP, PCSide);
126: PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP, PCSide *);
127: PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP, PetscReal, PetscReal, PetscReal, PetscInt);
128: PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP, PetscReal *, PetscReal *, PetscReal *, PetscInt *);
129: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP, PetscBool);
130: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP, PetscBool *);
131: PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP, PetscBool);
132: PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP, PetscBool *);
133: PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP, PetscBool);
134: PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP, PetscBool);
135: PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP, PetscBool *);
136: PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP, PetscBool);
137: PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP, PetscBool *);
138: PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP, Vec *);
139: PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP, Vec *);
140: PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP, PetscReal *);
141: PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP, PetscInt *);
142: PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP, PetscInt *);
143: PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP, PetscInt, Vec **, PetscInt, Vec **);
144: PETSC_DEPRECATED_FUNCTION("Use KSPCreateVecs() (since version 3.6)") static inline PetscErrorCode KSPGetVecs(KSP ksp, PetscInt n, Vec **x, PetscInt m, Vec **y) {
145:   return KSPCreateVecs(ksp, n, x, m, y);
146: }

148: PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *);
149: PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *);

151: PETSC_EXTERN PetscErrorCode KSPSetPC(KSP, PC);
152: PETSC_EXTERN PetscErrorCode KSPGetPC(KSP, PC *);

154: PETSC_EXTERN PetscErrorCode KSPMonitor(KSP, PetscInt, PetscReal);
155: PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void **));
156: PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
157: PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP, void *);
158: PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP, const PetscReal *[], PetscInt *);
159: PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP, PetscReal[], PetscInt, PetscBool);
160: PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP, const PetscReal *[], PetscInt *);
161: PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP, PetscReal[], PetscInt, PetscBool);

163: PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP, Vec, Vec *);
164: PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP, Vec, Vec, Vec *);
165: PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
166: PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP, PetscInt);

168: PETSC_EXTERN PetscErrorCode    PCKSPGetKSP(PC, KSP *);
169: PETSC_EXTERN PetscErrorCode    PCKSPSetKSP(PC, KSP);
170: PETSC_EXTERN PetscErrorCode    PCBJacobiGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
171: PETSC_EXTERN PetscErrorCode    PCASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
172: PETSC_EXTERN PetscErrorCode    PCGASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
173: PETSC_EXTERN PetscErrorCode    PCFieldSplitGetSubKSP(PC, PetscInt *, KSP *[]);
174: PETSC_EXTERN PetscErrorCode    PCFieldSplitSchurGetSubKSP(PC, PetscInt *, KSP *[]);
175: PETSC_EXTERN PetscErrorCode    PCMGGetSmoother(PC, PetscInt, KSP *);
176: PETSC_EXTERN PetscErrorCode    PCMGGetSmootherDown(PC, PetscInt, KSP *);
177: PETSC_EXTERN PetscErrorCode    PCMGGetSmootherUp(PC, PetscInt, KSP *);
178: PETSC_EXTERN PetscErrorCode    PCMGGetCoarseSolve(PC, KSP *);
179: PETSC_EXTERN PetscErrorCode    PCGalerkinGetKSP(PC, KSP *);
180: PETSC_EXTERN PetscErrorCode    PCDeflationGetCoarseKSP(PC, KSP *);
181: /*
182:   PCMGCoarseList contains the list of coarse space constructor currently registered
183:   These are added with PCMGRegisterCoarseSpaceConstructor()
184: */
185: PETSC_EXTERN PetscFunctionList PCMGCoarseList;
186: PETSC_EXTERN PetscErrorCode    PCMGRegisterCoarseSpaceConstructor(const char[], PetscErrorCode (*)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *));
187: PETSC_EXTERN PetscErrorCode    PCMGGetCoarseSpaceConstructor(const char[], PetscErrorCode (**)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *));

189: PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP, Vec, Vec *);
190: PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP, Vec, Vec, Vec *);

192: PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP, PetscReal);
193: PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP, PetscBool);
194: PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP, PetscReal, PetscReal);
195: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP, PetscReal, PetscReal, PetscReal, PetscReal);
196: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP, PetscBool);
197: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP, KSP *);
198: PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP, PetscReal *, PetscReal *);
199: PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP, PetscInt, PetscReal[], PetscReal[], PetscInt *);
200: PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP, PetscInt, PetscReal[], PetscReal[]);
201: PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP, PetscBool, PetscBool, PetscInt *, Vec[], PetscReal[], PetscReal[]);

203: /*E

205:   KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods

207:   Values:
208: $ `KSP_FCD_TRUNC_TYPE_STANDARD` - uses all (up to mmax) stored directions
209: $ `KSP_FCD_TRUNC_TYPE_NOTAY` - uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..

211:    Level: intermediate

213: .seealso `:` `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`, `KSPFCGSetTruncationType()`, `KSPFCGGetTruncationType()`
214: E*/
215: typedef enum {
216:   KSP_FCD_TRUNC_TYPE_STANDARD,
217:   KSP_FCD_TRUNC_TYPE_NOTAY
218: } KSPFCDTruncationType;
219: PETSC_EXTERN const char *const KSPFCDTruncationTypes[];

221: PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP, PetscInt);
222: PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP, PetscInt *);
223: PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP, PetscInt);
224: PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP, PetscInt *);
225: PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP, KSPFCDTruncationType);
226: PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP, KSPFCDTruncationType *);

228: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP, PetscInt);
229: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP, PetscInt *);
230: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP, PetscInt);
231: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP, PetscInt *);
232: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP, KSPFCDTruncationType);
233: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP, KSPFCDTruncationType *);

235: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP, PetscInt);
236: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP, PetscInt *);
237: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP, PetscInt);
238: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP, PetscInt *);
239: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP, KSPFCDTruncationType);
240: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP, KSPFCDTruncationType *);
241: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP, PetscBool);
242: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP, PetscBool *);

244: PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
245: PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt *);
246: PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP, PetscReal);
247: PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP, PetscReal);

249: PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
250: PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP, PetscErrorCode (*)(KSP, PetscInt));
251: PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP, PetscErrorCode (**)(KSP, PetscInt));
252: PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP, PetscInt);
253: PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP, PetscInt);

255: PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP, PetscInt);
256: PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);

258: PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP, PetscScalar);

260: PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP, PetscInt);
261: PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP, PetscInt *);
262: PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *));

264: PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP, PC *);
265: PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP, PC);
266: PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP, KSP *);
267: PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP, Mat);

269: PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationMat(KSP, Mat);
270: PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationMat(KSP, Mat *);
271: #if PetscDefined(HAVE_HPDDM)
272: PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMSetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMSetDeflationSpace(KSP ksp, Mat U) {
273:   return KSPHPDDMSetDeflationMat(ksp, U);
274: }
275: PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMGetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMGetDeflationSpace(KSP ksp, Mat *U) {
276:   return KSPHPDDMGetDeflationMat(ksp, U);
277: }
278: #endif
279: PETSC_DEPRECATED_FUNCTION("Use KSPMatSolve() (since version 3.14)") static inline PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X) {
280:   return KSPMatSolve(ksp, B, X);
281: }
282: /*E
283:     KSPHPDDMType - Type of Krylov method used by `KSPHPDDM`

285:     Level: intermediate

287:     Values:
288: $   `KSP_HPDDM_TYPE_GMRES` (default)
289: $   `KSP_HPDDM_TYPE_BGMRES`
290: $   `KSP_HPDDM_TYPE_CG`
291: $   `KSP_HPDDM_TYPE_BCG`
292: $   `KSP_HPDDM_TYPE_GCRODR`
293: $   `KSP_HPDDM_TYPE_BGCRODR`
294: $   `KSP_HPDDM_TYPE_BFBCG`
295: $   `KSP_HPDDM_TYPE_PREONLY`

297: .seealso: `KSPHPDDM`, `KSPHPDDMSetType()`
298: E*/
299: typedef enum {
300:   KSP_HPDDM_TYPE_GMRES   = 0,
301:   KSP_HPDDM_TYPE_BGMRES  = 1,
302:   KSP_HPDDM_TYPE_CG      = 2,
303:   KSP_HPDDM_TYPE_BCG     = 3,
304:   KSP_HPDDM_TYPE_GCRODR  = 4,
305:   KSP_HPDDM_TYPE_BGCRODR = 5,
306:   KSP_HPDDM_TYPE_BFBCG   = 6,
307:   KSP_HPDDM_TYPE_PREONLY = 7
308: } KSPHPDDMType;
309: PETSC_EXTERN const char *const KSPHPDDMTypes[];
310: /*E
311:     KSPHPDDMPrecision - Precision of Krylov bases used by `KSPHPDDM`

313:     Level: intermediate

315:     Values:
316: $   `KSP_HPDDM_PRECISION_HALF`
317: $   `KSP_HPDDM_PRECISION_SINGLE` (default when PETSc is configured --with-precision=single)
318: $   `KSP_HPDDM_PRECISION_DOUBLE` (default when PETSc is configured --with-precision=double)
319: $   `KSP_HPDDM_PRECISION_QUADRUPLE` (default when PETSc is configured --with-precision=__float128)

321: .seealso: `KSPHPDDM`
322: E*/
323: typedef enum {
324:   KSP_HPDDM_PRECISION_HALF      = 0,
325:   KSP_HPDDM_PRECISION_SINGLE    = 1,
326:   KSP_HPDDM_PRECISION_DOUBLE    = 2,
327:   KSP_HPDDM_PRECISION_QUADRUPLE = 3
328: } KSPHPDDMPrecision;
329: PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP, KSPHPDDMType);
330: PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP, KSPHPDDMType *);

332: /*E
333:     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.

335:    Level: advanced

337: .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`,
338:           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`

340: E*/
341: typedef enum {
342:   KSP_GMRES_CGS_REFINE_NEVER,
343:   KSP_GMRES_CGS_REFINE_IFNEEDED,
344:   KSP_GMRES_CGS_REFINE_ALWAYS
345: } KSPGMRESCGSRefinementType;
346: PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
347: /*MC
348:     KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process

350:    Level: advanced

352:    Note:
353:    Possibly unstable, but the fastest to compute

355: .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`,
356:           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
357:           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
358: M*/

360: /*MC
361:     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
362:           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
363:           poor orthogonality.

365:    Level: advanced

367:    Note: This is slower than `KSP_GMRES_CGS_REFINE_NEVER` because it requires an extra norm computation to
368:      estimate the orthogonality but is more stable.

370: .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`,
371:           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
372:           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
373: M*/

375: /*MC
376:     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.

378:    Level: advanced

380:    Notes:
381:    This is roughly twice the cost of `KSP_GMRES_CGS_REFINE_NEVER` because it performs the process twice
382:    but it saves the extra norm calculation needed by `KSP_GMRES_CGS_REFINE_IFNEEDED`.

384:    You should only use this if you absolutely know that the iterative refinement is needed.

386: .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`,
387:           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
388:           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
389: M*/

391: PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP, KSPGMRESCGSRefinementType);
392: PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP, KSPGMRESCGSRefinementType *);

394: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP, PetscInt, PetscInt, PetscReal, void *);
395: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP, PetscInt, PetscInt, PetscReal, void *);
396: PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *));

398: PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP, PetscReal);
399: PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP, PetscReal *);
400: PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP, PetscReal *);

402: PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP, PetscReal);
403: PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP, PetscBool);
404: PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP, PetscInt);
405: PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP, PetscBool);

407: PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
408: PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);
409: PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));

411: PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP, const char[], const char[], void *);
412: PETSC_EXTERN PetscErrorCode KSPMonitorLGCreate(MPI_Comm, const char[], const char[], const char[], PetscInt, const char *[], int, int, int, int, PetscDrawLG *);
413: PETSC_EXTERN PetscErrorCode KSPMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
414: PETSC_EXTERN PetscErrorCode KSPMonitorResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
415: PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
416: PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
417: PETSC_EXTERN PetscErrorCode KSPMonitorResidualShort(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
418: PETSC_EXTERN PetscErrorCode KSPMonitorResidualRange(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
419: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
420: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
421: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
422: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
423: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMax(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
424: PETSC_EXTERN PetscErrorCode KSPMonitorError(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
425: PETSC_EXTERN PetscErrorCode KSPMonitorErrorDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
426: PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
427: PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
428: PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
429: PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
430: PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
431: PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
432: PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
433: PETSC_EXTERN PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
434: PETSC_DEPRECATED_FUNCTION("Use KSPMonitorResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {
435:   return KSPMonitorResidual(ksp, n, rnorm, vf);
436: }
437: PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {
438:   return KSPMonitorTrueResidual(ksp, n, rnorm, vf);
439: }
440: PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidualMax() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {
441:   return KSPMonitorTrueResidualMax(ksp, n, rnorm, vf);
442: }

444: PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP, PetscInt, PetscReal, void *);
445: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp, PetscInt its, PetscReal fnorm, void *dummy);
446: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
447: PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP, PetscInt, PetscReal, void *);
448: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP, void **);
449: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void **);

451: PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP, Vec, Vec);
452: PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP, Vec, Vec, Vec, Vec, Vec);

454: PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP, Mat, Mat);
455: PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP, Mat *, Mat *);
456: PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP, PetscBool *, PetscBool *);
457: PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP, const char[]);
458: PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP, const char[]);
459: PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP, const char *[]);

461: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP, PetscBool);
462: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP, PetscBool *);
463: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP, PetscBool);
464: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP, PetscBool *);

466: PETSC_EXTERN PetscErrorCode KSPView(KSP, PetscViewer);
467: PETSC_EXTERN PetscErrorCode KSPLoad(KSP, PetscViewer);
468: PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP, PetscObject, const char[]);
469: PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP, PetscViewer);
470: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP, PetscErrorCode (*)(KSP, void *), void *vctx, PetscErrorCode (*)(void **));
471: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP);
472: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP);
473: PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP, PetscViewer);

475: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonView() (since version 3.14)") static inline PetscErrorCode KSPReasonView(KSP ksp, PetscViewer v) {
476:   return KSPConvergedReasonView(ksp, v);
477: }
478: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonViewFromOptions() (since version 3.14)") static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp) {
479:   return KSPConvergedReasonViewFromOptions(ksp);
480: }

482: #define KSP_FILE_CLASSID 1211223

484: PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP, PetscBool);
485: PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP, PetscBool);
486: PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP, Vec *);
487: PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP, PetscReal *, PetscReal *);
488: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
489: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
490: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);

492: PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC, KSP *);
493: PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC, KSP *);
494: PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC, KSP *);

496: /*E
497:     KSPNormType - Norm that is passed in the Krylov convergence
498:        test routines.

500:    Level: advanced

502:    Each solver only supports a subset of these and some may support different ones
503:    depending on left or right preconditioning, see `KSPSetPCSide()`

505:    Developer Note:
506:     this must match petsc/finclude/petscksp.h

508: .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`,
509:           `KSPSetConvergenceTest()`, `KSPSetPCSide()`
510: E*/
511: typedef enum {
512:   KSP_NORM_DEFAULT          = -1,
513:   KSP_NORM_NONE             = 0,
514:   KSP_NORM_PRECONDITIONED   = 1,
515:   KSP_NORM_UNPRECONDITIONED = 2,
516:   KSP_NORM_NATURAL          = 3
517: } KSPNormType;
518: #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
519: PETSC_EXTERN const char *const *const KSPNormTypes;

521: /*MC
522:     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
523:           possibly save some computation but means the convergence test cannot
524:           be based on a norm of a residual etc.

526:    Level: advanced

528:     Note: Some Krylov methods need to compute a residual norm (such as `KPSGMRES`) and then this option is ignored

530: .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`
531: M*/

533: /*MC
534:     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
535:        convergence test routine.

537:    Level: advanced

539: .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
540: M*/

542: /*MC
543:     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
544:        convergence test routine.

546:    Level: advanced

548: .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
549: M*/

551: /*MC
552:     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
553:        convergence test routine. This is only supported by  `KSPCG`, `KSPCR`, `KSPCGNE`, `KSPCGS`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`

555:    Level: advanced

557: .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSPSetConvergenceTest()`
558: M*/

560: PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP, KSPNormType);
561: PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP, KSPNormType *);
562: PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp, KSPNormType, PCSide, PetscInt);
563: PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP, PetscInt);
564: PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP, PetscBool);

566: #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM("Use KSP_DIVERGED_PC_FAILED (since version 3.11)")
567: /*E
568:     KSPConvergedReason - reason a Krylov method was said to have converged or diverged

570:    Level: beginner

572:    Notes:
573:     See `KSPGetConvergedReason()` for explanation of each value

575:    Developer Notes:
576:    This must match petsc/finclude/petscksp.h

578:    The string versions of these are `KSPConvergedReasons`; if you change
579:    any of the values here also change them that array of names.

581: .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()`
582: E*/
583: typedef enum { /* converged */
584:   KSP_CONVERGED_RTOL_NORMAL              = 1,
585:   KSP_CONVERGED_ATOL_NORMAL              = 9,
586:   KSP_CONVERGED_RTOL                     = 2,
587:   KSP_CONVERGED_ATOL                     = 3,
588:   KSP_CONVERGED_ITS                      = 4,
589:   KSP_CONVERGED_CG_NEG_CURVE             = 5,
590:   KSP_CONVERGED_CG_CONSTRAINED           = 6,
591:   KSP_CONVERGED_STEP_LENGTH              = 7,
592:   KSP_CONVERGED_HAPPY_BREAKDOWN          = 8,
593:   /* diverged */
594:   KSP_DIVERGED_NULL                      = -2,
595:   KSP_DIVERGED_ITS                       = -3,
596:   KSP_DIVERGED_DTOL                      = -4,
597:   KSP_DIVERGED_BREAKDOWN                 = -5,
598:   KSP_DIVERGED_BREAKDOWN_BICG            = -6,
599:   KSP_DIVERGED_NONSYMMETRIC              = -7,
600:   KSP_DIVERGED_INDEFINITE_PC             = -8,
601:   KSP_DIVERGED_NANORINF                  = -9,
602:   KSP_DIVERGED_INDEFINITE_MAT            = -10,
603:   KSP_DIVERGED_PC_FAILED                 = -11,
604:   KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11,

606:   KSP_CONVERGED_ITERATING = 0
607: } KSPConvergedReason;
608: PETSC_EXTERN const char *const *KSPConvergedReasons;

610: /*MC
611:      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if `KSPConvergedDefaultSetUIRNorm()` was called

613:    Level: beginner

615:    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
616:    for left preconditioning it is the 2-norm of the preconditioned residual, and the
617:    2-norm of the residual for right preconditioning

619:    See also `KSP_CONVERGED_ATOL` which may apply before this tolerance.

621: .seealso: `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`

623: M*/

625: /*MC
626:      KSP_CONVERGED_ATOL - norm(r) <= atol

628:    Level: beginner

630:    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
631:        for left preconditioning it is the 2-norm of the preconditioned residual, and the
632:        2-norm of the residual for right preconditioning

634:    See also `KSP_CONVERGED_RTOL` which may apply before this tolerance.

636: .seealso: `KSP_CONVERGED_RTOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`

638: M*/

640: /*MC
641:      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)

643:    Level: beginner

645:    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
646:        for left preconditioning it is the 2-norm of the preconditioned residual, and the
647:        2-norm of the residual for right preconditioning

649:    Level: beginner

651: .seealso: `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_RTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`

653: M*/

655: /*MC
656:      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
657:       reached

659:    Level: beginner

661: .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`

663: M*/

665: /*MC
666:      KSP_CONVERGED_ITS - Used by the `KSPPREONLY` solver after the single iteration of
667:            the preconditioner is applied. Also used when the `KSPConvergedSkip()` convergence
668:            test routine is set in KSP.

670:    Level: beginner

672: .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`

674: M*/

676: /*MC
677:      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
678:           method could not continue to enlarge the Krylov space. Could be due to a singular matrix or
679:           preconditioner. In `KSPHPDDM`, this is also returned when some search directions within a block
680:           are colinear.

682:    Level: beginner

684: .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`

686: M*/

688: /*MC
689:      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the `KSPBICG` method was detected so the
690:           method could not continue to enlarge the Krylov space.

692:    Level: beginner

694: .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`

696: M*/

698: /*MC
699:      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
700:         symmetric and this Krylov method (`KSPCG`, `KSPMINRES`, `KSPCR`) requires symmetry

702:    Level: beginner

704: .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`

706: M*/

708: /*MC
709:      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
710:         positive and negative eigenvalues) and this Krylov method (`KSPCG`) requires it to
711:         be positive definite

713:    Level: beginner

715:      Note:
716:     This can happen with the `PCICC` preconditioner, use -pc_factor_shift_positive_definite to force
717:   the `PCICC` preconditioner to generate a positive definite preconditioner

719: .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`

721: M*/

723: /*MC
724:      KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a
725:      zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner
726:      such as `PCFIELDSPLIT`.

728:    Level: beginner

730:     Notes:
731:     Run with -ksp_error_if_not_converged to stop the program when the error is detected and print an error message with details.

733: .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`

735: M*/

737: /*MC
738:      KSP_CONVERGED_ITERATING - This flag is returned if you call `KSPGetConvergedReason()`
739:         while the `KSPSolve()` is still running.

741:    Level: beginner

743: .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`

745: M*/

747: PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void *, PetscErrorCode (*)(void *));
748: PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *));
749: PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *));
750: PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP, void *);
751: PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
752: PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
753: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *);
754: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **);
755: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
756: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
757: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP, PetscBool);
758: PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
759: PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP, KSPConvergedReason *);
760: PETSC_EXTERN PetscErrorCode KSPGetConvergedReasonString(KSP, const char **);
761: PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP, PetscReal *, PetscReal *, PetscReal *, PetscReal *);

763: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") static inline void KSPDefaultConverged(void) { /* never called */
764: }
765: #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
766: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") static inline void KSPDefaultConvergedDestroy(void) { /* never called */
767: }
768: #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
769: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") static inline void KSPDefaultConvergedCreate(void) { /* never called */
770: }
771: #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
772: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUIRNorm(void) { /* never called */
773: }
774: #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
775: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */
776: }
777: #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
778: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") static inline void KSPSkipConverged(void) { /* never called */
779: }
780: #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)

782: PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP, MatType, Mat *);
783: PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") static inline PetscErrorCode KSPComputeExplicitOperator(KSP A, Mat *B) {
784:   return KSPComputeOperator(A, NULL, B);
785: }

787: /*E
788:     KSPCGType - Determines what type of CG to use

790:    Level: beginner

792: .seealso: `KSPCGSetType()`
793: E*/
794: typedef enum {
795:   KSP_CG_SYMMETRIC = 0,
796:   KSP_CG_HERMITIAN = 1
797: } KSPCGType;
798: PETSC_EXTERN const char *const KSPCGTypes[];

800: PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP, KSPCGType);
801: PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP, PetscBool);

803: PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP, PetscReal);
804: PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP, PetscReal *);
805: PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP, PetscReal *);

807: PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP, PetscReal *);
808: PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP, PetscReal *);
809: PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp, PetscReal *x) {
810:   return KSPGLTRGetMinEig(ksp, x);
811: }
812: PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetLambda(KSP ksp, PetscReal *x) {
813:   return KSPGLTRGetLambda(ksp, x);
814: }

816: PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP, const char[]);
817: PETSC_EXTERN PetscErrorCode KSPPythonGetType(KSP, const char *[]);

819: PETSC_EXTERN PetscErrorCode PCSetPreSolve(PC, PetscErrorCode (*)(PC, KSP));
820: PETSC_EXTERN PetscErrorCode PCPreSolve(PC, KSP);
821: PETSC_EXTERN PetscErrorCode PCPostSolve(PC, KSP);

823: #include <petscdrawtypes.h>
824: PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP, PetscInt, PetscReal, void *);

826: PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec));
827: PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec));

829: /*S
830:      KSPGuess - Abstract PETSc object that manages all initial guess methods in Krylov methods.

832:    Level: beginner

834: .seealso: `KSPCreate()`, `KSPSetGuessType()`, `KSPGuessType`
835: S*/
836: typedef struct _p_KSPGuess *KSPGuess;
837: /*J
838:     KSPGuessType - String with the name of a PETSc initial guess approach for Krylov methods.

840:    Level: beginner

842: .seealso: `KSPGuess`
843: J*/
844: typedef const char         *KSPGuessType;
845: #define KSPGUESSFISCHER "fischer"
846: #define KSPGUESSPOD     "pod"
847: PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[], PetscErrorCode (*)(KSPGuess));
848: PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP, KSPGuess);
849: PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP, KSPGuess *);
850: PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess, PetscViewer);
851: PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess *);
852: PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm, KSPGuess *);
853: PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess, KSPGuessType);
854: PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess, KSPGuessType *);
855: PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess, PetscReal);
856: PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
857: PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess, Vec, Vec);
858: PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess, Vec, Vec);
859: PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
860: PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess, PetscInt, PetscInt);
861: PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP, PetscInt, PetscInt);
862: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP, PetscBool);
863: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP, PetscBool *);

865: /*E
866:     MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines

868:     Level: intermediate

870: .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`, `MatCreateSchurComplementPmat()`
871: E*/
872: typedef enum {
873:   MAT_SCHUR_COMPLEMENT_AINV_DIAG,
874:   MAT_SCHUR_COMPLEMENT_AINV_LUMP,
875:   MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG,
876:   MAT_SCHUR_COMPLEMENT_AINV_FULL
877: } MatSchurComplementAinvType;
878: PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];

880: typedef enum {
881:   MAT_LMVM_SYMBROYDEN_SCALE_NONE     = 0,
882:   MAT_LMVM_SYMBROYDEN_SCALE_SCALAR   = 1,
883:   MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2,
884:   MAT_LMVM_SYMBROYDEN_SCALE_USER     = 3
885: } MatLMVMSymBroydenScaleType;
886: PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[];

888: PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat, Mat, Mat, Mat, Mat, Mat *);
889: PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat, KSP *);
890: PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat, KSP);
891: PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
892: PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
893: PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat, Mat *, Mat *, Mat *, Mat *, Mat *);
894: PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat, MatSchurComplementAinvType);
895: PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat, MatSchurComplementAinvType *);
896: PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat, MatReuse, Mat *);
897: PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat, Mat *);
898: PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatSchurComplementAinvType, MatReuse, Mat *);
899: PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat, Mat, Mat, Mat, MatSchurComplementAinvType, MatReuse, Mat *);

901: PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm, PetscInt, PetscInt, Mat *);
902: PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm, PetscInt, PetscInt, Mat *);
903: PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm, PetscInt, PetscInt, Mat *);
904: PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
905: PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
906: PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
907: PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
908: PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);

910: PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
911: PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool *);
912: PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
913: PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
914: PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
915: PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
916: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
917: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
918: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
919: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
920: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
921: PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
922: PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
923: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat *);
924: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC *);
925: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP *);
926: PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt);
927: PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt *);
928: PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt *);
929: PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar);
930: PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType);

932: PETSC_EXTERN PetscErrorCode KSPSetDM(KSP, DM);
933: PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP, PetscBool);
934: PETSC_EXTERN PetscErrorCode KSPGetDM(KSP, DM *);
935: PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP, void *);
936: PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP, void *);
937: PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP, PetscErrorCode (*func)(KSP, Vec, void *), void *);
938: PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP, PetscErrorCode (*)(KSP, Mat, Mat, void *), void *);
939: PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP, PetscErrorCode (*)(KSP, Vec, void *), void *);
940: PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM, PetscErrorCode (*)(KSP, Mat, Mat, void *), void *);
941: PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM, PetscErrorCode (**)(KSP, Mat, Mat, void *), void *);
942: PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM, PetscErrorCode (*)(KSP, Vec, void *), void *);
943: PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM, PetscErrorCode (**)(KSP, Vec, void *), void *);
944: PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM, PetscErrorCode (*)(KSP, Vec, void *), void *);
945: PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM, PetscErrorCode (**)(KSP, Vec, void *), void *);

947: PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM, Vec, Vec);
948: PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec, void (**)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), InsertMode, Vec);

950: PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *);
951: PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal);
952: #endif