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