Actual source code: petscksp.h

petsc-master 2020-07-10
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 application 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*);

161: PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec*);
162: PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec*);

164: PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
165: PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool );
166: PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal);
167: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal);
168: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP,PetscBool);
169: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*);
170: PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
171: PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt*);
172: PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]);
173: PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt*,Vec[],PetscReal[],PetscReal[]);

175: /*E

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

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

182:    Level: intermediate
183: .seealso : KSPFCG,KSPPIPEFCG,KSPPIPEGCR,KSPFCGSetTruncationType(),KSPFCGGetTruncationType()

185: E*/
186: typedef enum {KSP_FCD_TRUNC_TYPE_STANDARD,KSP_FCD_TRUNC_TYPE_NOTAY} KSPFCDTruncationType;
187: PETSC_EXTERN const char *const KSPFCDTruncationTypes[];

189: PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt);
190: PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*);
191: PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt);
192: PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*);
193: PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType);
194: PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*);

196: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt);
197: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*);
198: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt);
199: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*);
200: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType);
201: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*);

203: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt);
204: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*);
205: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt);
206: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*);
207: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType);
208: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*);
209: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool);
210: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*);

212: PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
213: PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*);
214: PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);

216: PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
217: PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
218: PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
219: PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
220: PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);

222: PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
223: PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);

225: PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP,PetscScalar);

227: PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
228: PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
229: PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));

231: PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP,PC*);
232: PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP,PC);
233: PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP,KSP*);
234: PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP,Mat);

236: PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationSpace(KSP,Mat);
237: PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationSpace(KSP,Mat*);
238: 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); }
239: /*E
240:     KSPHPDDMType - Type of Krylov method used by KSPHPDDM

242:     Level: intermediate

244:     Values:
245: +   KSP_HPDDM_TYPE_GMRES (default)
246: .   KSP_HPDDM_TYPE_BGMRES
247: .   KSP_HPDDM_TYPE_CG
248: .   KSP_HPDDM_TYPE_BCG
249: .   KSP_HPDDM_TYPE_GCRODR
250: .   KSP_HPDDM_TYPE_BGCRODR
251: .   KSP_HPDDM_TYPE_BFBCG
252: -   KSP_HPDDM_TYPE_PREONLY

254: .seealso: KSPHPDDM, KSPHPDDMSetType()
255: E*/
256: typedef enum {
257:   KSP_HPDDM_TYPE_GMRES = 0,
258:   KSP_HPDDM_TYPE_BGMRES = 1,
259:   KSP_HPDDM_TYPE_CG = 2,
260:   KSP_HPDDM_TYPE_BCG = 3,
261:   KSP_HPDDM_TYPE_GCRODR = 4,
262:   KSP_HPDDM_TYPE_BGCRODR = 5,
263:   KSP_HPDDM_TYPE_BFBCG = 6,
264:   KSP_HPDDM_TYPE_PREONLY = 7
265: } KSPHPDDMType;
266: PETSC_EXTERN const char *const KSPHPDDMTypes[];
267: PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP,KSPHPDDMType);
268: PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP,KSPHPDDMType*);

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

273:    Level: advanced

275: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
276:           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()

278: E*/
279: typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
280: PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
281: /*MC
282:     KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process

284:    Level: advanced

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

288: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
289:           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
290:           KSPGMRESModifiedGramSchmidtOrthogonalization()
291: M*/

293: /*MC
294:     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
295:           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
296:           poor orthogonality.

298:    Level: advanced

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

303: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
304:           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
305:           KSPGMRESModifiedGramSchmidtOrthogonalization()
306: M*/

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

311:    Level: advanced

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

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

318: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
319:           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
320:           KSPGMRESModifiedGramSchmidtOrthogonalization()
321: M*/

323: PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
324: PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);

326: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
327: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
328: PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));

330: PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
331: PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
332: PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);

334: PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
335: PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool );
336: PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
337: PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool);

339: PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
340: PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);
341: PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));

343: PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
344: PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
345: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
346: PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
347: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
348: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
349: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
350: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
351: PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
352: PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
353: PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*);
354: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**);
355: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**);
356: PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void*);
357: PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP,const char[],const char[],const char [],PetscErrorCode (*)(KSP,PetscInt,PetscReal,PetscViewerAndFormat*));

359: PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
360: PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);

362: PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat);
363: PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*);
364: PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool*,PetscBool*);
365: PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
366: PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
367: PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);

369: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool );
370: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool*);
371: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool );
372: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool*);

374: PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
375: PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);
376: PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP,PetscObject,const char[]);
377: PETSC_EXTERN PetscErrorCode KSPReasonView(KSP,PetscViewer);
378: PETSC_EXTERN PetscErrorCode KSPReasonViewFromOptions(KSP);

380: #define KSP_FILE_CLASSID 1211223

382: PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP,PetscBool);
383: PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP,PetscBool);
384: PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
385: PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP,PetscReal*,PetscReal*);

387: PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
388: PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
389: PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*);

391: /*E
392:     KSPNormType - Norm that is passed in the Krylov convergence
393:        test routines.

395:    Level: advanced

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

400:    Notes:
401:     this must match petsc/finclude/petscksp.h

403: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
404:           KSPSetConvergenceTest(), KSPSetPCSide()
405: E*/
406: typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
407: #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
408: PETSC_EXTERN const char *const*const KSPNormTypes;

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

415:    Level: advanced

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

419: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
420: M*/

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

426:    Level: advanced

428: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
429: M*/

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

435:    Level: advanced

437: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
438: M*/

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

444:    Level: advanced

446: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
447: M*/

449: PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
450: PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
451: PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
452: PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
453: PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool);

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

459:    Level: beginner

461:    Notes:
462:     See KSPGetConvergedReason() for explanation of each value

464:    Developer Notes:
465:     this must match petsc/finclude/petscksp.h

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

470: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
471: E*/
472: typedef enum {/* converged */
473:               KSP_CONVERGED_RTOL_NORMAL        =  1,
474:               KSP_CONVERGED_ATOL_NORMAL        =  9,
475:               KSP_CONVERGED_RTOL               =  2,
476:               KSP_CONVERGED_ATOL               =  3,
477:               KSP_CONVERGED_ITS                =  4,
478:               KSP_CONVERGED_CG_NEG_CURVE       =  5,
479:               KSP_CONVERGED_CG_CONSTRAINED     =  6,
480:               KSP_CONVERGED_STEP_LENGTH        =  7,
481:               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
482:               /* diverged */
483:               KSP_DIVERGED_NULL                = -2,
484:               KSP_DIVERGED_ITS                 = -3,
485:               KSP_DIVERGED_DTOL                = -4,
486:               KSP_DIVERGED_BREAKDOWN           = -5,
487:               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
488:               KSP_DIVERGED_NONSYMMETRIC        = -7,
489:               KSP_DIVERGED_INDEFINITE_PC       = -8,
490:               KSP_DIVERGED_NANORINF            = -9,
491:               KSP_DIVERGED_INDEFINITE_MAT      = -10,
492:               KSP_DIVERGED_PC_FAILED           = -11,
493:               KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED  = -11,

495:               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
496: PETSC_EXTERN const char *const*KSPConvergedReasons;

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

501:    Level: beginner

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

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

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

511: M*/

513: /*MC
514:      KSP_CONVERGED_ATOL - norm(r) <= atol

516:    Level: beginner

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

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

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

526: M*/

528: /*MC
529:      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)

531:    Level: beginner

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

537:    Level: beginner

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

541: M*/

543: /*MC
544:      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
545:       reached

547:    Level: beginner

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

551: M*/

553: /*MC
554:      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
555:            the preconditioner is applied. Also used when the KSPConvergedSkip() convergence
556:            test routine is set in KSP.

558:    Level: beginner

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

562: M*/

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

570:    Level: beginner

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

574: M*/

576: /*MC
577:      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
578:           method could not continue to enlarge the Krylov space.

580:    Level: beginner

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

584: M*/

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

590:    Level: beginner

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

594: M*/

596: /*MC
597:      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
598:         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
599:         be positive definite

601:    Level: beginner

603:      Notes:
604:     This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
605:   the PCICC preconditioner to generate a positive definite preconditioner

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

609: M*/

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

616:    Level: beginner

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


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

624: M*/

626: /*MC
627:      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
628:         while the KSPSolve() is still running.

630:    Level: beginner

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

634: M*/

636: PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void*,PetscErrorCode (*)(void*));
637: PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*));
638: PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*));
639: PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void**);
640: PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
641: PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
642: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void*);
643: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void**);
644: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
645: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
646: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP,PetscBool);
647: PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
648: PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason*);

650: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ }
651: #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
652: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ }
653: #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
654: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ }
655: #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
656: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ }
657: #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
658: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ }
659: #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
660: PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ }
661: #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)

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

666: /*E
667:     KSPCGType - Determines what type of CG to use

669:    Level: beginner

671: .seealso: KSPCGSetType()
672: E*/
673: typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
674: PETSC_EXTERN const char *const KSPCGTypes[];

676: PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
677: PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool );

679: PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP,PetscReal);
680: PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP,PetscReal*);
681: PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP,PetscReal*);

683: PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal*);
684: PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal*);
685: PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp,PetscReal *x) {return KSPGLTRGetMinEig(ksp,x);}
686: PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetLambda(KSP ksp,PetscReal *x) {return KSPGLTRGetLambda(ksp,x);}

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

690: PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
691: PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);

693:  #include <petscdrawtypes.h>
694: PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
695: PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*);
696: PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
697: PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
698: PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);

700: PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
701: PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));

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

706:    Level: beginner

708: .seealso:  KSPCreate(), KSPSetGuessType(), KSPGuessType
709: S*/
710: typedef struct _p_KSPGuess* KSPGuess;
711: /*J
712:     KSPGuessType - String with the name of a PETSc initial guess for Krylov methods.

714:    Level: beginner

716: .seealso: KSPGuess
717: J*/
718: typedef const char* KSPGuessType;
719: #define KSPGUESSFISCHER "fischer"
720: #define KSPGUESSPOD     "pod"
721: PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[],PetscErrorCode (*)(KSPGuess));
722: PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP,KSPGuess);
723: PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP,KSPGuess*);
724: PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess,PetscViewer);
725: PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess*);
726: PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm,KSPGuess*);
727: PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess,KSPGuessType);
728: PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess,KSPGuessType*);
729: PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
730: PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess,Vec,Vec);
731: PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess,Vec,Vec);
732: PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
733: PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess,PetscInt,PetscInt);
734: PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
735: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool);
736: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool*);

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

741:     Level: intermediate

743: .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat()
744: E*/
745: typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG} MatSchurComplementAinvType;
746: PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];

748: typedef enum {MAT_LMVM_SYMBROYDEN_SCALE_NONE       = 0,
749:               MAT_LMVM_SYMBROYDEN_SCALE_SCALAR     = 1,
750:               MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL   = 2,
751:               MAT_LMVM_SYMBROYDEN_SCALE_USER       = 3} MatLMVMSymBroydenScaleType;
752: PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[];

754: PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
755: PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
756: PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
757: PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
758: PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
759: PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
760: PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType);
761: PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*);
762: PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*);
763: PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*);
764: PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatSchurComplementAinvType,MatReuse,Mat*);
765: PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*);

767: PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm,PetscInt,PetscInt,Mat*);
768: PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm,PetscInt,PetscInt,Mat*);
769: PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm,PetscInt,PetscInt,Mat*);
770: PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
771: PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
772: PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
773: PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);
774: PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm,PetscInt,PetscInt,Mat*);

776: PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
777: PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool*);
778: PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
779: PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
780: PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
781: PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
782: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
783: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
784: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
785: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
786: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
787: PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
788: PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
789: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat*);
790: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC*);
791: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP*);
792: PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt);
793: PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt*);
794: PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt*);
795: PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar);
796: PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType);

798: PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
799: PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool );
800: PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
801: PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
802: PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
803: PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void*);
804: PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
805: PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
806: PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
807: PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*);
808: PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
809: PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
810: PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
811: PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);

813: PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec);
814: PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec,
815:                                            void (**)(PetscInt, PetscInt, PetscInt,
816:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
817:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
818:                                                      PetscReal, const PetscReal [], PetscInt, const PetscScalar[], PetscScalar []), InsertMode, Vec);
819: #endif