Actual source code: petscksp.h

petsc-3.4.5 2014-06-29
  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*);

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

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

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

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

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

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

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

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

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

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

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

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

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

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

172:    Level: advanced

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

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

183:    Level: advanced

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

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

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

197:    Level: advanced

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

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

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

210:    Level: advanced

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

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

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

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

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

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

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

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

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

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

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

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

273: PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
274: PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);

276: #define KSP_FILE_CLASSID 1211223

278: PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec);
279: PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);

281: PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
282: PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);

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

288:    Level: advanced

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

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

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

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

307:    Level: advanced

309:     Note: Some Krylov methods need to compute a residual norm and then this option is ignored

311: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
312: M*/

314: /*MC
315:     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the
316:        convergence test routine.

318:    Level: advanced

320: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
321: M*/

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

327:    Level: advanced

329: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
330: M*/

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

336:    Level: advanced

338: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
339: M*/

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

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

351:    Level: beginner

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

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

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

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

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

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

389:    Level: beginner

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

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

397: M*/

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

402:    Level: beginner

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

408:    Level: beginner

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

412: M*/

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

417:    Level: beginner

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

423:    Level: beginner

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

427: M*/

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

433:    Level: beginner

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

437: M*/

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


445:    Level: beginner


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

450: M*/

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

457:    Level: beginner

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

461: M*/

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


468:    Level: beginner


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

473: M*/

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

479:    Level: beginner

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

483: M*/

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

490:    Level: beginner

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

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

497: M*/

499: /*MC
500:      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
501:         while the KSPSolve() is still running.

503:    Level: beginner

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

507: M*/

509: PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*));
510: PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **);
511: PETSC_EXTERN PetscErrorCode KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
512: PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
513: PETSC_EXTERN PetscErrorCode KSPDefaultConvergedDestroy(void *);
514: PETSC_EXTERN PetscErrorCode KSPDefaultConvergedCreate(void **);
515: PETSC_EXTERN PetscErrorCode KSPDefaultConvergedSetUIRNorm(KSP);
516: PETSC_EXTERN PetscErrorCode KSPDefaultConvergedSetUMIRNorm(KSP);
517: PETSC_EXTERN PetscErrorCode KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
518: PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *);

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

522: /*E
523:     KSPCGType - Determines what type of CG to use

525:    Level: beginner

527: .seealso: KSPCGSetType()
528: E*/
529: typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
530: PETSC_EXTERN const char *const KSPCGTypes[];

532: PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
533: PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool );

535: PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal);
536: PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *);
537: PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *);

539: PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal);
540: PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *);
541: PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *);

543: PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal);
544: PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *);
545: PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *);
546: PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *);
547: PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *);

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

551: PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
552: PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);

554: #include <petscdrawtypes.h>
555: PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
556: PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*);
557: PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormDestroy(PetscDrawLG*);
558: PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
559: PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
560: PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG*);
561: PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);

563: PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
564: PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));

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

569: PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*);
570: PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*);
571: PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess);
572: PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec);
573: PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec);
574: PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess);

576: PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
577: PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess);
578: PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*);

580: PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
581: PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
582: PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
583: PETSC_EXTERN PetscErrorCode MatSchurComplementSet(Mat,Mat,Mat,Mat,Mat,Mat);
584: PETSC_EXTERN PetscErrorCode MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,Mat,MatStructure);
585: PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
586: PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatReuse,Mat *);

588: PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
589: PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool );
590: PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
591: PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
592: PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
593: PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *);
594: PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*);
595: PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
596: PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*);
597: PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,MatStructure*,void*),void*);
598: PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
599: PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
600: PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
601: PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);

603: #endif