xref: /petsc/include/petscts.h (revision 0b99f51463372431c8db4bb246d48fd57aaec3e2)
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 TSADAPTCFL   "cfl"
273 
274 /*MC
275    TSAdaptRegisterDynamic - adds a TSAdapt implementation
276 
277    Synopsis:
278    PetscErrorCode TSAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
279 
280    Not Collective
281 
282    Input Parameters:
283 +  name_scheme - name of user-defined adaptivity scheme
284 .  path - path (either absolute or relative) the library containing this scheme
285 .  name_create - name of routine to create method context
286 -  routine_create - routine to create method context
287 
288    Notes:
289    TSAdaptRegisterDynamic() may be called multiple times to add several user-defined families.
290 
291    If dynamic libraries are used, then the fourth input argument (routine_create)
292    is ignored.
293 
294    Sample usage:
295 .vb
296    TSAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
297                             "MySchemeCreate",MySchemeCreate);
298 .ve
299 
300    Then, your scheme can be chosen with the procedural interface via
301 $     TSAdaptSetType(ts,"my_scheme")
302    or at runtime via the option
303 $     -ts_adapt_type my_scheme
304 
305    Level: advanced
306 
307    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
308           and others of the form ${any_environmental_variable} occuring in pathname will be
309           replaced with appropriate values.
310 
311 .keywords: TSAdapt, register
312 
313 .seealso: TSAdaptRegisterAll()
314 M*/
315 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
316 #  define TSAdaptRegisterDynamic(a,b,c,d)  TSAdaptRegister(a,b,c,0)
317 #else
318 #  define TSAdaptRegisterDynamic(a,b,c,d)  TSAdaptRegister(a,b,c,d)
319 #endif
320 
321 extern PetscErrorCode TSGetAdapt(TS,TSAdapt*);
322 extern PetscErrorCode TSAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSAdapt));
323 extern PetscErrorCode TSAdaptRegisterAll(const char[]);
324 extern PetscErrorCode TSAdaptRegisterDestroy(void);
325 extern PetscErrorCode TSAdaptInitializePackage(const char[]);
326 extern PetscErrorCode TSAdaptFinalizePackage(void);
327 extern PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*);
328 extern PetscErrorCode TSAdaptSetType(TSAdapt,const TSAdaptType);
329 extern PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]);
330 extern PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
331 extern PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool);
332 extern PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**);
333 extern PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*);
334 extern PetscErrorCode TSAdaptView(TSAdapt,PetscViewer);
335 extern PetscErrorCode TSAdaptSetFromOptions(TSAdapt);
336 extern PetscErrorCode TSAdaptDestroy(TSAdapt*);
337 extern PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool);
338 extern PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal);
339 
340 /*S
341    TSGLAdapt - Abstract object that manages time-step adaptivity
342 
343    Level: beginner
344 
345    Developer Notes:
346    This functionality should be replaced by the TSAdapt.
347 
348 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType
349 S*/
350 typedef struct _p_TSGLAdapt *TSGLAdapt;
351 
352 /*J
353     TSGLAdaptType - String with the name of TSGLAdapt scheme or the creation function
354        with an optional dynamic library name, for example
355        http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate()
356 
357    Level: beginner
358 
359 .seealso: TSGLAdaptSetType(), TS
360 J*/
361 #define TSGLAdaptType  char*
362 #define TSGLADAPT_NONE "none"
363 #define TSGLADAPT_SIZE "size"
364 #define TSGLADAPT_BOTH "both"
365 
366 /*MC
367    TSGLAdaptRegisterDynamic - adds a TSGLAdapt implementation
368 
369    Synopsis:
370    PetscErrorCode TSGLAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
371 
372    Not Collective
373 
374    Input Parameters:
375 +  name_scheme - name of user-defined adaptivity scheme
376 .  path - path (either absolute or relative) the library containing this scheme
377 .  name_create - name of routine to create method context
378 -  routine_create - routine to create method context
379 
380    Notes:
381    TSGLAdaptRegisterDynamic() may be called multiple times to add several user-defined families.
382 
383    If dynamic libraries are used, then the fourth input argument (routine_create)
384    is ignored.
385 
386    Sample usage:
387 .vb
388    TSGLAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
389                             "MySchemeCreate",MySchemeCreate);
390 .ve
391 
392    Then, your scheme can be chosen with the procedural interface via
393 $     TSGLAdaptSetType(ts,"my_scheme")
394    or at runtime via the option
395 $     -ts_adapt_type my_scheme
396 
397    Level: advanced
398 
399    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
400           and others of the form ${any_environmental_variable} occuring in pathname will be
401           replaced with appropriate values.
402 
403 .keywords: TSGLAdapt, register
404 
405 .seealso: TSGLAdaptRegisterAll()
406 M*/
407 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
408 #  define TSGLAdaptRegisterDynamic(a,b,c,d)  TSGLAdaptRegister(a,b,c,0)
409 #else
410 #  define TSGLAdaptRegisterDynamic(a,b,c,d)  TSGLAdaptRegister(a,b,c,d)
411 #endif
412 
413 extern PetscErrorCode  TSGLAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSGLAdapt));
414 extern PetscErrorCode  TSGLAdaptRegisterAll(const char[]);
415 extern PetscErrorCode  TSGLAdaptRegisterDestroy(void);
416 extern PetscErrorCode  TSGLAdaptInitializePackage(const char[]);
417 extern PetscErrorCode  TSGLAdaptFinalizePackage(void);
418 extern PetscErrorCode  TSGLAdaptCreate(MPI_Comm,TSGLAdapt*);
419 extern PetscErrorCode  TSGLAdaptSetType(TSGLAdapt,const TSGLAdaptType);
420 extern PetscErrorCode  TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]);
421 extern PetscErrorCode  TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
422 extern PetscErrorCode  TSGLAdaptView(TSGLAdapt,PetscViewer);
423 extern PetscErrorCode  TSGLAdaptSetFromOptions(TSGLAdapt);
424 extern PetscErrorCode  TSGLAdaptDestroy(TSGLAdapt*);
425 
426 /*J
427     TSGLAcceptType - String with the name of TSGLAccept scheme or the function
428        with an optional dynamic library name, for example
429        http://www.mcs.anl.gov/petsc/lib.a:mytsglaccept()
430 
431    Level: beginner
432 
433 .seealso: TSGLSetAcceptType(), TS
434 J*/
435 #define TSGLAcceptType  char*
436 #define TSGLACCEPT_ALWAYS "always"
437 
438 typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *);
439 extern PetscErrorCode  TSGLAcceptRegister(const char[],const char[],const char[],TSGLAcceptFunction);
440 
441 /*MC
442    TSGLAcceptRegisterDynamic - adds a TSGL acceptance scheme
443 
444    Synopsis:
445    PetscErrorCode TSGLAcceptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
446 
447    Not Collective
448 
449    Input Parameters:
450 +  name_scheme - name of user-defined acceptance scheme
451 .  path - path (either absolute or relative) the library containing this scheme
452 .  name_create - name of routine to create method context
453 -  routine_create - routine to create method context
454 
455    Notes:
456    TSGLAcceptRegisterDynamic() may be called multiple times to add several user-defined families.
457 
458    If dynamic libraries are used, then the fourth input argument (routine_create)
459    is ignored.
460 
461    Sample usage:
462 .vb
463    TSGLAcceptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
464                              "MySchemeCreate",MySchemeCreate);
465 .ve
466 
467    Then, your scheme can be chosen with the procedural interface via
468 $     TSGLSetAcceptType(ts,"my_scheme")
469    or at runtime via the option
470 $     -ts_gl_accept_type my_scheme
471 
472    Level: advanced
473 
474    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
475           and others of the form ${any_environmental_variable} occuring in pathname will be
476           replaced with appropriate values.
477 
478 .keywords: TSGL, TSGLAcceptType, register
479 
480 .seealso: TSGLRegisterAll()
481 M*/
482 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
483 #  define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,0)
484 #else
485 #  define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,d)
486 #endif
487 
488 /*J
489   TSGLType - family of time integration method within the General Linear class
490 
491   Level: beginner
492 
493 .seealso: TSGLSetType(), TSGLRegister()
494 J*/
495 #define TSGLType char*
496 #define TSGL_IRKS   "irks"
497 
498 /*MC
499    TSGLRegisterDynamic - adds a TSGL implementation
500 
501    Synopsis:
502    PetscErrorCode TSGLRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
503 
504    Not Collective
505 
506    Input Parameters:
507 +  name_scheme - name of user-defined general linear scheme
508 .  path - path (either absolute or relative) the library containing this scheme
509 .  name_create - name of routine to create method context
510 -  routine_create - routine to create method context
511 
512    Notes:
513    TSGLRegisterDynamic() may be called multiple times to add several user-defined families.
514 
515    If dynamic libraries are used, then the fourth input argument (routine_create)
516    is ignored.
517 
518    Sample usage:
519 .vb
520    TSGLRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
521                        "MySchemeCreate",MySchemeCreate);
522 .ve
523 
524    Then, your scheme can be chosen with the procedural interface via
525 $     TSGLSetType(ts,"my_scheme")
526    or at runtime via the option
527 $     -ts_gl_type my_scheme
528 
529    Level: advanced
530 
531    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
532           and others of the form ${any_environmental_variable} occuring in pathname will be
533           replaced with appropriate values.
534 
535 .keywords: TSGL, register
536 
537 .seealso: TSGLRegisterAll()
538 M*/
539 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
540 #  define TSGLRegisterDynamic(a,b,c,d)       TSGLRegister(a,b,c,0)
541 #else
542 #  define TSGLRegisterDynamic(a,b,c,d)       TSGLRegister(a,b,c,d)
543 #endif
544 
545 extern PetscErrorCode  TSGLRegister(const char[],const char[],const char[],PetscErrorCode(*)(TS));
546 extern PetscErrorCode  TSGLRegisterAll(const char[]);
547 extern PetscErrorCode  TSGLRegisterDestroy(void);
548 extern PetscErrorCode  TSGLInitializePackage(const char[]);
549 extern PetscErrorCode  TSGLFinalizePackage(void);
550 extern PetscErrorCode  TSGLSetType(TS,const TSGLType);
551 extern PetscErrorCode  TSGLGetAdapt(TS,TSGLAdapt*);
552 extern PetscErrorCode  TSGLSetAcceptType(TS,const TSGLAcceptType);
553 
554 #define TSARKIMEXType char*
555 #define TSARKIMEX2D     "2d"
556 #define TSARKIMEX2E     "2e"
557 #define TSARKIMEXPRSSP2 "prssp2"
558 #define TSARKIMEX3      "3"
559 #define TSARKIMEXBPR3   "bpr3"
560 #define TSARKIMEXARS443 "ars443"
561 #define TSARKIMEX4      "4"
562 #define TSARKIMEX5      "5"
563 extern PetscErrorCode TSARKIMEXGetType(TS ts,const TSARKIMEXType*);
564 extern PetscErrorCode TSARKIMEXSetType(TS ts,const TSARKIMEXType);
565 extern PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool);
566 extern PetscErrorCode TSARKIMEXRegister(const TSARKIMEXType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[],const PetscReal[]);
567 extern PetscErrorCode TSARKIMEXFinalizePackage(void);
568 extern PetscErrorCode TSARKIMEXInitializePackage(const char path[]);
569 extern PetscErrorCode TSARKIMEXRegisterDestroy(void);
570 extern PetscErrorCode TSARKIMEXRegisterAll(void);
571 
572 #define TSRosWType char*
573 #define TSROSW2M          "2m"
574 #define TSROSW2P          "2p"
575 #define TSROSWRA3PW       "ra3pw"
576 #define TSROSWRA34PW2     "ra34pw2"
577 #define TSROSWRODAS3      "rodas3"
578 #define TSROSWSANDU3      "sandu3"
579 #define TSROSWASSP3P3S1C  "assp3p3s1c"
580 #define TSROSWLASSP3P4S2C "lassp3p4s2c"
581 #define TSROSWLLSSP3P3S2C "llssp3p3s2c"
582 #define TSROSWARK3        "ark3"
583 
584 extern PetscErrorCode TSRosWGetType(TS ts,const TSRosWType*);
585 extern PetscErrorCode TSRosWSetType(TS ts,const TSRosWType);
586 extern PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool);
587 extern PetscErrorCode TSRosWRegister(const TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[]);
588 extern PetscErrorCode TSRosWFinalizePackage(void);
589 extern PetscErrorCode TSRosWInitializePackage(const char path[]);
590 extern PetscErrorCode TSRosWRegisterDestroy(void);
591 extern PetscErrorCode TSRosWRegisterAll(void);
592 
593 /*
594        PETSc interface to Sundials
595 */
596 #ifdef PETSC_HAVE_SUNDIALS
597 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
598 extern const char *TSSundialsLmmTypes[];
599 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
600 extern const char *TSSundialsGramSchmidtTypes[];
601 extern PetscErrorCode   TSSundialsSetType(TS,TSSundialsLmmType);
602 extern PetscErrorCode   TSSundialsGetPC(TS,PC*);
603 extern PetscErrorCode   TSSundialsSetTolerance(TS,PetscReal,PetscReal);
604 extern PetscErrorCode   TSSundialsSetMinTimeStep(TS,PetscReal);
605 extern PetscErrorCode   TSSundialsSetMaxTimeStep(TS,PetscReal);
606 extern PetscErrorCode   TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
607 extern PetscErrorCode   TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
608 extern PetscErrorCode   TSSundialsSetGMRESRestart(TS,PetscInt);
609 extern PetscErrorCode   TSSundialsSetLinearTolerance(TS,PetscReal);
610 extern PetscErrorCode   TSSundialsMonitorInternalSteps(TS,PetscBool );
611 extern PetscErrorCode   TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
612 extern PetscErrorCode   TSSundialsSetMaxl(TS,PetscInt);
613 #endif
614 
615 extern PetscErrorCode   TSRKSetTolerance(TS,PetscReal);
616 
617 extern PetscErrorCode  TSThetaSetTheta(TS,PetscReal);
618 extern PetscErrorCode  TSThetaGetTheta(TS,PetscReal*);
619 extern PetscErrorCode  TSThetaGetEndpoint(TS,PetscBool*);
620 extern PetscErrorCode  TSThetaSetEndpoint(TS,PetscBool);
621 
622 extern PetscErrorCode  TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*);
623 extern PetscErrorCode  TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*);
624 extern PetscErrorCode  TSAlphaSetRadius(TS,PetscReal);
625 extern PetscErrorCode  TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal);
626 extern PetscErrorCode  TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*);
627 
628 extern PetscErrorCode  TSSetDM(TS,DM);
629 extern PetscErrorCode  TSGetDM(TS,DM*);
630 
631 extern PetscErrorCode  SNESTSFormFunction(SNES,Vec,Vec,void*);
632 extern PetscErrorCode  SNESTSFormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
633 
634 PETSC_EXTERN_CXX_END
635 #endif
636