xref: /petsc/include/petscksp.h (revision be1d678a52e6eff2808b2fa31ae986cdbf03c9fe)
1 /*
2    Defines the interface functions for the Krylov subspace accelerators.
3 */
4 #ifndef __PETSCKSP_H
5 #define __PETSCKSP_H
6 #include "petscpc.h"
7 PETSC_EXTERN_CXX_BEGIN
8 
9 EXTERN PetscErrorCode KSPInitializePackage(const char[]);
10 
11 /*S
12      KSP - Abstract PETSc object that manages all Krylov methods
13 
14    Level: beginner
15 
16   Concepts: Krylov methods
17 
18 .seealso:  KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP
19 S*/
20 typedef struct _p_KSP*     KSP;
21 
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()
26 
27    Level: beginner
28 
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 KSPTCQMR      "tcqmr"
37 #define KSPBCGS       "bcgs"
38 #define KSPBCGSL      "bcgsl"
39 #define KSPCGS        "cgs"
40 #define KSPTFQMR      "tfqmr"
41 #define KSPCR         "cr"
42 #define KSPLSQR       "lsqr"
43 #define KSPPREONLY    "preonly"
44 #define KSPQCG        "qcg"
45 #define KSPBICG       "bicg"
46 #define KSPFGMRES     "fgmres"
47 #define KSPMINRES     "minres"
48 #define KSPSYMMLQ     "symmlq"
49 #define KSPLGMRES     "lgmres"
50 #define KSPType char*
51 
52 /* Logging support */
53 extern PetscCookie KSP_COOKIE;
54 extern PetscEvent    KSP_GMRESOrthogonalization, KSP_SetUp, KSP_Solve;
55 
56 EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP *);
57 EXTERN PetscErrorCode KSPSetType(KSP,const KSPType);
58 EXTERN PetscErrorCode KSPSetUp(KSP);
59 EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
60 EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec);
61 EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec);
62 EXTERN PetscErrorCode KSPDestroy(KSP);
63 
64 extern PetscFList KSPList;
65 EXTERN PetscErrorCode KSPRegisterAll(const char[]);
66 EXTERN PetscErrorCode KSPRegisterDestroy(void);
67 
68 EXTERN PetscErrorCode KSPRegister(const char[],const char[],const char[],PetscErrorCode (*)(KSP));
69 
70 /*MC
71    KSPRegisterDynamic - Adds a method to the Krylov subspace solver package.
72 
73    Synopsis:
74    PetscErrorCode KSPRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(KSP))
75 
76    Not Collective
77 
78    Input Parameters:
79 +  name_solver - name of a new user-defined solver
80 .  path - path (either absolute or relative) the library containing this solver
81 .  name_create - name of routine to create method context
82 -  routine_create - routine to create method context
83 
84    Notes:
85    KSPRegisterDynamic() may be called multiple times to add several user-defined solvers.
86 
87    If dynamic libraries are used, then the fourth input argument (routine_create)
88    is ignored.
89 
90    Sample usage:
91 .vb
92    KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
93                "MySolverCreate",MySolverCreate);
94 .ve
95 
96    Then, your solver can be chosen with the procedural interface via
97 $     KSPSetType(ksp,"my_solver")
98    or at runtime via the option
99 $     -ksp_type my_solver
100 
101    Level: advanced
102 
103    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
104           and others of the form ${any_environmental_variable} occuring in pathname will be
105           replaced with appropriate values.
106          If your function is not being put into a shared library then use KSPRegister() instead
107 
108 .keywords: KSP, register
109 
110 .seealso: KSPRegisterAll(), KSPRegisterDestroy()
111 
112 M*/
113 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
114 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0)
115 #else
116 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d)
117 #endif
118 
119 EXTERN PetscErrorCode KSPGetType(KSP,KSPType *);
120 EXTERN PetscErrorCode KSPSetPreconditionerSide(KSP,PCSide);
121 EXTERN PetscErrorCode KSPGetPreconditionerSide(KSP,PCSide*);
122 EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
123 EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
124 EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscTruth);
125 EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscTruth *);
126 EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscTruth);
127 EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscTruth*);
128 EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscTruth);
129 EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscTruth);
130 EXTERN PetscErrorCode KSPGetRhs(KSP,Vec *);
131 EXTERN PetscErrorCode KSPGetSolution(KSP,Vec *);
132 EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*);
133 EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*);
134 EXTERN PetscErrorCode KSPSetNullSpace(KSP,MatNullSpace);
135 EXTERN PetscErrorCode KSPGetNullSpace(KSP,MatNullSpace*);
136 
137 EXTERN PetscErrorCode KSPSetPC(KSP,PC);
138 EXTERN PetscErrorCode KSPGetPC(KSP,PC*);
139 
140 EXTERN PetscErrorCode KSPSetMonitor(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void*));
141 EXTERN PetscErrorCode KSPClearMonitor(KSP);
142 EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **);
143 EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *);
144 EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscTruth);
145 
146 /* not sure where to put this */
147 EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*);
148 EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
149 EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
150 EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
151 
152 EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *);
153 EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *);
154 
155 EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
156 EXTERN PetscErrorCode KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal);
157 EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
158 EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *);
159 EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*);
160 
161 EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
162 EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);
163 
164 EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
165 EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
166 EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
167 EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
168 
169 EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
170 EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
171 
172 /*E
173     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
174 
175    Level: advanced
176 
177 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
178           KSPGMRESSetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
179 
180 E*/
181 typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
182 
183 /*M
184     KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process
185 
186    Level: advanced
187 
188    Note: Possible unstable, but the fastest to compute
189 
190 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
191           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
192           KSPGMRESModifiedGramSchmidtOrthogonalization()
193 M*/
194 
195 /*M
196     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
197           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
198           poor orthogonality.
199 
200    Level: advanced
201 
202    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
203      estimate the orthogonality but is more stable.
204 
205 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
206           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
207           KSPGMRESModifiedGramSchmidtOrthogonalization()
208 M*/
209 
210 /*M
211     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
212 
213    Level: advanced
214 
215    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
216      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
217 
218         You should only use this if you absolutely know that the iterative refinement is needed.
219 
220 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
221           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
222           KSPGMRESModifiedGramSchmidtOrthogonalization()
223 M*/
224 
225 EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
226 
227 EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
228 EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
229 EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
230 
231 EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
232 EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
233 EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);
234 
235 EXTERN PetscErrorCode KSPSetFromOptions(KSP);
236 EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
237 
238 EXTERN PetscErrorCode KSPSingularValueMonitor(KSP,PetscInt,PetscReal,void *);
239 EXTERN PetscErrorCode KSPDefaultMonitor(KSP,PetscInt,PetscReal,void *);
240 EXTERN PetscErrorCode KSPTrueMonitor(KSP,PetscInt,PetscReal,void *);
241 EXTERN PetscErrorCode KSPDefaultSMonitor(KSP,PetscInt,PetscReal,void *);
242 EXTERN PetscErrorCode KSPVecViewMonitor(KSP,PetscInt,PetscReal,void *);
243 EXTERN PetscErrorCode KSPGMRESKrylovMonitor(KSP,PetscInt,PetscReal,void *);
244 
245 EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
246 EXTERN PetscErrorCode KSPDefaultBuildSolution(KSP,Vec,Vec*);
247 EXTERN PetscErrorCode KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
248 
249 EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat,MatStructure);
250 EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*,MatStructure*);
251 EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
252 EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
253 EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,char*[]);
254 
255 EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscTruth);
256 EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscTruth*);
257 EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscTruth);
258 EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscTruth*);
259 
260 EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
261 
262 /*E
263     KSPNormType - Norm that is passed in the Krylov convergence
264        test routines.
265 
266    Level: advanced
267 
268    Notes: this must match finclude/petscksp.h
269 
270 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
271           KSPSetConvergenceTest()
272 E*/
273 typedef enum {KSP_NO_NORM               = 0,
274               KSP_PRECONDITIONED_NORM   = 1,
275               KSP_UNPRECONDITIONED_NORM = 2,
276               KSP_NATURAL_NORM          = 3} KSPNormType;
277 
278 /*M
279     KSP_NO_NORM - Do not compute a norm during the Krylov process. This will
280           possibly save some computation but means the convergence test cannot
281           be based on a norm of a residual etc.
282 
283    Level: advanced
284 
285     Note: Some Krylov methods need to compute a residual norm and then this is ignored
286 
287 .seealso: KSPNormType, KSPSetNormType(), KSP_PRECONDITIONED_NORM, KSP_UNPRECONDITIONED_NORM, KSP_NATURAL_NORM
288 M*/
289 
290 /*M
291     KSP_PRECONDITIONED_NORM - Compute the norm of the preconditioned residual and pass that to the
292        convergence test routine.
293 
294    Level: advanced
295 
296 .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_UNPRECONDITIONED_NORM, KSP_NATURAL_NORM, KSPSetConvergenceTest()
297 M*/
298 
299 /*M
300     KSP_UNPRECONDITIONED_NORM - Compute the norm of the true residual (b - A*x) and pass that to the
301        convergence test routine.
302 
303    Level: advanced
304 
305 .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_PRECONDITIONED_NORM, KSP_NATURAL_NORM, KSPSetConvergenceTest()
306 M*/
307 
308 /*M
309     KSP_NATURAL_NORM - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
310        convergence test routine.
311 
312    Level: advanced
313 
314 .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_PRECONDITIONED_NORM, KSP_UNPRECONDITIONED_NORM, KSPSetConvergenceTest()
315 M*/
316 
317 EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
318 
319 /*E
320     KSPConvergedReason - reason a Krylov method was said to
321          have converged or diverged
322 
323    Level: beginner
324 
325    Notes: this must match finclude/petscksp.h
326 
327    Developer note: The string versions of these are in
328      src/ksp/ksp/interface/itfunc.c called convergedreasons.
329      If these enums are changed you must change those.
330 
331 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
332 E*/
333 typedef enum {/* converged */
334               KSP_CONVERGED_RTOL               =  2,
335               KSP_CONVERGED_ATOL               =  3,
336               KSP_CONVERGED_ITS                =  4,
337               KSP_CONVERGED_QCG_NEG_CURVE      =  5,
338               KSP_CONVERGED_QCG_CONSTRAINED    =  6,
339               KSP_CONVERGED_STEP_LENGTH        =  7,
340               /* diverged */
341               KSP_DIVERGED_NULL                = -2,
342               KSP_DIVERGED_ITS                 = -3,
343               KSP_DIVERGED_DTOL                = -4,
344               KSP_DIVERGED_BREAKDOWN           = -5,
345               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
346               KSP_DIVERGED_NONSYMMETRIC        = -7,
347               KSP_DIVERGED_INDEFINITE_PC       = -8,
348               KSP_DIVERGED_NAN                 = -9,
349               KSP_DIVERGED_INDEFINITE_MAT      = -10,
350 
351               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
352 
353 /*MC
354      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)
355 
356    Level: beginner
357 
358    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
359        for left preconditioning it is the 2-norm of the preconditioned residual, and the
360        2-norm of the residual for right preconditioning
361 
362 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
363 
364 M*/
365 
366 /*MC
367      KSP_CONVERGED_ATOL - norm(r) <= atol
368 
369    Level: beginner
370 
371    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
372        for left preconditioning it is the 2-norm of the preconditioned residual, and the
373        2-norm of the residual for right preconditioning
374 
375    Level: beginner
376 
377 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
378 
379 M*/
380 
381 /*MC
382      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
383 
384    Level: beginner
385 
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
389 
390    Level: beginner
391 
392 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
393 
394 M*/
395 
396 /*MC
397      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
398       reached
399 
400    Level: beginner
401 
402 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
403 
404 M*/
405 
406 /*MC
407      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of the
408            preconditioner is applied.
409 
410 
411    Level: beginner
412 
413 
414 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
415 
416 M*/
417 
418 /*MC
419      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
420           method could not continue to enlarge the Krylov space.
421 
422    Level: beginner
423 
424 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
425 
426 M*/
427 
428 /*MC
429      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
430           method could not continue to enlarge the Krylov space.
431 
432 
433    Level: beginner
434 
435 
436 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
437 
438 M*/
439 
440 /*MC
441      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
442         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
443 
444    Level: beginner
445 
446 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
447 
448 M*/
449 
450 /*MC
451      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
452         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
453         be positive definite
454 
455    Level: beginner
456 
457      Notes: This can happen with the PCICC preconditioner, use -pc_icc_shift to force
458   the PCICC preconditioner to generate a positive definite preconditioner
459 
460 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
461 
462 M*/
463 
464 /*MC
465      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
466         while the KSPSolve() is still running.
467 
468    Level: beginner
469 
470 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
471 
472 M*/
473 
474 EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *);
475 EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **);
476 EXTERN PetscErrorCode KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
477 EXTERN PetscErrorCode KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
478 EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *);
479 
480 EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *);
481 
482 /*E
483     KSPCGType - Determines what type of CG to use
484 
485    Level: beginner
486 
487 .seealso: KSPCGSetType()
488 E*/
489 typedef enum {KSP_CG_SYMMETRIC=1,KSP_CG_HERMITIAN=2} KSPCGType;
490 
491 EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
492 
493 EXTERN PetscErrorCode PCPreSolve(PC,KSP);
494 EXTERN PetscErrorCode PCPostSolve(PC,KSP);
495 
496 EXTERN PetscErrorCode KSPLGMonitorCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
497 EXTERN PetscErrorCode KSPLGMonitor(KSP,PetscInt,PetscReal,void*);
498 EXTERN PetscErrorCode KSPLGMonitorDestroy(PetscDrawLG);
499 EXTERN PetscErrorCode KSPLGTrueMonitorCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
500 EXTERN PetscErrorCode KSPLGTrueMonitor(KSP,PetscInt,PetscReal,void*);
501 EXTERN PetscErrorCode KSPLGTrueMonitorDestroy(PetscDrawLG);
502 
503 PETSC_EXTERN_CXX_END
504 #endif
505