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