xref: /petsc/include/petscts.h (revision b90f4e8e2d1ae90880ffa0fd3113e4b748be9ce0)
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 extern PetscErrorCode   TSViewFromOptions(TS,const char[]);
207 
208 extern PetscErrorCode   TSSetApplicationContext(TS,void *);
209 extern PetscErrorCode   TSGetApplicationContext(TS,void *);
210 
211 extern PetscErrorCode   TSMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG *);
212 extern PetscErrorCode   TSMonitorLG(TS,PetscInt,PetscReal,Vec,void *);
213 extern PetscErrorCode   TSMonitorLGDestroy(PetscDrawLG*);
214 
215 /*S
216    TSGLAdapt - Abstract object that manages time-step adaptivity
217 
218    Level: beginner
219 
220 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType
221 S*/
222 typedef struct _p_TSGLAdapt *TSGLAdapt;
223 
224 /*E
225     TSGLAdaptType - String with the name of TSGLAdapt scheme or the creation function
226        with an optional dynamic library name, for example
227        http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate()
228 
229    Level: beginner
230 
231 .seealso: TSGLAdaptSetType(), TS
232 E*/
233 #define TSGLAdaptType  char*
234 #define TSGLADAPT_NONE "none"
235 #define TSGLADAPT_SIZE "size"
236 #define TSGLADAPT_BOTH "both"
237 
238 /*MC
239    TSGLAdaptRegisterDynamic - adds a TSGLAdapt implementation
240 
241    Synopsis:
242    PetscErrorCode TSGLAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
243 
244    Not Collective
245 
246    Input Parameters:
247 +  name_scheme - name of user-defined adaptivity scheme
248 .  path - path (either absolute or relative) the library containing this scheme
249 .  name_create - name of routine to create method context
250 -  routine_create - routine to create method context
251 
252    Notes:
253    TSGLAdaptRegisterDynamic() may be called multiple times to add several user-defined families.
254 
255    If dynamic libraries are used, then the fourth input argument (routine_create)
256    is ignored.
257 
258    Sample usage:
259 .vb
260    TSGLAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
261                             "MySchemeCreate",MySchemeCreate);
262 .ve
263 
264    Then, your scheme can be chosen with the procedural interface via
265 $     TSGLAdaptSetType(ts,"my_scheme")
266    or at runtime via the option
267 $     -ts_adapt_type my_scheme
268 
269    Level: advanced
270 
271    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
272           and others of the form ${any_environmental_variable} occuring in pathname will be
273           replaced with appropriate values.
274 
275 .keywords: TSGLAdapt, register
276 
277 .seealso: TSGLAdaptRegisterAll()
278 M*/
279 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
280 #  define TSGLAdaptRegisterDynamic(a,b,c,d)  TSGLAdaptRegister(a,b,c,0)
281 #else
282 #  define TSGLAdaptRegisterDynamic(a,b,c,d)  TSGLAdaptRegister(a,b,c,d)
283 #endif
284 
285 extern PetscErrorCode  TSGLAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSGLAdapt));
286 extern PetscErrorCode  TSGLAdaptRegisterAll(const char[]);
287 extern PetscErrorCode  TSGLAdaptRegisterDestroy(void);
288 extern PetscErrorCode  TSGLAdaptInitializePackage(const char[]);
289 extern PetscErrorCode  TSGLAdaptFinalizePackage(void);
290 extern PetscErrorCode  TSGLAdaptCreate(MPI_Comm,TSGLAdapt*);
291 extern PetscErrorCode  TSGLAdaptSetType(TSGLAdapt,const TSGLAdaptType);
292 extern PetscErrorCode  TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]);
293 extern PetscErrorCode  TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
294 extern PetscErrorCode  TSGLAdaptView(TSGLAdapt,PetscViewer);
295 extern PetscErrorCode  TSGLAdaptSetFromOptions(TSGLAdapt);
296 extern PetscErrorCode  TSGLAdaptDestroy(TSGLAdapt*);
297 
298 /*E
299     TSGLAcceptType - String with the name of TSGLAccept scheme or the function
300        with an optional dynamic library name, for example
301        http://www.mcs.anl.gov/petsc/lib.a:mytsglaccept()
302 
303    Level: beginner
304 
305 .seealso: TSGLSetAcceptType(), TS
306 E*/
307 #define TSGLAcceptType  char*
308 #define TSGLACCEPT_ALWAYS "always"
309 
310 typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *);
311 extern PetscErrorCode  TSGLAcceptRegister(const char[],const char[],const char[],TSGLAcceptFunction);
312 
313 /*MC
314    TSGLAcceptRegisterDynamic - adds a TSGL acceptance scheme
315 
316    Synopsis:
317    PetscErrorCode TSGLAcceptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
318 
319    Not Collective
320 
321    Input Parameters:
322 +  name_scheme - name of user-defined acceptance scheme
323 .  path - path (either absolute or relative) the library containing this scheme
324 .  name_create - name of routine to create method context
325 -  routine_create - routine to create method context
326 
327    Notes:
328    TSGLAcceptRegisterDynamic() may be called multiple times to add several user-defined families.
329 
330    If dynamic libraries are used, then the fourth input argument (routine_create)
331    is ignored.
332 
333    Sample usage:
334 .vb
335    TSGLAcceptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
336                              "MySchemeCreate",MySchemeCreate);
337 .ve
338 
339    Then, your scheme can be chosen with the procedural interface via
340 $     TSGLSetAcceptType(ts,"my_scheme")
341    or at runtime via the option
342 $     -ts_gl_accept_type my_scheme
343 
344    Level: advanced
345 
346    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
347           and others of the form ${any_environmental_variable} occuring in pathname will be
348           replaced with appropriate values.
349 
350 .keywords: TSGL, TSGLAcceptType, register
351 
352 .seealso: TSGLRegisterAll()
353 M*/
354 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
355 #  define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,0)
356 #else
357 #  define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,d)
358 #endif
359 
360 /*E
361   TSGLType - family of time integration method within the General Linear class
362 
363   Level: beginner
364 
365 .seealso: TSGLSetType(), TSGLRegister()
366 E*/
367 #define TSGLType char*
368 #define TSGL_IRKS   "irks"
369 
370 /*MC
371    TSGLRegisterDynamic - adds a TSGL implementation
372 
373    Synopsis:
374    PetscErrorCode TSGLRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
375 
376    Not Collective
377 
378    Input Parameters:
379 +  name_scheme - name of user-defined general linear scheme
380 .  path - path (either absolute or relative) the library containing this scheme
381 .  name_create - name of routine to create method context
382 -  routine_create - routine to create method context
383 
384    Notes:
385    TSGLRegisterDynamic() may be called multiple times to add several user-defined families.
386 
387    If dynamic libraries are used, then the fourth input argument (routine_create)
388    is ignored.
389 
390    Sample usage:
391 .vb
392    TSGLRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
393                        "MySchemeCreate",MySchemeCreate);
394 .ve
395 
396    Then, your scheme can be chosen with the procedural interface via
397 $     TSGLSetType(ts,"my_scheme")
398    or at runtime via the option
399 $     -ts_gl_type my_scheme
400 
401    Level: advanced
402 
403    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
404           and others of the form ${any_environmental_variable} occuring in pathname will be
405           replaced with appropriate values.
406 
407 .keywords: TSGL, register
408 
409 .seealso: TSGLRegisterAll()
410 M*/
411 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
412 #  define TSGLRegisterDynamic(a,b,c,d)       TSGLRegister(a,b,c,0)
413 #else
414 #  define TSGLRegisterDynamic(a,b,c,d)       TSGLRegister(a,b,c,d)
415 #endif
416 
417 extern PetscErrorCode  TSGLRegister(const char[],const char[],const char[],PetscErrorCode(*)(TS));
418 extern PetscErrorCode  TSGLRegisterAll(const char[]);
419 extern PetscErrorCode  TSGLRegisterDestroy(void);
420 extern PetscErrorCode  TSGLInitializePackage(const char[]);
421 extern PetscErrorCode  TSGLFinalizePackage(void);
422 extern PetscErrorCode  TSGLSetType(TS,const TSGLType);
423 extern PetscErrorCode  TSGLGetAdapt(TS,TSGLAdapt*);
424 extern PetscErrorCode  TSGLSetAcceptType(TS,const TSGLAcceptType);
425 
426 #define TSARKIMEXType char*
427 #define TSARKIMEX2D "2d"
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