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 the linear solves in PETSc (even those such as direct factorization-based solvers that
 15:          do not use Krylov accelerators).

 17:    Level: beginner

 19:    Notes:
 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:    Use `KSPSetType()` or the options database key `-ksp_type` to set the specific Krylov solver algorithm to use with a given `KSP` object

 25:    The `PC` object is used to control preconditioners in PETSc.

 27: .seealso: [](doc_linsolve), [](ch_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `SNES`, `TS`, `PC`, `KSP`, `KSPDestroy()`, `KSPCG`, `KSPGMRES`
 28: S*/
 29: typedef struct _p_KSP *KSP;

 31: /*J
 32:    KSPType - String with the name of a PETSc Krylov method. These are all the Krylov solvers that PETSc provides.

 34:    Level: beginner

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

 91: /* Logging support */
 92: PETSC_EXTERN PetscClassId KSP_CLASSID;
 93: PETSC_EXTERN PetscClassId KSPGUESS_CLASSID;
 94: PETSC_EXTERN PetscClassId DMKSP_CLASSID;

 96: PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm, KSP *);
 97: PETSC_EXTERN PetscErrorCode KSPSetType(KSP, KSPType);
 98: PETSC_EXTERN PetscErrorCode KSPGetType(KSP, KSPType *);
 99: PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
100: PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
101: PETSC_EXTERN PetscErrorCode KSPSolve(KSP, Vec, Vec);
102: PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP, Vec, Vec);
103: PETSC_EXTERN PetscErrorCode KSPSetUseExplicitTranspose(KSP, PetscBool);
104: PETSC_EXTERN PetscErrorCode KSPMatSolve(KSP, Mat, Mat);
105: PETSC_EXTERN PetscErrorCode KSPMatSolveTranspose(KSP, Mat, Mat);
106: PETSC_EXTERN PetscErrorCode KSPSetMatSolveBatchSize(KSP, PetscInt);
107: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPSetMatSolveBatchSize()", ) static inline PetscErrorCode KSPSetMatSolveBlockSize(KSP ksp, PetscInt n)
108: {
109:   return KSPSetMatSolveBatchSize(ksp, n);
110: }
111: PETSC_EXTERN PetscErrorCode KSPGetMatSolveBatchSize(KSP, PetscInt *);
112: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPGetMatSolveBatchSize()", ) static inline PetscErrorCode KSPGetMatSolveBlockSize(KSP ksp, PetscInt *n)
113: {
114:   return KSPGetMatSolveBatchSize(ksp, n);
115: }
116: PETSC_EXTERN PetscErrorCode KSPReset(KSP);
117: PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP);
118: PETSC_EXTERN PetscErrorCode KSPDestroy(KSP *);
119: PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP, PetscBool);
120: PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP, PetscBool *);
121: PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP, PetscBool);
122: PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP, PC, Vec);

124: PETSC_EXTERN PetscFunctionList KSPList;
125: PETSC_EXTERN PetscFunctionList KSPGuessList;
126: PETSC_EXTERN PetscFunctionList KSPMonitorList;
127: PETSC_EXTERN PetscFunctionList KSPMonitorCreateList;
128: PETSC_EXTERN PetscFunctionList KSPMonitorDestroyList;
129: PETSC_EXTERN PetscErrorCode    KSPRegister(const char[], PetscErrorCode (*)(KSP));
130: PETSC_EXTERN PetscErrorCode KSPMonitorRegister(const char[], PetscViewerType, PetscViewerFormat, PetscErrorCode (*)(KSP, PetscInt, PetscReal, PetscViewerAndFormat *), PetscErrorCode (*)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **), PetscErrorCode (*)(PetscViewerAndFormat **));

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

158: PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *);
159: PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *);

161: PETSC_EXTERN PetscErrorCode KSPSetPC(KSP, PC);
162: PETSC_EXTERN PetscErrorCode KSPGetPC(KSP, PC *);
163: PETSC_EXTERN PetscErrorCode KSPSetNestLevel(KSP, PetscInt);
164: PETSC_EXTERN PetscErrorCode KSPGetNestLevel(KSP, PetscInt *);

166: PETSC_EXTERN PetscErrorCode KSPMonitor(KSP, PetscInt, PetscReal);
167: PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void **));
168: PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
169: PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP, void *);
170: PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP, const PetscReal *[], PetscInt *);
171: PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP, PetscReal[], PetscCount, PetscBool);
172: PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP, const PetscReal *[], PetscInt *);
173: PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP, PetscReal[], PetscCount, PetscBool);

175: PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP, Vec, Vec *);
176: PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP, Vec, Vec, Vec *);
177: PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
178: PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP, PetscInt);

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

202: PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP, Vec, Vec *);
203: PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP, Vec, Vec, Vec *);

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

208:   Values:
209: + `KSP_CHEBYSHEV_FIRST`      - "classic" first-kind Chebyshev polynomial
210: . `KSP_CHEBYSHEV_FOURTH`     - fourth-kind Chebyshev polynomial
211: - `KSP_CHEBYSHEV_OPT_FOURTH` - optimized fourth-kind Chebyshev polynomial

213:   Level: intermediate

215: .seealso: [](ch_ksp), `KSPCHEBYSHEV`, `KSPChebyshevSetKind()`, `KSPChebyshevGetKind()`
216: E*/
217: typedef enum {
218:   KSP_CHEBYSHEV_FIRST,
219:   KSP_CHEBYSHEV_FOURTH,
220:   KSP_CHEBYSHEV_OPT_FOURTH
221: } KSPChebyshevKind;

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

236: /*E

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

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

244:    Level: intermediate

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

254: PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP, PetscInt);
255: PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP, PetscInt *);
256: PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP, PetscInt);
257: PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP, PetscInt *);
258: PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP, KSPFCDTruncationType);
259: PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP, KSPFCDTruncationType *);

261: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP, PetscInt);
262: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP, PetscInt *);
263: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP, PetscInt);
264: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP, PetscInt *);
265: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP, KSPFCDTruncationType);
266: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP, KSPFCDTruncationType *);

268: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP, PetscInt);
269: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP, PetscInt *);
270: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP, PetscInt);
271: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP, PetscInt *);
272: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP, KSPFCDTruncationType);
273: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP, KSPFCDTruncationType *);
274: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP, PetscBool);
275: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP, PetscBool *);
276: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *));

278: PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
279: PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt *);
280: PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP, PetscReal);
281: PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP, PetscReal);

283: PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
284: PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP, PetscErrorCode (*)(KSP, PetscInt));
285: PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP, PetscErrorCode (**)(KSP, PetscInt));
286: PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP, PetscInt);
287: PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP, PetscInt);

289: PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP, PetscInt);
290: PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);

292: PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP, PetscScalar);

294: PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP, PetscInt);
295: PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP, PetscInt *);
296: PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *));

298: PETSC_EXTERN PetscErrorCode KSPMINRESSetRadius(KSP, PetscReal);
299: PETSC_EXTERN PetscErrorCode KSPMINRESGetUseQLP(KSP, PetscBool *);
300: PETSC_EXTERN PetscErrorCode KSPMINRESSetUseQLP(KSP, PetscBool);

302: PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP, PC *);
303: PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP, PC);
304: PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP, KSP *);
305: PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP, Mat);

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

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

336:     Level: intermediate

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

352: /*E
353:     KSPHPDDMPrecision - Precision of Krylov bases used by `KSPHPDDM`

355:     Values:
356: +   `KSP_HPDDM_PRECISION_HALF`      - default when PETSc is configured `--with-precision=__fp16`
357: .   `KSP_HPDDM_PRECISION_SINGLE`    - default when PETSc is configured `--with-precision=single`
358: .   `KSP_HPDDM_PRECISION_DOUBLE`    - default when PETSc is configured `--with-precision=double`
359: -   `KSP_HPDDM_PRECISION_QUADRUPLE` - default when PETSc is configured `--with-precision=__float128`

361:     Level: intermediate

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

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

377:    Values:
378: +  `KSP_GMRES_CGS_REFINE_NEVER`    - one step of classical Gram-Schmidt
379: .  `KSP_GMRES_CGS_REFINE_IFNEEDED` - a second step is performed if the first step does not satisfy some criteria
380: -  `KSP_GMRES_CGS_REFINE_ALWAYS`   - always perform two steps

382:    Level: advanced

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

395: /*MC
396:    KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process

398:    Level: advanced

400:    Note:
401:    Possibly unstable, but the fastest to compute

403: .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
404:           `KSP`, `KSPGMRESGetOrthogonalization()`,
405:           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
406:           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
407: M*/

409: /*MC
410:     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
411:           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
412:           poor orthogonality.

414:    Level: advanced

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

420: .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
421:           `KSP`, `KSPGMRESGetOrthogonalization()`,
422:           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
423:           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
424: M*/

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

429:    Level: advanced

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

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

437: .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`,
438:           `KSP`, `KSPGMRESGetOrthogonalization()`,
439:           `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`,
440:           `KSPGMRESModifiedGramSchmidtOrthogonalization()`
441: M*/

443: PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP, KSPGMRESCGSRefinementType);
444: PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP, KSPGMRESCGSRefinementType *);

446: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP, PetscInt, PetscInt, PetscReal, void *);
447: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP, PetscInt, PetscInt, PetscReal, void *);
448: PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *));

450: PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP, PetscReal);
451: PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP, PetscReal *);
452: PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP, PetscReal *);

454: PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP, PetscReal);
455: PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP, PetscBool);
456: PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP, PetscInt);
457: PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP, PetscBool);

459: PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
460: PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);

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

498: PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP, PetscInt, PetscReal, void *);
499: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP, PetscInt, PetscReal, void *);
500: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **);
501: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceCreate(void *);
502: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceSetCoefficient(void *, PetscReal);
503: PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP, PetscInt, PetscReal, void *);
504: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP, void **);
505: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void **);

507: PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP, Vec, Vec);
508: PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP, Vec, Vec, Vec, Vec, Vec);

510: PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP, Mat, Mat);
511: PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP, Mat *, Mat *);
512: PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP, PetscBool *, PetscBool *);
513: PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP, const char[]);
514: PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP, const char[]);
515: PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP, const char *[]);

517: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP, PetscBool);
518: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP, PetscBool *);
519: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP, PetscBool);
520: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP, PetscBool *);

522: PETSC_EXTERN PetscErrorCode KSPView(KSP, PetscViewer);
523: PETSC_EXTERN PetscErrorCode KSPLoad(KSP, PetscViewer);
524: PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP, PetscObject, const char[]);
525: PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP, PetscViewer);
526: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP, PetscErrorCode (*)(KSP, void *), void *vctx, PetscErrorCode (*)(void **));
527: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP);
528: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP);
529: PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP, PetscViewer);

531: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonView()", ) static inline PetscErrorCode KSPReasonView(KSP ksp, PetscViewer v)
532: {
533:   return KSPConvergedReasonView(ksp, v);
534: }
535: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonViewFromOptions()", ) static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp)
536: {
537:   return KSPConvergedReasonViewFromOptions(ksp);
538: }

540: #define KSP_FILE_CLASSID 1211223

542: PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP, PetscBool);
543: PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP, PetscBool);
544: PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP, Vec *);
545: PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP, PetscReal *, PetscReal *);
546: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
547: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
548: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);

550: PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC, KSP *);
551: PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC, KSP *);
552: PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC, KSP *);
553: PETSC_EXTERN PetscErrorCode PCMPIGetKSP(PC, KSP *);

555: /*E
556:    KSPNormType - Norm calculated by the `KSP` and passed in the Krylov convergence
557:        test routines.

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

566:    Level: advanced

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

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

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

590:    Level: advanced

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

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

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

602:    Level: advanced

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

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

611:    Level: advanced

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

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

620:    Level: advanced

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

625: PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP, KSPNormType);
626: PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP, KSPNormType *);
627: PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp, KSPNormType, PCSide, PetscInt);
628: PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP, PetscInt);
629: PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP, PetscBool);

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

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

657:    Level: beginner

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

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

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

693:   KSP_CONVERGED_ITERATING = 0
694: } KSPConvergedReason;
695: PETSC_EXTERN const char *const *KSPConvergedReasons;

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

700:    Level: beginner

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

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

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

712: /*MC
713:    KSP_CONVERGED_ATOL - $||r|| \le atol$

715:    Level: beginner

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

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

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

727: /*MC
728:    KSP_DIVERGED_DTOL - $||r|| \ge dtol*||b||$

730:    Level: beginner

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

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

740: /*MC
741:    KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
742:    reached

744:    Level: beginner

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

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

754:    Level: beginner

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

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

765:    Level: beginner

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

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

774:    Level: beginner

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

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

783:    Level: beginner

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

788: /*MC
789:    KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
790:    positive and negative eigenvalues) and this Krylov method (`KSPCG`) requires it to
791:    be symmetric positive definite (SPD).

793:    Level: beginner

795:    Note:
796:    This can happen with the `PCICC` preconditioner, use the options database option `-pc_factor_shift_positive_definite` to force
797:    the `PCICC` preconditioner to generate a positive definite preconditioner

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

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

807:    Level: beginner

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

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

815: /*MC
816:    KSP_CONVERGED_ITERATING - This flag is returned if `KSPGetConvergedReason()` is called
817:    while `KSPSolve()` is still running.

819:    Level: beginner

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

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

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

867: PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP, MatType, Mat *);
868: PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPComputeOperator()", ) static inline PetscErrorCode KSPComputeExplicitOperator(KSP A, Mat *B)
869: {
870:   return KSPComputeOperator(A, PETSC_NULLPTR, B);
871: }

873: /*E
874:    KSPCGType - Determines what type of `KSPCG` to use

876:    Values:
877:  + `KSP_CG_SYMMETRIC` - the matrix is complex symmetric
878:  - `KSP_CG_HERMITIAN` - the matrix is complex Hermitian

880:    Level: beginner

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

890: PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP, KSPCGType);
891: PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP, PetscBool);

893: PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP, PetscReal);
894: PETSC_EXTERN PetscErrorCode KSPCGSetObjectiveTarget(KSP, PetscReal);
895: PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP, PetscReal *);
896: PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP, PetscReal *);

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

909: PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP, const char[]);
910: PETSC_EXTERN PetscErrorCode KSPPythonGetType(KSP, const char *[]);

912: PETSC_EXTERN PetscErrorCode PCSetPreSolve(PC, PetscErrorCode (*)(PC, KSP));
913: PETSC_EXTERN PetscErrorCode PCPreSolve(PC, KSP);
914: PETSC_EXTERN PetscErrorCode PCPostSolve(PC, KSP);

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

919: PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec));
920: PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec));

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

925:    Level: intermediate

927:    Note:
928:    These methods generate initial guesses based on a series of previous, related, linear solves. For example,
929:    in implicit time-stepping with `TS`.

931: .seealso: [](ch_ksp), `KSPCreate()`, `KSPGuessSetType()`, `KSPGuessType`
932: S*/
933: typedef struct _p_KSPGuess *KSPGuess;

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

938:    Values:
939:  + `KSPGUESSFISCHER` - methodology developed by Paul Fischer
940:  - `KSPGUESSPOD`     - methodology based on proper orthogonal decomposition (POD)

942:    Level: intermediate

944: .seealso: [](ch_ksp), `KSP`, `KSPGuess`
945: J*/
946: typedef const char *KSPGuessType;
947: #define KSPGUESSFISCHER "fischer"
948: #define KSPGUESSPOD     "pod"

950: PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[], PetscErrorCode (*)(KSPGuess));
951: PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP, KSPGuess);
952: PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP, KSPGuess *);
953: PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess, PetscViewer);
954: PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess *);
955: PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm, KSPGuess *);
956: PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess, KSPGuessType);
957: PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess, KSPGuessType *);
958: PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess, PetscReal);
959: PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
960: PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess, Vec, Vec);
961: PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess, Vec, Vec);
962: PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
963: PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess, PetscInt, PetscInt);
964: PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP, PetscInt, PetscInt);
965: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP, PetscBool);
966: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP, PetscBool *);

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

971:     Level: intermediate

973: .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`,
974:           `MatCreateSchurComplementPmat()`, `MatCreateSchurComplement()`
975: E*/
976: typedef enum {
977:   MAT_SCHUR_COMPLEMENT_AINV_DIAG,
978:   MAT_SCHUR_COMPLEMENT_AINV_LUMP,
979:   MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG,
980:   MAT_SCHUR_COMPLEMENT_AINV_FULL
981: } MatSchurComplementAinvType;
982: PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];

984: PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat, Mat, Mat, Mat, Mat, Mat *);
985: PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat, KSP *);
986: PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat, KSP);
987: PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
988: PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat);
989: PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat, Mat *, Mat *, Mat *, Mat *, Mat *);
990: PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat, MatSchurComplementAinvType);
991: PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat, MatSchurComplementAinvType *);
992: PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat, MatReuse, Mat *);
993: PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat, Mat *);
994: PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatSchurComplementAinvType, MatReuse, Mat *);
995: PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat, Mat, Mat, Mat, MatSchurComplementAinvType, MatReuse, Mat *);

997: PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm, PetscInt, PetscInt, Mat *);
998: PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm, PetscInt, PetscInt, Mat *);
999: PETSC_EXTERN PetscErrorCode MatCreateLMVMDBFGS(MPI_Comm, PetscInt, PetscInt, Mat *);
1000: PETSC_EXTERN PetscErrorCode MatCreateLMVMDDFP(MPI_Comm, PetscInt, PetscInt, Mat *);
1001: PETSC_EXTERN PetscErrorCode MatCreateLMVMDQN(MPI_Comm, PetscInt, PetscInt, Mat *);
1002: PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm, PetscInt, PetscInt, Mat *);
1003: PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1004: PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1005: PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1006: PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);
1007: PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm, PetscInt, PetscInt, Mat *);

1009: PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
1010: PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool *);
1011: PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
1012: PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
1013: PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
1014: PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
1015: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
1016: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
1017: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
1018: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
1019: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
1020: PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
1021: PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
1022: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat *);
1023: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC *);
1024: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP *);
1025: PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt);
1026: PETSC_EXTERN PetscErrorCode MatLMVMGetHistorySize(Mat, PetscInt *);
1027: PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt *);
1028: PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt *);
1029: PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar);

1031: /*E
1032:   MatLMVMSymBroydenScaleType - Scaling type for symmetric Broyden.

1034:   Values:
1035: + `MAT_LMVM_SYMBROYDEN_SCALE_NONE`     - No scaling
1036: . `MAT_LMVM_SYMBROYDEN_SCALE_SCALAR`   - scalar scaling
1037: . `MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL` - diagonal scaling
1038: - `MAT_LMVM_SYMBROYDEN_SCALE_USER`     - user-provided scale option

1040:   Level: intermediate

1042: .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMSymBroydenSetScaleType()`
1043: E*/
1044: typedef enum {
1045:   MAT_LMVM_SYMBROYDEN_SCALE_NONE     = 0,
1046:   MAT_LMVM_SYMBROYDEN_SCALE_SCALAR   = 1,
1047:   MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2,
1048:   MAT_LMVM_SYMBROYDEN_SCALE_USER     = 3
1049: } MatLMVMSymBroydenScaleType;
1050: PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[];

1052: PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType);

1054: /*E
1055:   MatLMVMDenseType - Memory storage strategy for dense variants of `MATLMVM`.

1057:   Values:
1058: + `MAT_LMVM_DENSE_REORDER` - reorders memory to minimize kernel launch
1059: - `MAT_LMVM_DENSE_INPLACE` - computes inplace to minimize memory movement

1061:   Level: intermediate

1063: .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMDenseSetType()`
1064: E*/
1065: typedef enum {
1066:   MAT_LMVM_DENSE_REORDER,
1067:   MAT_LMVM_DENSE_INPLACE
1068: } MatLMVMDenseType;
1069: PETSC_EXTERN const char *const MatLMVMDenseTypes[];

1071: PETSC_EXTERN PetscErrorCode MatLMVMDenseSetType(Mat, MatLMVMDenseType);

1073: PETSC_EXTERN PetscErrorCode KSPSetDM(KSP, DM);
1074: PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP, PetscBool);
1075: PETSC_EXTERN PetscErrorCode KSPGetDM(KSP, DM *);
1076: PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP, void *);
1077: PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP, void *);

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

1082:   Calling Sequence:
1083: + ksp  - `ksp` context
1084: . b    - output vector
1085: - ctx - [optional] user-defined function context

1087:   Level: beginner

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

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

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

1098:   Calling Sequence:
1099: + ksp - `KSP` context
1100: . A   - the operator that defines the linear system
1101: . P   - an operator from which to build the preconditioner (often the same as `A`)
1102: - ctx - [optional] user-defined function context

1104:   Level: beginner

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

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

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

1115:   Calling Sequence:
1116: + ksp  - `ksp` context
1117: . x    - output vector
1118: - ctx - [optional] user-defined function context

1120:   Level: beginner

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

1126: PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP, KSPComputeInitialGuessFn *, void *);
1127: PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM, KSPComputeOperatorsFn *, void *);
1128: PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM, KSPComputeOperatorsFn **, void *);
1129: PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM, KSPComputeRHSFn *, void *);
1130: PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM, KSPComputeRHSFn **, void *);
1131: PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM, KSPComputeInitialGuessFn *, void *);
1132: PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM, KSPComputeInitialGuessFn **, void *);

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

1138: PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *);
1139: PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal);

1141: PETSC_EXTERN PetscErrorCode PCBJKOKKOSSetKSP(PC, KSP);
1142: PETSC_EXTERN PetscErrorCode PCBJKOKKOSGetKSP(PC, KSP *);

1144: PETSC_EXTERN PetscErrorCode DMCopyDMKSP(DM, DM);