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