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