xref: /petsc/include/petscts.h (revision 4d7d938efd3480b430fca57a45eddde7b805a059)
1 /*
2    User interface for the timestepping package. This package
3    is for use in solving time-dependent PDEs.
4 */
5 #if !defined(__PETSCTS_H)
6 #define __PETSCTS_H
7 #include "petscsnes.h"
8 PETSC_EXTERN_CXX_BEGIN
9 
10 /*S
11      TS - Abstract PETSc object that manages all time-steppers (ODE integrators)
12 
13    Level: beginner
14 
15   Concepts: ODE solvers
16 
17 .seealso:  TSCreate(), TSSetType(), TSType, SNES, KSP, PC
18 S*/
19 typedef struct _p_TS* TS;
20 
21 /*E
22     TSType - String with the name of a PETSc TS method or the creation function
23        with an optional dynamic library name, for example
24        http://www.mcs.anl.gov/petsc/lib.a:mytscreate()
25 
26    Level: beginner
27 
28 .seealso: TSSetType(), TS
29 E*/
30 #define TSType char*
31 #define TSEULER           "euler"
32 #define TSBEULER          "beuler"
33 #define TSPSEUDO          "pseudo"
34 #define TSCN              "cn"
35 #define TSSUNDIALS        "sundials"
36 #define TSRK              "rk"
37 #define TSPYTHON          "python"
38 #define TSTHETA           "theta"
39 #define TSALPHA           "alpha"
40 #define TSGL              "gl"
41 #define TSSSP             "ssp"
42 #define TSARKIMEX         "arkimex"
43 
44 /*E
45     TSProblemType - Determines the type of problem this TS object is to be used to solve
46 
47    Level: beginner
48 
49 .seealso: TSCreate()
50 E*/
51 typedef enum {TS_LINEAR,TS_NONLINEAR} TSProblemType;
52 
53 typedef enum {
54   TS_CONVERGED_ITERATING      = 0,
55   TS_CONVERGED_TIME           = 1,
56   TS_CONVERGED_ITS            = 2,
57   TS_DIVERGED_NONLINEAR_SOLVE = -1,
58   TS_DIVERGED_STEP_REJECTED   = -2
59 } TSConvergedReason;
60 extern const char *const*TSConvergedReasons;
61 
62 /* Logging support */
63 extern PetscClassId  TS_CLASSID;
64 
65 extern PetscErrorCode   TSInitializePackage(const char[]);
66 
67 extern PetscErrorCode   TSCreate(MPI_Comm,TS*);
68 extern PetscErrorCode   TSDestroy(TS*);
69 
70 extern PetscErrorCode   TSSetProblemType(TS,TSProblemType);
71 extern PetscErrorCode   TSGetProblemType(TS,TSProblemType*);
72 extern PetscErrorCode   TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**));
73 extern PetscErrorCode   TSMonitorCancel(TS);
74 
75 extern PetscErrorCode   TSSetOptionsPrefix(TS,const char[]);
76 extern PetscErrorCode   TSAppendOptionsPrefix(TS,const char[]);
77 extern PetscErrorCode   TSGetOptionsPrefix(TS,const char *[]);
78 extern PetscErrorCode   TSSetFromOptions(TS);
79 extern PetscErrorCode   TSSetUp(TS);
80 extern PetscErrorCode   TSReset(TS);
81 
82 extern PetscErrorCode   TSSetSolution(TS,Vec);
83 extern PetscErrorCode   TSGetSolution(TS,Vec*);
84 
85 extern PetscErrorCode   TSSetDuration(TS,PetscInt,PetscReal);
86 extern PetscErrorCode   TSGetDuration(TS,PetscInt*,PetscReal*);
87 
88 extern PetscErrorCode   TSMonitorDefault(TS,PetscInt,PetscReal,Vec,void*);
89 extern PetscErrorCode   TSMonitorSolution(TS,PetscInt,PetscReal,Vec,void*);
90 extern PetscErrorCode   TSStep(TS);
91 extern PetscErrorCode   TSSolve(TS,Vec);
92 extern PetscErrorCode   TSGetConvergedReason(TS,TSConvergedReason*);
93 
94 extern PetscErrorCode   TSSetInitialTimeStep(TS,PetscReal,PetscReal);
95 extern PetscErrorCode   TSGetTimeStep(TS,PetscReal*);
96 extern PetscErrorCode   TSGetTime(TS,PetscReal*);
97 extern PetscErrorCode   TSSetTime(TS,PetscReal);
98 extern PetscErrorCode   TSGetTimeStepNumber(TS,PetscInt*);
99 extern PetscErrorCode   TSSetTimeStep(TS,PetscReal);
100 
101 typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*);
102 typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
103 extern PetscErrorCode   TSSetRHSFunction(TS,Vec,TSRHSFunction,void*);
104 extern PetscErrorCode   TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**);
105 extern PetscErrorCode   TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*);
106 extern PetscErrorCode   TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**);
107 
108 typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*);
109 typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*);
110 extern PetscErrorCode   TSSetIFunction(TS,Vec,TSIFunction,void*);
111 extern PetscErrorCode   TSGetIFunction(TS,Vec*,TSIFunction*,void**);
112 extern PetscErrorCode   TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*);
113 extern PetscErrorCode   TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**);
114 
115 extern PetscErrorCode   TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*);
116 extern PetscErrorCode   TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
117 extern PetscErrorCode   TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*);
118 extern PetscErrorCode   TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*);
119 extern PetscErrorCode   TSDefaultComputeJacobianColor(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
120 extern PetscErrorCode   TSDefaultComputeJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
121 
122 extern PetscErrorCode   TSSetPreStep(TS, PetscErrorCode (*)(TS));
123 extern PetscErrorCode   TSSetPostStep(TS, PetscErrorCode (*)(TS));
124 extern PetscErrorCode   TSPreStep(TS);
125 extern PetscErrorCode   TSPostStep(TS);
126 
127 extern PetscErrorCode   TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*);
128 extern PetscErrorCode   TSPseudoDefaultTimeStep(TS,PetscReal*,void*);
129 extern PetscErrorCode   TSPseudoComputeTimeStep(TS,PetscReal *);
130 
131 extern PetscErrorCode   TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*);
132 extern PetscErrorCode   TSPseudoDefaultVerifyTimeStep(TS,Vec,void*,PetscReal*,PetscBool *);
133 extern PetscErrorCode   TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *);
134 extern PetscErrorCode   TSPseudoSetTimeStepIncrement(TS,PetscReal);
135 extern PetscErrorCode   TSPseudoIncrementDtFromInitialDt(TS);
136 
137 extern PetscErrorCode   TSPythonSetType(TS,const char[]);
138 
139 extern PetscErrorCode   TSComputeRHSFunction(TS,PetscReal,Vec,Vec);
140 extern PetscErrorCode   TSComputeRHSJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*);
141 extern PetscErrorCode   TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool);
142 extern PetscErrorCode   TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,PetscBool);
143 
144 /* Dynamic creation and loading functions */
145 extern PetscFList TSList;
146 extern PetscBool  TSRegisterAllCalled;
147 extern PetscErrorCode   TSGetType(TS,const TSType*);
148 extern PetscErrorCode   TSSetType(TS,const TSType);
149 extern PetscErrorCode   TSRegister(const char[], const char[], const char[], PetscErrorCode (*)(TS));
150 extern PetscErrorCode   TSRegisterAll(const char[]);
151 extern PetscErrorCode   TSRegisterDestroy(void);
152 
153 /*MC
154   TSRegisterDynamic - Adds a creation method to the TS package.
155 
156   Synopsis:
157   PetscErrorCode TSRegisterDynamic(const char *name, const char *path, const char *func_name, PetscErrorCode (*create_func)(TS))
158 
159   Not Collective
160 
161   Input Parameters:
162 + name        - The name of a new user-defined creation routine
163 . path        - The path (either absolute or relative) of the library containing this routine
164 . func_name   - The name of the creation routine
165 - create_func - The creation routine itself
166 
167   Notes:
168   TSRegisterDynamic() may be called multiple times to add several user-defined tses.
169 
170   If dynamic libraries are used, then the fourth input argument (create_func) is ignored.
171 
172   Sample usage:
173 .vb
174   TSRegisterDynamic("my_ts", "/home/username/my_lib/lib/libO/solaris/libmy.a", "MyTSCreate", MyTSCreate);
175 .ve
176 
177   Then, your ts type can be chosen with the procedural interface via
178 .vb
179     TS ts;
180     TSCreate(MPI_Comm, &ts);
181     TSSetType(ts, "my_ts")
182 .ve
183   or at runtime via the option
184 .vb
185     -ts_type my_ts
186 .ve
187 
188   Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
189         If your function is not being put into a shared library then use TSRegister() instead
190 
191   Level: advanced
192 
193 .keywords: TS, register
194 .seealso: TSRegisterAll(), TSRegisterDestroy()
195 M*/
196 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
197 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,0)
198 #else
199 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,d)
200 #endif
201 
202 extern PetscErrorCode   TSGetSNES(TS,SNES*);
203 extern PetscErrorCode   TSGetKSP(TS,KSP*);
204 
205 extern PetscErrorCode   TSView(TS,PetscViewer);
206 
207 extern PetscErrorCode   TSSetApplicationContext(TS,void *);
208 extern PetscErrorCode   TSGetApplicationContext(TS,void *);
209 
210 extern PetscErrorCode   TSMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG *);
211 extern PetscErrorCode   TSMonitorLG(TS,PetscInt,PetscReal,Vec,void *);
212 extern PetscErrorCode   TSMonitorLGDestroy(PetscDrawLG*);
213 
214 /*S
215    TSGLAdapt - Abstract object that manages time-step adaptivity
216 
217    Level: beginner
218 
219 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType
220 S*/
221 typedef struct _p_TSGLAdapt *TSGLAdapt;
222 
223 /*E
224     TSGLAdaptType - String with the name of TSGLAdapt scheme or the creation function
225        with an optional dynamic library name, for example
226        http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate()
227 
228    Level: beginner
229 
230 .seealso: TSGLAdaptSetType(), TS
231 E*/
232 #define TSGLAdaptType  char*
233 #define TSGLADAPT_NONE "none"
234 #define TSGLADAPT_SIZE "size"
235 #define TSGLADAPT_BOTH "both"
236 
237 /*MC
238    TSGLAdaptRegisterDynamic - adds a TSGLAdapt implementation
239 
240    Synopsis:
241    PetscErrorCode TSGLAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
242 
243    Not Collective
244 
245    Input Parameters:
246 +  name_scheme - name of user-defined adaptivity scheme
247 .  path - path (either absolute or relative) the library containing this scheme
248 .  name_create - name of routine to create method context
249 -  routine_create - routine to create method context
250 
251    Notes:
252    TSGLAdaptRegisterDynamic() may be called multiple times to add several user-defined families.
253 
254    If dynamic libraries are used, then the fourth input argument (routine_create)
255    is ignored.
256 
257    Sample usage:
258 .vb
259    TSGLAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
260                             "MySchemeCreate",MySchemeCreate);
261 .ve
262 
263    Then, your scheme can be chosen with the procedural interface via
264 $     TSGLAdaptSetType(ts,"my_scheme")
265    or at runtime via the option
266 $     -ts_adapt_type my_scheme
267 
268    Level: advanced
269 
270    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
271           and others of the form ${any_environmental_variable} occuring in pathname will be
272           replaced with appropriate values.
273 
274 .keywords: TSGLAdapt, register
275 
276 .seealso: TSGLAdaptRegisterAll()
277 M*/
278 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
279 #  define TSGLAdaptRegisterDynamic(a,b,c,d)  TSGLAdaptRegister(a,b,c,0)
280 #else
281 #  define TSGLAdaptRegisterDynamic(a,b,c,d)  TSGLAdaptRegister(a,b,c,d)
282 #endif
283 
284 extern PetscErrorCode  TSGLAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSGLAdapt));
285 extern PetscErrorCode  TSGLAdaptRegisterAll(const char[]);
286 extern PetscErrorCode  TSGLAdaptRegisterDestroy(void);
287 extern PetscErrorCode  TSGLAdaptInitializePackage(const char[]);
288 extern PetscErrorCode  TSGLAdaptFinalizePackage(void);
289 extern PetscErrorCode  TSGLAdaptCreate(MPI_Comm,TSGLAdapt*);
290 extern PetscErrorCode  TSGLAdaptSetType(TSGLAdapt,const TSGLAdaptType);
291 extern PetscErrorCode  TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]);
292 extern PetscErrorCode  TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
293 extern PetscErrorCode  TSGLAdaptView(TSGLAdapt,PetscViewer);
294 extern PetscErrorCode  TSGLAdaptSetFromOptions(TSGLAdapt);
295 extern PetscErrorCode  TSGLAdaptDestroy(TSGLAdapt*);
296 
297 /*E
298     TSGLAcceptType - String with the name of TSGLAccept scheme or the function
299        with an optional dynamic library name, for example
300        http://www.mcs.anl.gov/petsc/lib.a:mytsglaccept()
301 
302    Level: beginner
303 
304 .seealso: TSGLSetAcceptType(), TS
305 E*/
306 #define TSGLAcceptType  char*
307 #define TSGLACCEPT_ALWAYS "always"
308 
309 typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *);
310 extern PetscErrorCode  TSGLAcceptRegister(const char[],const char[],const char[],TSGLAcceptFunction);
311 
312 /*MC
313    TSGLAcceptRegisterDynamic - adds a TSGL acceptance scheme
314 
315    Synopsis:
316    PetscErrorCode TSGLAcceptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
317 
318    Not Collective
319 
320    Input Parameters:
321 +  name_scheme - name of user-defined acceptance scheme
322 .  path - path (either absolute or relative) the library containing this scheme
323 .  name_create - name of routine to create method context
324 -  routine_create - routine to create method context
325 
326    Notes:
327    TSGLAcceptRegisterDynamic() may be called multiple times to add several user-defined families.
328 
329    If dynamic libraries are used, then the fourth input argument (routine_create)
330    is ignored.
331 
332    Sample usage:
333 .vb
334    TSGLAcceptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
335                              "MySchemeCreate",MySchemeCreate);
336 .ve
337 
338    Then, your scheme can be chosen with the procedural interface via
339 $     TSGLSetAcceptType(ts,"my_scheme")
340    or at runtime via the option
341 $     -ts_gl_accept_type my_scheme
342 
343    Level: advanced
344 
345    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
346           and others of the form ${any_environmental_variable} occuring in pathname will be
347           replaced with appropriate values.
348 
349 .keywords: TSGL, TSGLAcceptType, register
350 
351 .seealso: TSGLRegisterAll()
352 M*/
353 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
354 #  define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,0)
355 #else
356 #  define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,d)
357 #endif
358 
359 /*E
360   TSGLType - family of time integration method within the General Linear class
361 
362   Level: beginner
363 
364 .seealso: TSGLSetType(), TSGLRegister()
365 E*/
366 #define TSGLType char*
367 #define TSGL_IRKS   "irks"
368 
369 /*MC
370    TSGLRegisterDynamic - adds a TSGL implementation
371 
372    Synopsis:
373    PetscErrorCode TSGLRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
374 
375    Not Collective
376 
377    Input Parameters:
378 +  name_scheme - name of user-defined general linear scheme
379 .  path - path (either absolute or relative) the library containing this scheme
380 .  name_create - name of routine to create method context
381 -  routine_create - routine to create method context
382 
383    Notes:
384    TSGLRegisterDynamic() may be called multiple times to add several user-defined families.
385 
386    If dynamic libraries are used, then the fourth input argument (routine_create)
387    is ignored.
388 
389    Sample usage:
390 .vb
391    TSGLRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
392                        "MySchemeCreate",MySchemeCreate);
393 .ve
394 
395    Then, your scheme can be chosen with the procedural interface via
396 $     TSGLSetType(ts,"my_scheme")
397    or at runtime via the option
398 $     -ts_gl_type my_scheme
399 
400    Level: advanced
401 
402    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
403           and others of the form ${any_environmental_variable} occuring in pathname will be
404           replaced with appropriate values.
405 
406 .keywords: TSGL, register
407 
408 .seealso: TSGLRegisterAll()
409 M*/
410 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
411 #  define TSGLRegisterDynamic(a,b,c,d)       TSGLRegister(a,b,c,0)
412 #else
413 #  define TSGLRegisterDynamic(a,b,c,d)       TSGLRegister(a,b,c,d)
414 #endif
415 
416 extern PetscErrorCode  TSGLRegister(const char[],const char[],const char[],PetscErrorCode(*)(TS));
417 extern PetscErrorCode  TSGLRegisterAll(const char[]);
418 extern PetscErrorCode  TSGLRegisterDestroy(void);
419 extern PetscErrorCode  TSGLInitializePackage(const char[]);
420 extern PetscErrorCode  TSGLFinalizePackage(void);
421 extern PetscErrorCode  TSGLSetType(TS,const TSGLType);
422 extern PetscErrorCode  TSGLGetAdapt(TS,TSGLAdapt*);
423 extern PetscErrorCode  TSGLSetAcceptType(TS,const TSGLAcceptType);
424 
425 #define TSARKIMEXType char*
426 #define TSARKIMEX2D "2d"
427 #define TSARKIMEX2E "2e"
428 #define TSARKIMEX3  "3"
429 #define TSARKIMEX4  "4"
430 #define TSARKIMEX5  "5"
431 extern PetscErrorCode TSARKIMEXGetType(TS ts,const TSARKIMEXType*);
432 extern PetscErrorCode TSARKIMEXSetType(TS ts,const TSARKIMEXType);
433 extern PetscErrorCode TSARKIMEXRegister(const TSARKIMEXType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[]);
434 extern PetscErrorCode TSARKIMEXFinalizePackage(void);
435 extern PetscErrorCode TSARKIMEXInitializePackage(const char path[]);
436 extern PetscErrorCode TSARKIMEXRegisterDestroy(void);
437 extern PetscErrorCode TSARKIMEXRegisterAll(void);
438 
439 /*
440        PETSc interface to Sundials
441 */
442 #ifdef PETSC_HAVE_SUNDIALS
443 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
444 extern const char *TSSundialsLmmTypes[];
445 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
446 extern const char *TSSundialsGramSchmidtTypes[];
447 extern PetscErrorCode   TSSundialsSetType(TS,TSSundialsLmmType);
448 extern PetscErrorCode   TSSundialsGetPC(TS,PC*);
449 extern PetscErrorCode   TSSundialsSetTolerance(TS,PetscReal,PetscReal);
450 extern PetscErrorCode   TSSundialsSetMinTimeStep(TS,PetscReal);
451 extern PetscErrorCode   TSSundialsSetMaxTimeStep(TS,PetscReal);
452 extern PetscErrorCode   TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
453 extern PetscErrorCode   TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
454 extern PetscErrorCode   TSSundialsSetGMRESRestart(TS,PetscInt);
455 extern PetscErrorCode   TSSundialsSetLinearTolerance(TS,PetscReal);
456 extern PetscErrorCode   TSSundialsSetExactFinalTime(TS,PetscBool );
457 extern PetscErrorCode   TSSundialsMonitorInternalSteps(TS,PetscBool );
458 extern PetscErrorCode   TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
459 #endif
460 
461 extern PetscErrorCode   TSRKSetTolerance(TS,PetscReal);
462 
463 extern PetscErrorCode  TSThetaSetTheta(TS,PetscReal);
464 extern PetscErrorCode  TSThetaGetTheta(TS,PetscReal*);
465 
466 extern PetscErrorCode  TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*);
467 extern PetscErrorCode  TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*);
468 extern PetscErrorCode  TSAlphaSetRadius(TS,PetscReal);
469 extern PetscErrorCode  TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal);
470 extern PetscErrorCode  TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*);
471 
472 extern PetscErrorCode  TSSetDM(TS,DM);
473 extern PetscErrorCode  TSGetDM(TS,DM*);
474 
475 extern PetscErrorCode  SNESTSFormFunction(SNES,Vec,Vec,void*);
476 extern PetscErrorCode  SNESTSFormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
477 
478 PETSC_EXTERN_CXX_END
479 #endif
480