Actual source code: petscksp.h

  1: /*
  2:    Defines the interface functions for the Krylov subspace accelerators.
  3: */
  4: #ifndef __PETSCKSP_H
 6:  #include petscpc.h

  9: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPInitializePackage(const char[]);

 11: /*S
 12:      KSP - Abstract PETSc object that manages all Krylov methods

 14:    Level: beginner

 16:   Concepts: Krylov methods

 18: .seealso:  KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP
 19: S*/
 20: typedef struct _p_KSP*     KSP;

 22: /*E
 23:     KSPType - String with the name of a PETSc Krylov method or the creation function
 24:        with an optional dynamic library name, for example
 25:        http://www.mcs.anl.gov/petsc/lib.a:mykspcreate()

 27:    Level: beginner

 29: .seealso: KSPSetType(), KSP
 30: E*/
 31: #define KSPRICHARDSON "richardson"
 32: #define KSPCHEBYCHEV  "chebychev"
 33: #define KSPCG         "cg"
 34: #define   KSPCGNE       "cgne"
 35: #define KSPGMRES      "gmres"
 36: #define   KSPFGMRES     "fgmres" 
 37: #define   KSPLGMRES     "lgmres"
 38: #define KSPTCQMR      "tcqmr"
 39: #define KSPBCGS       "bcgs"
 40: #define KSPBCGSL      "bcgsl"
 41: #define KSPCGS        "cgs"
 42: #define KSPTFQMR      "tfqmr"
 43: #define KSPCR         "cr"
 44: #define KSPLSQR       "lsqr"
 45: #define KSPPREONLY    "preonly"
 46: #define KSPQCG        "qcg"
 47: #define KSPBICG       "bicg"
 48: #define KSPMINRES     "minres"
 49: #define KSPSYMMLQ     "symmlq"
 50: #define KSPType const char*

 52: /* Logging support */

 56: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCreate(MPI_Comm,KSP *);
 57: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetType(KSP,KSPType);
 58: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUp(KSP);
 59: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUpOnBlocks(KSP);
 60: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSolve(KSP,Vec,Vec);
 61: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSolveTranspose(KSP,Vec,Vec);
 62: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDestroy(KSP);

 65: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterAll(const char[]);
 66: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterDestroy(void);
 67: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegister(const char[],const char[],const char[],PetscErrorCode (*)(KSP));

 69: /*MC
 70:    KSPRegisterDynamic - Adds a method to the Krylov subspace solver package.

 72:    Synopsis:
 73:    PetscErrorCode KSPRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(KSP))

 75:    Not Collective

 77:    Input Parameters:
 78: +  name_solver - name of a new user-defined solver
 79: .  path - path (either absolute or relative) the library containing this solver
 80: .  name_create - name of routine to create method context
 81: -  routine_create - routine to create method context

 83:    Notes:
 84:    KSPRegisterDynamic() may be called multiple times to add several user-defined solvers.

 86:    If dynamic libraries are used, then the fourth input argument (routine_create)
 87:    is ignored.

 89:    Sample usage:
 90: .vb
 91:    KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
 92:                "MySolverCreate",MySolverCreate);
 93: .ve

 95:    Then, your solver can be chosen with the procedural interface via
 96: $     KSPSetType(ksp,"my_solver")
 97:    or at runtime via the option
 98: $     -ksp_type my_solver

100:    Level: advanced

102:    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
103:           and others of the form ${any_environmental_variable} occuring in pathname will be 
104:           replaced with appropriate values.
105:          If your function is not being put into a shared library then use KSPRegister() instead

107: .keywords: KSP, register

109: .seealso: KSPRegisterAll(), KSPRegisterDestroy()

111: M*/
112: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
113: #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0)
114: #else
115: #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d)
116: #endif

118: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetType(KSP,KSPType *);
119: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPreconditionerSide(KSP,PCSide);
120: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPreconditionerSide(KSP,PCSide*);
121: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
122: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
123: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetInitialGuessNonzero(KSP,PetscTruth);
124: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetInitialGuessNonzero(KSP,PetscTruth *);
125: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetInitialGuessKnoll(KSP,PetscTruth);
126: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetInitialGuessKnoll(KSP,PetscTruth*);
127: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetComputeEigenvalues(KSP,PetscTruth*);
128: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeEigenvalues(KSP,PetscTruth);
129: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetComputeSingularValues(KSP,PetscTruth*);
130: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeSingularValues(KSP,PetscTruth);
131: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetRhs(KSP,Vec *);
132: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetSolution(KSP,Vec *);
133: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualNorm(KSP,PetscReal*);
134: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetIterationNumber(KSP,PetscInt*);
135: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNullSpace(KSP,MatNullSpace);
136: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetNullSpace(KSP,MatNullSpace*);

138: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPC(KSP,PC);
139: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPC(KSP,PC*);

141: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetMonitor(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void*));
142: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPClearMonitor(KSP);
143: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetMonitorContext(KSP,void **);
144: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *);
145: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscTruth);

147: /* not sure where to put this */
148: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCKSPGetKSP(PC,KSP*);
149: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
150: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
151: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);

153: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildSolution(KSP,Vec,Vec *);
154: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildResidual(KSP,Vec,Vec,Vec *);

156: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRichardsonSetScale(KSP,PetscReal);
157: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal);
158: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
159: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *);
160: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*);

162: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetRestart(KSP, PetscInt);
163: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetHapTol(KSP,PetscReal);

165: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetPreAllocateVectors(KSP);
166: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
167: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
168: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);

170: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetAugDim(KSP,PetscInt);
171: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetConstant(KSP);

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

176:    Level: advanced

178: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
179:           KSPGMRESSetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()

181: E*/
182: typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
184: /*M
185:     KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process

187:    Level: advanced

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

191: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
192:           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
193:           KSPGMRESModifiedGramSchmidtOrthogonalization()
194: M*/

196: /*M
197:     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 
198:           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
199:           poor orthogonality.

201:    Level: advanced

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

206: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
207:           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
208:           KSPGMRESModifiedGramSchmidtOrthogonalization()
209: M*/

211: /*M
212:     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.

214:    Level: advanced

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

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

221: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
222:           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
223:           KSPGMRESModifiedGramSchmidtOrthogonalization()
224: M*/

226: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);

228: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
229: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
230: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));

232: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGSetTrustRegionRadius(KSP,PetscReal);
233: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetQuadratic(KSP,PetscReal*);
234: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetTrialStepNorm(KSP,PetscReal*);

236: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetXRes(KSP,PetscReal);
237: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetPol(KSP,PetscTruth);
238: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetEll(KSP,int);

240: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetFromOptions(KSP);
241: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAddOptionsChecker(PetscErrorCode (*)(KSP));

243: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSingularValueMonitor(KSP,PetscInt,PetscReal,void *);
244: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultMonitor(KSP,PetscInt,PetscReal,void *);
245: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPTrueMonitor(KSP,PetscInt,PetscReal,void *);
246: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultSMonitor(KSP,PetscInt,PetscReal,void *);
247: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPVecViewMonitor(KSP,PetscInt,PetscReal,void *);
248: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESKrylovMonitor(KSP,PetscInt,PetscReal,void *);

250: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPUnwindPreconditioner(KSP,Vec,Vec);
251: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildSolution(KSP,Vec,Vec*);
252: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
253: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);

255: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOperators(KSP,Mat,Mat,MatStructure);
256: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperators(KSP,Mat*,Mat*,MatStructure*);
257: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOptionsPrefix(KSP,const char[]);
258: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAppendOptionsPrefix(KSP,const char[]);
259: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOptionsPrefix(KSP,const char*[]);

261: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScale(KSP,PetscTruth);
262: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScale(KSP,PetscTruth*);
263: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScaleFix(KSP,PetscTruth);
264: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScaleFix(KSP,PetscTruth*);

266: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPView(KSP,PetscViewer);

268: /*E
269:     KSPNormType - Norm that is passed in the Krylov convergence
270:        test routines.

272:    Level: advanced

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

276: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
277:           KSPSetConvergenceTest()
278: E*/
279: typedef enum {KSP_NO_NORM               = 0,
280:               KSP_PRECONDITIONED_NORM   = 1,
281:               KSP_UNPRECONDITIONED_NORM = 2,
282:               KSP_NATURAL_NORM          = 3} KSPNormType;
284: /*M
285:     KSP_NO_NORM - Do not compute a norm during the Krylov process. This will 
286:           possibly save some computation but means the convergence test cannot
287:           be based on a norm of a residual etc.

289:    Level: advanced

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

293: .seealso: KSPNormType, KSPSetNormType(), KSP_PRECONDITIONED_NORM, KSP_UNPRECONDITIONED_NORM, KSP_NATURAL_NORM
294: M*/

296: /*M
297:     KSP_PRECONDITIONED_NORM - Compute the norm of the preconditioned residual and pass that to the 
298:        convergence test routine.

300:    Level: advanced

302: .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_UNPRECONDITIONED_NORM, KSP_NATURAL_NORM, KSPSetConvergenceTest()
303: M*/

305: /*M
306:     KSP_UNPRECONDITIONED_NORM - Compute the norm of the true residual (b - A*x) and pass that to the 
307:        convergence test routine.

309:    Level: advanced

311: .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_PRECONDITIONED_NORM, KSP_NATURAL_NORM, KSPSetConvergenceTest()
312: M*/

314: /*M
315:     KSP_NATURAL_NORM - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 
316:        convergence test routine.

318:    Level: advanced

320: .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_PRECONDITIONED_NORM, KSP_UNPRECONDITIONED_NORM, KSPSetConvergenceTest()
321: M*/

323: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNormType(KSP,KSPNormType);

325: /*E
326:     KSPConvergedReason - reason a Krylov method was said to 
327:          have converged or diverged

329:    Level: beginner

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

333:    Developer note: The string versions of these are in 
334:      src/ksp/ksp/interface/itfunc.c called convergedreasons.
335:      If these enums are changed you must change those.

337: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
338: E*/
339: typedef enum {/* converged */
340:               KSP_CONVERGED_RTOL               =  2,
341:               KSP_CONVERGED_ATOL               =  3,
342:               KSP_CONVERGED_ITS                =  4,
343:               KSP_CONVERGED_QCG_NEG_CURVE      =  5,
344:               KSP_CONVERGED_QCG_CONSTRAINED    =  6,
345:               KSP_CONVERGED_STEP_LENGTH        =  7,
346:               /* diverged */
347:               KSP_DIVERGED_NULL                = -2,
348:               KSP_DIVERGED_ITS                 = -3,
349:               KSP_DIVERGED_DTOL                = -4,
350:               KSP_DIVERGED_BREAKDOWN           = -5,
351:               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
352:               KSP_DIVERGED_NONSYMMETRIC        = -7,
353:               KSP_DIVERGED_INDEFINITE_PC       = -8,
354:               KSP_DIVERGED_NAN                 = -9,
355:               KSP_DIVERGED_INDEFINITE_MAT      = -10,
356: 
357:               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;

360: /*MC
361:      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)

363:    Level: beginner

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

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

371: M*/

373: /*MC
374:      KSP_CONVERGED_ATOL - norm(r) <= atol

376:    Level: beginner

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

382:    Level: beginner

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

386: M*/

388: /*MC
389:      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)

391:    Level: beginner

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

397:    Level: beginner

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

401: M*/

403: /*MC
404:      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 
405:       reached

407:    Level: beginner

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

411: M*/

413: /*MC
414:      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of the
415:            preconditioner is applied.


418:    Level: beginner


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

423: M*/

425: /*MC
426:      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
427:           method could not continue to enlarge the Krylov space.

429:    Level: beginner

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

433: M*/

435: /*MC
436:      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
437:           method could not continue to enlarge the Krylov space.


440:    Level: beginner


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

445: M*/

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

451:    Level: beginner

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

455: M*/

457: /*MC
458:      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
459:         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
460:         be positive definite

462:    Level: beginner

464:      Notes: This can happen with the PCICC preconditioner, use -pc_icc_shift to force 
465:   the PCICC preconditioner to generate a positive definite preconditioner

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

469: M*/

471: /*MC
472:      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
473:         while the KSPSolve() is still running.

475:    Level: beginner

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

479: M*/

481: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *);
482: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergenceContext(KSP,void **);
483: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
484: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
485: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergedReason(KSP,KSPConvergedReason *);

487: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExplicitOperator(KSP,Mat *);

489: /*E
490:     KSPCGType - Determines what type of CG to use

492:    Level: beginner

494: .seealso: KSPCGSetType()
495: E*/
496: typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;

499: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCGSetType(KSP,KSPCGType);

501: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPreSolve(PC,KSP);
502: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPostSolve(PC,KSP);

504: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMonitorCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
505: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMonitor(KSP,PetscInt,PetscReal,void*);
506: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMonitorDestroy(PetscDrawLG);
507: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGTrueMonitorCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
508: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGTrueMonitor(KSP,PetscInt,PetscReal,void*);
509: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGTrueMonitorDestroy(PetscDrawLG);

511: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPostPreSolve(PC,PetscErrorCode (*)(void*,KSP,Vec,Vec));
512: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPostSolve(PC,PetscErrorCode (*)(void*,KSP,Vec,Vec));

515: #endif