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