Actual source code: petscksp.h

petsc-master 2020-08-25
Report Typos and Errors
  1: /*
  2:    Defines the interface functions for the Krylov subspace accelerators.
  3: */
  4: #ifndef PETSCKSP_H
  5: #define PETSCKSP_H
  6:  #include <petscpc.h>

  8: PETSC_EXTERN PetscErrorCode KSPInitializePackage(void);

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

 14:    Level: beginner

 16:         Notes:
 17:     When a direct solver is used but no Krylov solver is used the KSP object is still used by with a
 18:        KSPType of KSPPREONLY (meaning Section 1.5 Writing Application Codes with PETSc of the preconditioner is only used as the linear solver).

 20: .seealso:  KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP, KSPDestroy(), KSPCG, KSPGMRES
 21: S*/
 22: typedef struct _p_KSP*     KSP;

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

 27:    Level: beginner

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

 81: /* Logging support */
 82: PETSC_EXTERN PetscClassId KSP_CLASSID;
 83: PETSC_EXTERN PetscClassId KSPGUESS_CLASSID;
 84: PETSC_EXTERN PetscClassId DMKSP_CLASSID;

 86: PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP*);
 87: PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType);
 88: PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType*);
 89: PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
 90: PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
 91: PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec);
 92: PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec);
 93: PETSC_EXTERN PetscErrorCode KSPMatSolve(KSP,Mat,Mat);
 94: PETSC_EXTERN PetscErrorCode KSPSetMatSolveBlockSize(KSP,PetscInt);
 95: PETSC_EXTERN PetscErrorCode KSPGetMatSolveBlockSize(KSP,PetscInt*);
 96: PETSC_EXTERN PetscErrorCode KSPReset(KSP);
 97: PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP);
 98: PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*);
 99: PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool);
100: PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP,PetscBool*);
101: PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool);
102: PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP,PC,Vec);

104: PETSC_EXTERN PetscFunctionList KSPList;
105: PETSC_EXTERN PetscFunctionList KSPGuessList;
106: PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP));

108: PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide);
109: PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*);
110: PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
111: PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
112: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool);
113: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool*);
114: PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool);
115: PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool*);
116: PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool);
117: PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool);
118: PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool*);
119: PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool);
120: PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool*);
121: PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec*);
122: PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec*);
123: PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*);
124: PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*);
125: PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*);
126: PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
127: PETSC_DEPRECATED_FUNCTION("Use KSPCreateVecs() (since version 3.6)") PETSC_STATIC_INLINE PetscErrorCode KSPGetVecs(KSP ksp,PetscInt n,Vec **x,PetscInt m,Vec **y) {return KSPCreateVecs(ksp,n,x,m,y);}

129: PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
130: PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);

132: PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC);
133: PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*);

135: PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal);
136: PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode (*)(void**));
137: PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
138: PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void**);
139: PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt*);
140: PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool);

142: PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*);
143: PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec*);
144: PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
145: PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt);

147: PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*);
148: PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC,KSP);
149: PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
150: PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
151: PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
152: PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
153: PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC,PetscInt*,KSP*[]);
154: PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*);
155: PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*);
156: PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*);
157: PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*);
158: PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP*);
159: PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC,KSP*);
160: /*
161:   PCMGCoarseList contains the list of coarse space constructor currently registered
162:   These are added with PCMGRegisterCoarseSpaceConstructor()
163: */
164: PETSC_EXTERN PetscFunctionList PCMGCoarseList;
165: PETSC_EXTERN PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char[], PetscErrorCode (*)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **));
166: PETSC_EXTERN PetscErrorCode PCMGGetCoarseSpaceConstructor(const char[], PetscErrorCode (**)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **));


169: PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec*);
170: PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec*);

172: PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
173: PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool);
174: PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal);
175: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal);
176: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP,PetscBool);
177: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*);
178: PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
179: PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt*);
180: PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]);
181: PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt*,Vec[],PetscReal[],PetscReal[]);

183: /*E

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

187:   KSP_FCD_TRUNC_TYPE_STANDARD uses all (up to mmax) stored directions
188:   KSP_FCD_TRUNC_TYPE_NOTAY uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..

190:    Level: intermediate
191: .seealso : KSPFCG,KSPPIPEFCG,KSPPIPEGCR,KSPFCGSetTruncationType(),KSPFCGGetTruncationType()

193: E*/
194: typedef enum {KSP_FCD_TRUNC_TYPE_STANDARD,KSP_FCD_TRUNC_TYPE_NOTAY} KSPFCDTruncationType;
195: PETSC_EXTERN const char *const KSPFCDTruncationTypes[];

197: PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt);
198: PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*);
199: PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt);
200: PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*);
201: PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType);
202: PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*);

204: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt);
205: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*);
206: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt);
207: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*);
208: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType);
209: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*);

211: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt);
212: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*);
213: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt);
214: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*);
215: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType);
216: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*);
217: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool);
218: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*);

220: PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
221: PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*);
222: PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);

224: PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
225: PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
226: PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
227: PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
228: PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);

230: PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
231: PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);

233: PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP,PetscScalar);

235: PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
236: PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
237: PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));

239: PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP,PC*);
240: PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP,PC);
241: PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP,KSP*);
242: PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP,Mat);

244: PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationSpace(KSP,Mat);
245: PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationSpace(KSP,Mat*);
246: PETSC_DEPRECATED_FUNCTION("Use KSPMatSolve() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X) { return KSPMatSolve(ksp, B, X); }
247: /*E
248:     KSPHPDDMType - Type of Krylov method used by KSPHPDDM

250:     Level: intermediate

252:     Values:
253: +   KSP_HPDDM_TYPE_GMRES (default)
254: .   KSP_HPDDM_TYPE_BGMRES
255: .   KSP_HPDDM_TYPE_CG
256: .   KSP_HPDDM_TYPE_BCG
257: .   KSP_HPDDM_TYPE_GCRODR
258: .   KSP_HPDDM_TYPE_BGCRODR
259: .   KSP_HPDDM_TYPE_BFBCG
260: -   KSP_HPDDM_TYPE_PREONLY

262: .seealso: KSPHPDDM, KSPHPDDMSetType()
263: E*/
264: typedef enum {
265:   KSP_HPDDM_TYPE_GMRES = 0,
266:   KSP_HPDDM_TYPE_BGMRES = 1,
267:   KSP_HPDDM_TYPE_CG = 2,
268:   KSP_HPDDM_TYPE_BCG = 3,
269:   KSP_HPDDM_TYPE_GCRODR = 4,
270:   KSP_HPDDM_TYPE_BGCRODR = 5,
271:   KSP_HPDDM_TYPE_BFBCG = 6,
272:   KSP_HPDDM_TYPE_PREONLY = 7
273: } KSPHPDDMType;
274: PETSC_EXTERN const char *const KSPHPDDMTypes[];
275: PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP,KSPHPDDMType);
276: PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP,KSPHPDDMType*);

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

281:    Level: advanced

283: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
284:           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()

286: E*/
287: typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
288: PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
289: /*MC
290:     KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process

292:    Level: advanced

294:    Note: Possible unstable, but the fastest to compute

296: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
297:           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
298:           KSPGMRESModifiedGramSchmidtOrthogonalization()
299: M*/

301: /*MC
302:     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
303:           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
304:           poor orthogonality.

306:    Level: advanced

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

311: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
312:           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
313:           KSPGMRESModifiedGramSchmidtOrthogonalization()
314: M*/

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

319:    Level: advanced

321:    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
322:      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.

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

326: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
327:           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
328:           KSPGMRESModifiedGramSchmidtOrthogonalization()
329: M*/

331: PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
332: PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);

334: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
335: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
336: PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));

338: PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
339: PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
340: PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);

342: PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
343: PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool);
344: PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
345: PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool);

347: PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
348: PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);
349: PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));

351: PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
352: PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
353: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
354: PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
355: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
356: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
357: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
358: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
359: PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
360: PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
361: PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*);
362: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**);
363: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**);
364: PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void*);
365: PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP,const char[],const char[],const char [],PetscErrorCode (*)(KSP,PetscInt,PetscReal,PetscViewerAndFormat*));

367: PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
368: PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);

370: PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat);
371: PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*);
372: PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool*,PetscBool*);
373: PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
374: PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
375: PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);

377: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool);
378: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool*);
379: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool);
380: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool*);

382: PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
383: PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);
384: PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP,PetscObject,const char[]);
385: PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP,PetscViewer,PetscViewerFormat);
386: PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP);

388: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonView() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode KSPReasonView(KSP ksp,PetscViewer v) {return KSPConvergedReasonView(ksp,v,PETSC_VIEWER_DEFAULT);}
389: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonViewFromOptions() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode KSPReasonViewFromOptions(KSP ksp) {return KSPConvergedReasonViewFromOptions(ksp);}

391: #define KSP_FILE_CLASSID 1211223

393: PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP,PetscBool);
394: PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP,PetscBool);
395: PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
396: PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP,PetscReal*,PetscReal*);

398: PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
399: PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
400: PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*);

402: /*E
403:     KSPNormType - Norm that is passed in the Krylov convergence
404:        test routines.

406:    Level: advanced

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

411:    Notes:
412:     this must match petsc/finclude/petscksp.h

414: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
415:           KSPSetConvergenceTest(), KSPSetPCSide()
416: E*/
417: typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
418: #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
419: PETSC_EXTERN const char *const*const KSPNormTypes;

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

426:    Level: advanced

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

430: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
431: M*/

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

437:    Level: advanced

439: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
440: M*/

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

446:    Level: advanced

448: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
449: M*/

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

455:    Level: advanced

457: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
458: M*/

460: PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
461: PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
462: PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
463: PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
464: PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool);

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

470:    Level: beginner

472:    Notes:
473:     See KSPGetConvergedReason() for explanation of each value

475:    Developer Notes:
476:     this must match petsc/finclude/petscksp.h

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

481: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances(), KSPConvergedReasonView()
482: E*/
483: typedef enum {/* converged */
484:               KSP_CONVERGED_RTOL_NORMAL        =  1,
485:               KSP_CONVERGED_ATOL_NORMAL        =  9,
486:               KSP_CONVERGED_RTOL               =  2,
487:               KSP_CONVERGED_ATOL               =  3,
488:               KSP_CONVERGED_ITS                =  4,
489:               KSP_CONVERGED_CG_NEG_CURVE       =  5,
490:               KSP_CONVERGED_CG_CONSTRAINED     =  6,
491:               KSP_CONVERGED_STEP_LENGTH        =  7,
492:               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
493:               /* diverged */
494:               KSP_DIVERGED_NULL                = -2,
495:               KSP_DIVERGED_ITS                 = -3,
496:               KSP_DIVERGED_DTOL                = -4,
497:               KSP_DIVERGED_BREAKDOWN           = -5,
498:               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
499:               KSP_DIVERGED_NONSYMMETRIC        = -7,
500:               KSP_DIVERGED_INDEFINITE_PC       = -8,
501:               KSP_DIVERGED_NANORINF            = -9,
502:               KSP_DIVERGED_INDEFINITE_MAT      = -10,
503:               KSP_DIVERGED_PC_FAILED           = -11,
504:               KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED  = -11,

506:               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
507: PETSC_EXTERN const char *const*KSPConvergedReasons;

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

512:    Level: beginner

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

518:    See also KSP_CONVERGED_ATOL which may apply before this tolerance.

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

522: M*/

524: /*MC
525:      KSP_CONVERGED_ATOL - norm(r) <= atol

527:    Level: beginner

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

533:    See also KSP_CONVERGED_RTOL which may apply before this tolerance.

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

537: M*/

539: /*MC
540:      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)

542:    Level: beginner

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

548:    Level: beginner

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

552: M*/

554: /*MC
555:      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
556:       reached

558:    Level: beginner

560: .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()

562: M*/

564: /*MC
565:      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
566:            the preconditioner is applied. Also used when the KSPConvergedSkip() convergence
567:            test routine is set in KSP.

569:    Level: beginner

571: .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()

573: M*/

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

581:    Level: beginner

583: .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()

585: M*/

587: /*MC
588:      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
589:           method could not continue to enlarge the Krylov space.

591:    Level: beginner

593: .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()

595: M*/

597: /*MC
598:      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
599:         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry

601:    Level: beginner

603: .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()

605: M*/

607: /*MC
608:      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
609:         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
610:         be positive definite

612:    Level: beginner

614:      Notes:
615:     This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
616:   the PCICC preconditioner to generate a positive definite preconditioner

618: .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()

620: M*/

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

627:    Level: beginner

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


633: .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()

635: M*/

637: /*MC
638:      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
639:         while the KSPSolve() is still running.

641:    Level: beginner

643: .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()

645: M*/

647: PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void*,PetscErrorCode (*)(void*));
648: PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*));
649: PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*));
650: PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void**);
651: PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
652: PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
653: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void*);
654: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void**);
655: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
656: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
657: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP,PetscBool);
658: PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
659: PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason*);

661: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ }
662: #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
663: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ }
664: #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
665: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ }
666: #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
667: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ }
668: #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
669: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ }
670: #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
671: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ }
672: #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)

674: PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP,MatType,Mat*);
675: PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPComputeExplicitOperator(KSP A,Mat* B) { return KSPComputeOperator(A,NULL,B); }

677: /*E
678:     KSPCGType - Determines what type of CG to use

680:    Level: beginner

682: .seealso: KSPCGSetType()
683: E*/
684: typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
685: PETSC_EXTERN const char *const KSPCGTypes[];

687: PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
688: PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool);

690: PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP,PetscReal);
691: PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP,PetscReal*);
692: PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP,PetscReal*);

694: PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal*);
695: PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal*);
696: PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp,PetscReal *x) {return KSPGLTRGetMinEig(ksp,x);}
697: PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetLambda(KSP ksp,PetscReal *x) {return KSPGLTRGetLambda(ksp,x);}

699: PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]);

701: PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
702: PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);

704:  #include <petscdrawtypes.h>
705: PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
706: PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*);
707: PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
708: PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
709: PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);

711: PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
712: PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));

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

717:    Level: beginner

719: .seealso:  KSPCreate(), KSPSetGuessType(), KSPGuessType
720: S*/
721: typedef struct _p_KSPGuess* KSPGuess;
722: /*J
723:     KSPGuessType - String with the name of a PETSc initial guess for Krylov methods.

725:    Level: beginner

727: .seealso: KSPGuess
728: J*/
729: typedef const char* KSPGuessType;
730: #define KSPGUESSFISCHER "fischer"
731: #define KSPGUESSPOD     "pod"
732: PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[],PetscErrorCode (*)(KSPGuess));
733: PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP,KSPGuess);
734: PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP,KSPGuess*);
735: PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess,PetscViewer);
736: PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess*);
737: PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm,KSPGuess*);
738: PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess,KSPGuessType);
739: PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess,KSPGuessType*);
740: PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
741: PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess,Vec,Vec);
742: PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess,Vec,Vec);
743: PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
744: PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess,PetscInt,PetscInt);
745: PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
746: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool);
747: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool*);

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

752:     Level: intermediate

754: .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat()
755: E*/
756: typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG} MatSchurComplementAinvType;
757: PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];

759: typedef enum {MAT_LMVM_SYMBROYDEN_SCALE_NONE       = 0,
760:               MAT_LMVM_SYMBROYDEN_SCALE_SCALAR     = 1,
761:               MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL   = 2,
762:               MAT_LMVM_SYMBROYDEN_SCALE_USER       = 3} MatLMVMSymBroydenScaleType;
763: PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[];

765: PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
766: PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
767: PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
768: PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
769: PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
770: PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
771: PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType);
772: PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*);
773: PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*);
774: PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*);
775: PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatSchurComplementAinvType,MatReuse,Mat*);
776: PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*);

778: PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm,PetscInt,PetscInt,Mat*);
779: PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm,PetscInt,PetscInt,Mat*);
780: PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm,PetscInt,PetscInt,Mat*);
781: PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
782: PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
783: PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
784: PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
785: PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);

787: PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
788: PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool*);
789: PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
790: PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
791: PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
792: PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
793: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
794: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
795: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
796: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
797: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
798: PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
799: PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
800: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat*);
801: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC*);
802: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP*);
803: PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt);
804: PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt*);
805: PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt*);
806: PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar);
807: PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType);

809: PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
810: PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool);
811: PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
812: PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
813: PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
814: PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void*);
815: PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
816: PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
817: PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
818: PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*);
819: PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
820: PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
821: PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
822: PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);

824: PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec);
825: PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec,
826:                                            void (**)(PetscInt, PetscInt, PetscInt,
827:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
828:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
829:                                                      PetscReal, const PetscReal [], PetscInt, const PetscScalar[], PetscScalar []), InsertMode, Vec);


832: PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, PetscInt, Vec[], Vec[], Mat *, void *);
833: PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, PetscInt, Vec[], Vec[], PetscReal);
834: #endif