xref: /petsc/include/petscts.h (revision a4e35b1925eceef64945ea472b84f2bf06a67b5e)
1 /*
2    User interface for the timestepping package. This package
3    is for use in solving time-dependent PDEs.
4 */
5 #ifndef PETSCTS_H
6 #define PETSCTS_H
7 
8 #include <petscsnes.h>
9 #include <petscconvest.h>
10 
11 /*I <petscts.h> I*/
12 
13 /* SUBMANSEC = TS */
14 
15 /*S
16      TS - Abstract PETSc object that manages all time-steppers (ODE integrators)
17 
18    Level: beginner
19 
20 .seealso: [](integrator_table), [](ch_ts), `TSCreate()`, `TSSetType()`, `TSType`, `SNES`, `KSP`, `PC`, `TSDestroy()`
21 S*/
22 typedef struct _p_TS *TS;
23 
24 /*J
25     TSType - String with the name of a PETSc `TS` method.
26 
27    Level: beginner
28 
29 .seealso: [](integrator_table), [](ch_ts), `TSSetType()`, `TS`, `TSRegister()`
30 J*/
31 typedef const char *TSType;
32 #define TSEULER           "euler"
33 #define TSBEULER          "beuler"
34 #define TSBASICSYMPLECTIC "basicsymplectic"
35 #define TSPSEUDO          "pseudo"
36 #define TSCN              "cn"
37 #define TSSUNDIALS        "sundials"
38 #define TSRK              "rk"
39 #define TSPYTHON          "python"
40 #define TSTHETA           "theta"
41 #define TSALPHA           "alpha"
42 #define TSALPHA2          "alpha2"
43 #define TSGLLE            "glle"
44 #define TSGLEE            "glee"
45 #define TSSSP             "ssp"
46 #define TSARKIMEX         "arkimex"
47 #define TSROSW            "rosw"
48 #define TSEIMEX           "eimex"
49 #define TSMIMEX           "mimex"
50 #define TSBDF             "bdf"
51 #define TSRADAU5          "radau5"
52 #define TSMPRK            "mprk"
53 #define TSDISCGRAD        "discgrad"
54 #define TSIRK             "irk"
55 
56 /*E
57     TSProblemType - Determines the type of problem this `TS` object is to be used to solve
58 
59    Values:
60  + `TS_LINEAR` - a linear ODE or DAE
61  - `TS_NONLINEAR` - a nonlinear ODE or DAE
62 
63    Level: beginner
64 
65 .seealso: [](ch_ts), `TS`, `TSCreate()`
66 E*/
67 typedef enum {
68   TS_LINEAR,
69   TS_NONLINEAR
70 } TSProblemType;
71 
72 /*E
73    TSEquationType - type of `TS` problem that is solved
74 
75    Level: beginner
76 
77    Values:
78 +  `TS_EQ_UNSPECIFIED` - (default)
79 .  `TS_EQ_EXPLICIT` - {ODE and DAE index 1, 2, 3, HI} F(t,U,U_t) := M(t) U_t - G(U,t) = 0
80 -  `TS_EQ_IMPLICIT` - {ODE and DAE index 1, 2, 3, HI} F(t,U,U_t) = 0
81 
82 .seealso: [](ch_ts), `TS`, `TSGetEquationType()`, `TSSetEquationType()`
83 E*/
84 typedef enum {
85   TS_EQ_UNSPECIFIED               = -1,
86   TS_EQ_EXPLICIT                  = 0,
87   TS_EQ_ODE_EXPLICIT              = 1,
88   TS_EQ_DAE_SEMI_EXPLICIT_INDEX1  = 100,
89   TS_EQ_DAE_SEMI_EXPLICIT_INDEX2  = 200,
90   TS_EQ_DAE_SEMI_EXPLICIT_INDEX3  = 300,
91   TS_EQ_DAE_SEMI_EXPLICIT_INDEXHI = 500,
92   TS_EQ_IMPLICIT                  = 1000,
93   TS_EQ_ODE_IMPLICIT              = 1001,
94   TS_EQ_DAE_IMPLICIT_INDEX1       = 1100,
95   TS_EQ_DAE_IMPLICIT_INDEX2       = 1200,
96   TS_EQ_DAE_IMPLICIT_INDEX3       = 1300,
97   TS_EQ_DAE_IMPLICIT_INDEXHI      = 1500
98 } TSEquationType;
99 PETSC_EXTERN const char *const *TSEquationTypes;
100 
101 /*E
102    TSConvergedReason - reason a `TS` method has converged (integrated to the requested time) or not
103 
104    Values:
105 +  `TS_CONVERGED_ITERATING`          - this only occurs if `TSGetConvergedReason()` is called during the `TSSolve()`
106 .  `TS_CONVERGED_TIME`               - the final time was reached
107 .  `TS_CONVERGED_ITS`                - the maximum number of iterations (time-steps) was reached prior to the final time
108 .  `TS_CONVERGED_USER`               - user requested termination
109 .  `TS_CONVERGED_EVENT`              - user requested termination on event detection
110 .  `TS_CONVERGED_PSEUDO_FATOL`       - stops when function norm decreased by a set amount, used only for `TSPSEUDO`
111 .  `TS_CONVERGED_PSEUDO_FRTOL`       - stops when function norm decreases below a set amount, used only for `TSPSEUDO`
112 .  `TS_DIVERGED_NONLINEAR_SOLVE`     - too many nonlinear solve failures have occurred
113 .  `TS_DIVERGED_STEP_REJECTED`       - too many steps were rejected
114 .  `TSFORWARD_DIVERGED_LINEAR_SOLVE` - tangent linear solve failed
115 -  `TSADJOINT_DIVERGED_LINEAR_SOLVE` - transposed linear solve failed
116 
117    Level: beginner
118 
119 .seealso: [](ch_ts), `TS`, `TSGetConvergedReason()`
120 E*/
121 typedef enum {
122   TS_CONVERGED_ITERATING          = 0,
123   TS_CONVERGED_TIME               = 1,
124   TS_CONVERGED_ITS                = 2,
125   TS_CONVERGED_USER               = 3,
126   TS_CONVERGED_EVENT              = 4,
127   TS_CONVERGED_PSEUDO_FATOL       = 5,
128   TS_CONVERGED_PSEUDO_FRTOL       = 6,
129   TS_DIVERGED_NONLINEAR_SOLVE     = -1,
130   TS_DIVERGED_STEP_REJECTED       = -2,
131   TSFORWARD_DIVERGED_LINEAR_SOLVE = -3,
132   TSADJOINT_DIVERGED_LINEAR_SOLVE = -4
133 } TSConvergedReason;
134 PETSC_EXTERN const char *const *TSConvergedReasons;
135 
136 /*MC
137    TS_CONVERGED_ITERATING - this only occurs if `TSGetConvergedReason()` is called during the `TSSolve()`
138 
139    Level: beginner
140 
141 .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`
142 M*/
143 
144 /*MC
145    TS_CONVERGED_TIME - the final time was reached
146 
147    Level: beginner
148 
149 .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxTime()`, `TSGetMaxTime()`, `TSGetSolveTime()`
150 M*/
151 
152 /*MC
153    TS_CONVERGED_ITS - the maximum number of iterations (time-steps) was reached prior to the final time
154 
155    Level: beginner
156 
157 .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxSteps()`, `TSGetMaxSteps()`
158 M*/
159 
160 /*MC
161    TS_CONVERGED_USER - user requested termination
162 
163    Level: beginner
164 
165 .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`
166 M*/
167 
168 /*MC
169    TS_CONVERGED_EVENT - user requested termination on event detection
170 
171    Level: beginner
172 
173 .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`
174 M*/
175 
176 /*MC
177    TS_CONVERGED_PSEUDO_FRTOL - stops when function norm decreased by a set amount, used only for `TSPSEUDO`
178 
179    Options Database Key:
180 .   -ts_pseudo_frtol <rtol> - use specified rtol
181 
182    Level: beginner
183 
184 .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`, `TS_CONVERGED_PSEUDO_FATOL`
185 M*/
186 
187 /*MC
188    TS_CONVERGED_PSEUDO_FATOL - stops when function norm decreases below a set amount, used only for `TSPSEUDO`
189 
190    Options Database Key:
191 .   -ts_pseudo_fatol <atol> - use specified atol
192 
193    Level: beginner
194 
195 .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`, `TS_CONVERGED_PSEUDO_FRTOL`
196 M*/
197 
198 /*MC
199    TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed
200 
201    Level: beginner
202 
203    Note:
204     See `TSSetMaxSNESFailures()` for how to allow more nonlinear solver failures.
205 
206 .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSGetSNES()`, `SNESGetConvergedReason()`, `TSSetMaxSNESFailures()`
207 M*/
208 
209 /*MC
210    TS_DIVERGED_STEP_REJECTED - too many steps were rejected
211 
212    Level: beginner
213 
214    Notes:
215     See `TSSetMaxStepRejections()` for how to allow more step rejections.
216 
217 .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxStepRejections()`
218 M*/
219 
220 /*E
221    TSExactFinalTimeOption - option for handling of final time step
222 
223    Values:
224 +  `TS_EXACTFINALTIME_STEPOVER`    - Don't do anything if requested final time is exceeded
225 .  `TS_EXACTFINALTIME_INTERPOLATE` - Interpolate back to final time
226 -  `TS_EXACTFINALTIME_MATCHSTEP` - Adapt final time step to match the final time requested
227 
228    Level: beginner
229 
230 .seealso: [](ch_ts), `TS`, `TSGetConvergedReason()`, `TSSetExactFinalTime()`, `TSGetExactFinalTime()`
231 E*/
232 typedef enum {
233   TS_EXACTFINALTIME_UNSPECIFIED = 0,
234   TS_EXACTFINALTIME_STEPOVER    = 1,
235   TS_EXACTFINALTIME_INTERPOLATE = 2,
236   TS_EXACTFINALTIME_MATCHSTEP   = 3
237 } TSExactFinalTimeOption;
238 PETSC_EXTERN const char *const TSExactFinalTimeOptions[];
239 
240 /* Logging support */
241 PETSC_EXTERN PetscClassId TS_CLASSID;
242 PETSC_EXTERN PetscClassId DMTS_CLASSID;
243 PETSC_EXTERN PetscClassId TSADAPT_CLASSID;
244 
245 PETSC_EXTERN PetscErrorCode TSInitializePackage(void);
246 PETSC_EXTERN PetscErrorCode TSFinalizePackage(void);
247 
248 PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm, TS *);
249 PETSC_EXTERN PetscErrorCode TSClone(TS, TS *);
250 PETSC_EXTERN PetscErrorCode TSDestroy(TS *);
251 
252 PETSC_EXTERN PetscErrorCode TSSetProblemType(TS, TSProblemType);
253 PETSC_EXTERN PetscErrorCode TSGetProblemType(TS, TSProblemType *);
254 PETSC_EXTERN PetscErrorCode TSMonitor(TS, PetscInt, PetscReal, Vec);
255 PETSC_EXTERN PetscErrorCode TSMonitorSet(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, void *), void *, PetscErrorCode (*)(void **));
256 PETSC_EXTERN PetscErrorCode TSMonitorSetFromOptions(TS, const char[], const char[], const char[], PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *), PetscErrorCode (*)(TS, PetscViewerAndFormat *));
257 PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS);
258 
259 PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS, const char[]);
260 PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS, const char[]);
261 PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS, const char *[]);
262 PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS);
263 PETSC_EXTERN PetscErrorCode TSSetUp(TS);
264 PETSC_EXTERN PetscErrorCode TSReset(TS);
265 
266 PETSC_EXTERN PetscErrorCode TSSetSolution(TS, Vec);
267 PETSC_EXTERN PetscErrorCode TSGetSolution(TS, Vec *);
268 
269 PETSC_EXTERN PetscErrorCode TS2SetSolution(TS, Vec, Vec);
270 PETSC_EXTERN PetscErrorCode TS2GetSolution(TS, Vec *, Vec *);
271 
272 PETSC_EXTERN PetscErrorCode TSGetSolutionComponents(TS, PetscInt *, Vec *);
273 PETSC_EXTERN PetscErrorCode TSGetAuxSolution(TS, Vec *);
274 PETSC_EXTERN PetscErrorCode TSGetTimeError(TS, PetscInt, Vec *);
275 PETSC_EXTERN PetscErrorCode TSSetTimeError(TS, Vec);
276 
277 PETSC_EXTERN PetscErrorCode TSSetRHSJacobianP(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), void *);
278 PETSC_EXTERN PetscErrorCode TSGetRHSJacobianP(TS, Mat *, PetscErrorCode (**)(TS, PetscReal, Vec, Mat, void *), void **);
279 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianP(TS, PetscReal, Vec, Mat);
280 PETSC_EXTERN PetscErrorCode TSSetIJacobianP(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Vec, PetscReal, Mat, void *), void *);
281 PETSC_EXTERN PetscErrorCode TSComputeIJacobianP(TS, PetscReal, Vec, Vec, PetscReal, Mat, PetscBool);
282 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS() then TSComputeRHSJacobianP()", ) PetscErrorCode TSComputeDRDPFunction(TS, PetscReal, Vec, Vec *);
283 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS() then TSComputeRHSJacobian()", ) PetscErrorCode TSComputeDRDUFunction(TS, PetscReal, Vec, Vec *);
284 PETSC_EXTERN PetscErrorCode TSSetIHessianProduct(TS, Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), void *);
285 PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUU(TS, PetscReal, Vec, Vec *, Vec, Vec *);
286 PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUP(TS, PetscReal, Vec, Vec *, Vec, Vec *);
287 PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPU(TS, PetscReal, Vec, Vec *, Vec, Vec *);
288 PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPP(TS, PetscReal, Vec, Vec *, Vec, Vec *);
289 PETSC_EXTERN PetscErrorCode TSSetRHSHessianProduct(TS, Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), void *);
290 PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS, PetscReal, Vec, Vec *, Vec, Vec *);
291 PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS, PetscReal, Vec, Vec *, Vec, Vec *);
292 PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS, PetscReal, Vec, Vec *, Vec, Vec *);
293 PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS, PetscReal, Vec, Vec *, Vec, Vec *);
294 PETSC_EXTERN PetscErrorCode TSSetCostHessianProducts(TS, PetscInt, Vec *, Vec *, Vec);
295 PETSC_EXTERN PetscErrorCode TSGetCostHessianProducts(TS, PetscInt *, Vec **, Vec **, Vec *);
296 PETSC_EXTERN PetscErrorCode TSComputeSNESJacobian(TS, Vec, Mat, Mat);
297 
298 /*S
299      TSTrajectory - Abstract PETSc object that stores the trajectory (solution of ODE/DAE at each time step)
300 
301    Level: advanced
302 
303 .seealso: [](ch_ts), `TS`, `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectorySetType()`, `TSTrajectoryDestroy()`, `TSTrajectoryReset()`
304 S*/
305 typedef struct _p_TSTrajectory *TSTrajectory;
306 
307 /*J
308     TSTrajectoryType - String with the name of a PETSc `TS` trajectory storage method
309 
310    Level: intermediate
311 
312 .seealso: [](ch_ts), `TS`, `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectoryDestroy()`
313 J*/
314 typedef const char *TSTrajectoryType;
315 #define TSTRAJECTORYBASIC         "basic"
316 #define TSTRAJECTORYSINGLEFILE    "singlefile"
317 #define TSTRAJECTORYMEMORY        "memory"
318 #define TSTRAJECTORYVISUALIZATION "visualization"
319 
320 PETSC_EXTERN PetscFunctionList TSTrajectoryList;
321 PETSC_EXTERN PetscClassId      TSTRAJECTORY_CLASSID;
322 PETSC_EXTERN PetscBool         TSTrajectoryRegisterAllCalled;
323 
324 PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS);
325 PETSC_EXTERN PetscErrorCode TSResetTrajectory(TS);
326 PETSC_EXTERN PetscErrorCode TSRemoveTrajectory(TS);
327 
328 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate(MPI_Comm, TSTrajectory *);
329 PETSC_EXTERN PetscErrorCode TSTrajectoryReset(TSTrajectory);
330 PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy(TSTrajectory *);
331 PETSC_EXTERN PetscErrorCode TSTrajectoryView(TSTrajectory, PetscViewer);
332 PETSC_EXTERN PetscErrorCode TSTrajectorySetType(TSTrajectory, TS, TSTrajectoryType);
333 PETSC_EXTERN PetscErrorCode TSTrajectoryGetType(TSTrajectory, TS, TSTrajectoryType *);
334 PETSC_EXTERN PetscErrorCode TSTrajectorySet(TSTrajectory, TS, PetscInt, PetscReal, Vec);
335 PETSC_EXTERN PetscErrorCode TSTrajectoryGet(TSTrajectory, TS, PetscInt, PetscReal *);
336 PETSC_EXTERN PetscErrorCode TSTrajectoryGetVecs(TSTrajectory, TS, PetscInt, PetscReal *, Vec, Vec);
337 PETSC_EXTERN PetscErrorCode TSTrajectoryGetUpdatedHistoryVecs(TSTrajectory, TS, PetscReal, Vec *, Vec *);
338 PETSC_EXTERN PetscErrorCode TSTrajectoryGetNumSteps(TSTrajectory, PetscInt *);
339 PETSC_EXTERN PetscErrorCode TSTrajectoryRestoreUpdatedHistoryVecs(TSTrajectory, Vec *, Vec *);
340 PETSC_EXTERN PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory, TS);
341 PETSC_EXTERN PetscErrorCode TSTrajectoryRegister(const char[], PetscErrorCode (*)(TSTrajectory, TS));
342 PETSC_EXTERN PetscErrorCode TSTrajectoryRegisterAll(void);
343 PETSC_EXTERN PetscErrorCode TSTrajectorySetUp(TSTrajectory, TS);
344 PETSC_EXTERN PetscErrorCode TSTrajectorySetUseHistory(TSTrajectory, PetscBool);
345 PETSC_EXTERN PetscErrorCode TSTrajectorySetMonitor(TSTrajectory, PetscBool);
346 PETSC_EXTERN PetscErrorCode TSTrajectorySetVariableNames(TSTrajectory, const char *const *);
347 PETSC_EXTERN PetscErrorCode TSTrajectorySetTransform(TSTrajectory, PetscErrorCode (*)(void *, Vec, Vec *), PetscErrorCode (*)(void *), void *);
348 PETSC_EXTERN PetscErrorCode TSTrajectorySetSolutionOnly(TSTrajectory, PetscBool);
349 PETSC_EXTERN PetscErrorCode TSTrajectoryGetSolutionOnly(TSTrajectory, PetscBool *);
350 PETSC_EXTERN PetscErrorCode TSTrajectorySetKeepFiles(TSTrajectory, PetscBool);
351 PETSC_EXTERN PetscErrorCode TSTrajectorySetDirname(TSTrajectory, const char[]);
352 PETSC_EXTERN PetscErrorCode TSTrajectorySetFiletemplate(TSTrajectory, const char[]);
353 PETSC_EXTERN PetscErrorCode TSGetTrajectory(TS, TSTrajectory *);
354 
355 PETSC_EXTERN PetscErrorCode TSSetCostGradients(TS, PetscInt, Vec *, Vec *);
356 PETSC_EXTERN PetscErrorCode TSGetCostGradients(TS, PetscInt *, Vec **, Vec **);
357 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 12, 0, "TSCreateQuadratureTS() and TSForwardSetSensitivities()", ) PetscErrorCode TSSetCostIntegrand(TS, PetscInt, Vec, PetscErrorCode (*)(TS, PetscReal, Vec, Vec, void *), PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, void *), PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, void *), PetscBool, void *);
358 PETSC_EXTERN PetscErrorCode TSGetCostIntegral(TS, Vec *);
359 PETSC_EXTERN PetscErrorCode TSComputeCostIntegrand(TS, PetscReal, Vec, Vec);
360 PETSC_EXTERN PetscErrorCode TSCreateQuadratureTS(TS, PetscBool, TS *);
361 PETSC_EXTERN PetscErrorCode TSGetQuadratureTS(TS, PetscBool *, TS *);
362 
363 PETSC_EXTERN PetscErrorCode TSAdjointSetFromOptions(TS, PetscOptionItems *);
364 PETSC_EXTERN PetscErrorCode TSAdjointMonitor(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *);
365 PETSC_EXTERN PetscErrorCode TSAdjointMonitorSet(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *), void *, PetscErrorCode (*)(void **));
366 PETSC_EXTERN PetscErrorCode TSAdjointMonitorCancel(TS);
367 PETSC_EXTERN PetscErrorCode TSAdjointMonitorSetFromOptions(TS, const char[], const char[], const char[], PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, PetscViewerAndFormat *), PetscErrorCode (*)(TS, PetscViewerAndFormat *));
368 
369 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSSetRHSJacobianP()", ) PetscErrorCode TSAdjointSetRHSJacobian(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), void *);
370 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSComputeRHSJacobianP()", ) PetscErrorCode TSAdjointComputeRHSJacobian(TS, PetscReal, Vec, Mat);
371 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS()", ) PetscErrorCode TSAdjointComputeDRDPFunction(TS, PetscReal, Vec, Vec *);
372 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS()", ) PetscErrorCode TSAdjointComputeDRDYFunction(TS, PetscReal, Vec, Vec *);
373 PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS);
374 PETSC_EXTERN PetscErrorCode TSAdjointSetSteps(TS, PetscInt);
375 
376 PETSC_EXTERN PetscErrorCode TSAdjointStep(TS);
377 PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS);
378 PETSC_EXTERN PetscErrorCode TSAdjointReset(TS);
379 PETSC_EXTERN PetscErrorCode TSAdjointCostIntegral(TS);
380 PETSC_EXTERN PetscErrorCode TSAdjointSetForward(TS, Mat);
381 PETSC_EXTERN PetscErrorCode TSAdjointResetForward(TS);
382 
383 PETSC_EXTERN PetscErrorCode TSForwardSetSensitivities(TS, PetscInt, Mat);
384 PETSC_EXTERN PetscErrorCode TSForwardGetSensitivities(TS, PetscInt *, Mat *);
385 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 12, 0, "TSCreateQuadratureTS()", ) PetscErrorCode TSForwardSetIntegralGradients(TS, PetscInt, Vec *);
386 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 12, 0, "TSForwardGetSensitivities()", ) PetscErrorCode TSForwardGetIntegralGradients(TS, PetscInt *, Vec **);
387 PETSC_EXTERN PetscErrorCode TSForwardSetUp(TS);
388 PETSC_EXTERN PetscErrorCode TSForwardReset(TS);
389 PETSC_EXTERN PetscErrorCode TSForwardCostIntegral(TS);
390 PETSC_EXTERN PetscErrorCode TSForwardStep(TS);
391 PETSC_EXTERN PetscErrorCode TSForwardSetInitialSensitivities(TS, Mat);
392 PETSC_EXTERN PetscErrorCode TSForwardGetStages(TS, PetscInt *, Mat *[]);
393 
394 PETSC_EXTERN PetscErrorCode TSSetMaxSteps(TS, PetscInt);
395 PETSC_EXTERN PetscErrorCode TSGetMaxSteps(TS, PetscInt *);
396 PETSC_EXTERN PetscErrorCode TSSetMaxTime(TS, PetscReal);
397 PETSC_EXTERN PetscErrorCode TSGetMaxTime(TS, PetscReal *);
398 PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS, TSExactFinalTimeOption);
399 PETSC_EXTERN PetscErrorCode TSGetExactFinalTime(TS, TSExactFinalTimeOption *);
400 PETSC_EXTERN PetscErrorCode TSSetTimeSpan(TS, PetscInt, PetscReal *);
401 PETSC_EXTERN PetscErrorCode TSGetTimeSpan(TS, PetscInt *, const PetscReal **);
402 PETSC_EXTERN PetscErrorCode TSGetTimeSpanSolutions(TS, PetscInt *, Vec **);
403 
404 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSSetTime()", ) PetscErrorCode TSSetInitialTimeStep(TS, PetscReal, PetscReal);
405 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSSetMax()", ) PetscErrorCode TSSetDuration(TS, PetscInt, PetscReal);
406 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetMax()", ) PetscErrorCode TSGetDuration(TS, PetscInt *, PetscReal *);
407 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetStepNumber()", ) PetscErrorCode TSGetTimeStepNumber(TS, PetscInt *);
408 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetStepNumber()", ) PetscErrorCode TSGetTotalSteps(TS, PetscInt *);
409 
410 PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
411 PETSC_EXTERN PetscErrorCode TSMonitorExtreme(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
412 
413 typedef struct _n_TSMonitorDrawCtx *TSMonitorDrawCtx;
414 PETSC_EXTERN PetscErrorCode         TSMonitorDrawCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorDrawCtx *);
415 PETSC_EXTERN PetscErrorCode         TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *);
416 PETSC_EXTERN PetscErrorCode         TSMonitorDrawSolution(TS, PetscInt, PetscReal, Vec, void *);
417 PETSC_EXTERN PetscErrorCode         TSMonitorDrawSolutionPhase(TS, PetscInt, PetscReal, Vec, void *);
418 PETSC_EXTERN PetscErrorCode         TSMonitorDrawError(TS, PetscInt, PetscReal, Vec, void *);
419 PETSC_EXTERN PetscErrorCode         TSMonitorDrawSolutionFunction(TS, PetscInt, PetscReal, Vec, void *);
420 
421 PETSC_EXTERN PetscErrorCode TSAdjointMonitorDefault(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, PetscViewerAndFormat *);
422 PETSC_EXTERN PetscErrorCode TSAdjointMonitorDrawSensi(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *);
423 
424 PETSC_EXTERN PetscErrorCode TSMonitorSolution(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
425 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS, PetscInt, PetscReal, Vec, void *);
426 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void *);
427 
428 PETSC_EXTERN PetscErrorCode TSStep(TS);
429 PETSC_EXTERN PetscErrorCode TSEvaluateWLTE(TS, NormType, PetscInt *, PetscReal *);
430 PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS, PetscInt, Vec, PetscBool *);
431 PETSC_EXTERN PetscErrorCode TSSolve(TS, Vec);
432 PETSC_EXTERN PetscErrorCode TSGetEquationType(TS, TSEquationType *);
433 PETSC_EXTERN PetscErrorCode TSSetEquationType(TS, TSEquationType);
434 PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS, TSConvergedReason *);
435 PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS, TSConvergedReason);
436 PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS, PetscReal *);
437 PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS, PetscInt *);
438 PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS, PetscInt *);
439 PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS, PetscInt *);
440 PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS, PetscInt);
441 PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS, PetscInt *);
442 PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS, PetscInt);
443 PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS, PetscBool);
444 PETSC_EXTERN PetscErrorCode TSRestartStep(TS);
445 PETSC_EXTERN PetscErrorCode TSRollBack(TS);
446 
447 PETSC_EXTERN PetscErrorCode TSGetStages(TS, PetscInt *, Vec *[]);
448 
449 PETSC_EXTERN PetscErrorCode TSGetTime(TS, PetscReal *);
450 PETSC_EXTERN PetscErrorCode TSSetTime(TS, PetscReal);
451 PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS, PetscReal *);
452 PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS, PetscReal *);
453 PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS, PetscReal);
454 PETSC_EXTERN PetscErrorCode TSGetStepNumber(TS, PetscInt *);
455 PETSC_EXTERN PetscErrorCode TSSetStepNumber(TS, PetscInt);
456 
457 /*S
458   TSRHSFunction - An `TS` right-hand-side evaluation function
459 
460   Calling Sequence:
461 + ts  - timestep context
462 . t   - current time
463 . u   - input vector
464 . F   - function vector
465 - ctx - [optional] user-defined function context
466 
467   Level: beginner
468 
469 .seealso: [](ch_ts), `TS`, `TSSetRHSFunction()`, `DMTSSetRHSFunction()`, `TSIFunction()`,
470 `TSIJacobian()`, `TSRHSJacobian()`
471 S*/
472 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS ts, PetscReal t, Vec u, Vec F, void *ctx);
473 
474 /*S
475   TSRHSJacobian - A `TS` right-hand-side Jacobian evaluation function
476 
477   Calling Sequence:
478 + ts   - the `TS` context obtained from `TSCreate()`
479 . t    - current time
480 . u    - input vector
481 . Amat - (approximate) Jacobian matrix
482 . Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat)
483 - ctx  - [optional] user-defined context for matrix evaluation routine
484 
485   Level: beginner
486 
487 .seealso: [](ch_ts), `TS`, `TSSetRHSJacobian()`, `DMTSSetRHSJacobian()`, `TSRHSFunction()`,
488 `TSIFunction()`, `TSIJacobian()`
489 S*/
490 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS ts, PetscReal t, Vec u, Mat Amat, Mat Pmat, void *ctx);
491 
492 /*S
493   TSRHSJacobianP - A function that computes the Jacobian of G w.r.t. the parameters P where U_t
494   = G(U,P,t), as well as the location to store the matrix.
495 
496   Calling Sequence:
497 + ts  - the `TS` context
498 . t   - current timestep
499 . U   - input vector (current ODE solution)
500 . A   - output matrix
501 - ctx - [optional] user-defined function context
502 
503   Level: beginner
504 
505 .seealso: [](ch_ts), `TS`, `TSSetRHSJacobianP()`, `TSGetRHSJacobianP()`
506 S*/
507 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobianP)(TS ts, PetscReal t, Vec U, Mat A, void *ctx);
508 
509 PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS, Vec, TSRHSFunction, void *);
510 PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS, Vec *, TSRHSFunction *, void **);
511 PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS, Mat, Mat, TSRHSJacobian, void *);
512 PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS, Mat *, Mat *, TSRHSJacobian *, void **);
513 PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS, PetscBool);
514 
515 /*S
516   TSSolutionFunction - A `TS` solution evaluation function
517 
518   Calling Sequence:
519 + ts  - timestep context
520 . t   - current time
521 . u   - output vector
522 - ctx - [optional] user-defined function context
523 
524   Level: advanced
525 
526 .seealso: [](ch_ts), `TS`, `TSSetSolutionFunction()`, `DMTSSetSolutionFunction()`
527 S*/
528 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS ts, PetscReal t, Vec u, void *ctx);
529 
530 PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS, TSSolutionFunction, void *);
531 
532 /*S
533   TSForcingFunction - A `TS` forcing function evaluation function
534 
535   Calling Sequence:
536 + ts  - timestep context
537 . t   - current time
538 . f   - output vector
539 - ctx - [optional] user-defined function context
540 
541   Level: advanced
542 
543 .seealso: [](ch_ts), `TS`, `TSSetForcingFunction()`, `DMTSSetForcingFunction()`
544 S*/
545 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSForcingFunction)(TS ts, PetscReal t, Vec f, void *ctx);
546 
547 PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS, TSForcingFunction, void *);
548 
549 /*S
550   TSIFunction - A `TS` implicit function evaluation function
551 
552   Calling Sequence:
553 + ts  - the `TS` context obtained from `TSCreate()`
554 . t   - time at step/stage being solved
555 . U   - state vector
556 . U_t - time derivative of state vector
557 . F   - function vector
558 - ctx - [optional] user-defined context for function
559 
560   Level: beginner
561 
562 .seealso: [](ch_ts), `TS`, `TSSetIFunction()`, `DMTSSetIFunction()`, `TSIJacobian()`, `TSRHSFunction()`, `TSRHSJacobian()`
563 S*/
564 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS ts, PetscReal t, Vec U, Vec U_tt, Vec F, void *ctx);
565 
566 /*S
567   TSIJacobian - A `TS` Jacobian evaluation function
568 
569   Calling Sequence:
570 + ts   - the `TS` context obtained from `TSCreate()`
571 . t    - time at step/stage being solved
572 . U    - state vector
573 . U_t  - time derivative of state vector
574 . a    - shift
575 . Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
576 . Pmat - matrix used for constructing preconditioner, usually the same as `Amat`
577 - ctx  - [optional] user-defined context for Jacobian evaluation routine
578 
579   Level: beginner
580 
581 .seealso: [](ch_ts), `TSSetIJacobian()`, `DMTSSetIJacobian()`, `TSIFunction()`, `TSRHSFunction()`, `TSRHSJacobian()`
582 S*/
583 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS ts, PetscReal t, Vec U, Vec U_t, PetscReal a, Mat Amat, Mat Pmat, void *ctx);
584 
585 PETSC_EXTERN PetscErrorCode TSSetIFunction(TS, Vec, TSIFunction, void *);
586 PETSC_EXTERN PetscErrorCode TSGetIFunction(TS, Vec *, TSIFunction *, void **);
587 PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS, Mat, Mat, TSIJacobian, void *);
588 PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS, Mat *, Mat *, TSIJacobian *, void **);
589 
590 /*S
591   TSI2Function - A `TS` implicit function evaluation function for 2nd order systems
592 
593   Calling Sequence:
594 + ts   - the `TS` context obtained from `TSCreate()`
595 . t    - time at step/stage being solved
596 . U    - state vector
597 . U_t  - time derivative of state vector
598 . U_tt - second time derivative of state vector
599 . F    - function vector
600 - ctx  - [optional] user-defined context for matrix evaluation routine (may be `NULL`)
601 
602   Level: advanced
603 
604 .seealso: [](ch_ts), `TS`, `TSSetI2Function()`, `DMTSSetI2Function()`
605 S*/
606 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Function)(TS ts, PetscReal t, Vec U, Vec U_t, Vec U_tt, Vec F, void *ctx);
607 
608 /*S
609   TSI2Jacobian - A `TS` implicit Jacobian evaluation function for 2nd order systems
610 
611   Calling Sequence:
612 + ts   - the `TS` context obtained from `TSCreate()`
613 . t    - time at step/stage being solved
614 . U    - state vector
615 . U_t  - time derivative of state vector
616 . U_tt - second time derivative of state vector
617 . v    - shift for U_t
618 . a    - shift for U_tt
619 . J    - Jacobian of G(U) = F(t,U,W+v*U,W'+a*U), equivalent to dF/dU + v*dF/dU_t  + a*dF/dU_tt
620 . jac  - matrix from which to construct the preconditioner, may be same as `J`
621 - ctx  - [optional] user-defined context for matrix evaluation routine
622 
623   Level: advanced
624 
625 .seealso: [](ch_ts), `TS`, `TSSetI2Jacobian()`, `DMTSSetI2Jacobian()`, `TSIFunction()`, `TSIJacobian()`, `TSRHSFunction()`, `TSRHSJacobian()`
626 S*/
627 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Jacobian)(TS ts, PetscReal t, Vec U, Vec U_t, Vec U_tt, PetscReal v, PetscReal a, Mat J, Mat Jac, void *ctx);
628 
629 PETSC_EXTERN PetscErrorCode TSSetI2Function(TS, Vec, TSI2Function, void *);
630 PETSC_EXTERN PetscErrorCode TSGetI2Function(TS, Vec *, TSI2Function *, void **);
631 PETSC_EXTERN PetscErrorCode TSSetI2Jacobian(TS, Mat, Mat, TSI2Jacobian, void *);
632 PETSC_EXTERN PetscErrorCode TSGetI2Jacobian(TS, Mat *, Mat *, TSI2Jacobian *, void **);
633 
634 PETSC_EXTERN PetscErrorCode TSRHSSplitSetIS(TS, const char[], IS);
635 PETSC_EXTERN PetscErrorCode TSRHSSplitGetIS(TS, const char[], IS *);
636 PETSC_EXTERN PetscErrorCode TSRHSSplitSetRHSFunction(TS, const char[], Vec, TSRHSFunction, void *);
637 PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTS(TS, const char[], TS *);
638 PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTSs(TS, PetscInt *, TS *[]);
639 PETSC_EXTERN PetscErrorCode TSSetUseSplitRHSFunction(TS, PetscBool);
640 PETSC_EXTERN PetscErrorCode TSGetUseSplitRHSFunction(TS, PetscBool *);
641 
642 PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS, PetscReal, Vec, Vec, void *);
643 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS, PetscReal, Vec, Mat, Mat, void *);
644 PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS, PetscReal, Vec, Vec, Vec, void *);
645 PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
646 PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS, PetscReal, Vec);
647 PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS, PetscReal, Vec);
648 PETSC_EXTERN PetscErrorCode TSComputeIJacobianDefaultColor(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
649 PETSC_EXTERN PetscErrorCode TSPruneIJacobianColor(TS, Mat, Mat);
650 
651 PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
652 PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS, PetscReal));
653 PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS, PetscReal, PetscInt, Vec *));
654 PETSC_EXTERN PetscErrorCode TSSetPostEvaluate(TS, PetscErrorCode (*)(TS));
655 PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
656 PETSC_EXTERN PetscErrorCode TSSetResize(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscBool *, void *), PetscErrorCode (*)(TS, PetscInt, Vec[], Vec[], void *), void *);
657 PETSC_EXTERN PetscErrorCode TSPreStep(TS);
658 PETSC_EXTERN PetscErrorCode TSPreStage(TS, PetscReal);
659 PETSC_EXTERN PetscErrorCode TSPostStage(TS, PetscReal, PetscInt, Vec *);
660 PETSC_EXTERN PetscErrorCode TSPostEvaluate(TS);
661 PETSC_EXTERN PetscErrorCode TSPostStep(TS);
662 PETSC_EXTERN PetscErrorCode TSResize(TS);
663 PETSC_EXTERN PetscErrorCode TSResizeRetrieveVec(TS, const char *, Vec *);
664 PETSC_EXTERN PetscErrorCode TSResizeRegisterVec(TS, const char *, Vec);
665 
666 PETSC_EXTERN PetscErrorCode TSInterpolate(TS, PetscReal, Vec);
667 PETSC_EXTERN PetscErrorCode TSSetTolerances(TS, PetscReal, Vec, PetscReal, Vec);
668 PETSC_EXTERN PetscErrorCode TSGetTolerances(TS, PetscReal *, Vec *, PetscReal *, Vec *);
669 PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *);
670 PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm(TS, Vec, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *);
671 PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS, PetscReal);
672 PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS, PetscReal *);
673 PETSC_EXTERN PetscErrorCode TSSetFunctionDomainError(TS, PetscErrorCode (*)(TS, PetscReal, Vec, PetscBool *));
674 PETSC_EXTERN PetscErrorCode TSFunctionDomainError(TS, PetscReal, Vec, PetscBool *);
675 
676 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS, PetscErrorCode (*)(TS, PetscReal *, void *), void *);
677 PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS, PetscReal *, void *);
678 PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS, PetscReal *);
679 PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS, PetscReal);
680 PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS, PetscErrorCode (*)(TS, Vec, void *, PetscReal *, PetscBool *), void *);
681 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS, Vec, void *, PetscReal *, PetscBool *);
682 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS, Vec, PetscReal *, PetscBool *);
683 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS, PetscReal);
684 PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);
685 
686 PETSC_EXTERN PetscErrorCode TSPythonSetType(TS, const char[]);
687 PETSC_EXTERN PetscErrorCode TSPythonGetType(TS, const char *[]);
688 
689 PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS, PetscReal, Vec, Vec);
690 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS, PetscReal, Vec, Mat, Mat);
691 PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS, PetscReal, Vec, Vec, Vec, PetscBool);
692 PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscBool);
693 PETSC_EXTERN PetscErrorCode TSComputeI2Function(TS, PetscReal, Vec, Vec, Vec, Vec);
694 PETSC_EXTERN PetscErrorCode TSComputeI2Jacobian(TS, PetscReal, Vec, Vec, Vec, PetscReal, PetscReal, Mat, Mat);
695 PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS, PetscReal, PetscReal, PetscReal *, PetscReal *);
696 
697 PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS, Vec, Vec);
698 
699 PETSC_EXTERN PetscErrorCode DMTSSetBoundaryLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *);
700 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM, TSRHSFunction, void *);
701 PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM, TSRHSFunction *, void **);
702 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionContextDestroy(DM, PetscErrorCode (*)(void *));
703 PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM, TSRHSJacobian, void *);
704 PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM, TSRHSJacobian *, void **);
705 PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobianContextDestroy(DM, PetscErrorCode (*)(void *));
706 PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM, TSIFunction, void *);
707 PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM, TSIFunction *, void **);
708 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionContextDestroy(DM, PetscErrorCode (*)(void *));
709 PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM, TSIJacobian, void *);
710 PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM, TSIJacobian *, void **);
711 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianContextDestroy(DM, PetscErrorCode (*)(void *));
712 PETSC_EXTERN PetscErrorCode DMTSSetI2Function(DM, TSI2Function, void *);
713 PETSC_EXTERN PetscErrorCode DMTSGetI2Function(DM, TSI2Function *, void **);
714 PETSC_EXTERN PetscErrorCode DMTSSetI2FunctionContextDestroy(DM, PetscErrorCode (*)(void *));
715 PETSC_EXTERN PetscErrorCode DMTSSetI2Jacobian(DM, TSI2Jacobian, void *);
716 PETSC_EXTERN PetscErrorCode DMTSGetI2Jacobian(DM, TSI2Jacobian *, void **);
717 PETSC_EXTERN PetscErrorCode DMTSSetI2JacobianContextDestroy(DM, PetscErrorCode (*)(void *));
718 
719 /*S
720   TSTransientVariable - A function to transform from state to transient variables
721 
722   Calling Sequence:
723 + ts  - timestep context
724 . p   - input vector (primitive form)
725 . c   - output vector, transient variables (conservative form)
726 - ctx - [optional] user-defined function context
727 
728   Level: advanced
729 
730 .seealso: [](ch_ts), `TS`, `TSSetTransientVariable()`, `DMTSSetTransientVariable()`
731 S*/
732 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSTransientVariable)(TS ts, Vec p, Vec c, void *ctx);
733 
734 PETSC_EXTERN PetscErrorCode TSSetTransientVariable(TS, TSTransientVariable, void *);
735 PETSC_EXTERN PetscErrorCode DMTSSetTransientVariable(DM, TSTransientVariable, void *);
736 PETSC_EXTERN PetscErrorCode DMTSGetTransientVariable(DM, TSTransientVariable *, void *);
737 PETSC_EXTERN PetscErrorCode TSComputeTransientVariable(TS, Vec, Vec);
738 PETSC_EXTERN PetscErrorCode TSHasTransientVariable(TS, PetscBool *);
739 
740 PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM, TSSolutionFunction, void *);
741 PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM, TSSolutionFunction *, void **);
742 PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM, TSForcingFunction, void *);
743 PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM, TSForcingFunction *, void **);
744 PETSC_EXTERN PetscErrorCode DMTSCheckResidual(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscReal *);
745 PETSC_EXTERN PetscErrorCode DMTSCheckJacobian(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscBool *, PetscReal *);
746 PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec);
747 
748 PETSC_EXTERN PetscErrorCode DMTSGetIFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, Vec, void *), void **);
749 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, Vec, void *), void *);
750 PETSC_EXTERN PetscErrorCode DMTSGetIJacobianLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void **);
751 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void *);
752 PETSC_EXTERN PetscErrorCode DMTSGetRHSFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, void *), void **);
753 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *);
754 PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrix(DM);
755 PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrixLumped(DM);
756 PETSC_EXTERN PetscErrorCode DMTSDestroyRHSMassMatrix(DM);
757 
758 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM, PetscErrorCode (*)(void *, PetscViewer), PetscErrorCode (*)(void **, PetscViewer));
759 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM, PetscErrorCode (*)(void *, PetscViewer), PetscErrorCode (*)(void **, PetscViewer));
760 
761 /*S
762   DMDATSRHSFunctionLocal - A local `TS` right hand side residual evaluation function for use with `DMDA`
763 
764   Calling Sequence:
765 + info - defines the subdomain to evaluate the residual on
766 . t    - time at which to evaluate residual
767 . x    - array of local state information
768 . f    - output array of local residual information
769 - ctx  - optional user context
770 
771   Level: beginner
772 
773 .seealso: `DMDA`, `DMDATSSetRHSFunctionLocal()`, `TSRHSFunction()`, `DMDATSRHSJacobianLocal()`, `DMDATSIJacobianLocal()`, `DMDATSIFunctionLocal()`
774 S*/
775 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo *info, PetscReal t, void *x, void *f, void *ctx);
776 
777 /*S
778   DMDATSRHSJacobianLocal - A local residual evaluation function for use with `DMDA`
779 
780   Calling Sequence:
781 + info - defines the subdomain to evaluate the residual on
782 . t    - time at which to evaluate residual
783 . x    - array of local state information
784 . J    - Jacobian matrix
785 . B    - matrix from which to construct the preconditioner; often same as `J`
786 - ctx  - optional context
787 
788   Level: beginner
789 
790 .seealso: `DMDA`, `DMDATSSetRHSJacobianLocal()`, `TSRHSJacobian()`, `DMDATSRHSFunctionLocal()`, `DMDATSIJacobianLocal()`, `DMDATSIFunctionLocal()`
791 S*/
792 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo *info, PetscReal t, void *x, Mat J, Mat B, void *ctx);
793 
794 /*S
795   DMDATSIFunctionLocal - A local residual evaluation function for use with `DMDA`
796 
797   Calling Sequence:
798 + info  - defines the subdomain to evaluate the residual on
799 . t     - time at which to evaluate residual
800 . x     - array of local state information
801 . xdot  - array of local time derivative information
802 . imode - output array of local function evaluation information
803 - ctx   - optional context
804 
805   Level: beginner
806 
807 .seealso: `DMDA`, `DMDATSSetIFunctionLocal()`, `DMDATSIJacobianLocal()`, `TSIFunction()`
808 S*/
809 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo *info, PetscReal t, void *x, void *xdot, void *imode, void *ctx);
810 
811 /*S
812   DMDATSIJacobianLocal - A local residual evaluation function for use with `DMDA`
813 
814   Calling Sequence:
815 + info  - defines the subdomain to evaluate the residual on
816 . t     - time at which to evaluate the jacobian
817 . x     - array of local state information
818 . xdot  - time derivative at this state
819 . shift - see `TSSetIJacobian()` for the meaning of this parameter
820 . J     - Jacobian matrix
821 . B     - matrix from which to construct the preconditioner; often same as `J`
822 - ctx   - optional context
823 
824   Level: beginner
825 
826 .seealso: `DMDA` `DMDATSSetIJacobianLocal()`, `TSIJacobian()`, `DMDATSIFunctionLocal()`, `DMDATSRHSFunctionLocal()`,  `DMDATSRHSJacob ianlocal()`
827 S*/
828 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo *info, PetscReal t, void *x, void *xdot, PetscReal shift, Mat J, Mat B, void *ctx);
829 
830 PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM, InsertMode, DMDATSRHSFunctionLocal, void *);
831 PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM, DMDATSRHSJacobianLocal, void *);
832 PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM, InsertMode, DMDATSIFunctionLocal, void *);
833 PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM, DMDATSIJacobianLocal, void *);
834 
835 typedef struct _n_TSMonitorLGCtx *TSMonitorLGCtx;
836 typedef struct {
837   Vec            ray;
838   VecScatter     scatter;
839   PetscViewer    viewer;
840   TSMonitorLGCtx lgctx;
841 } TSMonitorDMDARayCtx;
842 PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void **);
843 PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS, PetscInt, PetscReal, Vec, void *);
844 PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS, PetscInt, PetscReal, Vec, void *);
845 
846 /* Dynamic creation and loading functions */
847 PETSC_EXTERN PetscFunctionList TSList;
848 PETSC_EXTERN PetscErrorCode    TSGetType(TS, TSType *);
849 PETSC_EXTERN PetscErrorCode    TSSetType(TS, TSType);
850 PETSC_EXTERN PetscErrorCode    TSRegister(const char[], PetscErrorCode (*)(TS));
851 
852 PETSC_EXTERN PetscErrorCode TSGetSNES(TS, SNES *);
853 PETSC_EXTERN PetscErrorCode TSSetSNES(TS, SNES);
854 PETSC_EXTERN PetscErrorCode TSGetKSP(TS, KSP *);
855 
856 PETSC_EXTERN PetscErrorCode TSView(TS, PetscViewer);
857 PETSC_EXTERN PetscErrorCode TSLoad(TS, PetscViewer);
858 PETSC_EXTERN PetscErrorCode TSViewFromOptions(TS, PetscObject, const char[]);
859 PETSC_EXTERN PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory, PetscObject, const char[]);
860 
861 #define TS_FILE_CLASSID 1211225
862 
863 PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS, void *);
864 PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS, void *);
865 
866 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtx *);
867 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *);
868 PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS, PetscInt, PetscReal, Vec, void *);
869 PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS, PetscInt, PetscReal, Vec, void *);
870 PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS, const char *const *);
871 PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS, const char *const **);
872 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx, const char *const *);
873 PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS, const char *const *);
874 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx, const char *const *);
875 PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS, PetscErrorCode (*)(void *, Vec, Vec *), PetscErrorCode (*)(void *), void *);
876 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx, PetscErrorCode (*)(void *, Vec, Vec *), PetscErrorCode (*)(void *), void *);
877 PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS, PetscInt, PetscReal, Vec, void *);
878 PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS, PetscInt, PetscReal, Vec, void *);
879 PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS, PetscInt, PetscReal, Vec, void *);
880 PETSC_EXTERN PetscErrorCode TSMonitorError(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
881 PETSC_EXTERN PetscErrorCode TSDMSwarmMonitorMoments(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
882 
883 struct _n_TSMonitorLGCtxNetwork {
884   PetscInt     nlg;
885   PetscDrawLG *lg;
886   PetscBool    semilogy;
887   PetscInt     howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
888 };
889 typedef struct _n_TSMonitorLGCtxNetwork *TSMonitorLGCtxNetwork;
890 PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkDestroy(TSMonitorLGCtxNetwork *);
891 PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkCreate(TS, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtxNetwork *);
892 PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkSolution(TS, PetscInt, PetscReal, Vec, void *);
893 
894 typedef struct _n_TSMonitorEnvelopeCtx *TSMonitorEnvelopeCtx;
895 PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeCtxCreate(TS, TSMonitorEnvelopeCtx *);
896 PETSC_EXTERN PetscErrorCode             TSMonitorEnvelope(TS, PetscInt, PetscReal, Vec, void *);
897 PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeGetBounds(TS, Vec *, Vec *);
898 PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *);
899 
900 typedef struct _n_TSMonitorSPEigCtx *TSMonitorSPEigCtx;
901 PETSC_EXTERN PetscErrorCode          TSMonitorSPEigCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorSPEigCtx *);
902 PETSC_EXTERN PetscErrorCode          TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *);
903 PETSC_EXTERN PetscErrorCode          TSMonitorSPEig(TS, PetscInt, PetscReal, Vec, void *);
904 
905 typedef struct _n_TSMonitorSPCtx *TSMonitorSPCtx;
906 PETSC_EXTERN PetscErrorCode       TSMonitorSPCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscBool, PetscBool, TSMonitorSPCtx *);
907 PETSC_EXTERN PetscErrorCode       TSMonitorSPCtxDestroy(TSMonitorSPCtx *);
908 PETSC_EXTERN PetscErrorCode       TSMonitorSPSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
909 
910 typedef struct _n_TSMonitorHGCtx *TSMonitorHGCtx;
911 PETSC_EXTERN PetscErrorCode       TSMonitorHGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscInt, PetscBool, TSMonitorHGCtx *);
912 PETSC_EXTERN PetscErrorCode       TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
913 PETSC_EXTERN PetscErrorCode       TSMonitorHGCtxDestroy(TSMonitorHGCtx *);
914 PETSC_EXTERN PetscErrorCode       TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
915 
916 PETSC_EXTERN PetscErrorCode TSSetEventHandler(TS, PetscInt, PetscInt[], PetscBool[], PetscErrorCode (*)(TS, PetscReal, Vec, PetscScalar[], void *), PetscErrorCode (*)(TS, PetscInt, PetscInt[], PetscReal, Vec, PetscBool, void *), void *);
917 PETSC_EXTERN PetscErrorCode TSSetPostEventIntervalStep(TS, PetscReal);
918 PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS, PetscReal, PetscReal[]);
919 PETSC_EXTERN PetscErrorCode TSGetNumEvents(TS, PetscInt *);
920 
921 /*J
922    TSSSPType - string with the name of a `TSSSP` scheme.
923 
924    Level: beginner
925 
926 .seealso: [](ch_ts), `TSSSPSetType()`, `TS`, `TSSSP`
927 J*/
928 typedef const char *TSSSPType;
929 #define TSSSPRKS2  "rks2"
930 #define TSSSPRKS3  "rks3"
931 #define TSSSPRK104 "rk104"
932 
933 PETSC_EXTERN PetscErrorCode    TSSSPSetType(TS, TSSSPType);
934 PETSC_EXTERN PetscErrorCode    TSSSPGetType(TS, TSSSPType *);
935 PETSC_EXTERN PetscErrorCode    TSSSPSetNumStages(TS, PetscInt);
936 PETSC_EXTERN PetscErrorCode    TSSSPGetNumStages(TS, PetscInt *);
937 PETSC_EXTERN PetscErrorCode    TSSSPInitializePackage(void);
938 PETSC_EXTERN PetscErrorCode    TSSSPFinalizePackage(void);
939 PETSC_EXTERN PetscFunctionList TSSSPList;
940 
941 /*S
942    TSAdapt - Abstract object that manages time-step adaptivity
943 
944    Level: beginner
945 
946 .seealso: [](ch_ts), `TS`, `TSAdaptCreate()`, `TSAdaptType`
947 S*/
948 typedef struct _p_TSAdapt *TSAdapt;
949 
950 /*J
951     TSAdaptType - String with the name of `TSAdapt` scheme.
952 
953    Level: beginner
954 
955 .seealso: [](ch_ts), `TSAdaptSetType()`, `TS`, `TSAdapt`
956 J*/
957 typedef const char *TSAdaptType;
958 #define TSADAPTNONE    "none"
959 #define TSADAPTBASIC   "basic"
960 #define TSADAPTDSP     "dsp"
961 #define TSADAPTCFL     "cfl"
962 #define TSADAPTGLEE    "glee"
963 #define TSADAPTHISTORY "history"
964 
965 PETSC_EXTERN PetscErrorCode TSGetAdapt(TS, TSAdapt *);
966 PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[], PetscErrorCode (*)(TSAdapt));
967 PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void);
968 PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
969 PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm, TSAdapt *);
970 PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt, TSAdaptType);
971 PETSC_EXTERN PetscErrorCode TSAdaptGetType(TSAdapt, TSAdaptType *);
972 PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt, const char[]);
973 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
974 PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt, const char[], PetscInt, PetscInt, PetscReal, PetscReal, PetscBool);
975 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt, PetscInt *, const PetscInt **, const PetscInt **, const PetscReal **, const PetscReal **);
976 PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt, TS, PetscReal, PetscInt *, PetscReal *, PetscBool *);
977 PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt, TS, PetscReal, Vec, PetscBool *);
978 PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt, PetscViewer);
979 PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt, PetscViewer);
980 PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(TSAdapt, PetscOptionItems *);
981 PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt);
982 PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt *);
983 PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt, PetscBool);
984 PETSC_EXTERN PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt, PetscBool);
985 PETSC_EXTERN PetscErrorCode TSAdaptSetSafety(TSAdapt, PetscReal, PetscReal);
986 PETSC_EXTERN PetscErrorCode TSAdaptGetSafety(TSAdapt, PetscReal *, PetscReal *);
987 PETSC_EXTERN PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt, PetscReal);
988 PETSC_EXTERN PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt, PetscReal *);
989 PETSC_EXTERN PetscErrorCode TSAdaptSetClip(TSAdapt, PetscReal, PetscReal);
990 PETSC_EXTERN PetscErrorCode TSAdaptGetClip(TSAdapt, PetscReal *, PetscReal *);
991 PETSC_EXTERN PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt, PetscReal);
992 PETSC_EXTERN PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt, PetscReal *);
993 PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt, PetscReal, PetscReal);
994 PETSC_EXTERN PetscErrorCode TSAdaptGetStepLimits(TSAdapt, PetscReal *, PetscReal *);
995 PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt, PetscErrorCode (*)(TSAdapt, TS, PetscReal, Vec, PetscBool *));
996 PETSC_EXTERN PetscErrorCode TSAdaptHistorySetHistory(TSAdapt, PetscInt n, PetscReal hist[], PetscBool);
997 PETSC_EXTERN PetscErrorCode TSAdaptHistorySetTrajectory(TSAdapt, TSTrajectory, PetscBool);
998 PETSC_EXTERN PetscErrorCode TSAdaptHistoryGetStep(TSAdapt, PetscInt, PetscReal *, PetscReal *);
999 PETSC_EXTERN PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt, PetscInt);
1000 PETSC_EXTERN PetscErrorCode TSAdaptDSPSetFilter(TSAdapt, const char *);
1001 PETSC_EXTERN PetscErrorCode TSAdaptDSPSetPID(TSAdapt, PetscReal, PetscReal, PetscReal);
1002 
1003 /*S
1004    TSGLLEAdapt - Abstract object that manages time-step adaptivity for `TSGLLE`
1005 
1006    Level: beginner
1007 
1008    Developer Note:
1009    This functionality should be replaced by the `TSAdapt`.
1010 
1011 .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLEAdaptCreate()`, `TSGLLEAdaptType`
1012 S*/
1013 typedef struct _p_TSGLLEAdapt *TSGLLEAdapt;
1014 
1015 /*J
1016     TSGLLEAdaptType - String with the name of `TSGLLEAdapt` scheme
1017 
1018    Level: beginner
1019 
1020    Developer Note:
1021    This functionality should be replaced by the `TSAdaptType`.
1022 
1023 .seealso: [](ch_ts), `TSGLLEAdaptSetType()`, `TS`
1024 J*/
1025 typedef const char *TSGLLEAdaptType;
1026 #define TSGLLEADAPT_NONE "none"
1027 #define TSGLLEADAPT_SIZE "size"
1028 #define TSGLLEADAPT_BOTH "both"
1029 
1030 PETSC_EXTERN PetscErrorCode TSGLLEAdaptRegister(const char[], PetscErrorCode (*)(TSGLLEAdapt));
1031 PETSC_EXTERN PetscErrorCode TSGLLEAdaptInitializePackage(void);
1032 PETSC_EXTERN PetscErrorCode TSGLLEAdaptFinalizePackage(void);
1033 PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate(MPI_Comm, TSGLLEAdapt *);
1034 PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt, TSGLLEAdaptType);
1035 PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt, const char[]);
1036 PETSC_EXTERN PetscErrorCode TSGLLEAdaptChoose(TSGLLEAdapt, PetscInt, const PetscInt[], const PetscReal[], const PetscReal[], PetscInt, PetscReal, PetscReal, PetscInt *, PetscReal *, PetscBool *);
1037 PETSC_EXTERN PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt, PetscViewer);
1038 PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetFromOptions(TSGLLEAdapt, PetscOptionItems *);
1039 PETSC_EXTERN PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt *);
1040 
1041 /*J
1042     TSGLLEAcceptType - String with the name of `TSGLLEAccept` scheme
1043 
1044    Level: beginner
1045 
1046 .seealso: [](ch_ts), `TSGLLESetAcceptType()`, `TS`, `TSGLLEAccept`
1047 J*/
1048 typedef const char *TSGLLEAcceptType;
1049 #define TSGLLEACCEPT_ALWAYS "always"
1050 
1051 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLLEAcceptFunction)(TS, PetscReal, PetscReal, const PetscReal[], PetscBool *);
1052 PETSC_EXTERN PetscErrorCode TSGLLEAcceptRegister(const char[], TSGLLEAcceptFunction);
1053 
1054 /*J
1055   TSGLLEType - string with the name of a General Linear `TSGLLE` type
1056 
1057   Level: beginner
1058 
1059 .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLESetType()`, `TSGLLERegister()`, `TSGLLEAccept`
1060 J*/
1061 typedef const char *TSGLLEType;
1062 #define TSGLLE_IRKS "irks"
1063 
1064 PETSC_EXTERN PetscErrorCode TSGLLERegister(const char[], PetscErrorCode (*)(TS));
1065 PETSC_EXTERN PetscErrorCode TSGLLEInitializePackage(void);
1066 PETSC_EXTERN PetscErrorCode TSGLLEFinalizePackage(void);
1067 PETSC_EXTERN PetscErrorCode TSGLLESetType(TS, TSGLLEType);
1068 PETSC_EXTERN PetscErrorCode TSGLLEGetAdapt(TS, TSGLLEAdapt *);
1069 PETSC_EXTERN PetscErrorCode TSGLLESetAcceptType(TS, TSGLLEAcceptType);
1070 
1071 /*J
1072     TSEIMEXType - String with the name of an Extrapolated IMEX `TSEIMEX` type
1073 
1074    Level: beginner
1075 
1076 .seealso: [](ch_ts), `TSEIMEXSetType()`, `TS`, `TSEIMEX`, `TSEIMEXRegister()`
1077 J*/
1078 #define TSEIMEXType char *
1079 
1080 PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts, PetscInt);
1081 PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts, PetscInt, PetscInt);
1082 PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS, PetscBool);
1083 
1084 /*J
1085     TSRKType - String with the name of a Runge-Kutta `TSRK` type
1086 
1087    Level: beginner
1088 
1089 .seealso: [](ch_ts), `TS`, `TSRKSetType()`, `TS`, `TSRK`, `TSRKRegister()`
1090 J*/
1091 typedef const char *TSRKType;
1092 #define TSRK1FE "1fe"
1093 #define TSRK2A  "2a"
1094 #define TSRK2B  "2b"
1095 #define TSRK3   "3"
1096 #define TSRK3BS "3bs"
1097 #define TSRK4   "4"
1098 #define TSRK5F  "5f"
1099 #define TSRK5DP "5dp"
1100 #define TSRK5BS "5bs"
1101 #define TSRK6VR "6vr"
1102 #define TSRK7VR "7vr"
1103 #define TSRK8VR "8vr"
1104 
1105 PETSC_EXTERN PetscErrorCode TSRKGetOrder(TS, PetscInt *);
1106 PETSC_EXTERN PetscErrorCode TSRKGetType(TS, TSRKType *);
1107 PETSC_EXTERN PetscErrorCode TSRKSetType(TS, TSRKType);
1108 PETSC_EXTERN PetscErrorCode TSRKGetTableau(TS, PetscInt *, const PetscReal **, const PetscReal **, const PetscReal **, const PetscReal **, PetscInt *, const PetscReal **, PetscBool *);
1109 PETSC_EXTERN PetscErrorCode TSRKSetMultirate(TS, PetscBool);
1110 PETSC_EXTERN PetscErrorCode TSRKGetMultirate(TS, PetscBool *);
1111 PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
1112 PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void);
1113 PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void);
1114 PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void);
1115 
1116 /*J
1117    TSMPRKType - String with the name of a Partitioned Runge-Kutta `TSMPRK` type
1118 
1119    Level: beginner
1120 
1121 .seealso: [](ch_ts), `TSMPRKSetType()`, `TS`, `TSMPRK`, `TSMPRKRegister()`
1122 J*/
1123 typedef const char *TSMPRKType;
1124 #define TSMPRK2A22 "2a22"
1125 #define TSMPRK2A23 "2a23"
1126 #define TSMPRK2A32 "2a32"
1127 #define TSMPRK2A33 "2a33"
1128 #define TSMPRKP2   "p2"
1129 #define TSMPRKP3   "p3"
1130 
1131 PETSC_EXTERN PetscErrorCode TSMPRKGetType(TS ts, TSMPRKType *);
1132 PETSC_EXTERN PetscErrorCode TSMPRKSetType(TS ts, TSMPRKType);
1133 PETSC_EXTERN PetscErrorCode TSMPRKRegister(TSMPRKType, PetscInt, PetscInt, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscInt[], const PetscReal[], const PetscReal[], const PetscReal[], const PetscInt[], const PetscReal[], const PetscReal[], const PetscReal[]);
1134 PETSC_EXTERN PetscErrorCode TSMPRKInitializePackage(void);
1135 PETSC_EXTERN PetscErrorCode TSMPRKFinalizePackage(void);
1136 PETSC_EXTERN PetscErrorCode TSMPRKRegisterDestroy(void);
1137 
1138 /*J
1139     TSIRKType - String with the name of an implicit Runge-Kutta `TSIRK` type
1140 
1141    Level: beginner
1142 
1143 .seealso: [](ch_ts), `TSIRKSetType()`, `TS`, `TSIRK`, `TSIRKRegister()`
1144 J*/
1145 typedef const char *TSIRKType;
1146 #define TSIRKGAUSS "gauss"
1147 
1148 PETSC_EXTERN PetscErrorCode TSIRKGetType(TS, TSIRKType *);
1149 PETSC_EXTERN PetscErrorCode TSIRKSetType(TS, TSIRKType);
1150 PETSC_EXTERN PetscErrorCode TSIRKGetNumStages(TS, PetscInt *);
1151 PETSC_EXTERN PetscErrorCode TSIRKSetNumStages(TS, PetscInt);
1152 PETSC_EXTERN PetscErrorCode TSIRKRegister(const char[], PetscErrorCode (*function)(TS));
1153 PETSC_EXTERN PetscErrorCode TSIRKTableauCreate(TS, PetscInt, const PetscReal *, const PetscReal *, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, const PetscScalar *);
1154 PETSC_EXTERN PetscErrorCode TSIRKInitializePackage(void);
1155 PETSC_EXTERN PetscErrorCode TSIRKFinalizePackage(void);
1156 PETSC_EXTERN PetscErrorCode TSIRKRegisterDestroy(void);
1157 
1158 /*J
1159     TSGLEEType - String with the name of a General Linear with Error Estimation `TSGLEE` type
1160 
1161    Level: beginner
1162 
1163 .seealso: [](ch_ts), `TSGLEESetType()`, `TS`, `TSGLEE`, `TSGLEERegister()`
1164 J*/
1165 typedef const char *TSGLEEType;
1166 #define TSGLEEi1      "BE1"
1167 #define TSGLEE23      "23"
1168 #define TSGLEE24      "24"
1169 #define TSGLEE25I     "25i"
1170 #define TSGLEE35      "35"
1171 #define TSGLEEEXRK2A  "exrk2a"
1172 #define TSGLEERK32G1  "rk32g1"
1173 #define TSGLEERK285EX "rk285ex"
1174 
1175 /*J
1176     TSGLEEMode - String with the mode of error estimation for a General Linear with Error Estimation `TSGLEE` type
1177 
1178    Level: beginner
1179 
1180 .seealso: [](ch_ts), `TSGLEESetMode()`, `TS`, `TSGLEE`, `TSGLEERegister()`
1181 J*/
1182 PETSC_EXTERN PetscErrorCode TSGLEEGetType(TS ts, TSGLEEType *);
1183 PETSC_EXTERN PetscErrorCode TSGLEESetType(TS ts, TSGLEEType);
1184 PETSC_EXTERN PetscErrorCode TSGLEERegister(TSGLEEType, PetscInt, PetscInt, PetscInt, PetscReal, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
1185 PETSC_EXTERN PetscErrorCode TSGLEEFinalizePackage(void);
1186 PETSC_EXTERN PetscErrorCode TSGLEEInitializePackage(void);
1187 PETSC_EXTERN PetscErrorCode TSGLEERegisterDestroy(void);
1188 
1189 /*J
1190     TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX `TSARKIMEX` type
1191 
1192    Level: beginner
1193 
1194 .seealso: [](ch_ts), `TSARKIMEXSetType()`, `TS`, `TSARKIMEX`, `TSARKIMEXRegister()`
1195 J*/
1196 typedef const char *TSARKIMEXType;
1197 #define TSARKIMEX1BEE   "1bee"
1198 #define TSARKIMEXA2     "a2"
1199 #define TSARKIMEXL2     "l2"
1200 #define TSARKIMEXARS122 "ars122"
1201 #define TSARKIMEX2C     "2c"
1202 #define TSARKIMEX2D     "2d"
1203 #define TSARKIMEX2E     "2e"
1204 #define TSARKIMEXPRSSP2 "prssp2"
1205 #define TSARKIMEX3      "3"
1206 #define TSARKIMEXBPR3   "bpr3"
1207 #define TSARKIMEXARS443 "ars443"
1208 #define TSARKIMEX4      "4"
1209 #define TSARKIMEX5      "5"
1210 PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts, TSARKIMEXType *);
1211 PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts, TSARKIMEXType);
1212 PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS, PetscBool);
1213 PETSC_EXTERN PetscErrorCode TSARKIMEXGetFullyImplicit(TS, PetscBool *);
1214 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[]);
1215 PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void);
1216 PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
1217 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);
1218 
1219 /*J
1220     TSRosWType - String with the name of a Rosenbrock-W `TSROSW` type
1221 
1222    Level: beginner
1223 
1224 .seealso: [](ch_ts), `TSRosWSetType()`, `TS`, `TSROSW`, `TSRosWRegister()`
1225 J*/
1226 typedef const char *TSRosWType;
1227 #define TSROSW2M          "2m"
1228 #define TSROSW2P          "2p"
1229 #define TSROSWRA3PW       "ra3pw"
1230 #define TSROSWRA34PW2     "ra34pw2"
1231 #define TSROSWRODAS3      "rodas3"
1232 #define TSROSWSANDU3      "sandu3"
1233 #define TSROSWASSP3P3S1C  "assp3p3s1c"
1234 #define TSROSWLASSP3P4S2C "lassp3p4s2c"
1235 #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
1236 #define TSROSWARK3        "ark3"
1237 #define TSROSWTHETA1      "theta1"
1238 #define TSROSWTHETA2      "theta2"
1239 #define TSROSWGRK4T       "grk4t"
1240 #define TSROSWSHAMP4      "shamp4"
1241 #define TSROSWVELDD4      "veldd4"
1242 #define TSROSW4L          "4l"
1243 
1244 PETSC_EXTERN PetscErrorCode TSRosWGetType(TS, TSRosWType *);
1245 PETSC_EXTERN PetscErrorCode TSRosWSetType(TS, TSRosWType);
1246 PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS, PetscBool);
1247 PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
1248 PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal);
1249 PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void);
1250 PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
1251 PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);
1252 
1253 PETSC_EXTERN PetscErrorCode TSBDFSetOrder(TS, PetscInt);
1254 PETSC_EXTERN PetscErrorCode TSBDFGetOrder(TS, PetscInt *);
1255 
1256 /*J
1257   TSBasicSymplecticType - String with the name of a basic symplectic integration `TSBASICSYMPLECTIC` type
1258 
1259   Level: beginner
1260 
1261 .seealso: [](ch_ts), `TSBasicSymplecticSetType()`, `TS`, `TSBASICSYMPLECTIC`, `TSBasicSymplecticRegister()`
1262 J*/
1263 typedef const char *TSBasicSymplecticType;
1264 #define TSBASICSYMPLECTICSIEULER   "1"
1265 #define TSBASICSYMPLECTICVELVERLET "2"
1266 #define TSBASICSYMPLECTIC3         "3"
1267 #define TSBASICSYMPLECTIC4         "4"
1268 PETSC_EXTERN PetscErrorCode TSBasicSymplecticSetType(TS, TSBasicSymplecticType);
1269 PETSC_EXTERN PetscErrorCode TSBasicSymplecticGetType(TS, TSBasicSymplecticType *);
1270 PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegister(TSBasicSymplecticType, PetscInt, PetscInt, PetscReal[], PetscReal[]);
1271 PETSC_EXTERN PetscErrorCode TSBasicSymplecticInitializePackage(void);
1272 PETSC_EXTERN PetscErrorCode TSBasicSymplecticFinalizePackage(void);
1273 PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterDestroy(void);
1274 
1275 /*J
1276   TSDISCGRAD - The Discrete Gradient integrator is a timestepper for Hamiltonian systems designed to conserve the first integral (energy),
1277   but also has the property for some systems of monotonicity in a functional.
1278 
1279   Level: beginner
1280 
1281 .seealso: [](ch_ts), `TS`, TSDiscGradSetFormulation()`, `TSDiscGradGetFormulation()`
1282 J*/
1283 PETSC_EXTERN PetscErrorCode TSDiscGradSetFormulation(TS, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (*)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (*)(TS, PetscReal, Vec, Vec, void *), void *);
1284 PETSC_EXTERN PetscErrorCode TSDiscGradIsGonzalez(TS, PetscBool *);
1285 PETSC_EXTERN PetscErrorCode TSDiscGradUseGonzalez(TS, PetscBool);
1286 
1287 /*
1288        PETSc interface to Sundials
1289 */
1290 #ifdef PETSC_HAVE_SUNDIALS2
1291 typedef enum {
1292   SUNDIALS_ADAMS = 1,
1293   SUNDIALS_BDF   = 2
1294 } TSSundialsLmmType;
1295 PETSC_EXTERN const char *const TSSundialsLmmTypes[];
1296 typedef enum {
1297   SUNDIALS_MODIFIED_GS  = 1,
1298   SUNDIALS_CLASSICAL_GS = 2
1299 } TSSundialsGramSchmidtType;
1300 PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
1301 PETSC_EXTERN PetscErrorCode    TSSundialsSetType(TS, TSSundialsLmmType);
1302 PETSC_EXTERN PetscErrorCode    TSSundialsGetPC(TS, PC *);
1303 PETSC_EXTERN PetscErrorCode    TSSundialsSetTolerance(TS, PetscReal, PetscReal);
1304 PETSC_EXTERN PetscErrorCode    TSSundialsSetMinTimeStep(TS, PetscReal);
1305 PETSC_EXTERN PetscErrorCode    TSSundialsSetMaxTimeStep(TS, PetscReal);
1306 PETSC_EXTERN PetscErrorCode    TSSundialsGetIterations(TS, PetscInt *, PetscInt *);
1307 PETSC_EXTERN PetscErrorCode    TSSundialsSetGramSchmidtType(TS, TSSundialsGramSchmidtType);
1308 PETSC_EXTERN PetscErrorCode    TSSundialsSetGMRESRestart(TS, PetscInt);
1309 PETSC_EXTERN PetscErrorCode    TSSundialsSetLinearTolerance(TS, PetscReal);
1310 PETSC_EXTERN PetscErrorCode    TSSundialsMonitorInternalSteps(TS, PetscBool);
1311 PETSC_EXTERN PetscErrorCode    TSSundialsSetMaxl(TS, PetscInt);
1312 PETSC_EXTERN PetscErrorCode    TSSundialsSetMaxord(TS, PetscInt);
1313 PETSC_EXTERN PetscErrorCode    TSSundialsSetUseDense(TS, PetscBool);
1314 #endif
1315 
1316 PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS, PetscReal);
1317 PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS, PetscReal *);
1318 PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS, PetscBool *);
1319 PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS, PetscBool);
1320 
1321 PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS, PetscReal);
1322 PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS, PetscReal, PetscReal, PetscReal);
1323 PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS, PetscReal *, PetscReal *, PetscReal *);
1324 
1325 PETSC_EXTERN PetscErrorCode TSAlpha2SetRadius(TS, PetscReal);
1326 PETSC_EXTERN PetscErrorCode TSAlpha2SetParams(TS, PetscReal, PetscReal, PetscReal, PetscReal);
1327 PETSC_EXTERN PetscErrorCode TSAlpha2GetParams(TS, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
1328 
1329 PETSC_EXTERN PetscErrorCode TSSetDM(TS, DM);
1330 PETSC_EXTERN PetscErrorCode TSGetDM(TS, DM *);
1331 
1332 PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES, Vec, Vec, void *);
1333 PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES, Vec, Mat, Mat, void *);
1334 
1335 PETSC_EXTERN PetscErrorCode TSRHSJacobianTest(TS, PetscBool *);
1336 PETSC_EXTERN PetscErrorCode TSRHSJacobianTestTranspose(TS, PetscBool *);
1337 
1338 PETSC_EXTERN PetscErrorCode TSGetComputeInitialCondition(TS, PetscErrorCode (**)(TS, Vec));
1339 PETSC_EXTERN PetscErrorCode TSSetComputeInitialCondition(TS, PetscErrorCode (*)(TS, Vec));
1340 PETSC_EXTERN PetscErrorCode TSComputeInitialCondition(TS, Vec);
1341 PETSC_EXTERN PetscErrorCode TSGetComputeExactError(TS, PetscErrorCode (**)(TS, Vec, Vec));
1342 PETSC_EXTERN PetscErrorCode TSSetComputeExactError(TS, PetscErrorCode (*)(TS, Vec, Vec));
1343 PETSC_EXTERN PetscErrorCode TSComputeExactError(TS, Vec, Vec);
1344 PETSC_EXTERN PetscErrorCode PetscConvEstUseTS(PetscConvEst, PetscBool);
1345 
1346 PETSC_EXTERN PetscErrorCode TSSetMatStructure(TS, MatStructure);
1347 #endif
1348