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

 58: /* Logging support */

 61: EXTERN PetscErrorCode  KSPCreate(MPI_Comm,KSP *);
 62: EXTERN PetscErrorCode  KSPSetType(KSP,const KSPType);
 63: EXTERN PetscErrorCode  KSPSetUp(KSP);
 64: EXTERN PetscErrorCode  KSPSetUpOnBlocks(KSP);
 65: EXTERN PetscErrorCode  KSPSolve(KSP,Vec,Vec);
 66: EXTERN PetscErrorCode  KSPSolveTranspose(KSP,Vec,Vec);
 67: EXTERN PetscErrorCode  KSPDestroy(KSP);

 70: EXTERN PetscErrorCode  KSPRegisterAll(const char[]);
 71: EXTERN PetscErrorCode  KSPRegisterDestroy(void);
 72: EXTERN PetscErrorCode  KSPRegister(const char[],const char[],const char[],PetscErrorCode (*)(KSP));

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

 77:    Synopsis:
 78:    PetscErrorCode KSPRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(KSP))

 80:    Not Collective

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

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

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

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

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

105:    Level: advanced

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

112: .keywords: KSP, register

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

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

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

144: EXTERN PetscErrorCode  KSPSetPC(KSP,PC);
145: EXTERN PetscErrorCode  KSPGetPC(KSP,PC*);

147: EXTERN PetscErrorCode  KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void*));
148: EXTERN PetscErrorCode  KSPMonitorCancel(KSP);
149: EXTERN PetscErrorCode  KSPGetMonitorContext(KSP,void **);
150: EXTERN PetscErrorCode  KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *);
151: EXTERN PetscErrorCode  KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscTruth);

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

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

161: EXTERN PetscErrorCode  KSPBuildSolution(KSP,Vec,Vec *);
162: EXTERN PetscErrorCode  KSPBuildResidual(KSP,Vec,Vec,Vec *);

164: EXTERN PetscErrorCode  KSPRichardsonSetScale(KSP,PetscReal);
165: EXTERN PetscErrorCode  KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal);
166: EXTERN PetscErrorCode  KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
167: EXTERN PetscErrorCode  KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *);
168: EXTERN PetscErrorCode  KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*);

170: EXTERN PetscErrorCode  KSPGMRESSetRestart(KSP, PetscInt);
171: EXTERN PetscErrorCode  KSPGMRESSetHapTol(KSP,PetscReal);

173: EXTERN PetscErrorCode  KSPGMRESSetPreAllocateVectors(KSP);
174: EXTERN PetscErrorCode  KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
175: EXTERN PetscErrorCode  KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
176: EXTERN PetscErrorCode  KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);

178: EXTERN PetscErrorCode  KSPLGMRESSetAugDim(KSP,PetscInt);
179: EXTERN PetscErrorCode  KSPLGMRESSetConstant(KSP);

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

184:    Level: advanced

186: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
187:           KSPGMRESSetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()

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

195:    Level: advanced

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

199: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
200:           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
201:           KSPGMRESModifiedGramSchmidtOrthogonalization()
202: M*/

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

209:    Level: advanced

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

214: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
215:           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
216:           KSPGMRESModifiedGramSchmidtOrthogonalization()
217: M*/

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

222:    Level: advanced

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

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

229: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
230:           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
231:           KSPGMRESModifiedGramSchmidtOrthogonalization()
232: M*/

234: EXTERN PetscErrorCode  KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);

236: EXTERN PetscErrorCode  KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
237: EXTERN PetscErrorCode  KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
238: EXTERN PetscErrorCode  KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));

240: EXTERN PetscErrorCode  KSPQCGSetTrustRegionRadius(KSP,PetscReal);
241: EXTERN PetscErrorCode  KSPQCGGetQuadratic(KSP,PetscReal*);
242: EXTERN PetscErrorCode  KSPQCGGetTrialStepNorm(KSP,PetscReal*);

244: EXTERN PetscErrorCode  KSPBCGSLSetXRes(KSP,PetscReal);
245: EXTERN PetscErrorCode  KSPBCGSLSetPol(KSP,PetscTruth);
246: EXTERN PetscErrorCode  KSPBCGSLSetEll(KSP,int);

248: EXTERN PetscErrorCode  KSPSetFromOptions(KSP);
249: EXTERN PetscErrorCode  KSPAddOptionsChecker(PetscErrorCode (*)(KSP));

251: EXTERN PetscErrorCode  KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *);
252: EXTERN PetscErrorCode  KSPMonitorDefault(KSP,PetscInt,PetscReal,void *);
253: EXTERN PetscErrorCode  KSPMonitorRange(KSP,PetscInt,PetscReal,void *);
254: EXTERN PetscErrorCode  KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *);
255: EXTERN PetscErrorCode  KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *);
256: EXTERN PetscErrorCode  KSPMonitorSolution(KSP,PetscInt,PetscReal,void *);
257: EXTERN PetscErrorCode  KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *);

259: EXTERN PetscErrorCode  KSPUnwindPreconditioner(KSP,Vec,Vec);
260: EXTERN PetscErrorCode  KSPDefaultBuildSolution(KSP,Vec,Vec*);
261: EXTERN PetscErrorCode  KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
262: EXTERN PetscErrorCode  KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);

264: EXTERN PetscErrorCode  KSPSetOperators(KSP,Mat,Mat,MatStructure);
265: EXTERN PetscErrorCode  KSPGetOperators(KSP,Mat*,Mat*,MatStructure*);
266: EXTERN PetscErrorCode  KSPGetOperatorsSet(KSP,PetscTruth*,PetscTruth*);
267: EXTERN PetscErrorCode  KSPSetOptionsPrefix(KSP,const char[]);
268: EXTERN PetscErrorCode  KSPAppendOptionsPrefix(KSP,const char[]);
269: EXTERN PetscErrorCode  KSPGetOptionsPrefix(KSP,const char*[]);

271: EXTERN PetscErrorCode  KSPSetDiagonalScale(KSP,PetscTruth);
272: EXTERN PetscErrorCode  KSPGetDiagonalScale(KSP,PetscTruth*);
273: EXTERN PetscErrorCode  KSPSetDiagonalScaleFix(KSP,PetscTruth);
274: EXTERN PetscErrorCode  KSPGetDiagonalScaleFix(KSP,PetscTruth*);

276: EXTERN PetscErrorCode  KSPView(KSP,PetscViewer);

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

281: /*E
282:     KSPNormType - Norm that is passed in the Krylov convergence
283:        test routines.

285:    Level: advanced

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

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

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

302:    Level: advanced

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

306: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
307: M*/

309: /*MC
310:     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the 
311:        convergence test routine.

313:    Level: advanced

315: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
316: M*/

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

322:    Level: advanced

324: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
325: M*/

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

331:    Level: advanced

333: .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
334: M*/

336: EXTERN PetscErrorCode  KSPSetNormType(KSP,KSPNormType);
337: EXTERN PetscErrorCode  KSPGetNormType(KSP,KSPNormType*);
338: EXTERN PetscErrorCode  KSPSetCheckNormIteration(KSP,PetscInt);
339: EXTERN PetscErrorCode  KSPSetLagNorm(KSP,PetscTruth);

341: /*E
342:     KSPConvergedReason - reason a Krylov method was said to 
343:          have converged or diverged

345:    Level: beginner

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

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

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

354: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
355: E*/
356: typedef enum {/* converged */
357:               KSP_CONVERGED_RTOL               =  2,
358:               KSP_CONVERGED_ATOL               =  3,
359:               KSP_CONVERGED_ITS                =  4,
360:               KSP_CONVERGED_CG_NEG_CURVE       =  5,
361:               KSP_CONVERGED_CG_CONSTRAINED     =  6,
362:               KSP_CONVERGED_STEP_LENGTH        =  7,
363:               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
364:               /* diverged */
365:               KSP_DIVERGED_NULL                = -2,
366:               KSP_DIVERGED_ITS                 = -3,
367:               KSP_DIVERGED_DTOL                = -4,
368:               KSP_DIVERGED_BREAKDOWN           = -5,
369:               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
370:               KSP_DIVERGED_NONSYMMETRIC        = -7,
371:               KSP_DIVERGED_INDEFINITE_PC       = -8,
372:               KSP_DIVERGED_NAN                 = -9,
373:               KSP_DIVERGED_INDEFINITE_MAT      = -10,
374: 
375:               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;

378: /*MC
379:      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)

381:    Level: beginner

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

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

389: M*/

391: /*MC
392:      KSP_CONVERGED_ATOL - norm(r) <= atol

394:    Level: beginner

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

400:    Level: beginner

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

404: M*/

406: /*MC
407:      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)

409:    Level: beginner

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

415:    Level: beginner

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

419: M*/

421: /*MC
422:      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 
423:       reached

425:    Level: beginner

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

429: M*/

431: /*MC
432:      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 
433:            the preconditioner is applied. Also used when the KSPSkipConverged() convergence 
434:            test routine is set in KSP.


437:    Level: beginner


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

442: M*/

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

448:    Level: beginner

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

452: M*/

454: /*MC
455:      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
456:           method could not continue to enlarge the Krylov space.


459:    Level: beginner


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

464: M*/

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

470:    Level: beginner

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

474: M*/

476: /*MC
477:      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
478:         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
479:         be positive definite

481:    Level: beginner

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

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

488: M*/

490: /*MC
491:      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
492:         while the KSPSolve() is still running.

494:    Level: beginner

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

498: M*/

500: EXTERN PetscErrorCode  KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*));
501: EXTERN PetscErrorCode  KSPGetConvergenceContext(KSP,void **);
502: EXTERN PetscErrorCode  KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
503: EXTERN PetscErrorCode  KSPDefaultConvergedDestroy(void *);
504: EXTERN PetscErrorCode  KSPDefaultConvergedCreate(void **);
505: EXTERN PetscErrorCode  KSPDefaultConvergedSetUIRNorm(KSP);
506: EXTERN PetscErrorCode  KSPDefaultConvergedSetUMIRNorm(KSP);
507: EXTERN PetscErrorCode  KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
508: EXTERN PetscErrorCode  KSPGetConvergedReason(KSP,KSPConvergedReason *);

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

512: /*E
513:     KSPCGType - Determines what type of CG to use

515:    Level: beginner

517: .seealso: KSPCGSetType()
518: E*/
519: typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;

522: EXTERN PetscErrorCode  KSPCGSetType(KSP,KSPCGType);

524: EXTERN PetscErrorCode  KSPNASHSetRadius(KSP,PetscReal);
525: EXTERN PetscErrorCode  KSPNASHGetNormD(KSP,PetscReal *);
526: EXTERN PetscErrorCode  KSPNASHGetObjFcn(KSP,PetscReal *);

528: EXTERN PetscErrorCode  KSPSTCGSetRadius(KSP,PetscReal);
529: EXTERN PetscErrorCode  KSPSTCGGetNormD(KSP,PetscReal *);
530: EXTERN PetscErrorCode  KSPSTCGGetObjFcn(KSP,PetscReal *);

532: EXTERN PetscErrorCode  KSPGLTRSetRadius(KSP,PetscReal);
533: EXTERN PetscErrorCode  KSPGLTRGetNormD(KSP,PetscReal *);
534: EXTERN PetscErrorCode  KSPGLTRGetObjFcn(KSP,PetscReal *);
535: EXTERN PetscErrorCode  KSPGLTRGetMinEig(KSP,PetscReal *);
536: EXTERN PetscErrorCode  KSPGLTRGetLambda(KSP,PetscReal *);

538: EXTERN PetscErrorCode  KSPPythonSetType(KSP,const char[]);

540: EXTERN PetscErrorCode  PCPreSolve(PC,KSP);
541: EXTERN PetscErrorCode  PCPostSolve(PC,KSP);

543: EXTERN PetscErrorCode  KSPMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
544: EXTERN PetscErrorCode  KSPMonitorLG(KSP,PetscInt,PetscReal,void*);
545: EXTERN PetscErrorCode  KSPMonitorLGDestroy(PetscDrawLG);
546: EXTERN PetscErrorCode  KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
547: EXTERN PetscErrorCode  KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
548: EXTERN PetscErrorCode  KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG);
549: EXTERN PetscErrorCode  KSPMonitorLGRangeCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
550: EXTERN PetscErrorCode  KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
551: EXTERN PetscErrorCode  KSPMonitorLGRangeDestroy(PetscDrawLG);

553: EXTERN PetscErrorCode  PCShellSetPreSolve(PC,PetscErrorCode (*)(void*,KSP,Vec,Vec));
554: EXTERN PetscErrorCode  PCShellSetPostSolve(PC,PetscErrorCode (*)(void*,KSP,Vec,Vec));

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

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

561: EXTERN PetscErrorCode  KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*);
562: EXTERN PetscErrorCode  KSPFischerGuessDestroy(KSPFischerGuess);
563: EXTERN PetscErrorCode  KSPFischerGuessReset(KSPFischerGuess);
564: EXTERN PetscErrorCode  KSPFischerGuessUpdate(KSPFischerGuess,Vec);
565: EXTERN PetscErrorCode  KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec);
566: EXTERN PetscErrorCode  KSPFischerGuessSetFromOptions(KSPFischerGuess);

568: EXTERN PetscErrorCode  KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
569: EXTERN PetscErrorCode  KSPSetFischerGuess(KSP,KSPFischerGuess);
570: EXTERN PetscErrorCode  KSPGetFischerGuess(KSP,KSPFischerGuess*);

572: EXTERN PetscErrorCode  MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat*);
573: EXTERN PetscErrorCode  MatSchurComplementGetKSP(Mat,KSP*);
574: EXTERN PetscErrorCode  MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,MatStructure);
575: EXTERN PetscErrorCode  MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*);

578: #endif