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