xref: /petsc/include/petscts.h (revision 460c5f81cfdc5b4c06ddbb0232be5fd33ccb813f)
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, TSDestroy()
17 S*/
18 typedef struct _p_TS* TS;
19 
20 /*J
21     TSType - String with the name of a PETSc TS method.
22 
23    Level: beginner
24 
25 .seealso: TSSetType(), TS, TSRegister()
26 J*/
27 typedef const char* TSType;
28 #define TSEULER           "euler"
29 #define TSBEULER          "beuler"
30 #define TSPSEUDO          "pseudo"
31 #define TSCN              "cn"
32 #define TSSUNDIALS        "sundials"
33 #define TSRK              "rk"
34 #define TSPYTHON          "python"
35 #define TSTHETA           "theta"
36 #define TSALPHA           "alpha"
37 #define TSGL              "gl"
38 #define TSSSP             "ssp"
39 #define TSARKIMEX         "arkimex"
40 #define TSROSW            "rosw"
41 #define TSEIMEX           "eimex"
42 #define TSMIMEX           "mimex"
43 /*E
44     TSProblemType - Determines the type of problem this TS object is to be used to solve
45 
46    Level: beginner
47 
48 .seealso: TSCreate()
49 E*/
50 typedef enum {TS_LINEAR,TS_NONLINEAR} TSProblemType;
51 
52 /*E
53    TSEquationType - type of TS problem that is solved
54 
55    Level: beginner
56 
57    Developer Notes: this must match petsc/finclude/petscts.h
58 
59    Supported types are:
60      TS_EQ_UNSPECIFIED (default)
61      TS_EQ_EXPLICIT {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) := M(t) U_t - G(U,t) = 0
62      TS_EQ_IMPLICIT {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) = 0
63 
64 .seealso: TSGetEquationType(), TSSetEquationType()
65 E*/
66 typedef enum {
67   TS_EQ_UNSPECIFIED               = -1,
68   TS_EQ_EXPLICIT                  = 0,
69   TS_EQ_ODE_EXPLICIT              = 1,
70   TS_EQ_DAE_SEMI_EXPLICIT_INDEX1  = 100,
71   TS_EQ_DAE_SEMI_EXPLICIT_INDEX2  = 200,
72   TS_EQ_DAE_SEMI_EXPLICIT_INDEX3  = 300,
73   TS_EQ_DAE_SEMI_EXPLICIT_INDEXHI = 500,
74   TS_EQ_IMPLICIT                  = 1000,
75   TS_EQ_ODE_IMPLICIT              = 1001,
76   TS_EQ_DAE_IMPLICIT_INDEX1       = 1100,
77   TS_EQ_DAE_IMPLICIT_INDEX2       = 1200,
78   TS_EQ_DAE_IMPLICIT_INDEX3       = 1300,
79   TS_EQ_DAE_IMPLICIT_INDEXHI      = 1500
80 } TSEquationType;
81 PETSC_EXTERN const char *const*TSEquationTypes;
82 
83 /*E
84    TSConvergedReason - reason a TS method has converged or not
85 
86    Level: beginner
87 
88    Developer Notes: this must match petsc/finclude/petscts.h
89 
90    Each reason has its own manual page.
91 
92 .seealso: TSGetConvergedReason()
93 E*/
94 typedef enum {
95   TS_CONVERGED_ITERATING      = 0,
96   TS_CONVERGED_TIME           = 1,
97   TS_CONVERGED_ITS            = 2,
98   TS_CONVERGED_USER           = 3,
99   TS_CONVERGED_EVENT          = 4,
100   TS_DIVERGED_NONLINEAR_SOLVE = -1,
101   TS_DIVERGED_STEP_REJECTED   = -2
102 } TSConvergedReason;
103 PETSC_EXTERN const char *const*TSConvergedReasons;
104 
105 /*MC
106    TS_CONVERGED_ITERATING - this only occurs if TSGetConvergedReason() is called during the TSSolve()
107 
108    Level: beginner
109 
110 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt()
111 M*/
112 
113 /*MC
114    TS_CONVERGED_TIME - the final time was reached
115 
116    Level: beginner
117 
118 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetDuration(), TSGetSolveTime()
119 M*/
120 
121 /*MC
122    TS_CONVERGED_ITS - the maximum number of iterations (time-steps) was reached prior to the final time
123 
124    Level: beginner
125 
126 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetDuration()
127 M*/
128 
129 /*MC
130    TS_CONVERGED_USER - user requested termination
131 
132    Level: beginner
133 
134 .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TSSetDuration()
135 M*/
136 
137 /*MC
138    TS_CONVERGED_EVENT - user requested termination on event detection
139 
140    Level: beginner
141 
142 .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TSSetDuration()
143 M*/
144 
145 /*MC
146    TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed
147 
148    Level: beginner
149 
150    Notes: See TSSetMaxSNESFailures() for how to allow more nonlinear solver failures.
151 
152 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSGetSNES(), SNESGetConvergedReason(), TSSetMaxSNESFailures()
153 M*/
154 
155 /*MC
156    TS_DIVERGED_STEP_REJECTED - too many steps were rejected
157 
158    Level: beginner
159 
160    Notes: See TSSetMaxStepRejections() for how to allow more step rejections.
161 
162 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetMaxStepRejections()
163 M*/
164 
165 /*E
166    TSExactFinalTimeOption - option for handling of final time step
167 
168    Level: beginner
169 
170    Developer Notes: this must match petsc/finclude/petscts.h
171 
172 $  TS_EXACTFINALTIME_STEPOVER    - Don't do anything if final time is exceeded
173 $  TS_EXACTFINALTIME_INTERPOLATE - Interpolate back to final time
174 $  TS_EXACTFINALTIME_MATCHSTEP - Adapt final time step to match the final time
175 .seealso: TSGetConvergedReason(), TSSetExactFinalTime()
176 
177 E*/
178 typedef enum {TS_EXACTFINALTIME_STEPOVER=0,TS_EXACTFINALTIME_INTERPOLATE=1,TS_EXACTFINALTIME_MATCHSTEP=2} TSExactFinalTimeOption;
179 PETSC_EXTERN const char *const TSExactFinalTimeOptions[];
180 
181 
182 /* Logging support */
183 PETSC_EXTERN PetscClassId TS_CLASSID;
184 PETSC_EXTERN PetscClassId DMTS_CLASSID;
185 
186 PETSC_EXTERN PetscErrorCode TSInitializePackage(void);
187 
188 PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm,TS*);
189 PETSC_EXTERN PetscErrorCode TSDestroy(TS*);
190 
191 PETSC_EXTERN PetscErrorCode TSSetProblemType(TS,TSProblemType);
192 PETSC_EXTERN PetscErrorCode TSGetProblemType(TS,TSProblemType*);
193 PETSC_EXTERN PetscErrorCode TSMonitor(TS,PetscInt,PetscReal,Vec);
194 PETSC_EXTERN PetscErrorCode TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**));
195 PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS);
196 
197 PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS,const char[]);
198 PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS,const char[]);
199 PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS,const char *[]);
200 PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS);
201 PETSC_EXTERN PetscErrorCode TSSetUp(TS);
202 PETSC_EXTERN PetscErrorCode TSReset(TS);
203 
204 PETSC_EXTERN PetscErrorCode TSSetSolution(TS,Vec);
205 PETSC_EXTERN PetscErrorCode TSGetSolution(TS,Vec*);
206 
207 /*S
208      TSTrajectory - Abstract PETSc object that storing the trajectory (solution of ODE/ADE at each time step and stage)
209 
210    Level: advanced
211 
212   Concepts: ODE solvers, adjoint
213 
214 .seealso:  TSCreate(), TSSetType(), TSType, SNES, KSP, PC, TSDestroy()
215 S*/
216 typedef struct _p_TSTrajectory* TSTrajectory;
217 
218 /*J
219     TSTrajectoryType - String with the name of a PETSc TS trajectory storage method
220 
221    Level: intermediate
222 
223 .seealso: TSSetType(), TS, TSRegister(), TSTrajectoryCreate(), TSTrajectorySetType()
224 J*/
225 typedef const char* TSTrajectoryType;
226 #define TSTRAJECTORYBASIC      "basic"
227 #define TSTRAJECTORYSINGLEFILE "singlefile"
228 
229 PETSC_EXTERN PetscFunctionList TSTrajectoryList;
230 PETSC_EXTERN PetscClassId      TSTRAJECTORY_CLASSID;
231 PETSC_EXTERN PetscBool         TSTrajectoryRegisterAllCalled;
232 
233 PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS);
234 
235 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate(MPI_Comm,TSTrajectory*);
236 PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy(TSTrajectory*);
237 PETSC_EXTERN PetscErrorCode TSTrajectorySetType(TSTrajectory,const TSTrajectoryType);
238 PETSC_EXTERN PetscErrorCode TSTrajectorySet(TSTrajectory,TS,PetscInt,PetscReal,Vec);
239 PETSC_EXTERN PetscErrorCode TSTrajectoryGet(TSTrajectory,TS,PetscInt,PetscReal);
240 PETSC_EXTERN PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory);
241 PETSC_EXTERN PetscErrorCode TSTrajectoryRegisterAll(void);
242 
243 PETSC_EXTERN PetscErrorCode TSSetCostGradients(TS,PetscInt,Vec*,Vec*);
244 PETSC_EXTERN PetscErrorCode TSGetCostGradients(TS,PetscInt*,Vec**,Vec**);
245 PETSC_EXTERN PetscErrorCode TSSetCostIntegrand(TS,PetscInt,PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*),PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*),PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*),void*);
246 PETSC_EXTERN PetscErrorCode TSGetCostIntegral(TS,Vec*);
247 
248 PETSC_EXTERN PetscErrorCode TSAdjointSetRHSJacobian(TS,Mat,PetscErrorCode(*)(TS,PetscReal,Vec,Mat,void*),void*);
249 PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS);
250 PETSC_EXTERN PetscErrorCode TSAdjointSetSteps(TS,PetscInt);
251 
252 PETSC_EXTERN PetscErrorCode TSAdjointComputeRHSJacobian(TS,PetscReal,Vec,Mat);
253 PETSC_EXTERN PetscErrorCode TSAdjointStep(TS);
254 PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS);
255 PETSC_EXTERN PetscErrorCode TSAdjointComputeDRDPFunction(TS,PetscReal,Vec,Vec*);
256 PETSC_EXTERN PetscErrorCode TSAdjointComputeDRDYFunction(TS,PetscReal,Vec,Vec*);
257 PETSC_EXTERN PetscErrorCode TSAdjointComputeCostIntegrand(TS,PetscReal,Vec,Vec);
258 
259 PETSC_EXTERN PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal);
260 PETSC_EXTERN PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*);
261 PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS,TSExactFinalTimeOption);
262 
263 PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,void*);
264 
265 typedef struct _n_TSMonitorDrawCtx*  TSMonitorDrawCtx;
266 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorDrawCtx *);
267 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx*);
268 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS,PetscInt,PetscReal,Vec,void*);
269 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionPhase(TS,PetscInt,PetscReal,Vec,void*);
270 PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS,PetscInt,PetscReal,Vec,void*);
271 
272 PETSC_EXTERN PetscErrorCode TSMonitorSolutionBinary(TS,PetscInt,PetscReal,Vec,void*);
273 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*);
274 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void*);
275 
276 PETSC_EXTERN PetscErrorCode TSStep(TS);
277 PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS,PetscInt,Vec,PetscBool*);
278 PETSC_EXTERN PetscErrorCode TSSolve(TS,Vec);
279 PETSC_EXTERN PetscErrorCode TSGetEquationType(TS,TSEquationType*);
280 PETSC_EXTERN PetscErrorCode TSSetEquationType(TS,TSEquationType);
281 PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*);
282 PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS,TSConvergedReason);
283 PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS,PetscReal*);
284 PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS,PetscInt*);
285 PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS,PetscInt*);
286 PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS,PetscInt*);
287 PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS,PetscInt);
288 PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS,PetscInt*);
289 PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS,PetscInt);
290 PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS,PetscBool);
291 PETSC_EXTERN PetscErrorCode TSRollBack(TS);
292 PETSC_EXTERN PetscErrorCode TSGetTotalSteps(TS,PetscInt*);
293 
294 PETSC_EXTERN PetscErrorCode TSGetStages(TS,PetscInt*,Vec**);
295 
296 PETSC_EXTERN PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal);
297 PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS,PetscReal*);
298 PETSC_EXTERN PetscErrorCode TSGetTime(TS,PetscReal*);
299 PETSC_EXTERN PetscErrorCode TSSetTime(TS,PetscReal);
300 PETSC_EXTERN PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*);
301 PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS,PetscReal);
302 PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS,PetscReal*);
303 
304 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*);
305 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat,Mat,void*);
306 PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*);
307 PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**);
308 PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*);
309 PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**);
310 PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS,PetscBool);
311 
312 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS,PetscReal,Vec,void*);
313 PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS,TSSolutionFunction,void*);
314 PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS,PetscErrorCode (*)(TS,PetscReal,Vec,void*),void*);
315 
316 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*);
317 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
318 PETSC_EXTERN PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*);
319 PETSC_EXTERN PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**);
320 PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*);
321 PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**);
322 
323 PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*);
324 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat,Mat,void*);
325 PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*);
326 PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
327 PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS,PetscReal,Vec);
328 PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS,PetscReal,Vec);
329 PETSC_EXTERN PetscErrorCode TSComputeIJacobianDefaultColor(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
330 
331 PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
332 PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS,PetscReal));
333 PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS,PetscReal,PetscInt,Vec*));
334 PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
335 PETSC_EXTERN PetscErrorCode TSPreStep(TS);
336 PETSC_EXTERN PetscErrorCode TSPreStage(TS,PetscReal);
337 PETSC_EXTERN PetscErrorCode TSPostStage(TS,PetscReal,PetscInt,Vec*);
338 PETSC_EXTERN PetscErrorCode TSPostStep(TS);
339 PETSC_EXTERN PetscErrorCode TSSetRetainStages(TS,PetscBool);
340 PETSC_EXTERN PetscErrorCode TSInterpolate(TS,PetscReal,Vec);
341 PETSC_EXTERN PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec);
342 PETSC_EXTERN PetscErrorCode TSGetTolerances(TS,PetscReal*,Vec*,PetscReal*,Vec*);
343 PETSC_EXTERN PetscErrorCode TSSetDifferentialEquationsIS(TS,IS);
344 PETSC_EXTERN PetscErrorCode TSErrorWeightedNormInfinity(TS,Vec,Vec,PetscReal*);
345 PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm2(TS,Vec,Vec,PetscReal*);
346 PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS,Vec,Vec,NormType,PetscReal*);
347 PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal);
348 PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS,PetscReal*);
349 
350 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*);
351 PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS,PetscReal*,void*);
352 PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *);
353 PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal);
354 
355 PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*);
356 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS,Vec,void*,PetscReal*,PetscBool *);
357 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *);
358 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal);
359 PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);
360 
361 PETSC_EXTERN PetscErrorCode TSPythonSetType(TS,const char[]);
362 
363 PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec);
364 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat,Mat);
365 PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool);
366 PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,PetscBool);
367 PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS,PetscReal,PetscReal,PetscReal*,PetscReal*);
368 
369 PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS,Vec,Vec);
370 
371 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM,TSRHSFunction,void*);
372 PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM,TSRHSFunction*,void**);
373 PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM,TSRHSJacobian,void*);
374 PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM,TSRHSJacobian*,void**);
375 PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM,TSIFunction,void*);
376 PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM,TSIFunction*,void**);
377 PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM,TSIJacobian,void*);
378 PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM,TSIJacobian*,void**);
379 PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM,TSSolutionFunction,void*);
380 PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM,TSSolutionFunction*,void**);
381 PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM,PetscErrorCode (*)(TS,PetscReal,Vec,void*),void*);
382 PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM,PetscErrorCode (**)(TS,PetscReal,Vec,void*),void**);
383 PETSC_EXTERN PetscErrorCode DMTSGetMinRadius(DM,PetscReal*);
384 PETSC_EXTERN PetscErrorCode DMTSSetMinRadius(DM,PetscReal);
385 PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec, void (**)(const PetscReal[], PetscScalar *, void *), void **);
386 
387 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,Vec,void*),void*);
388 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*),void*);
389 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,void*),void*);
390 
391 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));
392 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));
393 
394 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*);
395 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*);
396 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*);
397 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*);
398 
399 PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*),void *);
400 PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*),void *);
401 PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*),void *);
402 PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*),void *);
403 
404 PETSC_EXTERN PetscErrorCode DMPlexTSGetGeometryFVM(DM, Vec *, Vec *, PetscReal *);
405 
406 typedef struct _n_TSMonitorLGCtx*  TSMonitorLGCtx;
407 typedef struct {
408   Vec            ray;
409   VecScatter     scatter;
410   PetscViewer    viewer;
411   TSMonitorLGCtx lgctx;
412 } TSMonitorDMDARayCtx;
413 PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void**);
414 PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS,PetscInt,PetscReal,Vec,void*);
415 PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS,PetscInt,PetscReal,Vec,void*);
416 
417 
418 /* Dynamic creation and loading functions */
419 PETSC_EXTERN PetscFunctionList TSList;
420 PETSC_EXTERN PetscErrorCode TSGetType(TS,TSType*);
421 PETSC_EXTERN PetscErrorCode TSSetType(TS,TSType);
422 PETSC_EXTERN PetscErrorCode TSRegister(const char[], PetscErrorCode (*)(TS));
423 
424 PETSC_EXTERN PetscErrorCode TSGetSNES(TS,SNES*);
425 PETSC_EXTERN PetscErrorCode TSSetSNES(TS,SNES);
426 PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*);
427 
428 PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer);
429 PETSC_EXTERN PetscErrorCode TSLoad(TS,PetscViewer);
430 PETSC_STATIC_INLINE PetscErrorCode TSViewFromOptions(TS A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);}
431 
432 
433 #define TS_FILE_CLASSID 1211225
434 
435 PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *);
436 PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *);
437 
438 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorLGCtx *);
439 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx*);
440 PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS,PetscInt,PetscReal,Vec,void *);
441 PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS,PetscInt,PetscReal,Vec,void *);
442 PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS,const char * const*);
443 PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS,const char *const **);
444 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx,const char * const *);
445 PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS,const char * const*);
446 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx,const char * const*);
447 PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void*);
448 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void *);
449 PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS,PetscInt,PetscReal,Vec,void *);
450 PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS,PetscInt,PetscReal,Vec,void *);
451 PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS,PetscInt,PetscReal,Vec,void *);
452 
453 typedef struct _n_TSMonitorEnvelopeCtx*  TSMonitorEnvelopeCtx;
454 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxCreate(TS,TSMonitorEnvelopeCtx*);
455 PETSC_EXTERN PetscErrorCode TSMonitorEnvelope(TS,PetscInt,PetscReal,Vec,void*);
456 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeGetBounds(TS,Vec*,Vec*);
457 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx*);
458 
459 typedef struct _n_TSMonitorSPEigCtx*  TSMonitorSPEigCtx;
460 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPEigCtx *);
461 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx*);
462 PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS,PetscInt,PetscReal,Vec,void *);
463 
464 PETSC_EXTERN PetscErrorCode TSSetEventMonitor(TS,PetscInt,PetscInt*,PetscBool*,PetscErrorCode (*)(TS,PetscReal,Vec,PetscScalar*,void*),PetscErrorCode (*)(TS,PetscInt,PetscInt[],PetscReal,Vec,PetscBool,void*),void*);
465 /*J
466    TSSSPType - string with the name of TSSSP scheme.
467 
468    Level: beginner
469 
470 .seealso: TSSSPSetType(), TS
471 J*/
472 typedef const char* TSSSPType;
473 #define TSSSPRKS2  "rks2"
474 #define TSSSPRKS3  "rks3"
475 #define TSSSPRK104 "rk104"
476 
477 PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,TSSSPType);
478 PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,TSSSPType*);
479 PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt);
480 PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*);
481 PETSC_EXTERN PetscErrorCode TSSSPFinalizePackage(void);
482 PETSC_EXTERN PetscErrorCode TSSSPInitializePackage(void);
483 
484 /*S
485    TSAdapt - Abstract object that manages time-step adaptivity
486 
487    Level: beginner
488 
489 .seealso: TS, TSAdaptCreate(), TSAdaptType
490 S*/
491 typedef struct _p_TSAdapt *TSAdapt;
492 
493 /*E
494     TSAdaptType - String with the name of TSAdapt scheme.
495 
496    Level: beginner
497 
498 .seealso: TSAdaptSetType(), TS
499 E*/
500 typedef const char *TSAdaptType;
501 #define TSADAPTBASIC "basic"
502 #define TSADAPTNONE  "none"
503 #define TSADAPTCFL   "cfl"
504 
505 PETSC_EXTERN PetscErrorCode TSGetAdapt(TS,TSAdapt*);
506 PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],PetscErrorCode (*)(TSAdapt));
507 PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void);
508 PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
509 PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*);
510 PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,TSAdaptType);
511 PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]);
512 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
513 PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool);
514 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**);
515 PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*);
516 PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscBool*);
517 PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer);
518 PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt,PetscViewer);
519 PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(PetscOptions*,TSAdapt);
520 PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt);
521 PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*);
522 PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool);
523 PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal);
524 PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt,PetscErrorCode(*)(TSAdapt,TS,PetscBool*));
525 
526 /*S
527    TSGLAdapt - Abstract object that manages time-step adaptivity
528 
529    Level: beginner
530 
531    Developer Notes:
532    This functionality should be replaced by the TSAdapt.
533 
534 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType
535 S*/
536 typedef struct _p_TSGLAdapt *TSGLAdapt;
537 
538 /*J
539     TSGLAdaptType - String with the name of TSGLAdapt scheme
540 
541    Level: beginner
542 
543 .seealso: TSGLAdaptSetType(), TS
544 J*/
545 typedef const char *TSGLAdaptType;
546 #define TSGLADAPT_NONE "none"
547 #define TSGLADAPT_SIZE "size"
548 #define TSGLADAPT_BOTH "both"
549 
550 PETSC_EXTERN PetscErrorCode TSGLAdaptRegister(const char[],PetscErrorCode (*)(TSGLAdapt));
551 PETSC_EXTERN PetscErrorCode TSGLAdaptInitializePackage(void);
552 PETSC_EXTERN PetscErrorCode TSGLAdaptFinalizePackage(void);
553 PETSC_EXTERN PetscErrorCode TSGLAdaptCreate(MPI_Comm,TSGLAdapt*);
554 PETSC_EXTERN PetscErrorCode TSGLAdaptSetType(TSGLAdapt,TSGLAdaptType);
555 PETSC_EXTERN PetscErrorCode TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]);
556 PETSC_EXTERN PetscErrorCode TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
557 PETSC_EXTERN PetscErrorCode TSGLAdaptView(TSGLAdapt,PetscViewer);
558 PETSC_EXTERN PetscErrorCode TSGLAdaptSetFromOptions(PetscOptions*,TSGLAdapt);
559 PETSC_EXTERN PetscErrorCode TSGLAdaptDestroy(TSGLAdapt*);
560 
561 /*J
562     TSGLAcceptType - String with the name of TSGLAccept scheme
563 
564    Level: beginner
565 
566 .seealso: TSGLSetAcceptType(), TS
567 J*/
568 typedef const char *TSGLAcceptType;
569 #define TSGLACCEPT_ALWAYS "always"
570 
571 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *);
572 PETSC_EXTERN PetscErrorCode TSGLAcceptRegister(const char[],TSGLAcceptFunction);
573 
574 /*J
575   TSGLType - family of time integration method within the General Linear class
576 
577   Level: beginner
578 
579 .seealso: TSGLSetType(), TSGLRegister()
580 J*/
581 typedef const char* TSGLType;
582 #define TSGL_IRKS   "irks"
583 
584 PETSC_EXTERN PetscErrorCode TSGLRegister(const char[],PetscErrorCode(*)(TS));
585 PETSC_EXTERN PetscErrorCode TSGLInitializePackage(void);
586 PETSC_EXTERN PetscErrorCode TSGLFinalizePackage(void);
587 PETSC_EXTERN PetscErrorCode TSGLSetType(TS,TSGLType);
588 PETSC_EXTERN PetscErrorCode TSGLGetAdapt(TS,TSGLAdapt*);
589 PETSC_EXTERN PetscErrorCode TSGLSetAcceptType(TS,TSGLAcceptType);
590 
591 /*J
592     TSEIMEXType - String with the name of an Extrapolated IMEX method.
593 
594    Level: beginner
595 
596 .seealso: TSEIMEXSetType(), TS, TSEIMEX, TSEIMEXRegister()
597 J*/
598 #define TSEIMEXType   char*
599 
600 PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts,PetscInt);
601 PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts,PetscInt,PetscInt);
602 PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS,PetscBool);
603 
604 /*J
605     TSRKType - String with the name of a Runge-Kutta method.
606 
607    Level: beginner
608 
609 .seealso: TSRKSetType(), TS, TSRK, TSRKRegister()
610 J*/
611 typedef const char* TSRKType;
612 #define TSRK1FE   "1fe"
613 #define TSRK2A    "2a"
614 #define TSRK3     "3"
615 #define TSRK3BS   "3bs"
616 #define TSRK4     "4"
617 #define TSRK5F    "5f"
618 #define TSRK5DP   "5dp"
619 PETSC_EXTERN PetscErrorCode TSRKGetType(TS ts,TSRKType*);
620 PETSC_EXTERN PetscErrorCode TSRKSetType(TS ts,TSRKType);
621 PETSC_EXTERN PetscErrorCode TSRKSetFullyImplicit(TS,PetscBool);
622 PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
623 PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void);
624 PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void);
625 PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void);
626 
627 /*J
628     TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX method.
629 
630    Level: beginner
631 
632 .seealso: TSARKIMEXSetType(), TS, TSARKIMEX, TSARKIMEXRegister()
633 J*/
634 typedef const char* TSARKIMEXType;
635 #define TSARKIMEX1BEE   "1bee"
636 #define TSARKIMEXA2     "a2"
637 #define TSARKIMEXL2     "l2"
638 #define TSARKIMEXARS122 "ars122"
639 #define TSARKIMEX2C     "2c"
640 #define TSARKIMEX2D     "2d"
641 #define TSARKIMEX2E     "2e"
642 #define TSARKIMEXPRSSP2 "prssp2"
643 #define TSARKIMEX3      "3"
644 #define TSARKIMEXBPR3   "bpr3"
645 #define TSARKIMEXARS443 "ars443"
646 #define TSARKIMEX4      "4"
647 #define TSARKIMEX5      "5"
648 PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType*);
649 PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType);
650 PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool);
651 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[]);
652 PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
653 PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void);
654 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);
655 
656 /*J
657     TSRosWType - String with the name of a Rosenbrock-W method.
658 
659    Level: beginner
660 
661 .seealso: TSRosWSetType(), TS, TSROSW, TSRosWRegister()
662 J*/
663 typedef const char* TSRosWType;
664 #define TSROSW2M          "2m"
665 #define TSROSW2P          "2p"
666 #define TSROSWRA3PW       "ra3pw"
667 #define TSROSWRA34PW2     "ra34pw2"
668 #define TSROSWRODAS3      "rodas3"
669 #define TSROSWSANDU3      "sandu3"
670 #define TSROSWASSP3P3S1C  "assp3p3s1c"
671 #define TSROSWLASSP3P4S2C "lassp3p4s2c"
672 #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
673 #define TSROSWARK3        "ark3"
674 #define TSROSWTHETA1      "theta1"
675 #define TSROSWTHETA2      "theta2"
676 #define TSROSWGRK4T       "grk4t"
677 #define TSROSWSHAMP4      "shamp4"
678 #define TSROSWVELDD4      "veldd4"
679 #define TSROSW4L          "4l"
680 
681 PETSC_EXTERN PetscErrorCode TSRosWGetType(TS ts,TSRosWType*);
682 PETSC_EXTERN PetscErrorCode TSRosWSetType(TS ts,TSRosWType);
683 PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool);
684 PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
685 PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
686 PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
687 PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void);
688 PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);
689 
690 /*
691        PETSc interface to Sundials
692 */
693 #ifdef PETSC_HAVE_SUNDIALS
694 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
695 PETSC_EXTERN const char *const TSSundialsLmmTypes[];
696 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
697 PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
698 PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType);
699 PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*);
700 PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal);
701 PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal);
702 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal);
703 PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
704 PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
705 PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt);
706 PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal);
707 PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool );
708 PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
709 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt);
710 #endif
711 
712 PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal);
713 PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*);
714 PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*);
715 PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool);
716 
717 PETSC_EXTERN PetscErrorCode TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*);
718 PETSC_EXTERN PetscErrorCode TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*);
719 PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal);
720 PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal);
721 PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*);
722 
723 PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM);
724 PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*);
725 
726 PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*);
727 PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat,Mat,void*);
728 
729 #endif
730