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