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