xref: /petsc/include/petscksp.h (revision 329f5518e9d4bb7ce96c0c5576cc53785c973973)
1 /* $Id: ksp.h,v 1.87 1999/12/08 22:17:32 balay Exp bsmith $ */
2 /*
3    Defines the interface functions for the Krylov subspace accelerators.
4 */
5 #ifndef __KSP_H
6 #define __KSP_H
7 #include "pc.h"
8 
9 #define KSP_COOKIE  PETSC_COOKIE+8
10 
11 typedef struct _p_KSP*     KSP;
12 
13 #define KSPRICHARDSON "richardson"
14 #define KSPCHEBYCHEV  "chebychev"
15 #define KSPCG         "cg"
16 #define KSPGMRES      "gmres"
17 #define KSPTCQMR      "tcqmr"
18 #define KSPBCGS       "bcgs"
19 #define KSPCGS        "cgs"
20 #define KSPTFQMR      "tfqmr"
21 #define KSPCR         "cr"
22 #define KSPLSQR       "lsqr"
23 #define KSPPREONLY    "preonly"
24 #define KSPQCG        "qcg"
25 #define KSPBICG       "bicg"
26 #define KSPFGMRES     "fgmres"
27 typedef char * KSPType;
28 
29 extern int KSPCreate(MPI_Comm,KSP *);
30 extern int KSPSetType(KSP,KSPType);
31 extern int KSPSetUp(KSP);
32 extern int KSPSolve(KSP,int *);
33 extern int KSPSolveTranspose(KSP,int *);
34 extern int KSPDestroy(KSP);
35 
36 extern FList KSPList;
37 extern int KSPRegisterAll(char *);
38 extern int KSPRegisterDestroy(void);
39 
40 extern int KSPRegister(char*,char*,char*,int(*)(KSP));
41 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
42 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0)
43 #else
44 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d)
45 #endif
46 
47 extern int KSPGetType(KSP,KSPType *);
48 extern int KSPSetPreconditionerSide(KSP,PCSide);
49 extern int KSPGetPreconditionerSide(KSP,PCSide*);
50 extern int KSPGetTolerances(KSP,double*,double*,double*,int*);
51 extern int KSPSetTolerances(KSP,double,double,double,int);
52 extern int KSPSetComputeResidual(KSP,PetscTruth);
53 extern int KSPSetUsePreconditionedResidual(KSP);
54 extern int KSPSetInitialGuessNonzero(KSP);
55 extern int KSPGetInitialGuessNonzero(KSP,PetscTruth *);
56 extern int KSPSetComputeEigenvalues(KSP);
57 extern int KSPSetComputeSingularValues(KSP);
58 extern int KSPSetRhs(KSP,Vec);
59 extern int KSPGetRhs(KSP,Vec *);
60 extern int KSPSetSolution(KSP,Vec);
61 extern int KSPGetSolution(KSP,Vec *);
62 extern int KSPGetResidualNorm(KSP,double*);
63 extern int KSPGetIterationNumber(KSP,int*);
64 
65 extern int KSPSetPC(KSP,PC);
66 extern int KSPGetPC(KSP,PC*);
67 
68 extern int KSPSetAvoidNorms(KSP);
69 
70 extern int KSPSetMonitor(KSP,int (*)(KSP,int,double,void*),void *,int (*)(void*));
71 extern int KSPClearMonitor(KSP);
72 extern int KSPGetMonitorContext(KSP,void **);
73 extern int KSPGetResidualHistory(KSP,double **,int *);
74 extern int KSPSetResidualHistory(KSP,double *,int,PetscTruth);
75 
76 
77 extern int KSPBuildSolution(KSP,Vec,Vec *);
78 extern int KSPBuildResidual(KSP,Vec,Vec,Vec *);
79 
80 extern int KSPRichardsonSetScale(KSP,double);
81 extern int KSPChebychevSetEigenvalues(KSP,double,double);
82 extern int KSPComputeExtremeSingularValues(KSP,double*,double*);
83 extern int KSPComputeEigenvalues(KSP,int,double*,double*,int *);
84 extern int KSPComputeEigenvaluesExplicitly(KSP,int,double*,double*);
85 
86 extern int KSPGMRESSetRestart(KSP,int);
87 extern int KSPGMRESSetPreAllocateVectors(KSP);
88 extern int KSPGMRESSetOrthogonalization(KSP,int (*)(KSP,int));
89 extern int KSPGMRESUnmodifiedGramSchmidtOrthogonalization(KSP,int);
90 extern int KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,int);
91 extern int KSPGMRESIROrthogonalization(KSP,int);
92 
93 extern int KSPFGMRESModifyPCNoChange(KSP,int,int,double,void*);
94 extern int KSPFGMRESModifyPCSLES(KSP,int,int,double,void*);
95 extern int KSPFGMRESSetModifyPC(KSP,int (*)(KSP,int,int,double,void*),void*,int(*)(void*));
96 
97 extern int KSPSetFromOptions(KSP);
98 extern int KSPSetTypeFromOptions(KSP);
99 extern int KSPAddOptionsChecker(int (*)(KSP));
100 
101 extern int KSPSingularValueMonitor(KSP,int,double,void *);
102 extern int KSPDefaultMonitor(KSP,int,double,void *);
103 extern int KSPTrueMonitor(KSP,int,double,void *);
104 extern int KSPDefaultSMonitor(KSP,int,double,void *);
105 extern int KSPVecViewMonitor(KSP,int,double,void *);
106 extern int KSPGMRESKrylovMonitor(KSP,int,double,void *);
107 
108 
109 extern int KSPResidual(KSP,Vec,Vec,Vec,Vec,Vec,Vec);
110 extern int KSPUnwindPreconditioner(KSP,Vec,Vec);
111 extern int KSPDefaultBuildSolution(KSP,Vec,Vec*);
112 extern int KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
113 
114 extern int KSPPrintHelp(KSP);
115 
116 extern int KSPSetOptionsPrefix(KSP,char*);
117 extern int KSPAppendOptionsPrefix(KSP,char*);
118 extern int KSPGetOptionsPrefix(KSP,char**);
119 
120 extern int KSPView(KSP,Viewer);
121 
122 typedef enum {/* converged */
123               KSP_CONVERGED_RTOL               =  2,
124               KSP_CONVERGED_ATOL               =  3,
125               KSP_CONVERGED_ITS                =  4,
126               KSP_CONVERGED_QCG_NEGATIVE_CURVE =  5,
127               KSP_CONVERGED_QCG_CONSTRAINED    =  6,
128               KSP_CONVERGED_STEP_LENGTH        =  7,
129               /* diverged */
130               KSP_DIVERGED_ITS                 = -3,
131               KSP_DIVERGED_DTOL                = -4,
132               KSP_DIVERGED_BREAKDOWN           = -5,
133 
134               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
135 
136 extern int KSPSetConvergenceTest(KSP,int (*)(KSP,int,double,KSPConvergedReason*,void*),void *);
137 extern int KSPGetConvergenceContext(KSP,void **);
138 extern int KSPDefaultConverged(KSP,int,double,KSPConvergedReason*,void *);
139 extern int KSPSkipConverged(KSP,int,double,KSPConvergedReason*,void *);
140 extern int KSPGetConvergedReason(KSP,KSPConvergedReason *);
141 
142 extern int KSPComputeExplicitOperator(KSP,Mat *);
143 
144 typedef enum {KSP_CG_SYMMETRIC=1,KSP_CG_HERMITIAN=2} KSPCGType;
145 extern int KSPCGSetType(KSP,KSPCGType);
146 
147 extern int PCPreSolve(PC,KSP);
148 extern int PCPostSolve(PC,KSP);
149 
150 extern int KSPLGMonitorCreate(char*,char*,int,int,int,int,DrawLG*);
151 extern int KSPLGMonitor(KSP,int,double,void*);
152 extern int KSPLGMonitorDestroy(DrawLG);
153 extern int KSPLGTrueMonitorCreate(MPI_Comm,char*,char*,int,int,int,int,DrawLG*);
154 extern int KSPLGTrueMonitor(KSP,int,double,void*);
155 extern int KSPLGTrueMonitorDestroy(DrawLG);
156 
157 #endif
158 
159 
160