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