Actual source code: petscksp.h

petsc-3.10.5 2019-03-28
Report Typos and Errors
  1: /*
  2:    Defines the interface functions for the Krylov subspace accelerators.
  3: */
  4: #ifndef __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:   Concepts: Krylov methods

 18:         Notes:
 19:     When a direct solver is used but no Krylov solver is used the KSP object is still used by with a
 20:        KSPType of KSPPREONLY (meaning application of the preconditioner is only used as the linear solver).

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

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

 29:    Level: beginner

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

 78: /* Logging support */
 79: PETSC_EXTERN PetscClassId KSP_CLASSID;
 80: PETSC_EXTERN PetscClassId KSPGUESS_CLASSID;
 81: PETSC_EXTERN PetscClassId DMKSP_CLASSID;

 83: PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP*);
 84: PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType);
 85: PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType*);
 86: PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
 87: PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
 88: PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec);
 89: PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec);
 90: PETSC_EXTERN PetscErrorCode KSPReset(KSP);
 91: PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*);
 92: PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool);
 93: PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool);

 95: PETSC_EXTERN PetscFunctionList KSPList;
 96: PETSC_EXTERN PetscFunctionList KSPGuessList;
 97: PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP));

 99: PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide);
100: PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*);
101: PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
102: PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
103: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool);
104: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool*);
105: PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool);
106: PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool*);
107: PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool);
108: PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool);
109: PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool*);
110: PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool);
111: PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool*);
112: PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec*);
113: PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec*);
114: PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*);
115: PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*);
116: PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*);
117: PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
118: PETSC_DEPRECATED("Use KSPCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode KSPGetVecs(KSP ksp,PetscInt n,Vec **x,PetscInt m,Vec **y) {return KSPCreateVecs(ksp,n,x,m,y);}

120: PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
121: PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);

123: PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC);
124: PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*);

126: PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal);
127: PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode (*)(void**));
128: PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
129: PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void**);
130: PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt*);
131: PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool );

133: PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*);
134: PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec*);
135: PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
136: PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt);

138: PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*);
139: PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC,KSP);
140: PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
141: PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
142: PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
143: PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
144: PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC,PetscInt*,KSP*[]);
145: PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*);
146: PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*);
147: PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*);
148: PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*);
149: PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP*);

151: PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec*);
152: PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec*);

154: PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
155: PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool );
156: PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal);
157: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal);
158: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP,PetscBool);
159: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*);
160: PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
161: PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt*);
162: PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]);
163: PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt*,Vec[],PetscReal[],PetscReal[]);

165: /*E

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

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

172:    Level: intermediate
173: .seealso : KSPFCG,KSPPIPEFCG,KSPPIPEGCR,KSPFCGSetTruncationType(),KSPFCGGetTruncationType()

175: E*/
176: typedef enum {KSP_FCD_TRUNC_TYPE_STANDARD,KSP_FCD_TRUNC_TYPE_NOTAY} KSPFCDTruncationType;
177: PETSC_EXTERN const char *const KSPFCDTruncationTypes[];

179: PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt);
180: PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*);
181: PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt);
182: PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*);
183: PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType);
184: PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*);

186: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt);
187: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*);
188: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt);
189: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*);
190: PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType);
191: PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*);

193: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt);
194: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*);
195: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt);
196: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*);
197: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType);
198: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*);
199: PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool);
200: PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*);

202: PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
203: PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*);
204: PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);

206: PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
207: PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
208: PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
209: PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
210: PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);

212: PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
213: PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);

215: PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP,PetscScalar);

217: PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
218: PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
219: PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));

221: PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP,PC*);
222: PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP,PC);
223: PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP,KSP*);
224: PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP,Mat);
225: /*E
226:     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.

228:    Level: advanced

230: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
231:           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()

233: E*/
234: typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
235: PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
236: /*MC
237:     KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process

239:    Level: advanced

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

243: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
244:           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
245:           KSPGMRESModifiedGramSchmidtOrthogonalization()
246: M*/

248: /*MC
249:     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
250:           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
251:           poor orthogonality.

253:    Level: advanced

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

258: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
259:           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
260:           KSPGMRESModifiedGramSchmidtOrthogonalization()
261: M*/

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

266:    Level: advanced

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

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

273: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
274:           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
275:           KSPGMRESModifiedGramSchmidtOrthogonalization()
276: M*/

278: PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
279: PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);

281: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
282: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
283: PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));

285: PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
286: PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
287: PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);

289: PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
290: PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool );
291: PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
292: PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool);

294: PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
295: PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);
296: PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));

298: PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
299: PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
300: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
301: PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
302: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
303: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
304: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
305: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
306: PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
307: PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
308: PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*);
309: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**);
310: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**);
311: PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void*);
312: PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP,const char[],const char[],const char [],PetscErrorCode (*)(KSP,PetscInt,PetscReal,PetscViewerAndFormat*));

314: PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
315: PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);

317: PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat);
318: PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*);
319: PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool*,PetscBool*);
320: PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
321: PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
322: PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);
323: PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt);
324: PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*);

326: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool );
327: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool*);
328: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool );
329: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool*);

331: PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
332: PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);
333: PETSC_STATIC_INLINE PetscErrorCode KSPViewFromOptions(KSP A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
334: PETSC_EXTERN PetscErrorCode KSPReasonView(KSP,PetscViewer);
335: PETSC_EXTERN PetscErrorCode KSPReasonViewFromOptions(KSP);

337: #define KSP_FILE_CLASSID 1211223

339: PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP,PetscBool);
340: PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP,PetscBool);
341: PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
342: PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP,PetscReal*,PetscReal*);

344: PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
345: PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
346: PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*);

348: /*E
349:     KSPNormType - Norm that is passed in the Krylov convergence
350:        test routines.

352:    Level: advanced

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

357:    Notes:
358:     this must match petsc/finclude/petscksp.h

360: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
361:           KSPSetConvergenceTest(), KSPSetPCSide()
362: E*/
363: typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
364: #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
365: PETSC_EXTERN const char *const*const KSPNormTypes;

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

372:    Level: advanced

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

376: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
377: M*/

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

383:    Level: advanced

385: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
386: M*/

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

392:    Level: advanced

394: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
395: M*/

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

401:    Level: advanced

403: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
404: M*/

406: PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
407: PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
408: PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
409: PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
410: PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool);

412: /*E
413:     KSPConvergedReason - reason a Krylov method was said to have converged or diverged

415:    Level: beginner

417:    Notes:
418:     See KSPGetConvergedReason() for explanation of each value

420:    Developer Notes:
421:     this must match petsc/finclude/petscksp.h

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

426: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
427: E*/
428: typedef enum {/* converged */
429:               KSP_CONVERGED_RTOL_NORMAL        =  1,
430:               KSP_CONVERGED_ATOL_NORMAL        =  9,
431:               KSP_CONVERGED_RTOL               =  2,
432:               KSP_CONVERGED_ATOL               =  3,
433:               KSP_CONVERGED_ITS                =  4,
434:               KSP_CONVERGED_CG_NEG_CURVE       =  5,
435:               KSP_CONVERGED_CG_CONSTRAINED     =  6,
436:               KSP_CONVERGED_STEP_LENGTH        =  7,
437:               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
438:               /* diverged */
439:               KSP_DIVERGED_NULL                = -2,
440:               KSP_DIVERGED_ITS                 = -3,
441:               KSP_DIVERGED_DTOL                = -4,
442:               KSP_DIVERGED_BREAKDOWN           = -5,
443:               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
444:               KSP_DIVERGED_NONSYMMETRIC        = -7,
445:               KSP_DIVERGED_INDEFINITE_PC       = -8,
446:               KSP_DIVERGED_NANORINF            = -9,
447:               KSP_DIVERGED_INDEFINITE_MAT      = -10,
448:               KSP_DIVERGED_PCSETUP_FAILED      = -11,

450:               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
451: PETSC_EXTERN const char *const*KSPConvergedReasons;

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

456:    Level: beginner

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

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

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

466: M*/

468: /*MC
469:      KSP_CONVERGED_ATOL - norm(r) <= atol

471:    Level: beginner

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

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

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

481: M*/

483: /*MC
484:      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)

486:    Level: beginner

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

492:    Level: beginner

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

496: M*/

498: /*MC
499:      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
500:       reached

502:    Level: beginner

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

506: M*/

508: /*MC
509:      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
510:            the preconditioner is applied. Also used when the KSPConvergedSkip() convergence
511:            test routine is set in KSP.

513:    Level: beginner

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

517: M*/

519: /*MC
520:      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
521:           method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or
522:           preconditioner.

524:    Level: beginner

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

528: M*/

530: /*MC
531:      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
532:           method could not continue to enlarge the Krylov space.

534:    Level: beginner

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

538: M*/

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

544:    Level: beginner

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

548: M*/

550: /*MC
551:      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
552:         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
553:         be positive definite

555:    Level: beginner

557:      Notes:
558:     This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
559:   the PCICC preconditioner to generate a positive definite preconditioner

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

563: M*/

565: /*MC
566:      KSP_DIVERGED_PCSETUP_FAILED - It was not possible to build the requested preconditioner. This is usually due to a 
567:      zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner
568:      such as PCFIELDSPLIT.

570:    Level: beginner

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


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

578: M*/

580: /*MC
581:      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
582:         while the KSPSolve() is still running.

584:    Level: beginner

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

588: M*/

590: PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void*,PetscErrorCode (*)(void*));
591: PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void**);
592: PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
593: PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
594: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void*);
595: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void**);
596: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
597: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
598: PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
599: PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason*);

601: PETSC_DEPRECATED("Use KSPConvergedDefault()") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ }
602: #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
603: PETSC_DEPRECATED("Use KSPConvergedDefaultDestroy()") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ }
604: #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
605: PETSC_DEPRECATED("Use KSPConvergedDefaultCreate()") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ }
606: #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
607: PETSC_DEPRECATED("Use KSPConvergedDefaultSetUIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ }
608: #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
609: PETSC_DEPRECATED("Use KSPConvergedDefaultSetUMIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ }
610: #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
611: PETSC_DEPRECATED("Use KSPConvergedSkip()") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ }
612: #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)

614: PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat*);

616: /*E
617:     KSPCGType - Determines what type of CG to use

619:    Level: beginner

621: .seealso: KSPCGSetType()
622: E*/
623: typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
624: PETSC_EXTERN const char *const KSPCGTypes[];

626: PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
627: PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool );

629: PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP,PetscReal);
630: PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP,PetscReal*);
631: PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP,PetscReal*);

633: PETSC_EXTERN PetscErrorCode KSPCGGLTRGetMinEig(KSP,PetscReal*);
634: PETSC_EXTERN PetscErrorCode KSPCGGLTRGetLambda(KSP,PetscReal*);

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

638: PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
639: PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);

641:  #include <petscdrawtypes.h>
642: PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
643: PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*);
644: PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
645: PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
646: PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);

648: PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
649: PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));

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

654:    Level: beginner

656:   Concepts: Krylov methods

658: .seealso:  KSPCreate(), KSPSetGuessType(), KSPGuessType
659: S*/
660: typedef struct _p_KSPGuess* KSPGuess;
661: /*J
662:     KSPGuessType - String with the name of a PETSc initial guess for Krylov methods.

664:    Level: beginner

666: .seealso: KSPGuess
667: J*/
668: typedef const char* KSPGuessType;
669: #define KSPGUESSFISCHER "fischer"
670: #define KSPGUESSPOD     "pod"
671: PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[],PetscErrorCode (*)(KSPGuess));
672: PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP,KSPGuess);
673: PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP,KSPGuess*);
674: PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess,PetscViewer);
675: PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess*);
676: PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm,KSPGuess*);
677: PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess,KSPGuessType);
678: PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess,KSPGuessType*);
679: PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
680: PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess,Vec,Vec);
681: PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess,Vec,Vec);
682: PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
683: PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess,PetscInt,PetscInt);
684: PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
685: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool);
686: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool*);

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

691:     Level: intermediate

693: .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat()
694: E*/
695: typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG} MatSchurComplementAinvType;
696: PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];

698: PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
699: PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
700: PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
701: PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
702: PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
703: PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
704: PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType);
705: PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*);
706: PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*);
707: PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*);
708: PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatSchurComplementAinvType,MatReuse,Mat*);
709: PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*);

711: PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm,PetscInt,PetscInt,Mat*);
712: PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm,PetscInt,PetscInt,Mat*);
713: PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm,PetscInt,PetscInt,Mat*);
714: PETSC_EXTERN PetscErrorCode MatCreateLMVMBrdn(MPI_Comm,PetscInt,PetscInt,Mat*);
715: PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBrdn(MPI_Comm,PetscInt,PetscInt,Mat*);
716: PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBrdn(MPI_Comm,PetscInt,PetscInt,Mat*);
717: PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBrdn(MPI_Comm,PetscInt,PetscInt,Mat*);
718: PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBrdn(MPI_Comm,PetscInt,PetscInt,Mat*);

720: PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
721: PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool*);
722: PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
723: PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
724: PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
725: PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
726: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
727: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
728: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
729: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
730: PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
731: PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
732: PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
733: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat*);
734: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC*);
735: PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP*);
736: PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt*);
737: PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt*);
738: PETSC_EXTERN PetscErrorCode MatSymBrdnSetDelta(Mat, PetscScalar);

740: PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
741: PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool );
742: PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
743: PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
744: PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
745: PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void*);
746: PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
747: PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
748: PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
749: PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*);
750: PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
751: PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
752: PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
753: PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);

755: PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec);
756: PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec,
757:                                            void (**)(PetscInt, PetscInt, PetscInt,
758:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
759:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
760:                                                      PetscReal, const PetscReal [], PetscInt, const PetscScalar[], PetscScalar []), InsertMode, Vec);
761: #endif