xref: /petsc/include/petscksp.h (revision dec1416f15364d8a66cef6f4b2a5a2aba5192d13)
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 
8 PETSC_EXTERN PetscErrorCode KSPInitializePackage(void);
9 
10 /*S
11      KSP - Abstract PETSc object that manages all Krylov methods. This is the object that manages the
12          linear solves in PETSc (even those such as direct solvers that do no use Krylov accelerators).
13 
14    Level: beginner
15 
16         Notes:
17     When a direct solver is used but no Krylov solver is used the KSP object is still used by with a
18        KSPType of KSPPREONLY (meaning application of the preconditioner is only used as the linear solver).
19 
20 .seealso:  KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP, KSPDestroy(), KSPCG, KSPGMRES
21 S*/
22 typedef struct _p_KSP*     KSP;
23 
24 /*J
25     KSPType - String with the name of a PETSc Krylov method.
26 
27    Level: beginner
28 
29 .seealso: KSPSetType(), KSP, KSPRegister(), KSPCreate(), KSPSetFromOptions()
30 J*/
31 typedef const char* KSPType;
32 #define KSPRICHARDSON "richardson"
33 #define KSPCHEBYSHEV  "chebyshev"
34 #define KSPCG         "cg"
35 #define KSPGROPPCG    "groppcg"
36 #define KSPPIPECG     "pipecg"
37 #define KSPPIPECGRR   "pipecgrr"
38 #define KSPPIPELCG     "pipelcg"
39 #define   KSPCGNE       "cgne"
40 #define   KSPNASH       "nash"
41 #define   KSPSTCG       "stcg"
42 #define   KSPGLTR       "gltr"
43 #define     KSPCGNASH  PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGNASH macro is deprecated use KSPNASH (since version 3.11)\"")  "nash"
44 #define     KSPCGSTCG  PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGSTCG macro is deprecated use KSPSTCG (since version 3.11)\"")  "stcg"
45 #define     KSPCGGLTR  PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGGLTR macro is deprecated use KSPSGLTR (since version 3.11)\"") "gltr"
46 #define KSPFCG        "fcg"
47 #define KSPPIPEFCG    "pipefcg"
48 #define KSPGMRES      "gmres"
49 #define KSPPIPEFGMRES "pipefgmres"
50 #define   KSPFGMRES     "fgmres"
51 #define   KSPLGMRES     "lgmres"
52 #define   KSPDGMRES     "dgmres"
53 #define   KSPPGMRES     "pgmres"
54 #define KSPTCQMR      "tcqmr"
55 #define KSPBCGS       "bcgs"
56 #define   KSPIBCGS      "ibcgs"
57 #define   KSPFBCGS      "fbcgs"
58 #define   KSPFBCGSR     "fbcgsr"
59 #define   KSPBCGSL      "bcgsl"
60 #define   KSPPIPEBCGS   "pipebcgs"
61 #define KSPCGS        "cgs"
62 #define KSPTFQMR      "tfqmr"
63 #define KSPCR         "cr"
64 #define KSPPIPECR     "pipecr"
65 #define KSPLSQR       "lsqr"
66 #define KSPPREONLY    "preonly"
67 #define KSPQCG        "qcg"
68 #define KSPBICG       "bicg"
69 #define KSPMINRES     "minres"
70 #define KSPSYMMLQ     "symmlq"
71 #define KSPLCD        "lcd"
72 #define KSPPYTHON     "python"
73 #define KSPGCR        "gcr"
74 #define KSPPIPEGCR    "pipegcr"
75 #define KSPTSIRM      "tsirm"
76 #define KSPCGLS       "cgls"
77 #define KSPFETIDP     "fetidp"
78 
79 /* Logging support */
80 PETSC_EXTERN PetscClassId KSP_CLASSID;
81 PETSC_EXTERN PetscClassId KSPGUESS_CLASSID;
82 PETSC_EXTERN PetscClassId DMKSP_CLASSID;
83 
84 PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP*);
85 PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType);
86 PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType*);
87 PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
88 PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
89 PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec);
90 PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec);
91 PETSC_EXTERN PetscErrorCode KSPReset(KSP);
92 PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP);
93 PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*);
94 PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool);
95 PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool);
96 PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP,PC,Vec);
97 
98 PETSC_EXTERN PetscFunctionList KSPList;
99 PETSC_EXTERN PetscFunctionList KSPGuessList;
100 PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP));
101 
102 PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide);
103 PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*);
104 PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
105 PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
106 PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool);
107 PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool*);
108 PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool);
109 PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool*);
110 PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool);
111 PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool);
112 PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool*);
113 PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool);
114 PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool*);
115 PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec*);
116 PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec*);
117 PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*);
118 PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*);
119 PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*);
120 PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
121 PETSC_DEPRECATED_FUNCTION("Use KSPCreateVecs() (since version 3.6)") PETSC_STATIC_INLINE PetscErrorCode KSPGetVecs(KSP ksp,PetscInt n,Vec **x,PetscInt m,Vec **y) {return KSPCreateVecs(ksp,n,x,m,y);}
122 
123 PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
124 PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
125 
126 PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC);
127 PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*);
128 
129 PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal);
130 PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode (*)(void**));
131 PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
132 PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void**);
133 PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt*);
134 PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool );
135 
136 PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*);
137 PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec*);
138 PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
139 PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt);
140 
141 PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*);
142 PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC,KSP);
143 PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
144 PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
145 PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
146 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
147 PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC,PetscInt*,KSP*[]);
148 PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*);
149 PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*);
150 PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*);
151 PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*);
152 PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP*);
153 
154 PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec*);
155 PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec*);
156 
157 PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
158 PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool );
159 PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal);
160 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal);
161 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP,PetscBool);
162 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*);
163 PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
164 PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt*);
165 PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]);
166 PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt*,Vec[],PetscReal[],PetscReal[]);
167 
168 /*E
169 
170   KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods
171 
172   KSP_FCD_TRUNC_TYPE_STANDARD uses all (up to mmax) stored directions
173   KSP_FCD_TRUNC_TYPE_NOTAY uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..
174 
175    Level: intermediate
176 .seealso : KSPFCG,KSPPIPEFCG,KSPPIPEGCR,KSPFCGSetTruncationType(),KSPFCGGetTruncationType()
177 
178 E*/
179 typedef enum {KSP_FCD_TRUNC_TYPE_STANDARD,KSP_FCD_TRUNC_TYPE_NOTAY} KSPFCDTruncationType;
180 PETSC_EXTERN const char *const KSPFCDTruncationTypes[];
181 
182 PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt);
183 PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*);
184 PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt);
185 PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*);
186 PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType);
187 PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*);
188 
189 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt);
190 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*);
191 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt);
192 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*);
193 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType);
194 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*);
195 
196 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt);
197 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*);
198 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt);
199 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*);
200 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType);
201 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*);
202 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool);
203 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*);
204 
205 PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
206 PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*);
207 PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);
208 
209 PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
210 PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
211 PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
212 PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
213 PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
214 
215 PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
216 PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
217 
218 PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP,PetscScalar);
219 
220 PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
221 PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
222 PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
223 
224 PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP,PC*);
225 PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP,PC);
226 PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP,KSP*);
227 PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP,Mat);
228 /*E
229     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
230 
231    Level: advanced
232 
233 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
234           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
235 
236 E*/
237 typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
238 PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
239 /*MC
240     KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
241 
242    Level: advanced
243 
244    Note: Possible unstable, but the fastest to compute
245 
246 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
247           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
248           KSPGMRESModifiedGramSchmidtOrthogonalization()
249 M*/
250 
251 /*MC
252     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
253           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
254           poor orthogonality.
255 
256    Level: advanced
257 
258    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
259      estimate the orthogonality but is more stable.
260 
261 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
262           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
263           KSPGMRESModifiedGramSchmidtOrthogonalization()
264 M*/
265 
266 /*MC
267     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
268 
269    Level: advanced
270 
271    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
272      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
273 
274         You should only use this if you absolutely know that the iterative refinement is needed.
275 
276 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
277           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
278           KSPGMRESModifiedGramSchmidtOrthogonalization()
279 M*/
280 
281 PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
282 PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);
283 
284 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
285 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
286 PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
287 
288 PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
289 PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
290 PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);
291 
292 PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
293 PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool );
294 PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
295 PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool);
296 
297 PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
298 PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP);
299 PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
300 
301 PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
302 PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
303 PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
304 PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
305 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
306 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
307 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
308 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
309 PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
310 PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,PetscViewerAndFormat*);
311 PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*);
312 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**);
313 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**);
314 PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void*);
315 PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP,const char[],const char[],const char [],PetscErrorCode (*)(KSP,PetscInt,PetscReal,PetscViewerAndFormat*));
316 
317 PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
318 PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
319 
320 PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat);
321 PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*);
322 PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool*,PetscBool*);
323 PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
324 PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
325 PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);
326 
327 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool );
328 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool*);
329 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool );
330 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool*);
331 
332 PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
333 PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);
334 PETSC_STATIC_INLINE PetscErrorCode KSPViewFromOptions(KSP A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
335 PETSC_EXTERN PetscErrorCode KSPReasonView(KSP,PetscViewer);
336 PETSC_EXTERN PetscErrorCode KSPReasonViewFromOptions(KSP);
337 
338 #define KSP_FILE_CLASSID 1211223
339 
340 PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP,PetscBool);
341 PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP,PetscBool);
342 PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
343 PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP,PetscReal*,PetscReal*);
344 
345 PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
346 PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
347 PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*);
348 
349 /*E
350     KSPNormType - Norm that is passed in the Krylov convergence
351        test routines.
352 
353    Level: advanced
354 
355    Each solver only supports a subset of these and some may support different ones
356    depending on left or right preconditioning, see KSPSetPCSide()
357 
358    Notes:
359     this must match petsc/finclude/petscksp.h
360 
361 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
362           KSPSetConvergenceTest(), KSPSetPCSide()
363 E*/
364 typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
365 #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
366 PETSC_EXTERN const char *const*const KSPNormTypes;
367 
368 /*MC
369     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
370           possibly save some computation but means the convergence test cannot
371           be based on a norm of a residual etc.
372 
373    Level: advanced
374 
375     Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored
376 
377 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
378 M*/
379 
380 /*MC
381     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
382        convergence test routine.
383 
384    Level: advanced
385 
386 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
387 M*/
388 
389 /*MC
390     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
391        convergence test routine.
392 
393    Level: advanced
394 
395 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
396 M*/
397 
398 /*MC
399     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
400        convergence test routine. This is only supported by  KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG, KSPPIPEFCG, KSPPIPEGCR
401 
402    Level: advanced
403 
404 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
405 M*/
406 
407 PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
408 PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
409 PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
410 PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
411 PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool);
412 
413 /*E
414     KSPConvergedReason - reason a Krylov method was said to have converged or diverged
415 
416    Level: beginner
417 
418    Notes:
419     See KSPGetConvergedReason() for explanation of each value
420 
421    Developer Notes:
422     this must match petsc/finclude/petscksp.h
423 
424       The string versions of these are KSPConvergedReasons; if you change
425       any of the values here also change them that array of names.
426 
427 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
428 E*/
429 #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM("Use KSP_DIVERGED_PC_FAILED (since version 3.11)")
430 typedef enum {/* converged */
431               KSP_CONVERGED_RTOL_NORMAL        =  1,
432               KSP_CONVERGED_ATOL_NORMAL        =  9,
433               KSP_CONVERGED_RTOL               =  2,
434               KSP_CONVERGED_ATOL               =  3,
435               KSP_CONVERGED_ITS                =  4,
436               KSP_CONVERGED_CG_NEG_CURVE       =  5,
437               KSP_CONVERGED_CG_CONSTRAINED     =  6,
438               KSP_CONVERGED_STEP_LENGTH        =  7,
439               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
440               /* diverged */
441               KSP_DIVERGED_NULL                = -2,
442               KSP_DIVERGED_ITS                 = -3,
443               KSP_DIVERGED_DTOL                = -4,
444               KSP_DIVERGED_BREAKDOWN           = -5,
445               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
446               KSP_DIVERGED_NONSYMMETRIC        = -7,
447               KSP_DIVERGED_INDEFINITE_PC       = -8,
448               KSP_DIVERGED_NANORINF            = -9,
449               KSP_DIVERGED_INDEFINITE_MAT      = -10,
450               KSP_DIVERGED_PC_FAILED           = -11,
451               KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED  = -11,
452 
453               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
454 PETSC_EXTERN const char *const*KSPConvergedReasons;
455 
456 /*MC
457      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if KSPConvergedDefaultSetUIRNorm() was called
458 
459    Level: beginner
460 
461    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
462        for left preconditioning it is the 2-norm of the preconditioned residual, and the
463        2-norm of the residual for right preconditioning
464 
465    See also KSP_CONVERGED_ATOL which may apply before this tolerance.
466 
467 .seealso:  KSP_CONVERGED_ATOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
468 
469 M*/
470 
471 /*MC
472      KSP_CONVERGED_ATOL - norm(r) <= atol
473 
474    Level: beginner
475 
476    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
477        for left preconditioning it is the 2-norm of the preconditioned residual, and the
478        2-norm of the residual for right preconditioning
479 
480    See also KSP_CONVERGED_RTOL which may apply before this tolerance.
481 
482 .seealso:  KSP_CONVERGED_RTOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
483 
484 M*/
485 
486 /*MC
487      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
488 
489    Level: beginner
490 
491    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
492        for left preconditioning it is the 2-norm of the preconditioned residual, and the
493        2-norm of the residual for right preconditioning
494 
495    Level: beginner
496 
497 .seealso:  KSP_CONVERGED_ATOL, KSP_DIVERGED_RTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
498 
499 M*/
500 
501 /*MC
502      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
503       reached
504 
505    Level: beginner
506 
507 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
508 
509 M*/
510 
511 /*MC
512      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
513            the preconditioner is applied. Also used when the KSPConvergedSkip() convergence
514            test routine is set in KSP.
515 
516    Level: beginner
517 
518 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
519 
520 M*/
521 
522 /*MC
523      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
524           method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or
525           preconditioner.
526 
527    Level: beginner
528 
529 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
530 
531 M*/
532 
533 /*MC
534      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
535           method could not continue to enlarge the Krylov space.
536 
537    Level: beginner
538 
539 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
540 
541 M*/
542 
543 /*MC
544      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
545         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
546 
547    Level: beginner
548 
549 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
550 
551 M*/
552 
553 /*MC
554      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
555         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
556         be positive definite
557 
558    Level: beginner
559 
560      Notes:
561     This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
562   the PCICC preconditioner to generate a positive definite preconditioner
563 
564 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
565 
566 M*/
567 
568 /*MC
569      KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a
570      zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner
571      such as PCFIELDSPLIT.
572 
573    Level: beginner
574 
575     Notes:
576     Run with -ksp_error_if_not_converged to stop the program when the error is detected and print an error message with details.
577 
578 
579 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
580 
581 M*/
582 
583 /*MC
584      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
585         while the KSPSolve() is still running.
586 
587    Level: beginner
588 
589 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
590 
591 M*/
592 
593 PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void*,PetscErrorCode (*)(void*));
594 PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*));
595 PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*));
596 PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void**);
597 PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
598 PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
599 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void*);
600 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void**);
601 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
602 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
603 PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
604 PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason*);
605 
606 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ }
607 #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
608 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ }
609 #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
610 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ }
611 #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
612 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ }
613 #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
614 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ }
615 #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
616 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ }
617 #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
618 
619 PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP,MatType,Mat*);
620 PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPComputeExplicitOperator(KSP A,Mat* B) { return KSPComputeOperator(A,NULL,B); }
621 
622 /*E
623     KSPCGType - Determines what type of CG to use
624 
625    Level: beginner
626 
627 .seealso: KSPCGSetType()
628 E*/
629 typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
630 PETSC_EXTERN const char *const KSPCGTypes[];
631 
632 PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
633 PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool );
634 
635 PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP,PetscReal);
636 PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP,PetscReal*);
637 PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP,PetscReal*);
638 
639 PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal*);
640 PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal*);
641 PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp,PetscReal *x) {return KSPGLTRGetMinEig(ksp,x);}
642 PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetLambda(KSP ksp,PetscReal *x) {return KSPGLTRGetLambda(ksp,x);}
643 
644 PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]);
645 
646 PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
647 PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);
648 
649 #include <petscdrawtypes.h>
650 PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
651 PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*);
652 PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
653 PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
654 PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
655 
656 PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
657 PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
658 
659 /*S
660      KSPGuess - Abstract PETSc object that manages all initial guess methods in Krylov methods.
661 
662    Level: beginner
663 
664 .seealso:  KSPCreate(), KSPSetGuessType(), KSPGuessType
665 S*/
666 typedef struct _p_KSPGuess* KSPGuess;
667 /*J
668     KSPGuessType - String with the name of a PETSc initial guess for Krylov methods.
669 
670    Level: beginner
671 
672 .seealso: KSPGuess
673 J*/
674 typedef const char* KSPGuessType;
675 #define KSPGUESSFISCHER "fischer"
676 #define KSPGUESSPOD     "pod"
677 PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[],PetscErrorCode (*)(KSPGuess));
678 PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP,KSPGuess);
679 PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP,KSPGuess*);
680 PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess,PetscViewer);
681 PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess*);
682 PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm,KSPGuess*);
683 PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess,KSPGuessType);
684 PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess,KSPGuessType*);
685 PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess);
686 PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess,Vec,Vec);
687 PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess,Vec,Vec);
688 PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess);
689 PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess,PetscInt,PetscInt);
690 PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
691 PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool);
692 PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool*);
693 
694 /*E
695     MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines
696 
697     Level: intermediate
698 
699 .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat()
700 E*/
701 typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG} MatSchurComplementAinvType;
702 PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
703 
704 PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
705 PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
706 PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
707 PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
708 PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
709 PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
710 PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType);
711 PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*);
712 PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*);
713 PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*);
714 PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatSchurComplementAinvType,MatReuse,Mat*);
715 PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*);
716 
717 PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm,PetscInt,PetscInt,Mat*);
718 PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm,PetscInt,PetscInt,Mat*);
719 PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm,PetscInt,PetscInt,Mat*);
720 PETSC_EXTERN PetscErrorCode MatCreateLMVMBrdn(MPI_Comm,PetscInt,PetscInt,Mat*);
721 PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBrdn(MPI_Comm,PetscInt,PetscInt,Mat*);
722 PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBrdn(MPI_Comm,PetscInt,PetscInt,Mat*);
723 PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBrdn(MPI_Comm,PetscInt,PetscInt,Mat*);
724 PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBrdn(MPI_Comm,PetscInt,PetscInt,Mat*);
725 
726 PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec);
727 PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool*);
728 PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec);
729 PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool);
730 PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat);
731 PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat);
732 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat);
733 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal);
734 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec);
735 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC);
736 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP);
737 PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec);
738 PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec);
739 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat*);
740 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC*);
741 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP*);
742 PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt*);
743 PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt*);
744 PETSC_EXTERN PetscErrorCode MatSymBrdnSetDelta(Mat, PetscScalar);
745 
746 PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
747 PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool );
748 PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
749 PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
750 PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
751 PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void*);
752 PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
753 PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
754 PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
755 PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*);
756 PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
757 PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
758 PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
759 PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
760 
761 PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec);
762 PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec,
763                                            void (**)(PetscInt, PetscInt, PetscInt,
764                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
765                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
766                                                      PetscReal, const PetscReal [], PetscInt, const PetscScalar[], PetscScalar []), InsertMode, Vec);
767 #endif
768