xref: /petsc/include/petscsnes.h (revision 0d6fbc72e08a63e91ce3cd64123f05c91ae55da2)
1 /*
2     User interface for the nonlinear solvers package.
3 */
4 #if !defined(__PETSCSNES_H)
5 #define __PETSCSNES_H
6 #include "petscksp.h"
7 PETSC_EXTERN_CXX_BEGIN
8 
9 /*S
10      SNES - Abstract PETSc object that manages all nonlinear solves
11 
12    Level: beginner
13 
14   Concepts: nonlinear solvers
15 
16 .seealso:  SNESCreate(), SNESSetType(), SNESType, TS, KSP, KSP, PC
17 S*/
18 typedef struct _p_SNES* SNES;
19 
20 /*E
21     SNESType - String with the name of a PETSc SNES method or the creation function
22        with an optional dynamic library name, for example
23        http://www.mcs.anl.gov/petsc/lib.a:mysnescreate()
24 
25    Level: beginner
26 
27 .seealso: SNESSetType(), SNES
28 E*/
29 #define SNESType char*
30 #define SNESLS      "ls"
31 #define SNESTR      "tr"
32 #define SNESPYTHON  "python"
33 #define SNESTEST    "test"
34 #define SNESPICARD  "picard"
35 #define SNESKSPONLY "ksponly"
36 #define SNESVI      "vi"
37 #define SNESNGMRES  "ngmres"
38 #define SNESSORQN   "sorqn"
39 #define SNESQN      "qn"
40 
41 /* Logging support */
42 extern PetscClassId  SNES_CLASSID;
43 
44 extern PetscErrorCode  SNESInitializePackage(const char[]);
45 
46 extern PetscErrorCode  SNESCreate(MPI_Comm,SNES*);
47 extern PetscErrorCode  SNESReset(SNES);
48 extern PetscErrorCode  SNESDestroy(SNES*);
49 extern PetscErrorCode  SNESSetType(SNES,const SNESType);
50 extern PetscErrorCode  SNESMonitor(SNES,PetscInt,PetscReal);
51 extern PetscErrorCode  SNESMonitorSet(SNES,PetscErrorCode(*)(SNES,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void**));
52 extern PetscErrorCode  SNESMonitorCancel(SNES);
53 extern PetscErrorCode  SNESSetConvergenceHistory(SNES,PetscReal[],PetscInt[],PetscInt,PetscBool );
54 extern PetscErrorCode  SNESGetConvergenceHistory(SNES,PetscReal*[],PetscInt *[],PetscInt *);
55 extern PetscErrorCode  SNESSetUp(SNES);
56 extern PetscErrorCode  SNESSolve(SNES,Vec,Vec);
57 extern PetscErrorCode  SNESSetErrorIfNotConverged(SNES,PetscBool );
58 extern PetscErrorCode  SNESGetErrorIfNotConverged(SNES,PetscBool  *);
59 
60 
61 extern PetscErrorCode  SNESAddOptionsChecker(PetscErrorCode (*)(SNES));
62 
63 extern PetscErrorCode  SNESSetUpdate(SNES, PetscErrorCode (*)(SNES, PetscInt));
64 extern PetscErrorCode  SNESDefaultUpdate(SNES, PetscInt);
65 
66 extern PetscFList SNESList;
67 extern PetscErrorCode  SNESRegisterDestroy(void);
68 extern PetscErrorCode  SNESRegisterAll(const char[]);
69 
70 extern PetscErrorCode  SNESRegister(const char[],const char[],const char[],PetscErrorCode (*)(SNES));
71 
72 /*MC
73    SNESRegisterDynamic - Adds a method to the nonlinear solver package.
74 
75    Synopsis:
76    PetscErrorCode SNESRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(SNES))
77 
78    Not collective
79 
80    Input Parameters:
81 +  name_solver - name of a new user-defined solver
82 .  path - path (either absolute or relative) the library containing this solver
83 .  name_create - name of routine to create method context
84 -  routine_create - routine to create method context
85 
86    Notes:
87    SNESRegisterDynamic() may be called multiple times to add several user-defined solvers.
88 
89    If dynamic libraries are used, then the fourth input argument (routine_create)
90    is ignored.
91 
92    Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
93    and others of the form ${any_environmental_variable} occuring in pathname will be
94    replaced with appropriate values.
95 
96    Sample usage:
97 .vb
98    SNESRegisterDynamic("my_solver",/home/username/my_lib/lib/libg/solaris/mylib.a,
99                 "MySolverCreate",MySolverCreate);
100 .ve
101 
102    Then, your solver can be chosen with the procedural interface via
103 $     SNESSetType(snes,"my_solver")
104    or at runtime via the option
105 $     -snes_type my_solver
106 
107    Level: advanced
108 
109     Note: If your function is not being put into a shared library then use SNESRegister() instead
110 
111 .keywords: SNES, nonlinear, register
112 
113 .seealso: SNESRegisterAll(), SNESRegisterDestroy()
114 M*/
115 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
116 #define SNESRegisterDynamic(a,b,c,d) SNESRegister(a,b,c,0)
117 #else
118 #define SNESRegisterDynamic(a,b,c,d) SNESRegister(a,b,c,d)
119 #endif
120 
121 extern PetscErrorCode  SNESGetKSP(SNES,KSP*);
122 extern PetscErrorCode  SNESSetKSP(SNES,KSP);
123 extern PetscErrorCode  SNESGetSolution(SNES,Vec*);
124 extern PetscErrorCode  SNESGetSolutionUpdate(SNES,Vec*);
125 extern PetscErrorCode  SNESGetRhs(SNES,Vec*);
126 extern PetscErrorCode  SNESView(SNES,PetscViewer);
127 
128 extern PetscErrorCode  SNESSetOptionsPrefix(SNES,const char[]);
129 extern PetscErrorCode  SNESAppendOptionsPrefix(SNES,const char[]);
130 extern PetscErrorCode  SNESGetOptionsPrefix(SNES,const char*[]);
131 extern PetscErrorCode  SNESSetFromOptions(SNES);
132 extern PetscErrorCode  SNESDefaultGetWork(SNES,PetscInt);
133 
134 extern PetscErrorCode  MatCreateSNESMF(SNES,Mat*);
135 extern PetscErrorCode  MatMFFDComputeJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
136 
137 extern PetscErrorCode  MatDAADSetSNES(Mat,SNES);
138 
139 extern PetscErrorCode  SNESGetType(SNES,const SNESType*);
140 extern PetscErrorCode  SNESMonitorDefault(SNES,PetscInt,PetscReal,void *);
141 extern PetscErrorCode  SNESMonitorRange(SNES,PetscInt,PetscReal,void *);
142 extern PetscErrorCode  SNESMonitorRatio(SNES,PetscInt,PetscReal,void *);
143 extern PetscErrorCode  SNESMonitorSetRatio(SNES,PetscViewer);
144 extern PetscErrorCode  SNESMonitorSolution(SNES,PetscInt,PetscReal,void *);
145 extern PetscErrorCode  SNESMonitorResidual(SNES,PetscInt,PetscReal,void *);
146 extern PetscErrorCode  SNESMonitorSolutionUpdate(SNES,PetscInt,PetscReal,void *);
147 extern PetscErrorCode  SNESMonitorDefaultShort(SNES,PetscInt,PetscReal,void *);
148 extern PetscErrorCode  SNESSetTolerances(SNES,PetscReal,PetscReal,PetscReal,PetscInt,PetscInt);
149 extern PetscErrorCode  SNESGetTolerances(SNES,PetscReal*,PetscReal*,PetscReal*,PetscInt*,PetscInt*);
150 extern PetscErrorCode  SNESSetTrustRegionTolerance(SNES,PetscReal);
151 extern PetscErrorCode  SNESGetFunctionNorm(SNES,PetscReal*);
152 extern PetscErrorCode  SNESGetIterationNumber(SNES,PetscInt*);
153 
154 extern PetscErrorCode  SNESGetNonlinearStepFailures(SNES,PetscInt*);
155 extern PetscErrorCode  SNESSetMaxNonlinearStepFailures(SNES,PetscInt);
156 extern PetscErrorCode  SNESGetMaxNonlinearStepFailures(SNES,PetscInt*);
157 extern PetscErrorCode  SNESGetNumberFunctionEvals(SNES,PetscInt*);
158 
159 extern PetscErrorCode  SNESSetLagPreconditioner(SNES,PetscInt);
160 extern PetscErrorCode  SNESGetLagPreconditioner(SNES,PetscInt*);
161 extern PetscErrorCode  SNESSetLagJacobian(SNES,PetscInt);
162 extern PetscErrorCode  SNESGetLagJacobian(SNES,PetscInt*);
163 extern PetscErrorCode  SNESSetGridSequence(SNES,PetscInt);
164 
165 extern PetscErrorCode  SNESGetLinearSolveIterations(SNES,PetscInt*);
166 extern PetscErrorCode  SNESGetLinearSolveFailures(SNES,PetscInt*);
167 extern PetscErrorCode  SNESSetMaxLinearSolveFailures(SNES,PetscInt);
168 extern PetscErrorCode  SNESGetMaxLinearSolveFailures(SNES,PetscInt*);
169 
170 extern PetscErrorCode  SNESKSPSetUseEW(SNES,PetscBool );
171 extern PetscErrorCode  SNESKSPGetUseEW(SNES,PetscBool *);
172 extern PetscErrorCode  SNESKSPSetParametersEW(SNES,PetscInt,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
173 extern PetscErrorCode  SNESKSPGetParametersEW(SNES,PetscInt*,PetscReal*,PetscReal*,PetscReal*,PetscReal*,PetscReal*,PetscReal*);
174 
175 extern PetscErrorCode  SNESMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
176 extern PetscErrorCode  SNESMonitorLG(SNES,PetscInt,PetscReal,void*);
177 extern PetscErrorCode  SNESMonitorLGDestroy(PetscDrawLG*);
178 extern PetscErrorCode  SNESMonitorLGRangeCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
179 extern PetscErrorCode  SNESMonitorLGRange(SNES,PetscInt,PetscReal,void*);
180 extern PetscErrorCode  SNESMonitorLGRangeDestroy(PetscDrawLG*);
181 
182 extern PetscErrorCode  SNESSetApplicationContext(SNES,void *);
183 extern PetscErrorCode  SNESGetApplicationContext(SNES,void *);
184 extern PetscErrorCode  SNESSetComputeApplicationContext(SNES,PetscErrorCode (*)(SNES,void**),PetscErrorCode (*)(void**));
185 
186 extern PetscErrorCode  SNESPythonSetType(SNES,const char[]);
187 
188 extern PetscErrorCode  SNESSetFunctionDomainError(SNES);
189 /*E
190     SNESConvergedReason - reason a SNES method was said to
191          have converged or diverged
192 
193    Level: beginner
194 
195    The two most common reasons for divergence are
196 $   1) an incorrectly coded or computed Jacobian or
197 $   2) failure or lack of convergence in the linear system (in this case we recommend
198 $      testing with -pc_type lu to eliminate the linear solver as the cause of the problem).
199 
200    Diverged Reasons:
201 .    SNES_DIVERGED_LOCAL_MIN - this can only occur when using the line-search variant of SNES.
202        The line search wants to minimize Q(alpha) = 1/2 || F(x + alpha s) ||^2_2  this occurs
203        at Q'(alpha) = s^T F'(x+alpha s)^T F(x+alpha s) = 0. If s is the Newton direction - F'(x)^(-1)F(x) then
204        you get Q'(alpha) = -F(x)^T F'(x)^(-1)^T F'(x+alpha s)F(x+alpha s); when alpha = 0
205        Q'(0) = - ||F(x)||^2_2 which is always NEGATIVE if F'(x) is invertible. This means the Newton
206        direction is a descent direction and the line search should succeed if alpha is small enough.
207 
208        If F'(x) is NOT invertible AND F'(x)^T F(x) = 0 then Q'(0) = 0 and the Newton direction
209        is NOT a descent direction so the line search will fail. All one can do at this point
210        is change the initial guess and try again.
211 
212        An alternative explanation: Newton's method can be regarded as replacing the function with
213        its linear approximation and minimizing the 2-norm of that. That is F(x+s) approx F(x) + F'(x)s
214        so we minimize || F(x) + F'(x) s ||^2_2; do this using Least Squares. If F'(x) is invertible then
215        s = - F'(x)^(-1)F(x) otherwise F'(x)^T F'(x) s = -F'(x)^T F(x). If F'(x)^T F(x) is NOT zero then there
216        exists a nontrival (that is F'(x)s != 0) solution to the equation and this direction is
217        s = - [F'(x)^T F'(x)]^(-1) F'(x)^T F(x) so Q'(0) = - F(x)^T F'(x) [F'(x)^T F'(x)]^(-T) F'(x)^T F(x)
218        = - (F'(x)^T F(x)) [F'(x)^T F'(x)]^(-T) (F'(x)^T F(x)). Since we are assuming (F'(x)^T F(x)) != 0
219        and F'(x)^T F'(x) has no negative eigenvalues Q'(0) < 0 so s is a descent direction and the line
220        search should succeed for small enough alpha.
221 
222        Note that this RARELY happens in practice. Far more likely the linear system is not being solved
223        (well enough?) or the Jacobian is wrong.
224 
225    SNES_DIVERGED_MAX_IT means that the solver reached the maximum number of iterations without satisfying any
226    convergence criteria. SNES_CONVERGED_ITS means that SNESSkipConverged() was chosen as the convergence test;
227    thus the usual convergence criteria have not been checked and may or may not be satisfied.
228 
229    Developer Notes: this must match finclude/petscsnes.h
230 
231        The string versions of these are in SNESConvergedReason, if you change any value here you must
232      also adjust that array.
233 
234    Each reason has its own manual page.
235 
236 .seealso: SNESSolve(), SNESGetConvergedReason(), KSPConvergedReason, SNESSetConvergenceTest()
237 E*/
238 typedef enum {/* converged */
239               SNES_CONVERGED_FNORM_ABS         =  2, /* ||F|| < atol */
240               SNES_CONVERGED_FNORM_RELATIVE    =  3, /* ||F|| < rtol*||F_initial|| */
241               SNES_CONVERGED_PNORM_RELATIVE    =  4, /* Newton computed step size small; || delta x || < stol */
242               SNES_CONVERGED_ITS               =  5, /* maximum iterations reached */
243               SNES_CONVERGED_TR_DELTA            =  7,
244               /* diverged */
245               SNES_DIVERGED_FUNCTION_DOMAIN     = -1, /* the new x location passed the function is not in the domain of F */
246               SNES_DIVERGED_FUNCTION_COUNT      = -2,
247               SNES_DIVERGED_LINEAR_SOLVE        = -3, /* the linear solve failed */
248               SNES_DIVERGED_FNORM_NAN           = -4,
249               SNES_DIVERGED_MAX_IT              = -5,
250               SNES_DIVERGED_LINE_SEARCH         = -6, /* the line search failed */
251               SNES_DIVERGED_LOCAL_MIN           = -8, /* || J^T b || is small, implies converged to local minimum of F() */
252               SNES_CONVERGED_ITERATING          =  0} SNESConvergedReason;
253 extern const char *const*SNESConvergedReasons;
254 
255 /*MC
256      SNES_CONVERGED_FNORM_ABS - 2-norm(F) <= abstol
257 
258    Level: beginner
259 
260 .seealso:  SNESSolve(), SNESGetConvergedReason(), SNESConvergedReason, SNESSetTolerances()
261 
262 M*/
263 
264 /*MC
265      SNES_CONVERGED_FNORM_RELATIVE - 2-norm(F) <= rtol*2-norm(F(x_0)) where x_0 is the initial guess
266 
267    Level: beginner
268 
269 .seealso:  SNESSolve(), SNESGetConvergedReason(), SNESConvergedReason, SNESSetTolerances()
270 
271 M*/
272 
273 /*MC
274      SNES_CONVERGED_PNORM_RELATIVE - The 2-norm of the last step <= stol * 2-norm(x) where x is the current
275           solution and stol is the 4th argument to SNESSetTolerances()
276 
277    Level: beginner
278 
279 .seealso:  SNESSolve(), SNESGetConvergedReason(), SNESConvergedReason, SNESSetTolerances()
280 
281 M*/
282 
283 /*MC
284      SNES_DIVERGED_FUNCTION_COUNT - The user provided function has been called more times then the final
285          argument to SNESSetTolerances()
286 
287    Level: beginner
288 
289 .seealso:  SNESSolve(), SNESGetConvergedReason(), SNESConvergedReason, SNESSetTolerances()
290 
291 M*/
292 
293 /*MC
294      SNES_DIVERGED_FNORM_NAN - the 2-norm of the current function evaluation is not-a-number (NaN), this
295       is usually caused by a division of 0 by 0.
296 
297    Level: beginner
298 
299 .seealso:  SNESSolve(), SNESGetConvergedReason(), SNESConvergedReason, SNESSetTolerances()
300 
301 M*/
302 
303 /*MC
304      SNES_DIVERGED_MAX_IT - SNESSolve() has reached the maximum number of iterations requested
305 
306    Level: beginner
307 
308 .seealso:  SNESSolve(), SNESGetConvergedReason(), SNESConvergedReason, SNESSetTolerances()
309 
310 M*/
311 
312 /*MC
313      SNES_DIVERGED_LINE_SEARCH - The line search has failed. This only occurs for a SNESType of SNESLS
314 
315    Level: beginner
316 
317 .seealso:  SNESSolve(), SNESGetConvergedReason(), SNESConvergedReason, SNESSetTolerances()
318 
319 M*/
320 
321 /*MC
322      SNES_DIVERGED_LOCAL_MIN - the algorithm seems to have stagnated at a local minimum that is not zero.
323         See the manual page for SNESConvergedReason for more details
324 
325    Level: beginner
326 
327 .seealso:  SNESSolve(), SNESGetConvergedReason(), SNESConvergedReason, SNESSetTolerances()
328 
329 M*/
330 
331 /*MC
332      SNES_CONERGED_ITERATING - this only occurs if SNESGetConvergedReason() is called during the SNESSolve()
333 
334    Level: beginner
335 
336 .seealso:  SNESSolve(), SNESGetConvergedReason(), SNESConvergedReason, SNESSetTolerances()
337 
338 M*/
339 
340 extern PetscErrorCode  SNESSetConvergenceTest(SNES,PetscErrorCode (*)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void*,PetscErrorCode (*)(void*));
341 extern PetscErrorCode  SNESDefaultConverged(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*);
342 extern PetscErrorCode  SNESSkipConverged(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*);
343 extern PetscErrorCode  SNESGetConvergedReason(SNES,SNESConvergedReason*);
344 
345 extern PetscErrorCode  SNESDAFormFunction(SNES,Vec,Vec,void*);
346 extern PetscErrorCode  SNESDAComputeJacobianWithAdic(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
347 extern PetscErrorCode  SNESDAComputeJacobianWithAdifor(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
348 extern PetscErrorCode  SNESDAComputeJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
349 
350 extern PetscErrorCode  SNESMeshFormFunction(SNES,Vec,Vec,void*);
351 extern PetscErrorCode SNESMeshFormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
352 
353 /* --------- Solving systems of nonlinear equations --------------- */
354 typedef PetscErrorCode (*SNESFunction)(SNES,Vec,Vec,void*);
355 typedef PetscErrorCode (*SNESJacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
356 extern PetscErrorCode  SNESSetFunction(SNES,Vec,SNESFunction,void*);
357 extern PetscErrorCode  SNESGetFunction(SNES,Vec*,SNESFunction*,void**);
358 extern PetscErrorCode  SNESComputeFunction(SNES,Vec,Vec);
359 extern PetscErrorCode  SNESSetJacobian(SNES,Mat,Mat,SNESJacobian,void*);
360 extern PetscErrorCode  SNESGetJacobian(SNES,Mat*,Mat*,SNESJacobian*,void**);
361 extern PetscErrorCode  SNESDefaultComputeJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
362 extern PetscErrorCode  SNESDefaultComputeJacobianColor(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
363 extern PetscErrorCode  SNESSetComputeInitialGuess(SNES,PetscErrorCode (*)(SNES,Vec,void*),void*);
364 
365 /* --------- Routines specifically for line search methods --------------- */
366 /*E
367     SNESLineSearchType - type of line search used in Newton's method as well as VI solvers and Picard solvers
368 
369     Level: beginner
370 
371 .seealso: SNESSetFromOptions(), SNESLineSearchSet()
372 E*/
373 typedef enum {SNES_LS_BASIC, SNES_LS_BASIC_NONORMS, SNES_LS_QUADRATIC, SNES_LS_CUBIC} SNESLineSearchType;
374 extern const char *const SNESLineSearchTypes[];
375 extern const char *SNESLineSearchTypeName(SNESLineSearchType); /* Does bounds checking, use this for viewing */
376 
377 extern PetscErrorCode  SNESLineSearchSet(SNES,PetscErrorCode(*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscBool *),void*);
378 extern PetscErrorCode  SNESLineSearchNo(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscBool *);
379 extern PetscErrorCode  SNESLineSearchNoNorms(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscBool *);
380 extern PetscErrorCode  SNESLineSearchCubic(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscBool *);
381 extern PetscErrorCode  SNESLineSearchQuadratic(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscBool *);
382 
383 extern PetscErrorCode  SNESLineSearchSetPostCheck(SNES,PetscErrorCode(*)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *),void*);
384 extern PetscErrorCode  SNESLineSearchSetPreCheck(SNES,PetscErrorCode(*)(SNES,Vec,Vec,void*,PetscBool *),void*);
385 extern PetscErrorCode  SNESLineSearchSetParams(SNES,PetscReal,PetscReal,PetscReal);
386 extern PetscErrorCode  SNESLineSearchGetParams(SNES,PetscReal*,PetscReal*,PetscReal*);
387 extern PetscErrorCode  SNESLineSearchSetMonitor(SNES,PetscBool );
388 
389 /* Routines for VI solver */
390 extern PetscErrorCode  SNESVISetVariableBounds(SNES,Vec,Vec);
391 extern PetscErrorCode  SNESVISetComputeVariableBounds(SNES, PetscErrorCode (*)(SNES,Vec,Vec));
392 extern PetscErrorCode  SNESVIGetInactiveSet(SNES,IS*);
393 extern PetscErrorCode  SNESVISetRedundancyCheck(SNES,PetscErrorCode(*)(SNES,IS,IS*,void*),void*);
394 #define SNES_VI_INF   1.0e20
395 #define SNES_VI_NINF -1.0e20
396 
397 extern PetscErrorCode  SNESTestLocalMin(SNES);
398 
399 /* Should this routine be private? */
400 extern PetscErrorCode  SNESComputeJacobian(SNES,Vec,Mat*,Mat*,MatStructure*);
401 
402 extern PetscErrorCode SNESSetDM(SNES,DM);
403 extern PetscErrorCode SNESGetDM(SNES,DM*);
404 extern PetscErrorCode SNESSetPC(SNES,SNES);
405 extern PetscErrorCode SNESGetPC(SNES,SNES*);
406 
407 /* Routines for Multiblock solver */
408 PetscErrorCode SNESMultiblockSetFields(SNES, const char [], PetscInt, const PetscInt *);
409 PetscErrorCode SNESMultiblockSetIS(SNES, const char [], IS);
410 PetscErrorCode SNESMultiblockSetBlockSize(SNES, PetscInt);
411 PetscErrorCode SNESMultiblockSetType(SNES, PCCompositeType);
412 
413 PETSC_EXTERN_CXX_END
414 #endif
415