xref: /petsc/include/petscts.h (revision b0f836d7bc5c231aef07e9463d99244d64e4598f)
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 /*J
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 J*/
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 #define TSROSW            "rosw"
44 
45 /*E
46     TSProblemType - Determines the type of problem this TS object is to be used to solve
47 
48    Level: beginner
49 
50 .seealso: TSCreate()
51 E*/
52 typedef enum {TS_LINEAR,TS_NONLINEAR} TSProblemType;
53 
54 typedef enum {
55   TS_CONVERGED_ITERATING      = 0,
56   TS_CONVERGED_TIME           = 1,
57   TS_CONVERGED_ITS            = 2,
58   TS_DIVERGED_NONLINEAR_SOLVE = -1,
59   TS_DIVERGED_STEP_REJECTED   = -2
60 } TSConvergedReason;
61 extern const char *const*TSConvergedReasons;
62 
63 /* Logging support */
64 extern PetscClassId  TS_CLASSID;
65 
66 extern PetscErrorCode   TSInitializePackage(const char[]);
67 
68 extern PetscErrorCode   TSCreate(MPI_Comm,TS*);
69 extern PetscErrorCode   TSDestroy(TS*);
70 
71 extern PetscErrorCode   TSSetProblemType(TS,TSProblemType);
72 extern PetscErrorCode   TSGetProblemType(TS,TSProblemType*);
73 extern PetscErrorCode   TSMonitor(TS,PetscInt,PetscReal,Vec);
74 extern PetscErrorCode   TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**));
75 extern PetscErrorCode   TSMonitorCancel(TS);
76 
77 extern PetscErrorCode   TSSetOptionsPrefix(TS,const char[]);
78 extern PetscErrorCode   TSAppendOptionsPrefix(TS,const char[]);
79 extern PetscErrorCode   TSGetOptionsPrefix(TS,const char *[]);
80 extern PetscErrorCode   TSSetFromOptions(TS);
81 extern PetscErrorCode   TSSetUp(TS);
82 extern PetscErrorCode   TSReset(TS);
83 
84 extern PetscErrorCode   TSSetSolution(TS,Vec);
85 extern PetscErrorCode   TSGetSolution(TS,Vec*);
86 
87 extern PetscErrorCode   TSSetDuration(TS,PetscInt,PetscReal);
88 extern PetscErrorCode   TSGetDuration(TS,PetscInt*,PetscReal*);
89 extern PetscErrorCode   TSSetExactFinalTime(TS,PetscBool);
90 
91 extern PetscErrorCode   TSMonitorDefault(TS,PetscInt,PetscReal,Vec,void*);
92 extern PetscErrorCode   TSMonitorSolution(TS,PetscInt,PetscReal,Vec,void*);
93 extern PetscErrorCode   TSMonitorSolutionCreate(TS,PetscViewer,void**);
94 extern PetscErrorCode   TSMonitorSolutionDestroy(void**);
95 extern PetscErrorCode   TSMonitorSolutionBinary(TS,PetscInt,PetscReal,Vec,void*);
96 extern PetscErrorCode   TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*);
97 extern PetscErrorCode   TSMonitorSolutionVTKDestroy(void*);
98 
99 extern PetscErrorCode   TSStep(TS);
100 extern PetscErrorCode   TSEvaluateStep(TS,PetscInt,Vec,PetscBool*);
101 extern PetscErrorCode   TSSolve(TS,Vec,PetscReal*);
102 extern PetscErrorCode   TSGetConvergedReason(TS,TSConvergedReason*);
103 extern PetscErrorCode   TSGetNonlinearSolveIterations(TS,PetscInt*);
104 extern PetscErrorCode   TSGetLinearSolveIterations(TS,PetscInt*);
105 
106 extern PetscErrorCode   TSSetInitialTimeStep(TS,PetscReal,PetscReal);
107 extern PetscErrorCode   TSGetTimeStep(TS,PetscReal*);
108 extern PetscErrorCode   TSGetTime(TS,PetscReal*);
109 extern PetscErrorCode   TSSetTime(TS,PetscReal);
110 extern PetscErrorCode   TSGetTimeStepNumber(TS,PetscInt*);
111 extern PetscErrorCode   TSSetTimeStep(TS,PetscReal);
112 
113 typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*);
114 typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
115 extern PetscErrorCode   TSSetRHSFunction(TS,Vec,TSRHSFunction,void*);
116 extern PetscErrorCode   TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**);
117 extern PetscErrorCode   TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*);
118 extern PetscErrorCode   TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**);
119 
120 typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*);
121 typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*);
122 extern PetscErrorCode   TSSetIFunction(TS,Vec,TSIFunction,void*);
123 extern PetscErrorCode   TSGetIFunction(TS,Vec*,TSIFunction*,void**);
124 extern PetscErrorCode   TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*);
125 extern PetscErrorCode   TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**);
126 
127 extern PetscErrorCode   TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*);
128 extern PetscErrorCode   TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
129 extern PetscErrorCode   TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*);
130 extern PetscErrorCode   TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*);
131 extern PetscErrorCode   TSDefaultComputeJacobianColor(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
132 extern PetscErrorCode   TSDefaultComputeJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
133 
134 extern PetscErrorCode   TSSetPreStep(TS, PetscErrorCode (*)(TS));
135 extern PetscErrorCode   TSSetPostStep(TS, PetscErrorCode (*)(TS));
136 extern PetscErrorCode   TSPreStep(TS);
137 extern PetscErrorCode   TSPostStep(TS);
138 extern PetscErrorCode   TSSetRetainStages(TS,PetscBool);
139 extern PetscErrorCode   TSInterpolate(TS,PetscReal,Vec);
140 extern PetscErrorCode   TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec);
141 extern PetscErrorCode   TSErrorNormWRMS(TS,Vec,PetscReal*);
142 extern PetscErrorCode   TSSetCFLTimeLocal(TS,PetscReal);
143 extern PetscErrorCode   TSGetCFLTime(TS,PetscReal*);
144 
145 extern PetscErrorCode   TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*);
146 extern PetscErrorCode   TSPseudoDefaultTimeStep(TS,PetscReal*,void*);
147 extern PetscErrorCode   TSPseudoComputeTimeStep(TS,PetscReal *);
148 extern PetscErrorCode   TSPseudoSetMaxTimeStep(TS,PetscReal);
149 
150 extern PetscErrorCode   TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*);
151 extern PetscErrorCode   TSPseudoDefaultVerifyTimeStep(TS,Vec,void*,PetscReal*,PetscBool *);
152 extern PetscErrorCode   TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *);
153 extern PetscErrorCode   TSPseudoSetTimeStepIncrement(TS,PetscReal);
154 extern PetscErrorCode   TSPseudoIncrementDtFromInitialDt(TS);
155 
156 extern PetscErrorCode   TSPythonSetType(TS,const char[]);
157 
158 extern PetscErrorCode   TSComputeRHSFunction(TS,PetscReal,Vec,Vec);
159 extern PetscErrorCode   TSComputeRHSJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*);
160 extern PetscErrorCode   TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool);
161 extern PetscErrorCode   TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,PetscBool);
162 
163 extern PetscErrorCode   TSVISetVariableBounds(TS,Vec,Vec);
164 
165 /* Dynamic creation and loading functions */
166 extern PetscFList TSList;
167 extern PetscBool  TSRegisterAllCalled;
168 extern PetscErrorCode   TSGetType(TS,const TSType*);
169 extern PetscErrorCode   TSSetType(TS,const TSType);
170 extern PetscErrorCode   TSRegister(const char[], const char[], const char[], PetscErrorCode (*)(TS));
171 extern PetscErrorCode   TSRegisterAll(const char[]);
172 extern PetscErrorCode   TSRegisterDestroy(void);
173 
174 /*MC
175   TSRegisterDynamic - Adds a creation method to the TS package.
176 
177   Synopsis:
178   PetscErrorCode TSRegisterDynamic(const char *name, const char *path, const char *func_name, PetscErrorCode (*create_func)(TS))
179 
180   Not Collective
181 
182   Input Parameters:
183 + name        - The name of a new user-defined creation routine
184 . path        - The path (either absolute or relative) of the library containing this routine
185 . func_name   - The name of the creation routine
186 - create_func - The creation routine itself
187 
188   Notes:
189   TSRegisterDynamic() may be called multiple times to add several user-defined tses.
190 
191   If dynamic libraries are used, then the fourth input argument (create_func) is ignored.
192 
193   Sample usage:
194 .vb
195   TSRegisterDynamic("my_ts", "/home/username/my_lib/lib/libO/solaris/libmy.a", "MyTSCreate", MyTSCreate);
196 .ve
197 
198   Then, your ts type can be chosen with the procedural interface via
199 .vb
200     TS ts;
201     TSCreate(MPI_Comm, &ts);
202     TSSetType(ts, "my_ts")
203 .ve
204   or at runtime via the option
205 .vb
206     -ts_type my_ts
207 .ve
208 
209   Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
210         If your function is not being put into a shared library then use TSRegister() instead
211 
212   Level: advanced
213 
214 .keywords: TS, register
215 .seealso: TSRegisterAll(), TSRegisterDestroy()
216 M*/
217 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
218 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,0)
219 #else
220 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,d)
221 #endif
222 
223 extern PetscErrorCode   TSGetSNES(TS,SNES*);
224 extern PetscErrorCode   TSGetKSP(TS,KSP*);
225 
226 extern PetscErrorCode   TSView(TS,PetscViewer);
227 
228 extern PetscErrorCode   TSSetApplicationContext(TS,void *);
229 extern PetscErrorCode   TSGetApplicationContext(TS,void *);
230 
231 extern PetscErrorCode   TSMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG *);
232 extern PetscErrorCode   TSMonitorLG(TS,PetscInt,PetscReal,Vec,void *);
233 extern PetscErrorCode   TSMonitorLGDestroy(PetscDrawLG*);
234 
235 /*J
236    TSSSPType - string with the name of TSSSP scheme.
237 
238    Level: beginner
239 
240 .seealso: TSSSPSetType(), TS
241 J*/
242 #define TSSSPType char*
243 #define TSSSPRKS2  "rks2"
244 #define TSSSPRKS3  "rks3"
245 #define TSSSPRK104 "rk104"
246 
247 extern PetscErrorCode TSSSPSetType(TS,const TSSSPType);
248 extern PetscErrorCode TSSSPGetType(TS,const TSSSPType*);
249 extern PetscErrorCode TSSSPSetNumStages(TS,PetscInt);
250 extern PetscErrorCode TSSSPGetNumStages(TS,PetscInt*);
251 
252 /*S
253    TSAdapt - Abstract object that manages time-step adaptivity
254 
255    Level: beginner
256 
257 .seealso: TS, TSAdaptCreate(), TSAdaptType
258 S*/
259 typedef struct _p_TSAdapt *TSAdapt;
260 
261 /*E
262     TSAdaptType - String with the name of TSAdapt scheme or the creation function
263        with an optional dynamic library name, for example
264        http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate()
265 
266    Level: beginner
267 
268 .seealso: TSAdaptSetType(), TS
269 E*/
270 #define TSAdaptType  char*
271 #define TSADAPTBASIC "basic"
272 #define TSADAPTNONE  "none"
273 #define TSADAPTCFL   "cfl"
274 
275 /*MC
276    TSAdaptRegisterDynamic - adds a TSAdapt implementation
277 
278    Synopsis:
279    PetscErrorCode TSAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
280 
281    Not Collective
282 
283    Input Parameters:
284 +  name_scheme - name of user-defined adaptivity scheme
285 .  path - path (either absolute or relative) the library containing this scheme
286 .  name_create - name of routine to create method context
287 -  routine_create - routine to create method context
288 
289    Notes:
290    TSAdaptRegisterDynamic() may be called multiple times to add several user-defined families.
291 
292    If dynamic libraries are used, then the fourth input argument (routine_create)
293    is ignored.
294 
295    Sample usage:
296 .vb
297    TSAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
298                             "MySchemeCreate",MySchemeCreate);
299 .ve
300 
301    Then, your scheme can be chosen with the procedural interface via
302 $     TSAdaptSetType(ts,"my_scheme")
303    or at runtime via the option
304 $     -ts_adapt_type my_scheme
305 
306    Level: advanced
307 
308    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
309           and others of the form ${any_environmental_variable} occuring in pathname will be
310           replaced with appropriate values.
311 
312 .keywords: TSAdapt, register
313 
314 .seealso: TSAdaptRegisterAll()
315 M*/
316 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
317 #  define TSAdaptRegisterDynamic(a,b,c,d)  TSAdaptRegister(a,b,c,0)
318 #else
319 #  define TSAdaptRegisterDynamic(a,b,c,d)  TSAdaptRegister(a,b,c,d)
320 #endif
321 
322 extern PetscErrorCode TSGetAdapt(TS,TSAdapt*);
323 extern PetscErrorCode TSAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSAdapt));
324 extern PetscErrorCode TSAdaptRegisterAll(const char[]);
325 extern PetscErrorCode TSAdaptRegisterDestroy(void);
326 extern PetscErrorCode TSAdaptInitializePackage(const char[]);
327 extern PetscErrorCode TSAdaptFinalizePackage(void);
328 extern PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*);
329 extern PetscErrorCode TSAdaptSetType(TSAdapt,const TSAdaptType);
330 extern PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]);
331 extern PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
332 extern PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool);
333 extern PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**);
334 extern PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*);
335 extern PetscErrorCode TSAdaptView(TSAdapt,PetscViewer);
336 extern PetscErrorCode TSAdaptSetFromOptions(TSAdapt);
337 extern PetscErrorCode TSAdaptDestroy(TSAdapt*);
338 extern PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool);
339 extern PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal);
340 
341 /*S
342    TSGLAdapt - Abstract object that manages time-step adaptivity
343 
344    Level: beginner
345 
346    Developer Notes:
347    This functionality should be replaced by the TSAdapt.
348 
349 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType
350 S*/
351 typedef struct _p_TSGLAdapt *TSGLAdapt;
352 
353 /*J
354     TSGLAdaptType - String with the name of TSGLAdapt scheme or the creation function
355        with an optional dynamic library name, for example
356        http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate()
357 
358    Level: beginner
359 
360 .seealso: TSGLAdaptSetType(), TS
361 J*/
362 #define TSGLAdaptType  char*
363 #define TSGLADAPT_NONE "none"
364 #define TSGLADAPT_SIZE "size"
365 #define TSGLADAPT_BOTH "both"
366 
367 /*MC
368    TSGLAdaptRegisterDynamic - adds a TSGLAdapt implementation
369 
370    Synopsis:
371    PetscErrorCode TSGLAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
372 
373    Not Collective
374 
375    Input Parameters:
376 +  name_scheme - name of user-defined adaptivity scheme
377 .  path - path (either absolute or relative) the library containing this scheme
378 .  name_create - name of routine to create method context
379 -  routine_create - routine to create method context
380 
381    Notes:
382    TSGLAdaptRegisterDynamic() may be called multiple times to add several user-defined families.
383 
384    If dynamic libraries are used, then the fourth input argument (routine_create)
385    is ignored.
386 
387    Sample usage:
388 .vb
389    TSGLAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
390                             "MySchemeCreate",MySchemeCreate);
391 .ve
392 
393    Then, your scheme can be chosen with the procedural interface via
394 $     TSGLAdaptSetType(ts,"my_scheme")
395    or at runtime via the option
396 $     -ts_adapt_type my_scheme
397 
398    Level: advanced
399 
400    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
401           and others of the form ${any_environmental_variable} occuring in pathname will be
402           replaced with appropriate values.
403 
404 .keywords: TSGLAdapt, register
405 
406 .seealso: TSGLAdaptRegisterAll()
407 M*/
408 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
409 #  define TSGLAdaptRegisterDynamic(a,b,c,d)  TSGLAdaptRegister(a,b,c,0)
410 #else
411 #  define TSGLAdaptRegisterDynamic(a,b,c,d)  TSGLAdaptRegister(a,b,c,d)
412 #endif
413 
414 extern PetscErrorCode  TSGLAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSGLAdapt));
415 extern PetscErrorCode  TSGLAdaptRegisterAll(const char[]);
416 extern PetscErrorCode  TSGLAdaptRegisterDestroy(void);
417 extern PetscErrorCode  TSGLAdaptInitializePackage(const char[]);
418 extern PetscErrorCode  TSGLAdaptFinalizePackage(void);
419 extern PetscErrorCode  TSGLAdaptCreate(MPI_Comm,TSGLAdapt*);
420 extern PetscErrorCode  TSGLAdaptSetType(TSGLAdapt,const TSGLAdaptType);
421 extern PetscErrorCode  TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]);
422 extern PetscErrorCode  TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
423 extern PetscErrorCode  TSGLAdaptView(TSGLAdapt,PetscViewer);
424 extern PetscErrorCode  TSGLAdaptSetFromOptions(TSGLAdapt);
425 extern PetscErrorCode  TSGLAdaptDestroy(TSGLAdapt*);
426 
427 /*J
428     TSGLAcceptType - String with the name of TSGLAccept scheme or the function
429        with an optional dynamic library name, for example
430        http://www.mcs.anl.gov/petsc/lib.a:mytsglaccept()
431 
432    Level: beginner
433 
434 .seealso: TSGLSetAcceptType(), TS
435 J*/
436 #define TSGLAcceptType  char*
437 #define TSGLACCEPT_ALWAYS "always"
438 
439 typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *);
440 extern PetscErrorCode  TSGLAcceptRegister(const char[],const char[],const char[],TSGLAcceptFunction);
441 
442 /*MC
443    TSGLAcceptRegisterDynamic - adds a TSGL acceptance scheme
444 
445    Synopsis:
446    PetscErrorCode TSGLAcceptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
447 
448    Not Collective
449 
450    Input Parameters:
451 +  name_scheme - name of user-defined acceptance scheme
452 .  path - path (either absolute or relative) the library containing this scheme
453 .  name_create - name of routine to create method context
454 -  routine_create - routine to create method context
455 
456    Notes:
457    TSGLAcceptRegisterDynamic() may be called multiple times to add several user-defined families.
458 
459    If dynamic libraries are used, then the fourth input argument (routine_create)
460    is ignored.
461 
462    Sample usage:
463 .vb
464    TSGLAcceptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
465                              "MySchemeCreate",MySchemeCreate);
466 .ve
467 
468    Then, your scheme can be chosen with the procedural interface via
469 $     TSGLSetAcceptType(ts,"my_scheme")
470    or at runtime via the option
471 $     -ts_gl_accept_type my_scheme
472 
473    Level: advanced
474 
475    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
476           and others of the form ${any_environmental_variable} occuring in pathname will be
477           replaced with appropriate values.
478 
479 .keywords: TSGL, TSGLAcceptType, register
480 
481 .seealso: TSGLRegisterAll()
482 M*/
483 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
484 #  define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,0)
485 #else
486 #  define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,d)
487 #endif
488 
489 /*J
490   TSGLType - family of time integration method within the General Linear class
491 
492   Level: beginner
493 
494 .seealso: TSGLSetType(), TSGLRegister()
495 J*/
496 #define TSGLType char*
497 #define TSGL_IRKS   "irks"
498 
499 /*MC
500    TSGLRegisterDynamic - adds a TSGL implementation
501 
502    Synopsis:
503    PetscErrorCode TSGLRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
504 
505    Not Collective
506 
507    Input Parameters:
508 +  name_scheme - name of user-defined general linear scheme
509 .  path - path (either absolute or relative) the library containing this scheme
510 .  name_create - name of routine to create method context
511 -  routine_create - routine to create method context
512 
513    Notes:
514    TSGLRegisterDynamic() may be called multiple times to add several user-defined families.
515 
516    If dynamic libraries are used, then the fourth input argument (routine_create)
517    is ignored.
518 
519    Sample usage:
520 .vb
521    TSGLRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
522                        "MySchemeCreate",MySchemeCreate);
523 .ve
524 
525    Then, your scheme can be chosen with the procedural interface via
526 $     TSGLSetType(ts,"my_scheme")
527    or at runtime via the option
528 $     -ts_gl_type my_scheme
529 
530    Level: advanced
531 
532    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
533           and others of the form ${any_environmental_variable} occuring in pathname will be
534           replaced with appropriate values.
535 
536 .keywords: TSGL, register
537 
538 .seealso: TSGLRegisterAll()
539 M*/
540 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
541 #  define TSGLRegisterDynamic(a,b,c,d)       TSGLRegister(a,b,c,0)
542 #else
543 #  define TSGLRegisterDynamic(a,b,c,d)       TSGLRegister(a,b,c,d)
544 #endif
545 
546 extern PetscErrorCode  TSGLRegister(const char[],const char[],const char[],PetscErrorCode(*)(TS));
547 extern PetscErrorCode  TSGLRegisterAll(const char[]);
548 extern PetscErrorCode  TSGLRegisterDestroy(void);
549 extern PetscErrorCode  TSGLInitializePackage(const char[]);
550 extern PetscErrorCode  TSGLFinalizePackage(void);
551 extern PetscErrorCode  TSGLSetType(TS,const TSGLType);
552 extern PetscErrorCode  TSGLGetAdapt(TS,TSGLAdapt*);
553 extern PetscErrorCode  TSGLSetAcceptType(TS,const TSGLAcceptType);
554 
555 #define TSARKIMEXType char*
556 #define TSARKIMEX2D     "2d"
557 #define TSARKIMEX2E     "2e"
558 #define TSARKIMEXPRSSP2 "prssp2"
559 #define TSARKIMEX3      "3"
560 #define TSARKIMEXBPR3   "bpr3"
561 #define TSARKIMEXARS443 "ars443"
562 #define TSARKIMEX4      "4"
563 #define TSARKIMEX5      "5"
564 extern PetscErrorCode TSARKIMEXGetType(TS ts,const TSARKIMEXType*);
565 extern PetscErrorCode TSARKIMEXSetType(TS ts,const TSARKIMEXType);
566 extern PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool);
567 extern PetscErrorCode TSARKIMEXRegister(const TSARKIMEXType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[],const PetscReal[]);
568 extern PetscErrorCode TSARKIMEXFinalizePackage(void);
569 extern PetscErrorCode TSARKIMEXInitializePackage(const char path[]);
570 extern PetscErrorCode TSARKIMEXRegisterDestroy(void);
571 extern PetscErrorCode TSARKIMEXRegisterAll(void);
572 
573 #define TSRosWType char*
574 #define TSROSW2M          "2m"
575 #define TSROSW2P          "2p"
576 #define TSROSWRA3PW       "ra3pw"
577 #define TSROSWRA34PW2     "ra34pw2"
578 #define TSROSWRODAS3      "rodas3"
579 #define TSROSWSANDU3      "sandu3"
580 #define TSROSWASSP3P3S1C  "assp3p3s1c"
581 #define TSROSWLASSP3P4S2C "lassp3p4s2c"
582 #define TSROSWLLSSP3P3S2C "llssp3p3s2c"
583 #define TSROSWARK3        "ark3"
584 
585 extern PetscErrorCode TSRosWGetType(TS ts,const TSRosWType*);
586 extern PetscErrorCode TSRosWSetType(TS ts,const TSRosWType);
587 extern PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool);
588 extern PetscErrorCode TSRosWRegister(const TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[]);
589 extern PetscErrorCode TSRosWFinalizePackage(void);
590 extern PetscErrorCode TSRosWInitializePackage(const char path[]);
591 extern PetscErrorCode TSRosWRegisterDestroy(void);
592 extern PetscErrorCode TSRosWRegisterAll(void);
593 
594 /*
595        PETSc interface to Sundials
596 */
597 #ifdef PETSC_HAVE_SUNDIALS
598 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
599 extern const char *TSSundialsLmmTypes[];
600 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
601 extern const char *TSSundialsGramSchmidtTypes[];
602 extern PetscErrorCode   TSSundialsSetType(TS,TSSundialsLmmType);
603 extern PetscErrorCode   TSSundialsGetPC(TS,PC*);
604 extern PetscErrorCode   TSSundialsSetTolerance(TS,PetscReal,PetscReal);
605 extern PetscErrorCode   TSSundialsSetMinTimeStep(TS,PetscReal);
606 extern PetscErrorCode   TSSundialsSetMaxTimeStep(TS,PetscReal);
607 extern PetscErrorCode   TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
608 extern PetscErrorCode   TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
609 extern PetscErrorCode   TSSundialsSetGMRESRestart(TS,PetscInt);
610 extern PetscErrorCode   TSSundialsSetLinearTolerance(TS,PetscReal);
611 extern PetscErrorCode   TSSundialsMonitorInternalSteps(TS,PetscBool );
612 extern PetscErrorCode   TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
613 extern PetscErrorCode   TSSundialsSetMaxl(TS,PetscInt);
614 #endif
615 
616 extern PetscErrorCode   TSRKSetTolerance(TS,PetscReal);
617 
618 extern PetscErrorCode  TSThetaSetTheta(TS,PetscReal);
619 extern PetscErrorCode  TSThetaGetTheta(TS,PetscReal*);
620 extern PetscErrorCode  TSThetaGetEndpoint(TS,PetscBool*);
621 extern PetscErrorCode  TSThetaSetEndpoint(TS,PetscBool);
622 
623 extern PetscErrorCode  TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*);
624 extern PetscErrorCode  TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*);
625 extern PetscErrorCode  TSAlphaSetRadius(TS,PetscReal);
626 extern PetscErrorCode  TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal);
627 extern PetscErrorCode  TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*);
628 
629 extern PetscErrorCode  TSSetDM(TS,DM);
630 extern PetscErrorCode  TSGetDM(TS,DM*);
631 
632 extern PetscErrorCode  SNESTSFormFunction(SNES,Vec,Vec,void*);
633 extern PetscErrorCode  SNESTSFormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
634 
635 PETSC_EXTERN_CXX_END
636 #endif
637