xref: /petsc/include/petscts.h (revision 43f9ffa7fee3e506d2e7ed45ea09f47c334a6d12)
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 typedef const char* TSType;
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 (*TSSolutionFunction)(TS,PetscReal,Vec,void*);
176 PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS,TSSolutionFunction,void*);
177 
178 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*);
179 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*);
180 PETSC_EXTERN PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*);
181 PETSC_EXTERN PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**);
182 PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*);
183 PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**);
184 
185 PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*);
186 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
187 PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*);
188 PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*);
189 PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS,PetscReal,Vec);
190 
191 PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
192 PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS,PetscReal));
193 PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
194 PETSC_EXTERN PetscErrorCode TSPreStep(TS);
195 PETSC_EXTERN PetscErrorCode TSPreStage(TS,PetscReal);
196 PETSC_EXTERN PetscErrorCode TSPostStep(TS);
197 PETSC_EXTERN PetscErrorCode TSSetRetainStages(TS,PetscBool);
198 PETSC_EXTERN PetscErrorCode TSInterpolate(TS,PetscReal,Vec);
199 PETSC_EXTERN PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec);
200 PETSC_EXTERN PetscErrorCode TSGetTolerances(TS,PetscReal*,Vec*,PetscReal*,Vec*);
201 PETSC_EXTERN PetscErrorCode TSErrorNormWRMS(TS,Vec,PetscReal*);
202 PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal);
203 PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS,PetscReal*);
204 
205 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*);
206 PETSC_EXTERN PetscErrorCode TSPseudoDefaultTimeStep(TS,PetscReal*,void*);
207 PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *);
208 PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal);
209 
210 PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*);
211 PETSC_EXTERN PetscErrorCode TSPseudoDefaultVerifyTimeStep(TS,Vec,void*,PetscReal*,PetscBool *);
212 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *);
213 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal);
214 PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);
215 
216 PETSC_EXTERN PetscErrorCode TSPythonSetType(TS,const char[]);
217 
218 PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec);
219 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*);
220 PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool);
221 PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,PetscBool);
222 
223 PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS,Vec,Vec);
224 
225 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM,TSRHSFunction,void*);
226 PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM,TSRHSFunction*,void**);
227 PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM,TSRHSJacobian,void*);
228 PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM,TSRHSJacobian*,void**);
229 PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM,TSIFunction,void*);
230 PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM,TSIFunction*,void**);
231 PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM,TSIJacobian,void*);
232 PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM,TSIJacobian*,void**);
233 PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM,TSSolutionFunction,void*);
234 PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM,TSSolutionFunction*,void**);
235 
236 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*);
237 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,MatStructure*,void*);
238 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*);
239 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,MatStructure*,void*);
240 
241 PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*),void *);
242 PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,MatStructure*,void*),void *);
243 PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*),void *);
244 PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,MatStructure*,void*),void *);
245 
246 /* Dynamic creation and loading functions */
247 PETSC_EXTERN PetscFList TSList;
248 PETSC_EXTERN PetscBool TSRegisterAllCalled;
249 PETSC_EXTERN PetscErrorCode TSGetType(TS,TSType*);
250 PETSC_EXTERN PetscErrorCode TSSetType(TS,TSType);
251 PETSC_EXTERN PetscErrorCode TSRegister(const char[], const char[], const char[], PetscErrorCode (*)(TS));
252 PETSC_EXTERN PetscErrorCode TSRegisterAll(const char[]);
253 PETSC_EXTERN PetscErrorCode TSRegisterDestroy(void);
254 
255 /*MC
256   TSRegisterDynamic - Adds a creation method to the TS package.
257 
258   Synopsis:
259   PetscErrorCode TSRegisterDynamic(const char *name, const char *path, const char *func_name, PetscErrorCode (*create_func)(TS))
260 
261   Not Collective
262 
263   Input Parameters:
264 + name        - The name of a new user-defined creation routine
265 . path        - The path (either absolute or relative) of the library containing this routine
266 . func_name   - The name of the creation routine
267 - create_func - The creation routine itself
268 
269   Notes:
270   TSRegisterDynamic() may be called multiple times to add several user-defined tses.
271 
272   If dynamic libraries are used, then the fourth input argument (create_func) is ignored.
273 
274   Sample usage:
275 .vb
276   TSRegisterDynamic("my_ts", "/home/username/my_lib/lib/libO/solaris/libmy.a", "MyTSCreate", MyTSCreate);
277 .ve
278 
279   Then, your ts type can be chosen with the procedural interface via
280 .vb
281     TS ts;
282     TSCreate(MPI_Comm, &ts);
283     TSSetType(ts, "my_ts")
284 .ve
285   or at runtime via the option
286 .vb
287     -ts_type my_ts
288 .ve
289 
290   Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
291         If your function is not being put into a shared library then use TSRegister() instead
292 
293   Level: advanced
294 
295 .keywords: TS, register
296 .seealso: TSRegisterAll(), TSRegisterDestroy()
297 M*/
298 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
299 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,0)
300 #else
301 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,d)
302 #endif
303 
304 PETSC_EXTERN PetscErrorCode TSGetSNES(TS,SNES*);
305 PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*);
306 
307 PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer);
308 
309 PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *);
310 PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *);
311 
312 PETSC_EXTERN PetscErrorCode TSMonitorLGCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG *);
313 PETSC_EXTERN PetscErrorCode TSMonitorLG(TS,PetscInt,PetscReal,Vec,void *);
314 PETSC_EXTERN PetscErrorCode TSMonitorLGDestroy(PetscDrawLG*);
315 
316 PETSC_EXTERN PetscErrorCode TSMonitorSolutionODECreate(MPI_Comm,PetscInt,const char[],const char[],int,int,int,int,PetscDrawLG *);
317 PETSC_EXTERN PetscErrorCode TSMonitorSolutionODE(TS,PetscInt,PetscReal,Vec,void *);
318 PETSC_EXTERN PetscErrorCode TSMonitorSolutionODEDestroy(PetscDrawLG*);
319 
320 PETSC_EXTERN PetscErrorCode TSMonitorErrorODECreate(MPI_Comm,PetscInt,const char[],const char[],int,int,int,int,PetscDrawLG *);
321 PETSC_EXTERN PetscErrorCode TSMonitorErrorODE(TS,PetscInt,PetscReal,Vec,void *);
322 PETSC_EXTERN PetscErrorCode TSMonitorErrorODEDestroy(PetscDrawLG*);
323 
324 /*J
325    TSSSPType - string with the name of TSSSP scheme.
326 
327    Level: beginner
328 
329 .seealso: TSSSPSetType(), TS
330 J*/
331 typedef const char* TSSSPType;
332 #define TSSSPRKS2  "rks2"
333 #define TSSSPRKS3  "rks3"
334 #define TSSSPRK104 "rk104"
335 
336 PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,TSSSPType);
337 PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,TSSSPType*);
338 PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt);
339 PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*);
340 
341 /*S
342    TSAdapt - Abstract object that manages time-step adaptivity
343 
344    Level: beginner
345 
346 .seealso: TS, TSAdaptCreate(), TSAdaptType
347 S*/
348 typedef struct _p_TSAdapt *TSAdapt;
349 
350 /*E
351     TSAdaptType - String with the name of TSAdapt scheme or the creation function
352        with an optional dynamic library name, for example
353        http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate()
354 
355    Level: beginner
356 
357 .seealso: TSAdaptSetType(), TS
358 E*/
359 typedef const char *TSAdaptType;
360 #define TSADAPTBASIC "basic"
361 #define TSADAPTNONE  "none"
362 #define TSADAPTCFL   "cfl"
363 
364 /*MC
365    TSAdaptRegisterDynamic - adds a TSAdapt implementation
366 
367    Synopsis:
368    PetscErrorCode TSAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
369 
370    Not Collective
371 
372    Input Parameters:
373 +  name_scheme - name of user-defined adaptivity scheme
374 .  path - path (either absolute or relative) the library containing this scheme
375 .  name_create - name of routine to create method context
376 -  routine_create - routine to create method context
377 
378    Notes:
379    TSAdaptRegisterDynamic() may be called multiple times to add several user-defined families.
380 
381    If dynamic libraries are used, then the fourth input argument (routine_create)
382    is ignored.
383 
384    Sample usage:
385 .vb
386    TSAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
387                             "MySchemeCreate",MySchemeCreate);
388 .ve
389 
390    Then, your scheme can be chosen with the procedural interface via
391 $     TSAdaptSetType(ts,"my_scheme")
392    or at runtime via the option
393 $     -ts_adapt_type my_scheme
394 
395    Level: advanced
396 
397    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
398           and others of the form ${any_environmental_variable} occuring in pathname will be
399           replaced with appropriate values.
400 
401 .keywords: TSAdapt, register
402 
403 .seealso: TSAdaptRegisterAll()
404 M*/
405 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
406 #  define TSAdaptRegisterDynamic(a,b,c,d)  TSAdaptRegister(a,b,c,0)
407 #else
408 #  define TSAdaptRegisterDynamic(a,b,c,d)  TSAdaptRegister(a,b,c,d)
409 #endif
410 
411 PETSC_EXTERN PetscErrorCode TSGetAdapt(TS,TSAdapt*);
412 PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSAdapt));
413 PETSC_EXTERN PetscErrorCode TSAdaptRegisterAll(const char[]);
414 PETSC_EXTERN PetscErrorCode TSAdaptRegisterDestroy(void);
415 PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(const char[]);
416 PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
417 PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*);
418 PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,TSAdaptType);
419 PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]);
420 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
421 PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool);
422 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**);
423 PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*);
424 PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscBool*);
425 PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer);
426 PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(TSAdapt);
427 PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*);
428 PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool);
429 PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal);
430 PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt,PetscErrorCode(*)(TSAdapt,TS,PetscBool*));
431 
432 /*S
433    TSGLAdapt - Abstract object that manages time-step adaptivity
434 
435    Level: beginner
436 
437    Developer Notes:
438    This functionality should be replaced by the TSAdapt.
439 
440 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType
441 S*/
442 typedef struct _p_TSGLAdapt *TSGLAdapt;
443 
444 /*J
445     TSGLAdaptType - String with the name of TSGLAdapt scheme or the creation function
446        with an optional dynamic library name, for example
447        http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate()
448 
449    Level: beginner
450 
451 .seealso: TSGLAdaptSetType(), TS
452 J*/
453 typedef const char *TSGLAdaptType;
454 #define TSGLADAPT_NONE "none"
455 #define TSGLADAPT_SIZE "size"
456 #define TSGLADAPT_BOTH "both"
457 
458 /*MC
459    TSGLAdaptRegisterDynamic - adds a TSGLAdapt implementation
460 
461    Synopsis:
462    PetscErrorCode TSGLAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
463 
464    Not Collective
465 
466    Input Parameters:
467 +  name_scheme - name of user-defined adaptivity scheme
468 .  path - path (either absolute or relative) the library containing this scheme
469 .  name_create - name of routine to create method context
470 -  routine_create - routine to create method context
471 
472    Notes:
473    TSGLAdaptRegisterDynamic() may be called multiple times to add several user-defined families.
474 
475    If dynamic libraries are used, then the fourth input argument (routine_create)
476    is ignored.
477 
478    Sample usage:
479 .vb
480    TSGLAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
481                             "MySchemeCreate",MySchemeCreate);
482 .ve
483 
484    Then, your scheme can be chosen with the procedural interface via
485 $     TSGLAdaptSetType(ts,"my_scheme")
486    or at runtime via the option
487 $     -ts_adapt_type my_scheme
488 
489    Level: advanced
490 
491    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
492           and others of the form ${any_environmental_variable} occuring in pathname will be
493           replaced with appropriate values.
494 
495 .keywords: TSGLAdapt, register
496 
497 .seealso: TSGLAdaptRegisterAll()
498 M*/
499 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
500 #  define TSGLAdaptRegisterDynamic(a,b,c,d)  TSGLAdaptRegister(a,b,c,0)
501 #else
502 #  define TSGLAdaptRegisterDynamic(a,b,c,d)  TSGLAdaptRegister(a,b,c,d)
503 #endif
504 
505 PETSC_EXTERN PetscErrorCode TSGLAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSGLAdapt));
506 PETSC_EXTERN PetscErrorCode TSGLAdaptRegisterAll(const char[]);
507 PETSC_EXTERN PetscErrorCode TSGLAdaptRegisterDestroy(void);
508 PETSC_EXTERN PetscErrorCode TSGLAdaptInitializePackage(const char[]);
509 PETSC_EXTERN PetscErrorCode TSGLAdaptFinalizePackage(void);
510 PETSC_EXTERN PetscErrorCode TSGLAdaptCreate(MPI_Comm,TSGLAdapt*);
511 PETSC_EXTERN PetscErrorCode TSGLAdaptSetType(TSGLAdapt,TSGLAdaptType);
512 PETSC_EXTERN PetscErrorCode TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]);
513 PETSC_EXTERN PetscErrorCode TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
514 PETSC_EXTERN PetscErrorCode TSGLAdaptView(TSGLAdapt,PetscViewer);
515 PETSC_EXTERN PetscErrorCode TSGLAdaptSetFromOptions(TSGLAdapt);
516 PETSC_EXTERN PetscErrorCode TSGLAdaptDestroy(TSGLAdapt*);
517 
518 /*J
519     TSGLAcceptType - String with the name of TSGLAccept scheme or the function
520        with an optional dynamic library name, for example
521        http://www.mcs.anl.gov/petsc/lib.a:mytsglaccept()
522 
523    Level: beginner
524 
525 .seealso: TSGLSetAcceptType(), TS
526 J*/
527 typedef const char *TSGLAcceptType;
528 #define TSGLACCEPT_ALWAYS "always"
529 
530 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *);
531 PETSC_EXTERN PetscErrorCode TSGLAcceptRegister(const char[],const char[],const char[],TSGLAcceptFunction);
532 
533 /*MC
534    TSGLAcceptRegisterDynamic - adds a TSGL acceptance scheme
535 
536    Synopsis:
537    PetscErrorCode TSGLAcceptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
538 
539    Not Collective
540 
541    Input Parameters:
542 +  name_scheme - name of user-defined acceptance scheme
543 .  path - path (either absolute or relative) the library containing this scheme
544 .  name_create - name of routine to create method context
545 -  routine_create - routine to create method context
546 
547    Notes:
548    TSGLAcceptRegisterDynamic() may be called multiple times to add several user-defined families.
549 
550    If dynamic libraries are used, then the fourth input argument (routine_create)
551    is ignored.
552 
553    Sample usage:
554 .vb
555    TSGLAcceptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
556                              "MySchemeCreate",MySchemeCreate);
557 .ve
558 
559    Then, your scheme can be chosen with the procedural interface via
560 $     TSGLSetAcceptType(ts,"my_scheme")
561    or at runtime via the option
562 $     -ts_gl_accept_type my_scheme
563 
564    Level: advanced
565 
566    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
567           and others of the form ${any_environmental_variable} occuring in pathname will be
568           replaced with appropriate values.
569 
570 .keywords: TSGL, TSGLAcceptType, register
571 
572 .seealso: TSGLRegisterAll()
573 M*/
574 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
575 #  define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,0)
576 #else
577 #  define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,d)
578 #endif
579 
580 /*J
581   TSGLType - family of time integration method within the General Linear class
582 
583   Level: beginner
584 
585 .seealso: TSGLSetType(), TSGLRegister()
586 J*/
587 typedef const char* TSGLType;
588 #define TSGL_IRKS   "irks"
589 
590 /*MC
591    TSGLRegisterDynamic - adds a TSGL implementation
592 
593    Synopsis:
594    PetscErrorCode TSGLRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
595 
596    Not Collective
597 
598    Input Parameters:
599 +  name_scheme - name of user-defined general linear scheme
600 .  path - path (either absolute or relative) the library containing this scheme
601 .  name_create - name of routine to create method context
602 -  routine_create - routine to create method context
603 
604    Notes:
605    TSGLRegisterDynamic() may be called multiple times to add several user-defined families.
606 
607    If dynamic libraries are used, then the fourth input argument (routine_create)
608    is ignored.
609 
610    Sample usage:
611 .vb
612    TSGLRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
613                        "MySchemeCreate",MySchemeCreate);
614 .ve
615 
616    Then, your scheme can be chosen with the procedural interface via
617 $     TSGLSetType(ts,"my_scheme")
618    or at runtime via the option
619 $     -ts_gl_type my_scheme
620 
621    Level: advanced
622 
623    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
624           and others of the form ${any_environmental_variable} occuring in pathname will be
625           replaced with appropriate values.
626 
627 .keywords: TSGL, register
628 
629 .seealso: TSGLRegisterAll()
630 M*/
631 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
632 #  define TSGLRegisterDynamic(a,b,c,d)       TSGLRegister(a,b,c,0)
633 #else
634 #  define TSGLRegisterDynamic(a,b,c,d)       TSGLRegister(a,b,c,d)
635 #endif
636 
637 PETSC_EXTERN PetscErrorCode TSGLRegister(const char[],const char[],const char[],PetscErrorCode(*)(TS));
638 PETSC_EXTERN PetscErrorCode TSGLRegisterAll(const char[]);
639 PETSC_EXTERN PetscErrorCode TSGLRegisterDestroy(void);
640 PETSC_EXTERN PetscErrorCode TSGLInitializePackage(const char[]);
641 PETSC_EXTERN PetscErrorCode TSGLFinalizePackage(void);
642 PETSC_EXTERN PetscErrorCode TSGLSetType(TS,TSGLType);
643 PETSC_EXTERN PetscErrorCode TSGLGetAdapt(TS,TSGLAdapt*);
644 PETSC_EXTERN PetscErrorCode TSGLSetAcceptType(TS,TSGLAcceptType);
645 
646 /*J
647     TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX method.
648 
649    Level: beginner
650 
651 .seealso: TSARKIMEXSetType(), TS, TSARKIMEX, TSARKIMEXRegister()
652 J*/
653 typedef const char* TSARKIMEXType;
654 #define TSARKIMEXA2     "a2"
655 #define TSARKIMEXL2     "l2"
656 #define TSARKIMEXARS122 "ars122"
657 #define TSARKIMEX2C     "2c"
658 #define TSARKIMEX2D     "2d"
659 #define TSARKIMEX2E     "2e"
660 #define TSARKIMEXPRSSP2 "prssp2"
661 #define TSARKIMEX3      "3"
662 #define TSARKIMEXBPR3   "bpr3"
663 #define TSARKIMEXARS443 "ars443"
664 #define TSARKIMEX4      "4"
665 #define TSARKIMEX5      "5"
666 PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType*);
667 PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType);
668 PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool);
669 PETSC_EXTERN PetscErrorCode TSARKIMEXRegister(TSARKIMEXType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[],const PetscReal[]);
670 PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
671 PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(const char path[]);
672 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);
673 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterAll(void);
674 
675 /*J
676     TSRosWType - String with the name of a Rosenbrock-W method.
677 
678    Level: beginner
679 
680 .seealso: TSRosWSetType(), TS, TSROSW, TSRosWRegister()
681 J*/
682 typedef const char* TSRosWType;
683 #define TSROSW2M          "2m"
684 #define TSROSW2P          "2p"
685 #define TSROSWRA3PW       "ra3pw"
686 #define TSROSWRA34PW2     "ra34pw2"
687 #define TSROSWRODAS3      "rodas3"
688 #define TSROSWSANDU3      "sandu3"
689 #define TSROSWASSP3P3S1C  "assp3p3s1c"
690 #define TSROSWLASSP3P4S2C "lassp3p4s2c"
691 #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
692 #define TSROSWARK3        "ark3"
693 #define TSROSWTHETA1      "theta1"
694 #define TSROSWTHETA2      "theta2"
695 #define TSROSWGRK4T       "grk4t"
696 #define TSROSWSHAMP4      "shamp4"
697 #define TSROSWVELDD4      "veldd4"
698 #define TSROSW4L          "4l"
699 
700 PETSC_EXTERN PetscErrorCode TSRosWGetType(TS ts,TSRosWType*);
701 PETSC_EXTERN PetscErrorCode TSRosWSetType(TS ts,TSRosWType);
702 PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool);
703 PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
704 PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
705 PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
706 PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(const char path[]);
707 PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);
708 PETSC_EXTERN PetscErrorCode TSRosWRegisterAll(void);
709 
710 /*
711        PETSc interface to Sundials
712 */
713 #ifdef PETSC_HAVE_SUNDIALS
714 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
715 PETSC_EXTERN const char *const TSSundialsLmmTypes[];
716 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
717 PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
718 PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType);
719 PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*);
720 PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal);
721 PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal);
722 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal);
723 PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
724 PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
725 PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt);
726 PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal);
727 PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool );
728 PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
729 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt);
730 #endif
731 
732 PETSC_EXTERN PetscErrorCode TSRKSetTolerance(TS,PetscReal);
733 
734 PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal);
735 PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*);
736 PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*);
737 PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool);
738 
739 PETSC_EXTERN PetscErrorCode TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*);
740 PETSC_EXTERN PetscErrorCode TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*);
741 PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal);
742 PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal);
743 PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*);
744 
745 PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM);
746 PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*);
747 
748 PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*);
749 PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
750 
751 #endif
752