xref: /petsc/include/petscts.h (revision deb2cd25b4974e90da8065ee93cb8b2449f9ae06)
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 TSSetSNES(TS,SNES);
348 PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*);
349 
350 PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer);
351 PETSC_EXTERN PetscErrorCode TSLoad(TS,PetscViewer);
352 
353 #define TS_FILE_CLASSID 1211225
354 
355 PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *);
356 PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *);
357 
358 typedef struct _n_TSMonitorLGCtx*  TSMonitorLGCtx;
359 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorLGCtx *);
360 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx*);
361 PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS,PetscInt,PetscReal,Vec,void *);
362 PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS,PetscInt,PetscReal,Vec,void *);
363 PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS,PetscInt,PetscReal,Vec,void *);
364 PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS,PetscInt,PetscReal,Vec,void *);
365 PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS,PetscInt,PetscReal,Vec,void *);
366 
367 typedef struct _n_TSMonitorSPEigCtx*  TSMonitorSPEigCtx;
368 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPEigCtx *);
369 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx*);
370 PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS,PetscInt,PetscReal,Vec,void *);
371 
372 /*J
373    TSSSPType - string with the name of TSSSP scheme.
374 
375    Level: beginner
376 
377 .seealso: TSSSPSetType(), TS
378 J*/
379 typedef const char* TSSSPType;
380 #define TSSSPRKS2  "rks2"
381 #define TSSSPRKS3  "rks3"
382 #define TSSSPRK104 "rk104"
383 
384 PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,TSSSPType);
385 PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,TSSSPType*);
386 PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt);
387 PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*);
388 
389 /*S
390    TSAdapt - Abstract object that manages time-step adaptivity
391 
392    Level: beginner
393 
394 .seealso: TS, TSAdaptCreate(), TSAdaptType
395 S*/
396 typedef struct _p_TSAdapt *TSAdapt;
397 
398 /*E
399     TSAdaptType - String with the name of TSAdapt scheme or the creation function
400        with an optional dynamic library name, for example
401        http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate()
402 
403    Level: beginner
404 
405 .seealso: TSAdaptSetType(), TS
406 E*/
407 typedef const char *TSAdaptType;
408 #define TSADAPTBASIC "basic"
409 #define TSADAPTNONE  "none"
410 #define TSADAPTCFL   "cfl"
411 
412 /*MC
413    TSAdaptRegisterDynamic - adds a TSAdapt implementation
414 
415    Synopsis:
416    PetscErrorCode TSAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
417 
418    Not Collective
419 
420    Input Parameters:
421 +  name_scheme - name of user-defined adaptivity scheme
422 .  path - path (either absolute or relative) the library containing this scheme
423 .  name_create - name of routine to create method context
424 -  routine_create - routine to create method context
425 
426    Notes:
427    TSAdaptRegisterDynamic() may be called multiple times to add several user-defined families.
428 
429    If dynamic libraries are used, then the fourth input argument (routine_create)
430    is ignored.
431 
432    Sample usage:
433 .vb
434    TSAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
435                             "MySchemeCreate",MySchemeCreate);
436 .ve
437 
438    Then, your scheme can be chosen with the procedural interface via
439 $     TSAdaptSetType(ts,"my_scheme")
440    or at runtime via the option
441 $     -ts_adapt_type my_scheme
442 
443    Level: advanced
444 
445    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
446           and others of the form ${any_environmental_variable} occuring in pathname will be
447           replaced with appropriate values.
448 
449 .keywords: TSAdapt, register
450 
451 .seealso: TSAdaptRegisterAll()
452 M*/
453 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
454 #  define TSAdaptRegisterDynamic(a,b,c,d)  TSAdaptRegister(a,b,c,0)
455 #else
456 #  define TSAdaptRegisterDynamic(a,b,c,d)  TSAdaptRegister(a,b,c,d)
457 #endif
458 
459 PETSC_EXTERN PetscErrorCode TSGetTSAdapt(TS,TSAdapt*);
460 PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSAdapt));
461 PETSC_EXTERN PetscErrorCode TSAdaptRegisterAll(const char[]);
462 PETSC_EXTERN PetscErrorCode TSAdaptRegisterDestroy(void);
463 PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(const char[]);
464 PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
465 PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*);
466 PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,TSAdaptType);
467 PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]);
468 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
469 PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool);
470 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**);
471 PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*);
472 PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscBool*);
473 PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer);
474 PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt,PetscViewer);
475 PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(TSAdapt);
476 PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*);
477 PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool);
478 PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal);
479 PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt,PetscErrorCode(*)(TSAdapt,TS,PetscBool*));
480 
481 /*S
482    TSGLAdapt - Abstract object that manages time-step adaptivity
483 
484    Level: beginner
485 
486    Developer Notes:
487    This functionality should be replaced by the TSAdapt.
488 
489 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType
490 S*/
491 typedef struct _p_TSGLAdapt *TSGLAdapt;
492 
493 /*J
494     TSGLAdaptType - String with the name of TSGLAdapt scheme or the creation function
495        with an optional dynamic library name, for example
496        http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate()
497 
498    Level: beginner
499 
500 .seealso: TSGLAdaptSetType(), TS
501 J*/
502 typedef const char *TSGLAdaptType;
503 #define TSGLADAPT_NONE "none"
504 #define TSGLADAPT_SIZE "size"
505 #define TSGLADAPT_BOTH "both"
506 
507 /*MC
508    TSGLAdaptRegisterDynamic - adds a TSGLAdapt implementation
509 
510    Synopsis:
511    PetscErrorCode TSGLAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
512 
513    Not Collective
514 
515    Input Parameters:
516 +  name_scheme - name of user-defined adaptivity scheme
517 .  path - path (either absolute or relative) the library containing this scheme
518 .  name_create - name of routine to create method context
519 -  routine_create - routine to create method context
520 
521    Notes:
522    TSGLAdaptRegisterDynamic() may be called multiple times to add several user-defined families.
523 
524    If dynamic libraries are used, then the fourth input argument (routine_create)
525    is ignored.
526 
527    Sample usage:
528 .vb
529    TSGLAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
530                             "MySchemeCreate",MySchemeCreate);
531 .ve
532 
533    Then, your scheme can be chosen with the procedural interface via
534 $     TSGLAdaptSetType(ts,"my_scheme")
535    or at runtime via the option
536 $     -ts_adapt_type my_scheme
537 
538    Level: advanced
539 
540    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
541           and others of the form ${any_environmental_variable} occuring in pathname will be
542           replaced with appropriate values.
543 
544 .keywords: TSGLAdapt, register
545 
546 .seealso: TSGLAdaptRegisterAll()
547 M*/
548 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
549 #  define TSGLAdaptRegisterDynamic(a,b,c,d)  TSGLAdaptRegister(a,b,c,0)
550 #else
551 #  define TSGLAdaptRegisterDynamic(a,b,c,d)  TSGLAdaptRegister(a,b,c,d)
552 #endif
553 
554 PETSC_EXTERN PetscErrorCode TSGLAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSGLAdapt));
555 PETSC_EXTERN PetscErrorCode TSGLAdaptRegisterAll(const char[]);
556 PETSC_EXTERN PetscErrorCode TSGLAdaptRegisterDestroy(void);
557 PETSC_EXTERN PetscErrorCode TSGLAdaptInitializePackage(const char[]);
558 PETSC_EXTERN PetscErrorCode TSGLAdaptFinalizePackage(void);
559 PETSC_EXTERN PetscErrorCode TSGLAdaptCreate(MPI_Comm,TSGLAdapt*);
560 PETSC_EXTERN PetscErrorCode TSGLAdaptSetType(TSGLAdapt,TSGLAdaptType);
561 PETSC_EXTERN PetscErrorCode TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]);
562 PETSC_EXTERN PetscErrorCode TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
563 PETSC_EXTERN PetscErrorCode TSGLAdaptView(TSGLAdapt,PetscViewer);
564 PETSC_EXTERN PetscErrorCode TSGLAdaptSetFromOptions(TSGLAdapt);
565 PETSC_EXTERN PetscErrorCode TSGLAdaptDestroy(TSGLAdapt*);
566 
567 /*J
568     TSGLAcceptType - String with the name of TSGLAccept scheme or the function
569        with an optional dynamic library name, for example
570        http://www.mcs.anl.gov/petsc/lib.a:mytsglaccept()
571 
572    Level: beginner
573 
574 .seealso: TSGLSetAcceptType(), TS
575 J*/
576 typedef const char *TSGLAcceptType;
577 #define TSGLACCEPT_ALWAYS "always"
578 
579 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *);
580 PETSC_EXTERN PetscErrorCode TSGLAcceptRegister(const char[],const char[],const char[],TSGLAcceptFunction);
581 
582 /*MC
583    TSGLAcceptRegisterDynamic - adds a TSGL acceptance scheme
584 
585    Synopsis:
586    PetscErrorCode TSGLAcceptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
587 
588    Not Collective
589 
590    Input Parameters:
591 +  name_scheme - name of user-defined acceptance scheme
592 .  path - path (either absolute or relative) the library containing this scheme
593 .  name_create - name of routine to create method context
594 -  routine_create - routine to create method context
595 
596    Notes:
597    TSGLAcceptRegisterDynamic() may be called multiple times to add several user-defined families.
598 
599    If dynamic libraries are used, then the fourth input argument (routine_create)
600    is ignored.
601 
602    Sample usage:
603 .vb
604    TSGLAcceptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
605                              "MySchemeCreate",MySchemeCreate);
606 .ve
607 
608    Then, your scheme can be chosen with the procedural interface via
609 $     TSGLSetAcceptType(ts,"my_scheme")
610    or at runtime via the option
611 $     -ts_gl_accept_type my_scheme
612 
613    Level: advanced
614 
615    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
616           and others of the form ${any_environmental_variable} occuring in pathname will be
617           replaced with appropriate values.
618 
619 .keywords: TSGL, TSGLAcceptType, register
620 
621 .seealso: TSGLRegisterAll()
622 M*/
623 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
624 #  define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,0)
625 #else
626 #  define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,d)
627 #endif
628 
629 /*J
630   TSGLType - family of time integration method within the General Linear class
631 
632   Level: beginner
633 
634 .seealso: TSGLSetType(), TSGLRegister()
635 J*/
636 typedef const char* TSGLType;
637 #define TSGL_IRKS   "irks"
638 
639 /*MC
640    TSGLRegisterDynamic - adds a TSGL implementation
641 
642    Synopsis:
643    PetscErrorCode TSGLRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
644 
645    Not Collective
646 
647    Input Parameters:
648 +  name_scheme - name of user-defined general linear scheme
649 .  path - path (either absolute or relative) the library containing this scheme
650 .  name_create - name of routine to create method context
651 -  routine_create - routine to create method context
652 
653    Notes:
654    TSGLRegisterDynamic() may be called multiple times to add several user-defined families.
655 
656    If dynamic libraries are used, then the fourth input argument (routine_create)
657    is ignored.
658 
659    Sample usage:
660 .vb
661    TSGLRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
662                        "MySchemeCreate",MySchemeCreate);
663 .ve
664 
665    Then, your scheme can be chosen with the procedural interface via
666 $     TSGLSetType(ts,"my_scheme")
667    or at runtime via the option
668 $     -ts_gl_type my_scheme
669 
670    Level: advanced
671 
672    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
673           and others of the form ${any_environmental_variable} occuring in pathname will be
674           replaced with appropriate values.
675 
676 .keywords: TSGL, register
677 
678 .seealso: TSGLRegisterAll()
679 M*/
680 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
681 #  define TSGLRegisterDynamic(a,b,c,d)       TSGLRegister(a,b,c,0)
682 #else
683 #  define TSGLRegisterDynamic(a,b,c,d)       TSGLRegister(a,b,c,d)
684 #endif
685 
686 PETSC_EXTERN PetscErrorCode TSGLRegister(const char[],const char[],const char[],PetscErrorCode(*)(TS));
687 PETSC_EXTERN PetscErrorCode TSGLRegisterAll(const char[]);
688 PETSC_EXTERN PetscErrorCode TSGLRegisterDestroy(void);
689 PETSC_EXTERN PetscErrorCode TSGLInitializePackage(const char[]);
690 PETSC_EXTERN PetscErrorCode TSGLFinalizePackage(void);
691 PETSC_EXTERN PetscErrorCode TSGLSetType(TS,TSGLType);
692 PETSC_EXTERN PetscErrorCode TSGLGetAdapt(TS,TSGLAdapt*);
693 PETSC_EXTERN PetscErrorCode TSGLSetAcceptType(TS,TSGLAcceptType);
694 
695 /*J
696     TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX method.
697 
698    Level: beginner
699 
700 .seealso: TSARKIMEXSetType(), TS, TSARKIMEX, TSARKIMEXRegister()
701 J*/
702 typedef const char* TSARKIMEXType;
703 #define TSARKIMEXA2     "a2"
704 #define TSARKIMEXL2     "l2"
705 #define TSARKIMEXARS122 "ars122"
706 #define TSARKIMEX2C     "2c"
707 #define TSARKIMEX2D     "2d"
708 #define TSARKIMEX2E     "2e"
709 #define TSARKIMEXPRSSP2 "prssp2"
710 #define TSARKIMEX3      "3"
711 #define TSARKIMEXBPR3   "bpr3"
712 #define TSARKIMEXARS443 "ars443"
713 #define TSARKIMEX4      "4"
714 #define TSARKIMEX5      "5"
715 PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType*);
716 PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType);
717 PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool);
718 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[]);
719 PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
720 PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(const char path[]);
721 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);
722 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterAll(void);
723 
724 /*J
725     TSRosWType - String with the name of a Rosenbrock-W method.
726 
727    Level: beginner
728 
729 .seealso: TSRosWSetType(), TS, TSROSW, TSRosWRegister()
730 J*/
731 typedef const char* TSRosWType;
732 #define TSROSW2M          "2m"
733 #define TSROSW2P          "2p"
734 #define TSROSWRA3PW       "ra3pw"
735 #define TSROSWRA34PW2     "ra34pw2"
736 #define TSROSWRODAS3      "rodas3"
737 #define TSROSWSANDU3      "sandu3"
738 #define TSROSWASSP3P3S1C  "assp3p3s1c"
739 #define TSROSWLASSP3P4S2C "lassp3p4s2c"
740 #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
741 #define TSROSWARK3        "ark3"
742 #define TSROSWTHETA1      "theta1"
743 #define TSROSWTHETA2      "theta2"
744 #define TSROSWGRK4T       "grk4t"
745 #define TSROSWSHAMP4      "shamp4"
746 #define TSROSWVELDD4      "veldd4"
747 #define TSROSW4L          "4l"
748 
749 PETSC_EXTERN PetscErrorCode TSRosWGetType(TS ts,TSRosWType*);
750 PETSC_EXTERN PetscErrorCode TSRosWSetType(TS ts,TSRosWType);
751 PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool);
752 PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
753 PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
754 PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
755 PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(const char path[]);
756 PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);
757 PETSC_EXTERN PetscErrorCode TSRosWRegisterAll(void);
758 
759 /*
760        PETSc interface to Sundials
761 */
762 #ifdef PETSC_HAVE_SUNDIALS
763 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
764 PETSC_EXTERN const char *const TSSundialsLmmTypes[];
765 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
766 PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
767 PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType);
768 PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*);
769 PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal);
770 PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal);
771 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal);
772 PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
773 PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
774 PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt);
775 PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal);
776 PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool );
777 PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
778 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt);
779 #endif
780 
781 PETSC_EXTERN PetscErrorCode TSRKSetTolerance(TS,PetscReal);
782 
783 PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal);
784 PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*);
785 PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*);
786 PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool);
787 
788 PETSC_EXTERN PetscErrorCode TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*);
789 PETSC_EXTERN PetscErrorCode TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*);
790 PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal);
791 PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal);
792 PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*);
793 
794 PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM);
795 PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*);
796 
797 PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*);
798 PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
799 
800 #endif
801