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