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  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 KSPType const char*
 32: #define KSPRICHARDSON "richardson"
 33: #define KSPCHEBYCHEV  "chebychev"
 34: #define KSPCG         "cg"
 35: #define   KSPCGNE       "cgne"
 36: #define   KSPSTCG       "stcg"
 37: #define   KSPGLTR       "gltr"
 38: #define KSPGMRES      "gmres"
 39: #define   KSPFGMRES     "fgmres" 
 40: #define   KSPLGMRES     "lgmres"
 41: #define KSPTCQMR      "tcqmr"
 42: #define KSPBCGS       "bcgs"
 43: #define KSPBCGSL      "bcgsl"
 44: #define KSPCGS        "cgs"
 45: #define KSPTFQMR      "tfqmr"
 46: #define KSPCR         "cr"
 47: #define KSPLSQR       "lsqr"
 48: #define KSPPREONLY    "preonly"
 49: #define KSPQCG        "qcg"
 50: #define KSPBICG       "bicg"
 51: #define KSPMINRES     "minres"
 52: #define KSPSYMMLQ     "symmlq"
 53: #define KSPLCD        "lcd"

 55: /* Logging support */

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

 67: EXTERN PetscErrorCode  KSPRegisterAll(const char[]);
 68: EXTERN PetscErrorCode  KSPRegisterDestroy(void);
 69: EXTERN PetscErrorCode  KSPRegister(const char[],const char[],const char[],PetscErrorCode (*)(KSP));

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

 74:    Synopsis:
 75:    PetscErrorCode KSPRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(KSP))

 77:    Not Collective

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

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

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

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

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

102:    Level: advanced

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

109: .keywords: KSP, register

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

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

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

141: EXTERN PetscErrorCode  KSPSetPC(KSP,PC);
142: EXTERN PetscErrorCode  KSPGetPC(KSP,PC*);

144: EXTERN PetscErrorCode  KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void*));
145: EXTERN PetscErrorCode  KSPMonitorCancel(KSP);
146: EXTERN PetscErrorCode  KSPGetMonitorContext(KSP,void **);
147: EXTERN PetscErrorCode  KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *);
148: EXTERN PetscErrorCode  KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscTruth);

150: /* not sure where to put this */
151: EXTERN PetscErrorCode  PCKSPGetKSP(PC,KSP*);
152: EXTERN PetscErrorCode  PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
153: EXTERN PetscErrorCode  PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
154: EXTERN PetscErrorCode  PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);

156: EXTERN PetscErrorCode  PCGalerkinGetKSP(PC,KSP *);

158: EXTERN PetscErrorCode  KSPBuildSolution(KSP,Vec,Vec *);
159: EXTERN PetscErrorCode  KSPBuildResidual(KSP,Vec,Vec,Vec *);

161: EXTERN PetscErrorCode  KSPRichardsonSetScale(KSP,PetscReal);
162: EXTERN PetscErrorCode  KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal);
163: EXTERN PetscErrorCode  KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
164: EXTERN PetscErrorCode  KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *);
165: EXTERN PetscErrorCode  KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*);

167: EXTERN PetscErrorCode  KSPGMRESSetRestart(KSP, PetscInt);
168: EXTERN PetscErrorCode  KSPGMRESSetHapTol(KSP,PetscReal);

170: EXTERN PetscErrorCode  KSPGMRESSetPreAllocateVectors(KSP);
171: EXTERN PetscErrorCode  KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
172: EXTERN PetscErrorCode  KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
173: EXTERN PetscErrorCode  KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);

175: EXTERN PetscErrorCode  KSPLGMRESSetAugDim(KSP,PetscInt);
176: EXTERN PetscErrorCode  KSPLGMRESSetConstant(KSP);

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

181:    Level: advanced

183: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
184:           KSPGMRESSetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()

186: E*/
187: typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
189: /*MC
190:     KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process

192:    Level: advanced

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

196: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
197:           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
198:           KSPGMRESModifiedGramSchmidtOrthogonalization()
199: M*/

201: /*MC
202:     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 
203:           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
204:           poor orthogonality.

206:    Level: advanced

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

211: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
212:           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
213:           KSPGMRESModifiedGramSchmidtOrthogonalization()
214: M*/

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

219:    Level: advanced

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

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

226: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
227:           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
228:           KSPGMRESModifiedGramSchmidtOrthogonalization()
229: M*/

231: EXTERN PetscErrorCode  KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);

233: EXTERN PetscErrorCode  KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
234: EXTERN PetscErrorCode  KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
235: EXTERN PetscErrorCode  KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));

237: EXTERN PetscErrorCode  KSPQCGSetTrustRegionRadius(KSP,PetscReal);
238: EXTERN PetscErrorCode  KSPQCGGetQuadratic(KSP,PetscReal*);
239: EXTERN PetscErrorCode  KSPQCGGetTrialStepNorm(KSP,PetscReal*);

241: EXTERN PetscErrorCode  KSPBCGSLSetXRes(KSP,PetscReal);
242: EXTERN PetscErrorCode  KSPBCGSLSetPol(KSP,PetscTruth);
243: EXTERN PetscErrorCode  KSPBCGSLSetEll(KSP,int);

245: EXTERN PetscErrorCode  KSPSetFromOptions(KSP);
246: EXTERN PetscErrorCode  KSPAddOptionsChecker(PetscErrorCode (*)(KSP));

248: EXTERN PetscErrorCode  KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *);
249: EXTERN PetscErrorCode  KSPMonitorDefault(KSP,PetscInt,PetscReal,void *);
250: EXTERN PetscErrorCode  KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *);
251: EXTERN PetscErrorCode  KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *);
252: EXTERN PetscErrorCode  KSPMonitorSolution(KSP,PetscInt,PetscReal,void *);
253: EXTERN PetscErrorCode  KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *);

255: EXTERN PetscErrorCode  KSPUnwindPreconditioner(KSP,Vec,Vec);
256: EXTERN PetscErrorCode  KSPDefaultBuildSolution(KSP,Vec,Vec*);
257: EXTERN PetscErrorCode  KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
258: EXTERN PetscErrorCode  KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);

260: EXTERN PetscErrorCode  KSPSetOperators(KSP,Mat,Mat,MatStructure);
261: EXTERN PetscErrorCode  KSPGetOperators(KSP,Mat*,Mat*,MatStructure*);
262: EXTERN PetscErrorCode  KSPGetOperatorsSet(KSP,PetscTruth*,PetscTruth*);
263: EXTERN PetscErrorCode  KSPSetOptionsPrefix(KSP,const char[]);
264: EXTERN PetscErrorCode  KSPAppendOptionsPrefix(KSP,const char[]);
265: EXTERN PetscErrorCode  KSPGetOptionsPrefix(KSP,const char*[]);

267: EXTERN PetscErrorCode  KSPSetDiagonalScale(KSP,PetscTruth);
268: EXTERN PetscErrorCode  KSPGetDiagonalScale(KSP,PetscTruth*);
269: EXTERN PetscErrorCode  KSPSetDiagonalScaleFix(KSP,PetscTruth);
270: EXTERN PetscErrorCode  KSPGetDiagonalScaleFix(KSP,PetscTruth*);

272: EXTERN PetscErrorCode  KSPView(KSP,PetscViewer);

274: EXTERN PetscErrorCode  KSPLSQRSetStandardErrorVec(KSP,Vec);
275: EXTERN PetscErrorCode  KSPLSQRGetStandardErrorVec(KSP,Vec*);

277: /*E
278:     KSPNormType - Norm that is passed in the Krylov convergence
279:        test routines.

281:    Level: advanced

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

285: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
286:           KSPSetConvergenceTest()
287: E*/
288: typedef enum {KSP_NORM_NO = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
290: /*MC
291:     KSP_NORM_NO - Do not compute a norm during the Krylov process. This will 
292:           possibly save some computation but means the convergence test cannot
293:           be based on a norm of a residual etc.

295:    Level: advanced

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

299: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
300: M*/

302: /*MC
303:     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the 
304:        convergence test routine.

306:    Level: advanced

308: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
309: M*/

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

315:    Level: advanced

317: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
318: M*/

320: /*MC
321:     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 
322:        convergence test routine.

324:    Level: advanced

326: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
327: M*/

329: EXTERN PetscErrorCode  KSPSetNormType(KSP,KSPNormType);
330: EXTERN PetscErrorCode  KSPGetNormType(KSP,KSPNormType*);

332: /*E
333:     KSPConvergedReason - reason a Krylov method was said to 
334:          have converged or diverged

336:    Level: beginner

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

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

344: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
345: E*/
346: typedef enum {/* converged */
347:               KSP_CONVERGED_RTOL               =  2,
348:               KSP_CONVERGED_ATOL               =  3,
349:               KSP_CONVERGED_ITS                =  4,
350:               KSP_CONVERGED_CG_NEG_CURVE       =  5,
351:               KSP_CONVERGED_CG_CONSTRAINED     =  6,
352:               KSP_CONVERGED_STEP_LENGTH        =  7,
353:               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
354:               /* diverged */
355:               KSP_DIVERGED_NULL                = -2,
356:               KSP_DIVERGED_ITS                 = -3,
357:               KSP_DIVERGED_DTOL                = -4,
358:               KSP_DIVERGED_BREAKDOWN           = -5,
359:               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
360:               KSP_DIVERGED_NONSYMMETRIC        = -7,
361:               KSP_DIVERGED_INDEFINITE_PC       = -8,
362:               KSP_DIVERGED_NAN                 = -9,
363:               KSP_DIVERGED_INDEFINITE_MAT      = -10,
364: 
365:               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;

368: /*MC
369:      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)

371:    Level: beginner

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

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

379: M*/

381: /*MC
382:      KSP_CONVERGED_ATOL - norm(r) <= atol

384:    Level: beginner

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

390:    Level: beginner

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

394: M*/

396: /*MC
397:      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)

399:    Level: beginner

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

405:    Level: beginner

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

409: M*/

411: /*MC
412:      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 
413:       reached

415:    Level: beginner

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

419: M*/

421: /*MC
422:      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of the
423:            preconditioner is applied.


426:    Level: beginner


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

431: M*/

433: /*MC
434:      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
435:           method could not continue to enlarge the Krylov space.

437:    Level: beginner

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

441: M*/

443: /*MC
444:      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
445:           method could not continue to enlarge the Krylov space.


448:    Level: beginner


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

453: M*/

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

459:    Level: beginner

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

463: M*/

465: /*MC
466:      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
467:         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
468:         be positive definite

470:    Level: beginner

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

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

477: M*/

479: /*MC
480:      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
481:         while the KSPSolve() is still running.

483:    Level: beginner

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

487: M*/

489: EXTERN PetscErrorCode  KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *);
490: EXTERN PetscErrorCode  KSPGetConvergenceContext(KSP,void **);
491: EXTERN PetscErrorCode  KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
492: EXTERN PetscErrorCode  KSPDefaultConvergedSetUIRNorm(KSP);
493: EXTERN PetscErrorCode  KSPDefaultConvergedSetUMIRNorm(KSP);
494: EXTERN PetscErrorCode  KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
495: EXTERN PetscErrorCode  KSPGetConvergedReason(KSP,KSPConvergedReason *);

497: EXTERN PetscErrorCode  KSPComputeExplicitOperator(KSP,Mat *);

499: /*E
500:     KSPCGType - Determines what type of CG to use

502:    Level: beginner

504: .seealso: KSPCGSetType()
505: E*/
506: typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;

509: EXTERN PetscErrorCode  KSPCGSetType(KSP,KSPCGType);

511: EXTERN PetscErrorCode  KSPSTCGSetRadius(KSP,PetscReal);
512: EXTERN PetscErrorCode  KSPSTCGGetNormD(KSP,PetscReal *);
513: EXTERN PetscErrorCode  KSPSTCGGetObjFcn(KSP,PetscReal *);

515: EXTERN PetscErrorCode  KSPGLTRSetRadius(KSP,PetscReal);
516: EXTERN PetscErrorCode  KSPGLTRGetNormD(KSP,PetscReal *);
517: EXTERN PetscErrorCode  KSPGLTRGetObjFcn(KSP,PetscReal *);
518: EXTERN PetscErrorCode  KSPGLTRGetMinEig(KSP,PetscReal *);
519: EXTERN PetscErrorCode  KSPGLTRGetLambda(KSP,PetscReal *);

521: EXTERN PetscErrorCode  PCPreSolve(PC,KSP);
522: EXTERN PetscErrorCode  PCPostSolve(PC,KSP);

524: EXTERN PetscErrorCode  KSPMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
525: EXTERN PetscErrorCode  KSPMonitorLG(KSP,PetscInt,PetscReal,void*);
526: EXTERN PetscErrorCode  KSPMonitorLGDestroy(PetscDrawLG);
527: EXTERN PetscErrorCode  KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
528: EXTERN PetscErrorCode  KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
529: EXTERN PetscErrorCode  KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG);

531: EXTERN PetscErrorCode  PCShellSetPreSolve(PC,PetscErrorCode (*)(void*,KSP,Vec,Vec));
532: EXTERN PetscErrorCode  PCShellSetPostSolve(PC,PetscErrorCode (*)(void*,KSP,Vec,Vec));

534: EXTERN PetscErrorCode  MatMFFDKSPMonitor(KSP,PetscInt,PetscReal,void *);

537: #endif