Actual source code: petscksp.h

petsc-3.5.4 2015-05-23
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: When a direct solver is used but no Krylov solver is used the KSP object is still used by with a
 19:        KSPType of KSPPREONLY (meaning application of the preconditioner is only used as the linear solver).

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

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

 28:    Level: beginner

 30: .seealso: KSPSetType(), KSP, KSPRegister(), KSPCreate(), KSPSetFromOptions()
 31: J*/
 32: typedef const char* KSPType;
 33: #define KSPRICHARDSON "richardson"
 34: #define KSPCHEBYSHEV  "chebyshev"
 35: #define KSPCG         "cg"
 36: #define KSPGROPPCG    "groppcg"
 37: #define KSPPIPECG     "pipecg"
 38: #define   KSPCGNE       "cgne"
 39: #define   KSPNASH       "nash"
 40: #define   KSPSTCG       "stcg"
 41: #define   KSPGLTR       "gltr"
 42: #define KSPGMRES      "gmres"
 43: #define   KSPFGMRES     "fgmres"
 44: #define   KSPLGMRES     "lgmres"
 45: #define   KSPDGMRES     "dgmres"
 46: #define   KSPPGMRES     "pgmres"
 47: #define KSPTCQMR      "tcqmr"
 48: #define KSPBCGS       "bcgs"
 49: #define   KSPIBCGS      "ibcgs"
 50: #define   KSPFBCGS      "fbcgs"
 51: #define   KSPFBCGSR     "fbcgsr"
 52: #define   KSPBCGSL      "bcgsl"
 53: #define KSPCGS        "cgs"
 54: #define KSPTFQMR      "tfqmr"
 55: #define KSPCR         "cr"
 56: #define KSPPIPECR     "pipecr"
 57: #define KSPLSQR       "lsqr"
 58: #define KSPPREONLY    "preonly"
 59: #define KSPQCG        "qcg"
 60: #define KSPBICG       "bicg"
 61: #define KSPMINRES     "minres"
 62: #define KSPSYMMLQ     "symmlq"
 63: #define KSPLCD        "lcd"
 64: #define KSPPYTHON     "python"
 65: #define KSPGCR        "gcr"
 66: #define KSPSPECEST    "specest"

 68: /* Logging support */
 69: PETSC_EXTERN PetscClassId KSP_CLASSID;
 70: PETSC_EXTERN PetscClassId DMKSP_CLASSID;

 72: PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP *);
 73: PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType);
 74: PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
 75: PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
 76: PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec);
 77: PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec);
 78: PETSC_EXTERN PetscErrorCode KSPReset(KSP);
 79: PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*);
 80: PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool);

 82: PETSC_EXTERN PetscFunctionList KSPList;
 83: PETSC_EXTERN PetscBool         KSPRegisterAllCalled;
 84: PETSC_EXTERN PetscErrorCode KSPRegisterAll(void);
 85: PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP));
 86: PETSC_EXTERN PetscErrorCode KSPMatRegisterAll(void);

 88: PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType *);
 89: PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide);
 90: PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*);
 91: PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
 92: PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
 93: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool );
 94: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool  *);
 95: PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool );
 96: PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool *);
 97: PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool );
 98: PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool  *);
 99: PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool *);
100: PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool );
101: PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool *);
102: PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool );
103: PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec *);
104: PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec *);
105: PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*);
106: PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*);
107: PETSC_EXTERN PetscErrorCode KSPSetNullSpace(KSP,MatNullSpace);
108: PETSC_EXTERN PetscErrorCode KSPGetNullSpace(KSP,MatNullSpace*);
109: PETSC_EXTERN PetscErrorCode KSPGetVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);

111: PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
112: PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);

114: PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC);
115: PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*);

117: PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal);
118: PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void**));
119: PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
120: PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **);
121: PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *);
122: PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool );

124: PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*);
125: PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec *);
126: PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
127: PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt);

129: PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*);
130: PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
131: PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
132: PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
133: PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
134: PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*);
135: PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*);
136: PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*);
137: PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*);
138: PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP *);

140: PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *);
141: PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *);

143: PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
144: PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool );
145: PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal);
146: PETSC_EXTERN PetscErrorCode KSPChebyshevSetEstimateEigenvalues(KSP,PetscReal,PetscReal,PetscReal,PetscReal);
147: PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetRandom(KSP,PetscRandom);
148: PETSC_EXTERN PetscErrorCode KSPChebyshevSetNewMatrix(KSP);
149: PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
150: PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt *);
151: PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]);

153: PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
154: PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*);
155: PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);

157: PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
158: PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
159: PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
160: PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
161: PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);

163: PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
164: PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);

166: PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
167: PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
168: PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));

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

173:    Level: advanced

175: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
176:           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()

178: E*/
179: typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
180: PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
181: /*MC
182:     KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process

184:    Level: advanced

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

188: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
189:           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
190:           KSPGMRESModifiedGramSchmidtOrthogonalization()
191: M*/

193: /*MC
194:     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
195:           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
196:           poor orthogonality.

198:    Level: advanced

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

203: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
204:           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
205:           KSPGMRESModifiedGramSchmidtOrthogonalization()
206: M*/

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

211:    Level: advanced

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

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

218: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
219:           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
220:           KSPGMRESModifiedGramSchmidtOrthogonalization()
221: M*/

223: PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
224: PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);

226: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
227: PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
228: PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));

230: PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
231: PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
232: PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);

234: PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
235: PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool );
236: PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
237: PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool);

239: PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
240: PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));

242: PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *);
243: PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *);
244: PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *);
245: PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *);
246: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
247: PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
248: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *);
249: PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,void *);
250: PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *);
251: PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *);
252: PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*);
253: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**);
254: PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**);
255: PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *);

257: PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
258: PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);

260: PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat);
261: PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*);
262: PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *);
263: PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
264: PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
265: PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);
266: PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt);
267: PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*);

269: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool );
270: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *);
271: PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool );
272: PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *);

274: PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
275: PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);
276: PETSC_STATIC_INLINE PetscErrorCode KSPViewFromOptions(KSP A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);}

278: #define KSP_FILE_CLASSID 1211223

280: PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec);
281: PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);

283: PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
284: PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);

286: /*E
287:     KSPNormType - Norm that is passed in the Krylov convergence
288:        test routines.

290:    Level: advanced

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

295:    Notes: this must match finclude/petscksp.h

297: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
298:           KSPSetConvergenceTest(), KSPSetPCSide()
299: E*/
300: typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
301: #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
302: PETSC_EXTERN const char *const*const KSPNormTypes;

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

309:    Level: advanced

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

313: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
314: M*/

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

320:    Level: advanced

322: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
323: M*/

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

329:    Level: advanced

331: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
332: M*/

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

338:    Level: advanced

340: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
341: M*/

343: PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
344: PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
345: PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
346: PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
347: PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool);

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

352:    Level: beginner

354:    Notes: See KSPGetConvergedReason() for explanation of each value

356:    Developer notes: this must match finclude/petscksp.h

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

361: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
362: E*/
363: typedef enum {/* converged */
364:               KSP_CONVERGED_RTOL_NORMAL        =  1,
365:               KSP_CONVERGED_ATOL_NORMAL        =  9,
366:               KSP_CONVERGED_RTOL               =  2,
367:               KSP_CONVERGED_ATOL               =  3,
368:               KSP_CONVERGED_ITS                =  4,
369:               KSP_CONVERGED_CG_NEG_CURVE       =  5,
370:               KSP_CONVERGED_CG_CONSTRAINED     =  6,
371:               KSP_CONVERGED_STEP_LENGTH        =  7,
372:               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
373:               /* diverged */
374:               KSP_DIVERGED_NULL                = -2,
375:               KSP_DIVERGED_ITS                 = -3,
376:               KSP_DIVERGED_DTOL                = -4,
377:               KSP_DIVERGED_BREAKDOWN           = -5,
378:               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
379:               KSP_DIVERGED_NONSYMMETRIC        = -7,
380:               KSP_DIVERGED_INDEFINITE_PC       = -8,
381:               KSP_DIVERGED_NANORINF            = -9,
382:               KSP_DIVERGED_INDEFINITE_MAT      = -10,

384:               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
385: PETSC_EXTERN const char *const*KSPConvergedReasons;

387: /*MC
388:      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)

390:    Level: beginner

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

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

398: M*/

400: /*MC
401:      KSP_CONVERGED_ATOL - norm(r) <= atol

403:    Level: beginner

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

409:    Level: beginner

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

413: M*/

415: /*MC
416:      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)

418:    Level: beginner

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

424:    Level: beginner

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

428: M*/

430: /*MC
431:      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
432:       reached

434:    Level: beginner

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

438: M*/

440: /*MC
441:      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
442:            the preconditioner is applied. Also used when the KSPConvergedSkip() convergence
443:            test routine is set in KSP.

445:    Level: beginner

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

449: M*/

451: /*MC
452:      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
453:           method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or
454:           preconditioner.

456:    Level: beginner

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

460: M*/

462: /*MC
463:      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
464:           method could not continue to enlarge the Krylov space.

466:    Level: beginner

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

470: M*/

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

476:    Level: beginner

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

480: M*/

482: /*MC
483:      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
484:         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
485:         be positive definite

487:    Level: beginner

489:      Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
490:   the PCICC preconditioner to generate a positive definite preconditioner

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

494: M*/

496: /*MC
497:      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
498:         while the KSPSolve() is still running.

500:    Level: beginner

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

504: M*/

506: PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*));
507: PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **);
508: PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
509: PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
510: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *);
511: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **);
512: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
513: PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
514: PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
515: PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *);

517: PETSC_DEPRECATED("Use KSPConvergedDefault()") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ }
518: #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
519: PETSC_DEPRECATED("Use KSPConvergedDefaultDestroy()") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ }
520: #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
521: PETSC_DEPRECATED("Use KSPConvergedDefaultCreate()") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ }
522: #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
523: PETSC_DEPRECATED("Use KSPConvergedDefaultSetUIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ }
524: #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
525: PETSC_DEPRECATED("Use KSPConvergedDefaultSetUMIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ }
526: #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
527: PETSC_DEPRECATED("Use KSPConvergedSkip()") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ }
528: #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)

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

532: /*E
533:     KSPCGType - Determines what type of CG to use

535:    Level: beginner

537: .seealso: KSPCGSetType()
538: E*/
539: typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
540: PETSC_EXTERN const char *const KSPCGTypes[];

542: PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
543: PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool );

545: PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal);
546: PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *);
547: PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *);

549: PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal);
550: PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *);
551: PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *);

553: PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal);
554: PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *);
555: PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *);
556: PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *);
557: PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *);

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

561: PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
562: PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);

564: #include <petscdrawtypes.h>
565: PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
566: PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*);
567: PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormDestroy(PetscDrawLG*);
568: PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
569: PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
570: PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG*);
571: PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);

573: PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
574: PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));

576: /* see src/ksp/ksp/interface/iguess.c */
577: typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool  monitor;Mat mat; KSP ksp;}* KSPFischerGuess;

579: PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*);
580: PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*);
581: PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess);
582: PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec);
583: PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec);
584: PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess);

586: PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
587: PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess);
588: PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*);

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

593:     Level: intermediate

595: .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat()
596: E*/
597: typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP} MatSchurComplementAinvType;
598: PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];

600: PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
601: PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
602: PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
603: PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
604: PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
605: PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
606: PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType);
607: PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*);
608: PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*);
609: PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*);
610: PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatSchurComplementAinvType,MatReuse,Mat *);
611: PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*);

613: PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
614: PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool );
615: PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
616: PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
617: PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
618: PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *);
619: PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
620: PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
621: PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
622: PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*);
623: PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
624: PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
625: PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
626: PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);

628: #endif