Actual source code: petscksp.h

  1: /*
  2:    Defines the interface functions for the Krylov subspace accelerators.
  3: */
  4: #pragma once

  6: #include <petscpc.h>

  8: /* SUBMANSEC = KSP */

 10: PETSC_EXTERN PetscErrorCode KSPInitializePackage(void);
 11: PETSC_EXTERN PetscErrorCode KSPFinalizePackage(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 factorization-based 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` (or equivalently `KSPNONE`), meaning that only application of the preconditioner is used as the linear solver.

 23: .seealso: [](doc_linsolve), [](ch_ksp), `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: [](doc_linsolve), [](ch_ksp), `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(3, 11, 0, "KSPNASH", ) "nash"
 49: #define KSPCGSTCG     PETSC_DEPRECATED_MACRO(3, 11, 0, "KSPSTCG", ) "stcg"
 50: #define KSPCGGLTR     PETSC_DEPRECATED_MACRO(3, 11, 0, "KSPSGLTR", ) "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 KSPMatSolveTranspose(KSP, Mat, Mat);
102: PETSC_EXTERN PetscErrorCode KSPSetMatSolveBatchSize(KSP, PetscInt);
103: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPSetMatSolveBatchSize()", ) static inline PetscErrorCode KSPSetMatSolveBlockSize(KSP ksp, PetscInt n)
104: {
105:   return KSPSetMatSolveBatchSize(ksp, n);
106: }
107: PETSC_EXTERN PetscErrorCode KSPGetMatSolveBatchSize(KSP, PetscInt *);
108: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPGetMatSolveBatchSize()", ) static inline PetscErrorCode KSPGetMatSolveBlockSize(KSP ksp, PetscInt *n)
109: {
110:   return KSPGetMatSolveBatchSize(ksp, n);
111: }
112: PETSC_EXTERN PetscErrorCode KSPReset(KSP);
113: PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP);
114: PETSC_EXTERN PetscErrorCode KSPDestroy(KSP *);
115: PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP, PetscBool);
116: PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP, PetscBool *);
117: PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP, PetscBool);
118: PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP, PC, Vec);

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

128: PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP, PCSide);
129: PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP, PCSide *);
130: PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP, PetscReal, PetscReal, PetscReal, PetscInt);
131: PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP, PetscReal *, PetscReal *, PetscReal *, PetscInt *);
132: PETSC_EXTERN PetscErrorCode KSPSetMinimumIterations(KSP, PetscInt);
133: PETSC_EXTERN PetscErrorCode KSPGetMinimumIterations(KSP, PetscInt *);
134: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP, PetscBool);
135: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP, PetscBool *);
136: PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP, PetscBool);
137: PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP, PetscBool *);
138: PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP, PetscBool);
139: PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP, PetscBool);
140: PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP, PetscBool *);
141: PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP, PetscBool);
142: PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP, PetscBool *);
143: PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP, Vec *);
144: PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP, Vec *);
145: PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP, PetscReal *);
146: PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP, PetscInt *);
147: PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP, PetscInt *);
148: PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP, PetscInt, Vec **, PetscInt, Vec **);
149: PETSC_DEPRECATED_FUNCTION(3, 6, 0, "KSPCreateVecs()", ) static inline PetscErrorCode KSPGetVecs(KSP ksp, PetscInt n, Vec **x, PetscInt m, Vec **y)
150: {
151:   return KSPCreateVecs(ksp, n, x, m, y);
152: }

154: PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *);
155: PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *);

157: PETSC_EXTERN PetscErrorCode KSPSetPC(KSP, PC);
158: PETSC_EXTERN PetscErrorCode KSPGetPC(KSP, PC *);
159: PETSC_EXTERN PetscErrorCode KSPSetNestLevel(KSP, PetscInt);
160: PETSC_EXTERN PetscErrorCode KSPGetNestLevel(KSP, PetscInt *);

162: PETSC_EXTERN PetscErrorCode KSPMonitor(KSP, PetscInt, PetscReal);
163: PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void **));
164: PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
165: PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP, void *);
166: PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP, const PetscReal *[], PetscInt *);
167: PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP, PetscReal[], PetscCount, PetscBool);
168: PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP, const PetscReal *[], PetscInt *);
169: PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP, PetscReal[], PetscCount, PetscBool);

171: PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP, Vec, Vec *);
172: PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP, Vec, Vec, Vec *);
173: PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
174: PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP, PetscInt);

176: PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC, KSP *);
177: PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC, KSP);
178: PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
179: PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
180: PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]);
181: PETSC_EXTERN PetscErrorCode PCPatchGetSubKSP(PC, PetscInt *, KSP *[]);
182: PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC, PetscInt *, KSP *[]);
183: PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC, PetscInt *, KSP *[]);
184: PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC, PetscInt, KSP *);
185: PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC, PetscInt, KSP *);
186: PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC, PetscInt, KSP *);
187: PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC, KSP *);
188: PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC, KSP *);
189: PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC, KSP *);
190: /*
191:   PCMGCoarseList contains the list of coarse space constructor currently registered
192:   These are added with PCMGRegisterCoarseSpaceConstructor()
193: */
194: PETSC_EXTERN PetscFunctionList PCMGCoarseList;
195: PETSC_EXTERN PetscErrorCode    PCMGRegisterCoarseSpaceConstructor(const char[], PetscErrorCode (*)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *));
196: PETSC_EXTERN PetscErrorCode    PCMGGetCoarseSpaceConstructor(const char[], PetscErrorCode (**)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *));

198: PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP, Vec, Vec *);
199: PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP, Vec, Vec, Vec *);

201: /*E
202:   KSPChebyshevKind - Which kind of Chebyshev polynomial to use with `KSPCHEBYSHEV`

204:   Values:
205: + `KSP_CHEBYSHEV_FIRST`      - "classic" first-kind Chebyshev polynomial
206: . `KSP_CHEBYSHEV_FOURTH`     - fourth-kind Chebyshev polynomial
207: - `KSP_CHEBYSHEV_OPT_FOURTH` - optimized fourth-kind Chebyshev polynomial

209:   Level: intermediate

211: .seealso: [](ch_ksp), `KSPCHEBYSHEV`, `KSPChebyshevSetKind()`, `KSPChebyshevGetKind()`
212: E*/
213: typedef enum {
214:   KSP_CHEBYSHEV_FIRST,
215:   KSP_CHEBYSHEV_FOURTH,
216:   KSP_CHEBYSHEV_OPT_FOURTH
217: } KSPChebyshevKind;

219: PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP, PetscReal);
220: PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP, PetscBool);
221: PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP, PetscReal, PetscReal);
222: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP, PetscReal, PetscReal, PetscReal, PetscReal);
223: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP, PetscBool);
224: PETSC_EXTERN PetscErrorCode KSPChebyshevSetKind(KSP, KSPChebyshevKind);
225: PETSC_EXTERN PetscErrorCode KSPChebyshevGetKind(KSP, KSPChebyshevKind *);
226: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP, KSP *);
227: PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP, PetscReal *, PetscReal *);
228: PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP, PetscInt, PetscReal[], PetscReal[], PetscInt *);
229: PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP, PetscInt, PetscReal[], PetscReal[]);
230: PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP, PetscBool, PetscBool, PetscInt *, Vec[], PetscReal[], PetscReal[]);

232: /*E

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

236:   Values:
237: + `KSP_FCD_TRUNC_TYPE_STANDARD` - uses all (up to mmax) stored directions
238: - `KSP_FCD_TRUNC_TYPE_NOTAY`    - uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..

240:    Level: intermediate

242: .seealso: [](ch_ksp), `KSP`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`, `KSPFCGSetTruncationType()`, `KSPFCGGetTruncationType()`
243: E*/
244: typedef enum {
245:   KSP_FCD_TRUNC_TYPE_STANDARD,
246:   KSP_FCD_TRUNC_TYPE_NOTAY
247: } KSPFCDTruncationType;
248: PETSC_EXTERN const char *const KSPFCDTruncationTypes[];

250: PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP, PetscInt);
251: PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP, PetscInt *);
252: PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP, PetscInt);
253: PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP, PetscInt *);
254: PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP, KSPFCDTruncationType);
255: PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP, KSPFCDTruncationType *);

257: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP, PetscInt);
258: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP, PetscInt *);
259: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP, PetscInt);
260: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP, PetscInt *);
261: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP, KSPFCDTruncationType);
262: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP, KSPFCDTruncationType *);

264: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP, PetscInt);
265: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP, PetscInt *);
266: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP, PetscInt);
267: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP, PetscInt *);
268: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP, KSPFCDTruncationType);
269: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP, KSPFCDTruncationType *);
270: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP, PetscBool);
271: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP, PetscBool *);
272: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *));

274: PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
275: PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt *);
276: PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP, PetscReal);
277: PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP, PetscReal);

279: PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
280: PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP, PetscErrorCode (*)(KSP, PetscInt));
281: PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP, PetscErrorCode (**)(KSP, PetscInt));
282: PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP, PetscInt);
283: PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP, PetscInt);

285: PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP, PetscInt);
286: PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);

288: PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP, PetscScalar);

290: PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP, PetscInt);
291: PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP, PetscInt *);
292: PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *));

294: PETSC_EXTERN PetscErrorCode KSPMINRESSetRadius(KSP, PetscReal);
295: PETSC_EXTERN PetscErrorCode KSPMINRESGetUseQLP(KSP, PetscBool *);
296: PETSC_EXTERN PetscErrorCode KSPMINRESSetUseQLP(KSP, PetscBool);

298: PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP, PC *);
299: PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP, PC);
300: PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP, KSP *);
301: PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP, Mat);

303: PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationMat(KSP, Mat);
304: PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationMat(KSP, Mat *);
305: #if PetscDefined(HAVE_HPDDM)
306: PETSC_DEPRECATED_FUNCTION(3, 18, 0, "KSPHPDDMSetDeflationMat()", ) static inline PetscErrorCode KSPHPDDMSetDeflationSpace(KSP ksp, Mat U)
307: {
308:   return KSPHPDDMSetDeflationMat(ksp, U);
309: }
310: PETSC_DEPRECATED_FUNCTION(3, 18, 0, "KSPHPDDMGetDeflationMat()", ) static inline PetscErrorCode KSPHPDDMGetDeflationSpace(KSP ksp, Mat *U)
311: {
312:   return KSPHPDDMGetDeflationMat(ksp, U);
313: }
314: #endif
315: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPMatSolve()", ) static inline PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X)
316: {
317:   return KSPMatSolve(ksp, B, X);
318: }
319: /*E
320:     KSPHPDDMType - Type of Krylov method used by `KSPHPDDM`

322:     Values:
323: +   `KSP_HPDDM_TYPE_GMRES` (default) - Generalized Minimal Residual method
324: .   `KSP_HPDDM_TYPE_BGMRES`          - block GMRES
325: .   `KSP_HPDDM_TYPE_CG`              - Conjugate Gradient
326: .   `KSP_HPDDM_TYPE_BCG`             - block CG
327: .   `KSP_HPDDM_TYPE_GCRODR`          - Generalized Conjugate Residual method with inner Orthogonalization and Deflated Restarting
328: .   `KSP_HPDDM_TYPE_BGCRODR`         - block GCRODR
329: .   `KSP_HPDDM_TYPE_BFBCG`           - breakdown-free BCG
330: -   `KSP_HPDDM_TYPE_PREONLY`         - apply the preconditioner only

332:     Level: intermediate

334: .seealso: [](ch_ksp), `KSPHPDDM`, `KSPHPDDMSetType()`
335: E*/
336: typedef enum {
337:   KSP_HPDDM_TYPE_GMRES   = 0,
338:   KSP_HPDDM_TYPE_BGMRES  = 1,
339:   KSP_HPDDM_TYPE_CG      = 2,
340:   KSP_HPDDM_TYPE_BCG     = 3,
341:   KSP_HPDDM_TYPE_GCRODR  = 4,
342:   KSP_HPDDM_TYPE_BGCRODR = 5,
343:   KSP_HPDDM_TYPE_BFBCG   = 6,
344:   KSP_HPDDM_TYPE_PREONLY = 7
345: } KSPHPDDMType;
346: PETSC_EXTERN const char *const KSPHPDDMTypes[];

348: /*E
349:     KSPHPDDMPrecision - Precision of Krylov bases used by `KSPHPDDM`

351:     Values:
352: +   `KSP_HPDDM_PRECISION_HALF`      - default when PETSc is configured `--with-precision=__fp16`
353: .   `KSP_HPDDM_PRECISION_SINGLE`    - default when PETSc is configured `--with-precision=single`
354: .   `KSP_HPDDM_PRECISION_DOUBLE`    - default when PETSc is configured `--with-precision=double`
355: -   `KSP_HPDDM_PRECISION_QUADRUPLE` - default when PETSc is configured `--with-precision=__float128`

357:     Level: intermediate

359: .seealso: [](ch_ksp), `KSP`, `KSPHPDDM`
360: E*/
361: typedef enum {
362:   KSP_HPDDM_PRECISION_HALF      = 0,
363:   KSP_HPDDM_PRECISION_SINGLE    = 1,
364:   KSP_HPDDM_PRECISION_DOUBLE    = 2,
365:   KSP_HPDDM_PRECISION_QUADRUPLE = 3
366: } KSPHPDDMPrecision;
367: PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP, KSPHPDDMType);
368: PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP, KSPHPDDMType *);

370: /*E
371:    KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed in the GMRES solvers

373:    Values:
374: +  `KSP_GMRES_CGS_REFINE_NEVER`    - one step of classical Gram-Schmidt
375: .  `KSP_GMRES_CGS_REFINE_IFNEEDED` - a second step is performed if the first step does not satisfy some criteria
376: -  `KSP_GMRES_CGS_REFINE_ALWAYS`   - always perform two steps

378:    Level: advanced

380: .seealso: [](ch_ksp), `KSP`, `KSPGMRES`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
381:           `KSPGMRESGetOrthogonalization()`,
382:           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`
383: E*/
384: typedef enum {
385:   KSP_GMRES_CGS_REFINE_NEVER,
386:   KSP_GMRES_CGS_REFINE_IFNEEDED,
387:   KSP_GMRES_CGS_REFINE_ALWAYS
388: } KSPGMRESCGSRefinementType;
389: PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];

391: /*MC
392:    KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process

394:    Level: advanced

396:    Note:
397:    Possibly unstable, but the fastest to compute

399: .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
400:           `KSP`, `KSPGMRESGetOrthogonalization()`,
401:           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
402:           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
403: M*/

405: /*MC
406:     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
407:           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
408:           poor orthogonality.

410:    Level: advanced

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

416: .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
417:           `KSP`, `KSPGMRESGetOrthogonalization()`,
418:           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
419:           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
420: M*/

422: /*MC
423:    KSP_GMRES_CGS_REFINE_ALWAYS - Do two steps of the classical (unmodified) Gram-Schmidt process.

425:    Level: advanced

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

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

433: .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
434:           `KSP`, `KSPGMRESGetOrthogonalization()`,
435:           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
436:           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
437: M*/

439: PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP, KSPGMRESCGSRefinementType);
440: PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP, KSPGMRESCGSRefinementType *);

442: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP, PetscInt, PetscInt, PetscReal, void *);
443: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP, PetscInt, PetscInt, PetscReal, void *);
444: PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *));

446: PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP, PetscReal);
447: PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP, PetscReal *);
448: PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP, PetscReal *);

450: PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP, PetscReal);
451: PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP, PetscBool);
452: PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP, PetscInt);
453: PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP, PetscBool);

455: PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
456: PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);

458: PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP, const char[], const char[], void *);
459: PETSC_EXTERN PetscErrorCode KSPMonitorLGCreate(MPI_Comm, const char[], const char[], const char[], PetscInt, const char *[], int, int, int, int, PetscDrawLG *);
460: PETSC_EXTERN PetscErrorCode KSPMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
461: PETSC_EXTERN PetscErrorCode KSPMonitorResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
462: PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
463: PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
464: PETSC_EXTERN PetscErrorCode KSPMonitorResidualShort(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
465: PETSC_EXTERN PetscErrorCode KSPMonitorResidualRange(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
466: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
467: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
468: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
469: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
470: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMax(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
471: PETSC_EXTERN PetscErrorCode KSPMonitorError(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
472: PETSC_EXTERN PetscErrorCode KSPMonitorErrorDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
473: PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
474: PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
475: PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
476: PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
477: PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
478: PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
479: PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
480: PETSC_EXTERN PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
481: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorResidual()", ) static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
482: {
483:   return KSPMonitorResidual(ksp, n, rnorm, vf);
484: }
485: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorTrueResidual()", ) static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
486: {
487:   return KSPMonitorTrueResidual(ksp, n, rnorm, vf);
488: }
489: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorTrueResidualMax()", ) static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
490: {
491:   return KSPMonitorTrueResidualMax(ksp, n, rnorm, vf);
492: }

494: PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP, PetscInt, PetscReal, void *);
495: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP, PetscInt, PetscReal, void *);
496: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **);
497: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceCreate(void *);
498: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceSetCoefficient(void *, PetscReal);
499: PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP, PetscInt, PetscReal, void *);
500: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP, void **);
501: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void **);

503: PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP, Vec, Vec);
504: PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP, Vec, Vec, Vec, Vec, Vec);

506: PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP, Mat, Mat);
507: PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP, Mat *, Mat *);
508: PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP, PetscBool *, PetscBool *);
509: PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP, const char[]);
510: PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP, const char[]);
511: PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP, const char *[]);

513: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP, PetscBool);
514: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP, PetscBool *);
515: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP, PetscBool);
516: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP, PetscBool *);

518: PETSC_EXTERN PetscErrorCode KSPView(KSP, PetscViewer);
519: PETSC_EXTERN PetscErrorCode KSPLoad(KSP, PetscViewer);
520: PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP, PetscObject, const char[]);
521: PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP, PetscViewer);
522: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP, PetscErrorCode (*)(KSP, void *), void *vctx, PetscErrorCode (*)(void **));
523: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP);
524: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP);
525: PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP, PetscViewer);

527: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonView()", ) static inline PetscErrorCode KSPReasonView(KSP ksp, PetscViewer v)
528: {
529:   return KSPConvergedReasonView(ksp, v);
530: }
531: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonViewFromOptions()", ) static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp)
532: {
533:   return KSPConvergedReasonViewFromOptions(ksp);
534: }

536: #define KSP_FILE_CLASSID 1211223

538: PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP, PetscBool);
539: PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP, PetscBool);
540: PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP, Vec *);
541: PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP, PetscReal *, PetscReal *);
542: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
543: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
544: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);

546: PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC, KSP *);
547: PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC, KSP *);
548: PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC, KSP *);
549: PETSC_EXTERN PetscErrorCode PCMPIGetKSP(PC, KSP *);

551: /*E
552:    KSPNormType - Norm calculated by the `KSP` and passed in the Krylov convergence
553:        test routines.

555:    Values:
556: +  `KSP_NORM_DEFAULT`          - use the default for the current `KSPType`
557: .  `KSP_NORM_NONE`             - use no norm calculation
558: .  `KSP_NORM_PRECONDITIONED`   - use the preconditioned residual norm
559: .  `KSP_NORM_UNPRECONDITIONED` - use the unpreconditioned residual norm
560: -  `KSP_NORM_NATURAL`          - use the natural norm (the norm induced by the linear operator)

562:    Level: advanced

564:    Note:
565:    Each solver only supports a subset of these and some may support different ones
566:    depending on whether left or right preconditioning is used, see `KSPSetPCSide()`

568: .seealso: [](ch_ksp), `KSP`, `PCSide`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`,
569:           `KSPSetConvergenceTest()`, `KSPSetPCSide()`, `KSP_NORM_DEFAULT`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`
570: E*/
571: typedef enum {
572:   KSP_NORM_DEFAULT          = -1,
573:   KSP_NORM_NONE             = 0,
574:   KSP_NORM_PRECONDITIONED   = 1,
575:   KSP_NORM_UNPRECONDITIONED = 2,
576:   KSP_NORM_NATURAL          = 3
577: } KSPNormType;
578: #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
579: PETSC_EXTERN const char *const *const KSPNormTypes;

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

586:    Level: advanced

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

591: .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`
592: M*/

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

598:    Level: advanced

600: .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
601: M*/

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

607:    Level: advanced

609: .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()`
610: M*/

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

616:    Level: advanced

618: .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSPSetConvergenceTest()`
619: M*/

621: PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP, KSPNormType);
622: PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP, KSPNormType *);
623: PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp, KSPNormType, PCSide, PetscInt);
624: PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP, PetscInt);
625: PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP, PetscBool);

627: #define KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED   KSP_CONVERGED_CG_NEG_CURVE PETSC_DEPRECATED_ENUM(3, 19, 0, "KSP_CONVERGED_NEG_CURVE", )
628: #define KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED KSP_CONVERGED_CG_CONSTRAINED PETSC_DEPRECATED_ENUM(3, 19, 0, "KSP_CONVERGED_STEP_LENGTH", )
629: #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED  KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM(3, 11, 0, "KSP_DIVERGED_PC_FAILED", )
630: /*E
631:    KSPConvergedReason - reason a Krylov method was determined to have converged or diverged

633:    Values:
634: +  `KSP_CONVERGED_RTOL_NORMAL`     - requested decrease in the residual for the normal equations
635: .  `KSP_CONVERGED_ATOL_NORMAL`     - requested absolute value in the residual for the normal equations
636: .  `KSP_CONVERGED_RTOL`            - requested decrease in the residual
637: .  `KSP_CONVERGED_ATOL`            - requested absolute value in the residual
638: .  `KSP_CONVERGED_ITS`             - requested number of iterations
639: .  `KSP_CONVERGED_NEG_CURVE`       - see note below
640: .  `KSP_CONVERGED_STEP_LENGTH`     - see note below
641: .  `KSP_CONVERGED_HAPPY_BREAKDOWN` - happy breakdown (meaning early convergence of the `KSPType` occurred).
642: .  `KSP_DIVERGED_NULL`             - breakdown when solving the Hessenberg system within GMRES
643: .  `KSP_DIVERGED_ITS`              - requested number of iterations
644: .  `KSP_DIVERGED_DTOL`             - large increase in the residual norm
645: .  `KSP_DIVERGED_BREAKDOWN`        - breakdown in the Krylov method
646: .  `KSP_DIVERGED_BREAKDOWN_BICG`   - breakdown in the `KSPBGCS` Krylov method
647: .  `KSP_DIVERGED_NONSYMMETRIC`     - the operator or preonditioner was not symmetric for a `KSPType` that requires symmetry
648: .  `KSP_DIVERGED_INDEFINITE_PC`    - the preconditioner was indefinite for a `KSPType` that requires it be definite
649: .  `KSP_DIVERGED_NANORINF`         - a not a number of infinity was detected in a vector during the computation
650: .  `KSP_DIVERGED_INDEFINITE_MAT`   - the operator was indefinite for a `KSPType` that requires it be definite
651: -  `KSP_DIVERGED_PC_FAILED`        - the action of the preconditioner failed for some reason

653:    Level: beginner

655:    Note:
656:    The values `KSP_CONVERGED_NEG_CURVE`, and `KSP_CONVERGED_STEP_LENGTH` are returned only by `KSPCG`, `KSPMINRES` and by
657:    the special `KSPNASH`, `KSPSTCG`, and `KSPGLTR` solvers which are used by the `SNESNEWTONTR` (trust region) solver.

659:    Developer Note:
660:    The string versions of these are `KSPConvergedReasons`; if you change
661:    any of the values here also change them that array of names.

663: .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()`
664: E*/
665: typedef enum { /* converged */
666:   KSP_CONVERGED_RTOL_NORMAL               = 1,
667:   KSP_CONVERGED_ATOL_NORMAL               = 9,
668:   KSP_CONVERGED_RTOL                      = 2,
669:   KSP_CONVERGED_ATOL                      = 3,
670:   KSP_CONVERGED_ITS                       = 4,
671:   KSP_CONVERGED_NEG_CURVE                 = 5,
672:   KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED   = 5,
673:   KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED = 6,
674:   KSP_CONVERGED_STEP_LENGTH               = 6,
675:   KSP_CONVERGED_HAPPY_BREAKDOWN           = 7,
676:   /* diverged */
677:   KSP_DIVERGED_NULL                      = -2,
678:   KSP_DIVERGED_ITS                       = -3,
679:   KSP_DIVERGED_DTOL                      = -4,
680:   KSP_DIVERGED_BREAKDOWN                 = -5,
681:   KSP_DIVERGED_BREAKDOWN_BICG            = -6,
682:   KSP_DIVERGED_NONSYMMETRIC              = -7,
683:   KSP_DIVERGED_INDEFINITE_PC             = -8,
684:   KSP_DIVERGED_NANORINF                  = -9,
685:   KSP_DIVERGED_INDEFINITE_MAT            = -10,
686:   KSP_DIVERGED_PC_FAILED                 = -11,
687:   KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11,

689:   KSP_CONVERGED_ITERATING = 0
690: } KSPConvergedReason;
691: PETSC_EXTERN const char *const *KSPConvergedReasons;

693: /*MC
694:    KSP_CONVERGED_RTOL - $||r|| \le rtol*||b||$ or $rtol*||b - A*x_0||$ if `KSPConvergedDefaultSetUIRNorm()` was called

696:    Level: beginner

698:    Notes:
699:    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
700:    for left preconditioning it is the 2-norm of the preconditioned residual, and the
701:    2-norm of the residual for right preconditioning

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

705: .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
706: M*/

708: /*MC
709:    KSP_CONVERGED_ATOL - $||r|| \le atol$

711:    Level: beginner

713:    Notes:
714:    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
715:    for left preconditioning it is the 2-norm of the preconditioned residual, and the
716:    2-norm of the residual for right preconditioning

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

720: .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_RTOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
721: M*/

723: /*MC
724:    KSP_DIVERGED_DTOL - $||r|| \ge dtol*||b||$

726:    Level: beginner

728:    Note:
729:    See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default
730:    for left preconditioning it is the 2-norm of the preconditioned residual, and the
731:    2-norm of the residual for right preconditioning

733: .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_CONVERGED_RTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
734: M*/

736: /*MC
737:    KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
738:    reached

740:    Level: beginner

742: .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
743: M*/

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

750:    Level: beginner

752: .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
753: M*/

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

761:    Level: beginner

763: .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
764: M*/

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

770:    Level: beginner

772: .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
773: M*/

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

779:    Level: beginner

781: .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
782: M*/

784: /*MC
785:    KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
786:    positive and negative eigenvalues) and this Krylov method (`KSPCG`) requires it to
787:    be positive definite

789:    Level: beginner

791:    Note:
792:    This can happen with the `PCICC` preconditioner, use the options database option `-pc_factor_shift_positive_definite` to force
793:    the `PCICC` preconditioner to generate a positive definite preconditioner

795: .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
796: M*/

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

803:    Level: beginner

805:    Note:
806:    Run with `-ksp_error_if_not_converged` to stop the program when the error is detected and print an error message with details.

808: .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
809: M*/

811: /*MC
812:    KSP_CONVERGED_ITERATING - This flag is returned if `KSPGetConvergedReason()` is called
813:    while `KSPSolve()` is still running.

815:    Level: beginner

817: .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()`
818: M*/

820: PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void *, PetscErrorCode (*)(void *));
821: PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *));
822: PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *));
823: PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP, void *);
824: PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
825: PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
826: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *);
827: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **);
828: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
829: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
830: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP, PetscBool);
831: PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
832: PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP, KSPConvergedReason *);
833: PETSC_EXTERN PetscErrorCode KSPGetConvergedReasonString(KSP, const char **);
834: PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
835: PETSC_EXTERN PetscErrorCode KSPSetConvergedNegativeCurvature(KSP, PetscBool);
836: PETSC_EXTERN PetscErrorCode KSPGetConvergedNegativeCurvature(KSP, PetscBool *);

838: PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefault()", ) static inline void KSPDefaultConverged(void)
839: { /* never called */
840: }
841: #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
842: PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultDestroy()", ) static inline void KSPDefaultConvergedDestroy(void)
843: { /* never called */
844: }
845: #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
846: PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultCreate()", ) static inline void KSPDefaultConvergedCreate(void)
847: { /* never called */
848: }
849: #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
850: PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultSetUIRNorm()", ) static inline void KSPDefaultConvergedSetUIRNorm(void)
851: { /* never called */
852: }
853: #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
854: PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultSetUMIRNorm()", ) static inline void KSPDefaultConvergedSetUMIRNorm(void)
855: { /* never called */
856: }
857: #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
858: PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedSkip()", ) static inline void KSPSkipConverged(void)
859: { /* never called */
860: }
861: #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)

863: PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP, MatType, Mat *);
864: PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPComputeOperator()", ) static inline PetscErrorCode KSPComputeExplicitOperator(KSP A, Mat *B)
865: {
866:   return KSPComputeOperator(A, PETSC_NULLPTR, B);
867: }

869: /*E
870:    KSPCGType - Determines what type of `KSPCG` to use

872:    Values:
873:  + `KSP_CG_SYMMETRIC` - the matrix is complex symmetric
874:  - `KSP_CG_HERMITIAN` - the matrix is complex Hermitian

876:    Level: beginner

878: .seealso: [](ch_ksp), `KSPCG`, `KSP`, `KSPCGSetType()`
879: E*/
880: typedef enum {
881:   KSP_CG_SYMMETRIC = 0,
882:   KSP_CG_HERMITIAN = 1
883: } KSPCGType;
884: PETSC_EXTERN const char *const KSPCGTypes[];

886: PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP, KSPCGType);
887: PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP, PetscBool);

889: PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP, PetscReal);
890: PETSC_EXTERN PetscErrorCode KSPCGSetObjectiveTarget(KSP, PetscReal);
891: PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP, PetscReal *);
892: PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP, PetscReal *);

894: PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP, PetscReal *);
895: PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP, PetscReal *);
896: PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPGLTRGetMinEig()", ) static inline PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp, PetscReal *x)
897: {
898:   return KSPGLTRGetMinEig(ksp, x);
899: }
900: PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPGLTRGetLambda()", ) static inline PetscErrorCode KSPCGGLTRGetLambda(KSP ksp, PetscReal *x)
901: {
902:   return KSPGLTRGetLambda(ksp, x);
903: }

905: PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP, const char[]);
906: PETSC_EXTERN PetscErrorCode KSPPythonGetType(KSP, const char *[]);

908: PETSC_EXTERN PetscErrorCode PCSetPreSolve(PC, PetscErrorCode (*)(PC, KSP));
909: PETSC_EXTERN PetscErrorCode PCPreSolve(PC, KSP);
910: PETSC_EXTERN PetscErrorCode PCPostSolve(PC, KSP);

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

915: PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec));
916: PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec));

918: /*S
919:    KSPGuess - Abstract PETSc object that manages all initial guess generation methods for Krylov methods.

921:    Level: intermediate

923: .seealso: [](ch_ksp), `KSPCreate()`, `KSPGuessSetType()`, `KSPGuessType`
924: S*/
925: typedef struct _p_KSPGuess *KSPGuess;

927: /*J
928:    KSPGuessType - String with the name of a PETSc initial guess approach for Krylov methods.

930:    Values:
931:  + `KSPGUESSFISCHER` - methodology developed by Paul Fischer
932:  - `KSPGUESSPOD`     - methodology based on proper orthogonal decomposition

934:    Level: intermediate

936: .seealso: [](ch_ksp), `KSP`, `KSPGuess`
937: J*/
938: typedef const char *KSPGuessType;
939: #define KSPGUESSFISCHER "fischer"
940: #define KSPGUESSPOD     "pod"

942: PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[], PetscErrorCode (*)(KSPGuess));
943: PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP, KSPGuess);
944: PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP, KSPGuess *);
945: PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess, PetscViewer);
946: PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess *);
947: PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm, KSPGuess *);
948: PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess, KSPGuessType);
949: PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess, KSPGuessType *);
950: PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess, PetscReal);
951: PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
952: PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess, Vec, Vec);
953: PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess, Vec, Vec);
954: PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
955: PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess, PetscInt, PetscInt);
956: PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP, PetscInt, PetscInt);
957: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP, PetscBool);
958: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP, PetscBool *);

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

963:     Level: intermediate

965: .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`,
966:           `MatCreateSchurComplementPmat()`
967: E*/
968: typedef enum {
969:   MAT_SCHUR_COMPLEMENT_AINV_DIAG,
970:   MAT_SCHUR_COMPLEMENT_AINV_LUMP,
971:   MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG,
972:   MAT_SCHUR_COMPLEMENT_AINV_FULL
973: } MatSchurComplementAinvType;
974: PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];

976: PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat, Mat, Mat, Mat, Mat, Mat *);
977: PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat, KSP *);
978: PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat, KSP);
979: PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
980: PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
981: PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat, Mat *, Mat *, Mat *, Mat *, Mat *);
982: PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat, MatSchurComplementAinvType);
983: PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat, MatSchurComplementAinvType *);
984: PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat, MatReuse, Mat *);
985: PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat, Mat *);
986: PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatSchurComplementAinvType, MatReuse, Mat *);
987: PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat, Mat, Mat, Mat, MatSchurComplementAinvType, MatReuse, Mat *);

989: PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm, PetscInt, PetscInt, Mat *);
990: PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm, PetscInt, PetscInt, Mat *);
991: PETSC_EXTERN PetscErrorCode MatCreateLMVMDBFGS(MPI_Comm, PetscInt, PetscInt, Mat *);
992: PETSC_EXTERN PetscErrorCode MatCreateLMVMDDFP(MPI_Comm, PetscInt, PetscInt, Mat *);
993: PETSC_EXTERN PetscErrorCode MatCreateLMVMDQN(MPI_Comm, PetscInt, PetscInt, Mat *);
994: PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm, PetscInt, PetscInt, Mat *);
995: PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
996: PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
997: PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
998: PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
999: PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);

1001: PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
1002: PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool *);
1003: PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
1004: PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
1005: PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
1006: PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
1007: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
1008: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
1009: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
1010: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
1011: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
1012: PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
1013: PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
1014: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat *);
1015: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC *);
1016: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP *);
1017: PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt);
1018: PETSC_EXTERN PetscErrorCode MatLMVMGetHistorySize(Mat, PetscInt *);
1019: PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt *);
1020: PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt *);
1021: PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar);

1023: /*E
1024:   MatLMVMSymBroydenScaleType - Scaling type for symmetric Broyden.

1026:   Values:
1027: + `MAT_LMVM_SYMBROYDEN_SCALE_NONE`     - No scaling
1028: . `MAT_LMVM_SYMBROYDEN_SCALE_SCALAR`   - scalar scaling
1029: . `MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL` - diagonal scaling
1030: - `MAT_LMVM_SYMBROYDEN_SCALE_USER`     - user-provided scale option

1032:   Level: intermediate

1034: .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMSymBroydenSetScaleType()`
1035: E*/
1036: typedef enum {
1037:   MAT_LMVM_SYMBROYDEN_SCALE_NONE     = 0,
1038:   MAT_LMVM_SYMBROYDEN_SCALE_SCALAR   = 1,
1039:   MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2,
1040:   MAT_LMVM_SYMBROYDEN_SCALE_USER     = 3
1041: } MatLMVMSymBroydenScaleType;
1042: PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[];

1044: PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType);

1046: /*E
1047:   MatLMVMDenseType - Memory storage strategy for dense variants `MATLMVM`.

1049:   Values:
1050: + `MAT_LMVM_DENSE_REORDER` - reorders memory to minimize kernel launch
1051: - `MAT_LMVM_DENSE_INPLACE` - computes inplace to minimize memory movement

1053:   Level: intermediate

1055: .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMDenseSetType()`
1056: E*/
1057: typedef enum {
1058:   MAT_LMVM_DENSE_REORDER,
1059:   MAT_LMVM_DENSE_INPLACE
1060: } MatLMVMDenseType;
1061: PETSC_EXTERN const char *const MatLMVMDenseTypes[];

1063: PETSC_EXTERN PetscErrorCode MatLMVMDenseSetType(Mat, MatLMVMDenseType);

1065: PETSC_EXTERN PetscErrorCode KSPSetDM(KSP, DM);
1066: PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP, PetscBool);
1067: PETSC_EXTERN PetscErrorCode KSPGetDM(KSP, DM *);
1068: PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP, void *);
1069: PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP, void *);

1071: /*S
1072:   KSPComputeRHSFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeRHS()`

1074:   Calling Sequence:
1075: + ksp  - `ksp` context
1076: . b    - output vector
1077: - ctx - [optional] user-defined function context

1079:   Level: beginner

1081: .seealso: [](ch_snes), `KSP`, `KSPSetComputeRHS()`, `SNESGetFunction()`, `KSPComputeInitialGuessFn`, `KSPComputeOperatorsFn`
1082: S*/
1083: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(KSPComputeRHSFn)(KSP ksp, Vec b, void *ctx);

1085: PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP, KSPComputeRHSFn *, void *);

1087: /*S
1088:   KSPComputeOperatorsFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeOperators()`

1090:   Calling Sequence:
1091: + ksp - `KSP` context
1092: . A   - the operator that defines the linear system
1093: . P   - an operator from which to build the preconditioner (often the same as `A`)
1094: - ctx - [optional] user-defined function context

1096:   Level: beginner

1098: .seealso: [](ch_snes), `KSP`, `KSPSetComputeRHS()`, `SNESGetFunction()`, `KSPComputeRHSFn`, `KSPComputeInitialGuessFn`
1099: S*/
1100: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(KSPComputeOperatorsFn)(KSP ksp, Mat A, Mat P, void *ctx);

1102: PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP, KSPComputeOperatorsFn, void *);

1104: /*S
1105:   KSPComputeInitialGuessFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeInitialGuess()`

1107:   Calling Sequence:
1108: + ksp  - `ksp` context
1109: . x    - output vector
1110: - ctx - [optional] user-defined function context

1112:   Level: beginner

1114: .seealso: [](ch_snes), `KSP`, `KSPSetComputeInitialGuess()`, `SNESGetFunction()`, `KSPComputeRHSFn`, `KSPComputeOperatorsFn`
1115: S*/
1116: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(KSPComputeInitialGuessFn)(KSP ksp, Vec x, void *ctx);

1118: PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP, KSPComputeInitialGuessFn *, void *);
1119: PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM, KSPComputeOperatorsFn *, void *);
1120: PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM, KSPComputeOperatorsFn **, void *);
1121: PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM, KSPComputeRHSFn *, void *);
1122: PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM, KSPComputeRHSFn **, void *);
1123: PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM, KSPComputeInitialGuessFn *, void *);
1124: PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM, KSPComputeInitialGuessFn **, void *);

1126: PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM, Vec, Vec);
1127: 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);
1128: PETSC_EXTERN PetscErrorCode DMSwarmProjectFields(DM, DM, PetscInt, const char **, Vec[], ScatterMode mode);

1130: PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *);
1131: PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal);

1133: PETSC_EXTERN PetscErrorCode PCBJKOKKOSSetKSP(PC, KSP);
1134: PETSC_EXTERN PetscErrorCode PCBJKOKKOSGetKSP(PC, KSP *);

1136: PETSC_EXTERN PetscErrorCode DMCopyDMKSP(DM, DM);