xref: /petsc/include/petscts.h (revision 0e9bae810fdaeb60e2713eaa8ddb89f42e079fd1)
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 
9 /*S
10      TS - Abstract PETSc object that manages all time-steppers (ODE integrators)
11 
12    Level: beginner
13 
14   Concepts: ODE solvers
15 
16 .seealso:  TSCreate(), TSSetType(), TSType, SNES, KSP, PC
17 S*/
18 typedef struct _p_TS* TS;
19 
20 /*J
21     TSType - String with the name of a PETSc TS method or the creation function
22        with an optional dynamic library name, for example
23        http://www.mcs.anl.gov/petsc/lib.a:mytscreate()
24 
25    Level: beginner
26 
27 .seealso: TSSetType(), TS
28 J*/
29 #define TSType char*
30 #define TSEULER           "euler"
31 #define TSBEULER          "beuler"
32 #define TSPSEUDO          "pseudo"
33 #define TSCN              "cn"
34 #define TSSUNDIALS        "sundials"
35 #define TSRK              "rk"
36 #define TSPYTHON          "python"
37 #define TSTHETA           "theta"
38 #define TSALPHA           "alpha"
39 #define TSGL              "gl"
40 #define TSSSP             "ssp"
41 #define TSARKIMEX         "arkimex"
42 #define TSROSW            "rosw"
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 /*E
54    TSConvergedReason - reason a TS method has converged or not
55 
56    Level: beginner
57 
58    Developer Notes: this must match finclude/petscts.h
59 
60    Each reason has its own manual page.
61 
62 .seealso: TSGetConvergedReason()
63 E*/
64 typedef enum {
65   TS_CONVERGED_ITERATING      = 0,
66   TS_CONVERGED_TIME           = 1,
67   TS_CONVERGED_ITS            = 2,
68   TS_DIVERGED_NONLINEAR_SOLVE = -1,
69   TS_DIVERGED_STEP_REJECTED   = -2
70 } TSConvergedReason;
71 PETSC_EXTERN const char *const*TSConvergedReasons;
72 
73 /*MC
74    TS_CONVERGED_ITERATING - this only occurs if TSGetConvergedReason() is called during the TSSolve()
75 
76    Level: beginner
77 
78 .seealso: TSSolve(), TSConvergedReason(), TSGetAdapt()
79 M*/
80 
81 /*MC
82    TS_CONVERGED_TIME - the final time was reached
83 
84    Level: beginner
85 
86 .seealso: TSSolve(), TSConvergedReason(), TSGetAdapt(), TSSetDuration()
87 M*/
88 
89 /*MC
90    TS_CONVERGED_ITS - the maximum number of iterations was reached prior to the final time
91 
92    Level: beginner
93 
94 .seealso: TSSolve(), TSConvergedReason(), TSGetAdapt(), TSSetDuration()
95 M*/
96 
97 /*MC
98    TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed
99 
100    Level: beginner
101 
102 .seealso: TSSolve(), TSConvergedReason(), TSGetAdapt(), TSGetSNES(), SNESGetConvergedReason()
103 M*/
104 
105 /*MC
106    TS_DIVERGED_STEP_REJECTED - too many steps were rejected
107 
108    Level: beginner
109 
110 .seealso: TSSolve(), TSConvergedReason(), TSGetAdapt()
111 M*/
112 
113 /* Logging support */
114 PETSC_EXTERN PetscClassId TS_CLASSID;
115 
116 PETSC_EXTERN PetscErrorCode TSInitializePackage(const char[]);
117 
118 PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm,TS*);
119 PETSC_EXTERN PetscErrorCode TSDestroy(TS*);
120 
121 PETSC_EXTERN PetscErrorCode TSSetProblemType(TS,TSProblemType);
122 PETSC_EXTERN PetscErrorCode TSGetProblemType(TS,TSProblemType*);
123 PETSC_EXTERN PetscErrorCode TSMonitor(TS,PetscInt,PetscReal,Vec);
124 PETSC_EXTERN PetscErrorCode TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**));
125 PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS);
126 
127 PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS,const char[]);
128 PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS,const char[]);
129 PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS,const char *[]);
130 PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS);
131 PETSC_EXTERN PetscErrorCode TSSetUp(TS);
132 PETSC_EXTERN PetscErrorCode TSReset(TS);
133 
134 PETSC_EXTERN PetscErrorCode TSSetSolution(TS,Vec);
135 PETSC_EXTERN PetscErrorCode TSGetSolution(TS,Vec*);
136 
137 PETSC_EXTERN PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal);
138 PETSC_EXTERN PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*);
139 PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS,PetscBool);
140 
141 PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,void*);
142 PETSC_EXTERN PetscErrorCode TSMonitorSolution(TS,PetscInt,PetscReal,Vec,void*);
143 PETSC_EXTERN PetscErrorCode TSMonitorSolutionCreate(TS,PetscViewer,void**);
144 PETSC_EXTERN PetscErrorCode TSMonitorSolutionDestroy(void**);
145 PETSC_EXTERN PetscErrorCode TSMonitorSolutionBinary(TS,PetscInt,PetscReal,Vec,void*);
146 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*);
147 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void*);
148 
149 PETSC_EXTERN PetscErrorCode TSStep(TS);
150 PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS,PetscInt,Vec,PetscBool*);
151 PETSC_EXTERN PetscErrorCode TSSolve(TS,Vec,PetscReal*);
152 PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*);
153 PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS,PetscInt*);
154 PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS,PetscInt*);
155 PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS,PetscInt*);
156 PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS,PetscInt);
157 PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS,PetscInt*);
158 PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS,PetscInt);
159 PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS,PetscBool);
160 
161 PETSC_EXTERN PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal);
162 PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS,PetscReal*);
163 PETSC_EXTERN PetscErrorCode TSGetTime(TS,PetscReal*);
164 PETSC_EXTERN PetscErrorCode TSSetTime(TS,PetscReal);
165 PETSC_EXTERN PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*);
166 PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS,PetscReal);
167 
168 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*);
169 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
170 PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*);
171 PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**);
172 PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*);
173 PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**);
174 
175 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*);
176 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*);
177 PETSC_EXTERN PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*);
178 PETSC_EXTERN PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**);
179 PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*);
180 PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**);
181 
182 PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*);
183 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
184 PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*);
185 PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*);
186 
187 PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
188 PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS,PetscReal));
189 PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
190 PETSC_EXTERN PetscErrorCode TSPreStep(TS);
191 PETSC_EXTERN PetscErrorCode TSPreStage(TS,PetscReal);
192 PETSC_EXTERN PetscErrorCode TSPostStep(TS);
193 PETSC_EXTERN PetscErrorCode TSSetRetainStages(TS,PetscBool);
194 PETSC_EXTERN PetscErrorCode TSInterpolate(TS,PetscReal,Vec);
195 PETSC_EXTERN PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec);
196 PETSC_EXTERN PetscErrorCode TSGetTolerances(TS,PetscReal*,Vec*,PetscReal*,Vec*);
197 PETSC_EXTERN PetscErrorCode TSErrorNormWRMS(TS,Vec,PetscReal*);
198 PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal);
199 PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS,PetscReal*);
200 
201 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*);
202 PETSC_EXTERN PetscErrorCode TSPseudoDefaultTimeStep(TS,PetscReal*,void*);
203 PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *);
204 PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal);
205 
206 PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*);
207 PETSC_EXTERN PetscErrorCode TSPseudoDefaultVerifyTimeStep(TS,Vec,void*,PetscReal*,PetscBool *);
208 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *);
209 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal);
210 PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);
211 
212 PETSC_EXTERN PetscErrorCode TSPythonSetType(TS,const char[]);
213 
214 PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec);
215 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*);
216 PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool);
217 PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,PetscBool);
218 
219 PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS,Vec,Vec);
220 
221 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM,TSRHSFunction,void*);
222 PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM,TSRHSFunction*,void**);
223 PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM,TSRHSJacobian,void*);
224 PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM,TSRHSJacobian*,void**);
225 PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM,TSIFunction,void*);
226 PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM,TSIFunction*,void**);
227 PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM,TSIJacobian,void*);
228 PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM,TSIJacobian*,void**);
229 
230 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*);
231 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,MatStructure*,void*);
232 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*);
233 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,MatStructure*,void*);
234 
235 PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*),void *);
236 PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,MatStructure*,void*),void *);
237 PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*),void *);
238 PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,MatStructure*,void*),void *);
239 
240 /* Dynamic creation and loading functions */
241 PETSC_EXTERN PetscFList TSList;
242 PETSC_EXTERN PetscBool TSRegisterAllCalled;
243 PETSC_EXTERN PetscErrorCode TSGetType(TS,const TSType*);
244 PETSC_EXTERN PetscErrorCode TSSetType(TS,const TSType);
245 PETSC_EXTERN PetscErrorCode TSRegister(const char[], const char[], const char[], PetscErrorCode (*)(TS));
246 PETSC_EXTERN PetscErrorCode TSRegisterAll(const char[]);
247 PETSC_EXTERN PetscErrorCode TSRegisterDestroy(void);
248 
249 /*MC
250   TSRegisterDynamic - Adds a creation method to the TS package.
251 
252   Synopsis:
253   PetscErrorCode TSRegisterDynamic(const char *name, const char *path, const char *func_name, PetscErrorCode (*create_func)(TS))
254 
255   Not Collective
256 
257   Input Parameters:
258 + name        - The name of a new user-defined creation routine
259 . path        - The path (either absolute or relative) of the library containing this routine
260 . func_name   - The name of the creation routine
261 - create_func - The creation routine itself
262 
263   Notes:
264   TSRegisterDynamic() may be called multiple times to add several user-defined tses.
265 
266   If dynamic libraries are used, then the fourth input argument (create_func) is ignored.
267 
268   Sample usage:
269 .vb
270   TSRegisterDynamic("my_ts", "/home/username/my_lib/lib/libO/solaris/libmy.a", "MyTSCreate", MyTSCreate);
271 .ve
272 
273   Then, your ts type can be chosen with the procedural interface via
274 .vb
275     TS ts;
276     TSCreate(MPI_Comm, &ts);
277     TSSetType(ts, "my_ts")
278 .ve
279   or at runtime via the option
280 .vb
281     -ts_type my_ts
282 .ve
283 
284   Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
285         If your function is not being put into a shared library then use TSRegister() instead
286 
287   Level: advanced
288 
289 .keywords: TS, register
290 .seealso: TSRegisterAll(), TSRegisterDestroy()
291 M*/
292 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
293 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,0)
294 #else
295 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,d)
296 #endif
297 
298 PETSC_EXTERN PetscErrorCode TSGetSNES(TS,SNES*);
299 PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*);
300 
301 PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer);
302 
303 PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *);
304 PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *);
305 
306 PETSC_EXTERN PetscErrorCode TSMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG *);
307 PETSC_EXTERN PetscErrorCode TSMonitorLG(TS,PetscInt,PetscReal,Vec,void *);
308 PETSC_EXTERN PetscErrorCode TSMonitorLGDestroy(PetscDrawLG*);
309 
310 /*J
311    TSSSPType - string with the name of TSSSP scheme.
312 
313    Level: beginner
314 
315 .seealso: TSSSPSetType(), TS
316 J*/
317 #define TSSSPType char*
318 #define TSSSPRKS2  "rks2"
319 #define TSSSPRKS3  "rks3"
320 #define TSSSPRK104 "rk104"
321 
322 PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,const TSSSPType);
323 PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,const TSSSPType*);
324 PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt);
325 PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*);
326 
327 /*S
328    TSAdapt - Abstract object that manages time-step adaptivity
329 
330    Level: beginner
331 
332 .seealso: TS, TSAdaptCreate(), TSAdaptType
333 S*/
334 typedef struct _p_TSAdapt *TSAdapt;
335 
336 /*E
337     TSAdaptType - String with the name of TSAdapt scheme or the creation function
338        with an optional dynamic library name, for example
339        http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate()
340 
341    Level: beginner
342 
343 .seealso: TSAdaptSetType(), TS
344 E*/
345 #define TSAdaptType  char*
346 #define TSADAPTBASIC "basic"
347 #define TSADAPTNONE  "none"
348 #define TSADAPTCFL   "cfl"
349 
350 /*MC
351    TSAdaptRegisterDynamic - adds a TSAdapt implementation
352 
353    Synopsis:
354    PetscErrorCode TSAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
355 
356    Not Collective
357 
358    Input Parameters:
359 +  name_scheme - name of user-defined adaptivity scheme
360 .  path - path (either absolute or relative) the library containing this scheme
361 .  name_create - name of routine to create method context
362 -  routine_create - routine to create method context
363 
364    Notes:
365    TSAdaptRegisterDynamic() may be called multiple times to add several user-defined families.
366 
367    If dynamic libraries are used, then the fourth input argument (routine_create)
368    is ignored.
369 
370    Sample usage:
371 .vb
372    TSAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
373                             "MySchemeCreate",MySchemeCreate);
374 .ve
375 
376    Then, your scheme can be chosen with the procedural interface via
377 $     TSAdaptSetType(ts,"my_scheme")
378    or at runtime via the option
379 $     -ts_adapt_type my_scheme
380 
381    Level: advanced
382 
383    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
384           and others of the form ${any_environmental_variable} occuring in pathname will be
385           replaced with appropriate values.
386 
387 .keywords: TSAdapt, register
388 
389 .seealso: TSAdaptRegisterAll()
390 M*/
391 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
392 #  define TSAdaptRegisterDynamic(a,b,c,d)  TSAdaptRegister(a,b,c,0)
393 #else
394 #  define TSAdaptRegisterDynamic(a,b,c,d)  TSAdaptRegister(a,b,c,d)
395 #endif
396 
397 PETSC_EXTERN PetscErrorCode TSGetAdapt(TS,TSAdapt*);
398 PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSAdapt));
399 PETSC_EXTERN PetscErrorCode TSAdaptRegisterAll(const char[]);
400 PETSC_EXTERN PetscErrorCode TSAdaptRegisterDestroy(void);
401 PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(const char[]);
402 PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
403 PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*);
404 PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,const TSAdaptType);
405 PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]);
406 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
407 PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool);
408 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**);
409 PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*);
410 PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscBool*);
411 PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer);
412 PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(TSAdapt);
413 PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*);
414 PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool);
415 PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal);
416 PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt,PetscErrorCode(*)(TSAdapt,TS,PetscBool*));
417 
418 /*S
419    TSGLAdapt - Abstract object that manages time-step adaptivity
420 
421    Level: beginner
422 
423    Developer Notes:
424    This functionality should be replaced by the TSAdapt.
425 
426 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType
427 S*/
428 typedef struct _p_TSGLAdapt *TSGLAdapt;
429 
430 /*J
431     TSGLAdaptType - String with the name of TSGLAdapt scheme or the creation function
432        with an optional dynamic library name, for example
433        http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate()
434 
435    Level: beginner
436 
437 .seealso: TSGLAdaptSetType(), TS
438 J*/
439 #define TSGLAdaptType  char*
440 #define TSGLADAPT_NONE "none"
441 #define TSGLADAPT_SIZE "size"
442 #define TSGLADAPT_BOTH "both"
443 
444 /*MC
445    TSGLAdaptRegisterDynamic - adds a TSGLAdapt implementation
446 
447    Synopsis:
448    PetscErrorCode TSGLAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
449 
450    Not Collective
451 
452    Input Parameters:
453 +  name_scheme - name of user-defined adaptivity scheme
454 .  path - path (either absolute or relative) the library containing this scheme
455 .  name_create - name of routine to create method context
456 -  routine_create - routine to create method context
457 
458    Notes:
459    TSGLAdaptRegisterDynamic() may be called multiple times to add several user-defined families.
460 
461    If dynamic libraries are used, then the fourth input argument (routine_create)
462    is ignored.
463 
464    Sample usage:
465 .vb
466    TSGLAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
467                             "MySchemeCreate",MySchemeCreate);
468 .ve
469 
470    Then, your scheme can be chosen with the procedural interface via
471 $     TSGLAdaptSetType(ts,"my_scheme")
472    or at runtime via the option
473 $     -ts_adapt_type my_scheme
474 
475    Level: advanced
476 
477    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
478           and others of the form ${any_environmental_variable} occuring in pathname will be
479           replaced with appropriate values.
480 
481 .keywords: TSGLAdapt, register
482 
483 .seealso: TSGLAdaptRegisterAll()
484 M*/
485 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
486 #  define TSGLAdaptRegisterDynamic(a,b,c,d)  TSGLAdaptRegister(a,b,c,0)
487 #else
488 #  define TSGLAdaptRegisterDynamic(a,b,c,d)  TSGLAdaptRegister(a,b,c,d)
489 #endif
490 
491 PETSC_EXTERN PetscErrorCode TSGLAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSGLAdapt));
492 PETSC_EXTERN PetscErrorCode TSGLAdaptRegisterAll(const char[]);
493 PETSC_EXTERN PetscErrorCode TSGLAdaptRegisterDestroy(void);
494 PETSC_EXTERN PetscErrorCode TSGLAdaptInitializePackage(const char[]);
495 PETSC_EXTERN PetscErrorCode TSGLAdaptFinalizePackage(void);
496 PETSC_EXTERN PetscErrorCode TSGLAdaptCreate(MPI_Comm,TSGLAdapt*);
497 PETSC_EXTERN PetscErrorCode TSGLAdaptSetType(TSGLAdapt,const TSGLAdaptType);
498 PETSC_EXTERN PetscErrorCode TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]);
499 PETSC_EXTERN PetscErrorCode TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
500 PETSC_EXTERN PetscErrorCode TSGLAdaptView(TSGLAdapt,PetscViewer);
501 PETSC_EXTERN PetscErrorCode TSGLAdaptSetFromOptions(TSGLAdapt);
502 PETSC_EXTERN PetscErrorCode TSGLAdaptDestroy(TSGLAdapt*);
503 
504 /*J
505     TSGLAcceptType - String with the name of TSGLAccept scheme or the function
506        with an optional dynamic library name, for example
507        http://www.mcs.anl.gov/petsc/lib.a:mytsglaccept()
508 
509    Level: beginner
510 
511 .seealso: TSGLSetAcceptType(), TS
512 J*/
513 #define TSGLAcceptType  char*
514 #define TSGLACCEPT_ALWAYS "always"
515 
516 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *);
517 PETSC_EXTERN PetscErrorCode TSGLAcceptRegister(const char[],const char[],const char[],TSGLAcceptFunction);
518 
519 /*MC
520    TSGLAcceptRegisterDynamic - adds a TSGL acceptance scheme
521 
522    Synopsis:
523    PetscErrorCode TSGLAcceptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
524 
525    Not Collective
526 
527    Input Parameters:
528 +  name_scheme - name of user-defined acceptance scheme
529 .  path - path (either absolute or relative) the library containing this scheme
530 .  name_create - name of routine to create method context
531 -  routine_create - routine to create method context
532 
533    Notes:
534    TSGLAcceptRegisterDynamic() may be called multiple times to add several user-defined families.
535 
536    If dynamic libraries are used, then the fourth input argument (routine_create)
537    is ignored.
538 
539    Sample usage:
540 .vb
541    TSGLAcceptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
542                              "MySchemeCreate",MySchemeCreate);
543 .ve
544 
545    Then, your scheme can be chosen with the procedural interface via
546 $     TSGLSetAcceptType(ts,"my_scheme")
547    or at runtime via the option
548 $     -ts_gl_accept_type my_scheme
549 
550    Level: advanced
551 
552    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
553           and others of the form ${any_environmental_variable} occuring in pathname will be
554           replaced with appropriate values.
555 
556 .keywords: TSGL, TSGLAcceptType, register
557 
558 .seealso: TSGLRegisterAll()
559 M*/
560 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
561 #  define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,0)
562 #else
563 #  define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,d)
564 #endif
565 
566 /*J
567   TSGLType - family of time integration method within the General Linear class
568 
569   Level: beginner
570 
571 .seealso: TSGLSetType(), TSGLRegister()
572 J*/
573 #define TSGLType char*
574 #define TSGL_IRKS   "irks"
575 
576 /*MC
577    TSGLRegisterDynamic - adds a TSGL implementation
578 
579    Synopsis:
580    PetscErrorCode TSGLRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
581 
582    Not Collective
583 
584    Input Parameters:
585 +  name_scheme - name of user-defined general linear scheme
586 .  path - path (either absolute or relative) the library containing this scheme
587 .  name_create - name of routine to create method context
588 -  routine_create - routine to create method context
589 
590    Notes:
591    TSGLRegisterDynamic() may be called multiple times to add several user-defined families.
592 
593    If dynamic libraries are used, then the fourth input argument (routine_create)
594    is ignored.
595 
596    Sample usage:
597 .vb
598    TSGLRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
599                        "MySchemeCreate",MySchemeCreate);
600 .ve
601 
602    Then, your scheme can be chosen with the procedural interface via
603 $     TSGLSetType(ts,"my_scheme")
604    or at runtime via the option
605 $     -ts_gl_type my_scheme
606 
607    Level: advanced
608 
609    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
610           and others of the form ${any_environmental_variable} occuring in pathname will be
611           replaced with appropriate values.
612 
613 .keywords: TSGL, register
614 
615 .seealso: TSGLRegisterAll()
616 M*/
617 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
618 #  define TSGLRegisterDynamic(a,b,c,d)       TSGLRegister(a,b,c,0)
619 #else
620 #  define TSGLRegisterDynamic(a,b,c,d)       TSGLRegister(a,b,c,d)
621 #endif
622 
623 PETSC_EXTERN PetscErrorCode TSGLRegister(const char[],const char[],const char[],PetscErrorCode(*)(TS));
624 PETSC_EXTERN PetscErrorCode TSGLRegisterAll(const char[]);
625 PETSC_EXTERN PetscErrorCode TSGLRegisterDestroy(void);
626 PETSC_EXTERN PetscErrorCode TSGLInitializePackage(const char[]);
627 PETSC_EXTERN PetscErrorCode TSGLFinalizePackage(void);
628 PETSC_EXTERN PetscErrorCode TSGLSetType(TS,const TSGLType);
629 PETSC_EXTERN PetscErrorCode TSGLGetAdapt(TS,TSGLAdapt*);
630 PETSC_EXTERN PetscErrorCode TSGLSetAcceptType(TS,const TSGLAcceptType);
631 
632 /*J
633     TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX method.
634 
635    Level: beginner
636 
637 .seealso: TSARKIMEXSetType(), TS, TSARKIMEX, TSARKIMEXRegister()
638 J*/
639 #define TSARKIMEXType char*
640 #define TSARKIMEXA2     "a2"
641 #define TSARKIMEXL2     "l2"
642 #define TSARKIMEXARS122 "ars122"
643 #define TSARKIMEX2C     "2c"
644 #define TSARKIMEX2D     "2d"
645 #define TSARKIMEX2E     "2e"
646 #define TSARKIMEXPRSSP2 "prssp2"
647 #define TSARKIMEX3      "3"
648 #define TSARKIMEXBPR3   "bpr3"
649 #define TSARKIMEXARS443 "ars443"
650 #define TSARKIMEX4      "4"
651 #define TSARKIMEX5      "5"
652 PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,const TSARKIMEXType*);
653 PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,const TSARKIMEXType);
654 PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool);
655 PETSC_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[]);
656 PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
657 PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(const char path[]);
658 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);
659 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterAll(void);
660 
661 /*J
662     TSRosWType - String with the name of a Rosenbrock-W method.
663 
664    Level: beginner
665 
666 .seealso: TSRosWSetType(), TS, TSROSW, TSRosWRegister()
667 J*/
668 #define TSRosWType char*
669 #define TSROSW2M          "2m"
670 #define TSROSW2P          "2p"
671 #define TSROSWRA3PW       "ra3pw"
672 #define TSROSWRA34PW2     "ra34pw2"
673 #define TSROSWRODAS3      "rodas3"
674 #define TSROSWSANDU3      "sandu3"
675 #define TSROSWASSP3P3S1C  "assp3p3s1c"
676 #define TSROSWLASSP3P4S2C "lassp3p4s2c"
677 #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
678 #define TSROSWARK3        "ark3"
679 #define TSROSWTHETA1      "theta1"
680 #define TSROSWTHETA2      "theta2"
681 #define TSROSWGRK4T       "grk4t"
682 #define TSROSWSHAMP4      "shamp4"
683 #define TSROSWVELDD4      "veldd4"
684 #define TSROSW4L          "4l"
685 
686 PETSC_EXTERN PetscErrorCode TSRosWGetType(TS ts,const TSRosWType*);
687 PETSC_EXTERN PetscErrorCode TSRosWSetType(TS ts,const TSRosWType);
688 PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool);
689 PETSC_EXTERN PetscErrorCode TSRosWRegister(const TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
690 PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(const TSRosWType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
691 PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
692 PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(const char path[]);
693 PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);
694 PETSC_EXTERN PetscErrorCode TSRosWRegisterAll(void);
695 
696 /*
697        PETSc interface to Sundials
698 */
699 #ifdef PETSC_HAVE_SUNDIALS
700 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
701 PETSC_EXTERN const char *const TSSundialsLmmTypes[];
702 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
703 PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
704 PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType);
705 PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*);
706 PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal);
707 PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal);
708 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal);
709 PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
710 PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
711 PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt);
712 PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal);
713 PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool );
714 PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
715 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt);
716 #endif
717 
718 PETSC_EXTERN PetscErrorCode TSRKSetTolerance(TS,PetscReal);
719 
720 PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal);
721 PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*);
722 PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*);
723 PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool);
724 
725 PETSC_EXTERN PetscErrorCode TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*);
726 PETSC_EXTERN PetscErrorCode TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*);
727 PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal);
728 PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal);
729 PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*);
730 
731 PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM);
732 PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*);
733 
734 PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*);
735 PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
736 
737 #endif
738