Actual source code: petscts.h

  1: /*
  2:    User interface for the timestepping package. This package
  3:    is for use in solving time-dependent PDEs.
  4: */
  5: #pragma once

  7: #include <petscsnes.h>
  8: #include <petscconvest.h>

 10: /*I <petscts.h> I*/

 12: /* SUBMANSEC = TS */

 14: /*S
 15:    TS - Abstract PETSc object that manages all time-steppers (ODE integrators)

 17:    Level: beginner

 19: .seealso: [](integrator_table), [](ch_ts), `TSCreate()`, `TSSetType()`, `TSType`, `SNES`, `KSP`, `PC`, `TSDestroy()`
 20: S*/
 21: typedef struct _p_TS *TS;

 23: /*J
 24:    TSType - String with the name of a PETSc `TS` method.

 26:    Level: beginner

 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"

 56: /*E
 57:    TSProblemType - Determines the type of problem this `TS` object is to be used to solve

 59:    Values:
 60:  + `TS_LINEAR`    - a linear ODE or DAE
 61:  - `TS_NONLINEAR` - a nonlinear ODE or DAE

 63:    Level: beginner

 65: .seealso: [](ch_ts), `TS`, `TSCreate()`
 66: E*/
 67: typedef enum {
 68:   TS_LINEAR,
 69:   TS_NONLINEAR
 70: } TSProblemType;

 72: /*E
 73:    TSEquationType - type of `TS` problem that is solved

 75:    Values:
 76: +  `TS_EQ_UNSPECIFIED` - (default)
 77: .  `TS_EQ_EXPLICIT`    - {ODE and DAE index 1, 2, 3, HI} F(t,U,U_t) := M(t) U_t - G(U,t) = 0
 78: -  `TS_EQ_IMPLICIT`    - {ODE and DAE index 1, 2, 3, HI} F(t,U,U_t) = 0

 80:    Level: beginner

 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;

101: /*E
102:    TSConvergedReason - reason a `TS` method has converged (integrated to the requested time) or not

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

117:    Level: beginner

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;

136: /*MC
137:    TS_CONVERGED_ITERATING - this only occurs if `TSGetConvergedReason()` is called during the `TSSolve()`

139:    Level: beginner

141: .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`
142: M*/

144: /*MC
145:    TS_CONVERGED_TIME - the final time was reached

147:    Level: beginner

149: .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxTime()`, `TSGetMaxTime()`, `TSGetSolveTime()`
150: M*/

152: /*MC
153:    TS_CONVERGED_ITS - the maximum number of iterations (time-steps) was reached prior to the final time

155:    Level: beginner

157: .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxSteps()`, `TSGetMaxSteps()`
158: M*/

160: /*MC
161:    TS_CONVERGED_USER - user requested termination

163:    Level: beginner

165: .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`
166: M*/

168: /*MC
169:    TS_CONVERGED_EVENT - user requested termination on event detection

171:    Level: beginner

173: .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`
174: M*/

176: /*MC
177:    TS_CONVERGED_PSEUDO_FRTOL - stops when function norm decreased by a set amount, used only for `TSPSEUDO`

179:    Options Database Key:
180: .   -ts_pseudo_frtol <rtol> - use specified rtol

182:    Level: beginner

184: .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`, `TS_CONVERGED_PSEUDO_FATOL`
185: M*/

187: /*MC
188:    TS_CONVERGED_PSEUDO_FATOL - stops when function norm decreases below a set amount, used only for `TSPSEUDO`

190:    Options Database Key:
191: .   -ts_pseudo_fatol <atol> - use specified atol

193:    Level: beginner

195: .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`, `TS_CONVERGED_PSEUDO_FRTOL`
196: M*/

198: /*MC
199:    TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed

201:    Level: beginner

203:    Note:
204:    See `TSSetMaxSNESFailures()` for how to allow more nonlinear solver failures.

206: .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSGetSNES()`, `SNESGetConvergedReason()`, `TSSetMaxSNESFailures()`
207: M*/

209: /*MC
210:    TS_DIVERGED_STEP_REJECTED - too many steps were rejected

212:    Level: beginner

214:    Notes:
215:    See `TSSetMaxStepRejections()` for how to allow more step rejections.

217: .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxStepRejections()`
218: M*/

220: /*E
221:    TSExactFinalTimeOption - option for handling of final time step

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

228:    Level: beginner

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[];

240: /* Logging support */
241: PETSC_EXTERN PetscClassId TS_CLASSID;
242: PETSC_EXTERN PetscClassId DMTS_CLASSID;
243: PETSC_EXTERN PetscClassId TSADAPT_CLASSID;

245: PETSC_EXTERN PetscErrorCode TSInitializePackage(void);
246: PETSC_EXTERN PetscErrorCode TSFinalizePackage(void);

248: PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm, TS *);
249: PETSC_EXTERN PetscErrorCode TSClone(TS, TS *);
250: PETSC_EXTERN PetscErrorCode TSDestroy(TS *);

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);

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);

266: PETSC_EXTERN PetscErrorCode TSSetSolution(TS, Vec);
267: PETSC_EXTERN PetscErrorCode TSGetSolution(TS, Vec *);

269: PETSC_EXTERN PetscErrorCode TS2SetSolution(TS, Vec, Vec);
270: PETSC_EXTERN PetscErrorCode TS2GetSolution(TS, Vec *, Vec *);

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);

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);

298: /*S
299:    TSTrajectory - Abstract PETSc object that stores the trajectory (solution of ODE/DAE at each time step)

301:    Level: advanced

303: .seealso: [](ch_ts), `TS`, `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectorySetType()`, `TSTrajectoryDestroy()`, `TSTrajectoryReset()`
304: S*/
305: typedef struct _p_TSTrajectory *TSTrajectory;

307: /*J
308:    TSTrajectoryType - String with the name of a PETSc `TS` trajectory storage method

310:    Level: intermediate

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"

320: PETSC_EXTERN PetscFunctionList TSTrajectoryList;
321: PETSC_EXTERN PetscClassId      TSTRAJECTORY_CLASSID;
322: PETSC_EXTERN PetscBool         TSTrajectoryRegisterAllCalled;

324: PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS);
325: PETSC_EXTERN PetscErrorCode TSResetTrajectory(TS);
326: PETSC_EXTERN PetscErrorCode TSRemoveTrajectory(TS);

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 *);

355: typedef enum {
356:   TJ_REVOLVE,
357:   TJ_CAMS,
358:   TJ_PETSC
359: } TSTrajectoryMemoryType;
360: PETSC_EXTERN const char *const TSTrajectoryMemoryTypes[];

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);

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 *);

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 *));

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);

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);

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 *[]);

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 **);

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 *);

423: PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
424: PETSC_EXTERN PetscErrorCode TSMonitorExtreme(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);

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 *);

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 *);

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 *);

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: PETSC_EXTERN PetscErrorCode TSGetStepRollBack(TS, PetscBool *);

461: PETSC_EXTERN PetscErrorCode TSGetStages(TS, PetscInt *, Vec *[]);

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);

471: /*S
472:   TSRHSFunctionFn - A prototype of a `TS` right-hand-side evaluation function that would be passed to `TSSetRHSFunction()`

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

481:   Level: beginner

483:   Note:
484:   The deprecated `TSRHSFunction` still works as a replacement for `TSRHSFunctionFn` *.

486: .seealso: [](ch_ts), `TS`, `TSSetRHSFunction()`, `DMTSSetRHSFunction()`, `TSIFunctionFn`,
487: `TSIJacobianFn`, `TSRHSJacobianFn`
488: S*/
489: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSRHSFunctionFn)(TS ts, PetscReal t, Vec u, Vec F, void *ctx);

491: PETSC_EXTERN_TYPEDEF typedef TSRHSFunctionFn *TSRHSFunction;

493: /*S
494:   TSRHSJacobianFn - A prototype of a `TS` right-hand-side Jacobian evaluation function that would be passed to `TSSetRHSJacobian()`

496:   Calling Sequence:
497: + ts   - the `TS` context obtained from `TSCreate()`
498: . t    - current time
499: . u    - input vector
500: . Amat - (approximate) Jacobian matrix
501: . Pmat - matrix from which preconditioner is to be constructed (usually the same as `Amat`)
502: - ctx  - [optional] user-defined context for matrix evaluation routine

504:   Level: beginner

506:   Note:
507:   The deprecated `TSRHSJacobian` still works as a replacement for `TSRHSJacobianFn` *.

509: .seealso: [](ch_ts), `TS`, `TSSetRHSJacobian()`, `DMTSSetRHSJacobian()`, `TSRHSFunctionFn`,
510: `TSIFunctionFn`, `TSIJacobianFn`
511: S*/
512: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSRHSJacobianFn)(TS ts, PetscReal t, Vec u, Mat Amat, Mat Pmat, void *ctx);

514: PETSC_EXTERN_TYPEDEF typedef TSRHSJacobianFn *TSRHSJacobian;

516: /*S
517:   TSRHSJacobianPFn - A prototype of a function that computes the Jacobian of G w.r.t. the parameters P where
518:   U_t = G(U,P,t), as well as the location to store the matrix that would be passed to `TSSetRHSJacobianP()`

520:   Calling Sequence:
521: + ts  - the `TS` context
522: . t   - current timestep
523: . U   - input vector (current ODE solution)
524: . A   - output matrix
525: - ctx - [optional] user-defined function context

527:   Level: beginner

529:   Note:
530:   The deprecated `TSRHSJacobianP` still works as a replacement for `TSRHSJacobianPFn` *.

532: .seealso: [](ch_ts), `TS`, `TSSetRHSJacobianP()`, `TSGetRHSJacobianP()`
533: S*/
534: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSRHSJacobianPFn)(TS ts, PetscReal t, Vec U, Mat A, void *ctx);

536: PETSC_EXTERN_TYPEDEF typedef TSRHSJacobianPFn *TSRHSJacobianP;

538: PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS, Vec, TSRHSFunctionFn *, void *);
539: PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS, Vec *, TSRHSFunctionFn **, void **);
540: PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS, Mat, Mat, TSRHSJacobianFn *, void *);
541: PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS, Mat *, Mat *, TSRHSJacobianFn **, void **);
542: PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS, PetscBool);

544: /*S
545:   TSSolutionFn - A prototype of a `TS` solution evaluation function that would be passed to `TSSetSolutionFunction()`

547:   Calling Sequence:
548: + ts  - timestep context
549: . t   - current time
550: . u   - output vector
551: - ctx - [optional] user-defined function context

553:   Level: advanced

555:   Note:
556:   The deprecated `TSSolutionFunction` still works as a replacement for `TSSolutionFn` *.

558: .seealso: [](ch_ts), `TS`, `TSSetSolutionFunction()`, `DMTSSetSolutionFunction()`
559: S*/
560: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSSolutionFn)(TS ts, PetscReal t, Vec u, void *ctx);

562: PETSC_EXTERN_TYPEDEF typedef TSSolutionFn *TSSolutionFunction;

564: PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS, TSSolutionFn *, void *);

566: /*S
567:   TSForcingFn - A prototype of a `TS` forcing function evaluation function that would be passed to `TSSetForcingFunction()`

569:   Calling Sequence:
570: + ts  - timestep context
571: . t   - current time
572: . f   - output vector
573: - ctx - [optional] user-defined function context

575:   Level: advanced

577:   Note:
578:   The deprecated `TSForcingFunction` still works as a replacement for `TSForcingFn` *.

580: .seealso: [](ch_ts), `TS`, `TSSetForcingFunction()`, `DMTSSetForcingFunction()`
581: S*/
582: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSForcingFn)(TS ts, PetscReal t, Vec f, void *ctx);

584: PETSC_EXTERN_TYPEDEF typedef TSForcingFn *TSForcingFunction;

586: PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS, TSForcingFn *, void *);

588: /*S
589:   TSIFunctionFn - A prototype of a `TS` implicit function evaluation function that would be passed to `TSSetIFunction()

591:   Calling Sequence:
592: + ts  - the `TS` context obtained from `TSCreate()`
593: . t   - time at step/stage being solved
594: . U   - state vector
595: . U_t - time derivative of state vector
596: . F   - function vector
597: - ctx - [optional] user-defined context for function

599:   Level: beginner

601:   Note:
602:   The deprecated `TSIFunction` still works as a replacement for `TSIFunctionFn` *.

604: .seealso: [](ch_ts), `TS`, `TSSetIFunction()`, `DMTSSetIFunction()`, `TSIJacobianFn`, `TSRHSFunctionFn`, `TSRHSJacobianFn`
605: S*/
606: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSIFunctionFn)(TS ts, PetscReal t, Vec U, Vec U_t, Vec F, void *ctx);

608: PETSC_EXTERN_TYPEDEF typedef TSIFunctionFn *TSIFunction;

610: /*S
611:   TSIJacobianFn - A prototype of a `TS` Jacobian evaluation function that would be passed to `TSSetIJacobian()`

613:   Calling Sequence:
614: + ts   - the `TS` context obtained from `TSCreate()`
615: . t    - time at step/stage being solved
616: . U    - state vector
617: . U_t  - time derivative of state vector
618: . a    - shift
619: . Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
620: . Pmat - matrix used for constructing preconditioner, usually the same as `Amat`
621: - ctx  - [optional] user-defined context for Jacobian evaluation routine

623:   Level: beginner

625:   Note:
626:   The deprecated `TSIJacobian` still works as a replacement for `TSIJacobianFn` *.

628: .seealso: [](ch_ts), `TSSetIJacobian()`, `DMTSSetIJacobian()`, `TSIFunctionFn`, `TSRHSFunctionFn`, `TSRHSJacobianFn`
629: S*/
630: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSIJacobianFn)(TS ts, PetscReal t, Vec U, Vec U_t, PetscReal a, Mat Amat, Mat Pmat, void *ctx);

632: PETSC_EXTERN_TYPEDEF typedef TSIJacobianFn *TSIJacobian;

634: PETSC_EXTERN PetscErrorCode TSSetIFunction(TS, Vec, TSIFunctionFn *, void *);
635: PETSC_EXTERN PetscErrorCode TSGetIFunction(TS, Vec *, TSIFunctionFn **, void **);
636: PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS, Mat, Mat, TSIJacobianFn *, void *);
637: PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS, Mat *, Mat *, TSIJacobianFn **, void **);

639: /*S
640:   TSI2FunctionFn - A prototype of a `TS` implicit function evaluation function for 2nd order systems that would be passed to `TSSetI2Function()`

642:   Calling Sequence:
643: + ts   - the `TS` context obtained from `TSCreate()`
644: . t    - time at step/stage being solved
645: . U    - state vector
646: . U_t  - time derivative of state vector
647: . U_tt - second time derivative of state vector
648: . F    - function vector
649: - ctx  - [optional] user-defined context for matrix evaluation routine (may be `NULL`)

651:   Level: advanced

653:   Note:
654:   The deprecated `TSI2Function` still works as a replacement for `TSI2FunctionFn` *.

656: .seealso: [](ch_ts), `TS`, `TSSetI2Function()`, `DMTSSetI2Function()`, `TSIFunctionFn`
657: S*/
658: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSI2FunctionFn)(TS ts, PetscReal t, Vec U, Vec U_t, Vec U_tt, Vec F, void *ctx);

660: PETSC_EXTERN_TYPEDEF typedef TSI2FunctionFn *TSI2Function;

662: /*S
663:   TSI2JacobianFn - A prototype of a `TS` implicit Jacobian evaluation function for 2nd order systems that would be passed to `TSSetI2Jacobian()`

665:   Calling Sequence:
666: + ts   - the `TS` context obtained from `TSCreate()`
667: . t    - time at step/stage being solved
668: . U    - state vector
669: . U_t  - time derivative of state vector
670: . U_tt - second time derivative of state vector
671: . v    - shift for U_t
672: . a    - shift for U_tt
673: . 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
674: . jac  - matrix from which to construct the preconditioner, may be same as `J`
675: - ctx  - [optional] user-defined context for matrix evaluation routine

677:   Level: advanced

679:   Note:
680:   The deprecated `TSI2Jacobian` still works as a replacement for `TSI2JacobianFn` *.

682: .seealso: [](ch_ts), `TS`, `TSSetI2Jacobian()`, `DMTSSetI2Jacobian()`, `TSIFunctionFn`, `TSIJacobianFn`, `TSRHSFunctionFn`, `TSRHSJacobianFn`
683: S*/
684: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSI2JacobianFn)(TS ts, PetscReal t, Vec U, Vec U_t, Vec U_tt, PetscReal v, PetscReal a, Mat J, Mat Jac, void *ctx);

686: PETSC_EXTERN_TYPEDEF typedef TSI2JacobianFn *TSI2Jacobian;

688: PETSC_EXTERN PetscErrorCode TSSetI2Function(TS, Vec, TSI2FunctionFn *, void *);
689: PETSC_EXTERN PetscErrorCode TSGetI2Function(TS, Vec *, TSI2FunctionFn **, void **);
690: PETSC_EXTERN PetscErrorCode TSSetI2Jacobian(TS, Mat, Mat, TSI2JacobianFn *, void *);
691: PETSC_EXTERN PetscErrorCode TSGetI2Jacobian(TS, Mat *, Mat *, TSI2JacobianFn **, void **);

693: PETSC_EXTERN PetscErrorCode TSRHSSplitSetIS(TS, const char[], IS);
694: PETSC_EXTERN PetscErrorCode TSRHSSplitGetIS(TS, const char[], IS *);
695: PETSC_EXTERN PetscErrorCode TSRHSSplitSetRHSFunction(TS, const char[], Vec, TSRHSFunctionFn *, void *);
696: PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTS(TS, const char[], TS *);
697: PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTSs(TS, PetscInt *, TS *[]);
698: PETSC_EXTERN PetscErrorCode TSSetUseSplitRHSFunction(TS, PetscBool);
699: PETSC_EXTERN PetscErrorCode TSGetUseSplitRHSFunction(TS, PetscBool *);

701: PETSC_EXTERN TSRHSFunctionFn TSComputeRHSFunctionLinear;
702: PETSC_EXTERN TSRHSJacobianFn TSComputeRHSJacobianConstant;
703: PETSC_EXTERN PetscErrorCode  TSComputeIFunctionLinear(TS, PetscReal, Vec, Vec, Vec, void *);
704: PETSC_EXTERN PetscErrorCode  TSComputeIJacobianConstant(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
705: PETSC_EXTERN PetscErrorCode  TSComputeSolutionFunction(TS, PetscReal, Vec);
706: PETSC_EXTERN PetscErrorCode  TSComputeForcingFunction(TS, PetscReal, Vec);
707: PETSC_EXTERN PetscErrorCode  TSComputeIJacobianDefaultColor(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
708: PETSC_EXTERN PetscErrorCode  TSPruneIJacobianColor(TS, Mat, Mat);

710: PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
711: PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS, PetscReal));
712: PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS, PetscReal, PetscInt, Vec *));
713: PETSC_EXTERN PetscErrorCode TSSetPostEvaluate(TS, PetscErrorCode (*)(TS));
714: PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
715: PETSC_EXTERN PetscErrorCode TSSetResize(TS, PetscBool, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscBool *, void *), PetscErrorCode (*)(TS, PetscInt, Vec[], Vec[], void *), void *);
716: PETSC_EXTERN PetscErrorCode TSPreStep(TS);
717: PETSC_EXTERN PetscErrorCode TSPreStage(TS, PetscReal);
718: PETSC_EXTERN PetscErrorCode TSPostStage(TS, PetscReal, PetscInt, Vec *);
719: PETSC_EXTERN PetscErrorCode TSPostEvaluate(TS);
720: PETSC_EXTERN PetscErrorCode TSPostStep(TS);
721: PETSC_EXTERN PetscErrorCode TSResize(TS);
722: PETSC_EXTERN PetscErrorCode TSResizeRetrieveVec(TS, const char *, Vec *);
723: PETSC_EXTERN PetscErrorCode TSResizeRegisterVec(TS, const char *, Vec);

725: PETSC_EXTERN PetscErrorCode TSInterpolate(TS, PetscReal, Vec);
726: PETSC_EXTERN PetscErrorCode TSSetTolerances(TS, PetscReal, Vec, PetscReal, Vec);
727: PETSC_EXTERN PetscErrorCode TSGetTolerances(TS, PetscReal *, Vec *, PetscReal *, Vec *);
728: PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *);
729: PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm(TS, Vec, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *);
730: PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS, PetscReal);
731: PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS, PetscReal *);
732: PETSC_EXTERN PetscErrorCode TSSetFunctionDomainError(TS, PetscErrorCode (*)(TS, PetscReal, Vec, PetscBool *));
733: PETSC_EXTERN PetscErrorCode TSFunctionDomainError(TS, PetscReal, Vec, PetscBool *);

735: PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS, PetscErrorCode (*)(TS, PetscReal *, void *), void *);
736: PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS, PetscReal *, void *);
737: PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS, PetscReal *);
738: PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS, PetscReal);
739: PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS, PetscErrorCode (*)(TS, Vec, void *, PetscReal *, PetscBool *), void *);
740: PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS, Vec, void *, PetscReal *, PetscBool *);
741: PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS, Vec, PetscReal *, PetscBool *);
742: PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS, PetscReal);
743: PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);

745: PETSC_EXTERN PetscErrorCode TSPythonSetType(TS, const char[]);
746: PETSC_EXTERN PetscErrorCode TSPythonGetType(TS, const char *[]);

748: PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS, PetscReal, Vec, Vec);
749: PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS, PetscReal, Vec, Mat, Mat);
750: PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS, PetscReal, Vec, Vec, Vec, PetscBool);
751: PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscBool);
752: PETSC_EXTERN PetscErrorCode TSComputeI2Function(TS, PetscReal, Vec, Vec, Vec, Vec);
753: PETSC_EXTERN PetscErrorCode TSComputeI2Jacobian(TS, PetscReal, Vec, Vec, Vec, PetscReal, PetscReal, Mat, Mat);
754: PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS, PetscReal, PetscReal, PetscReal *, PetscReal *);

756: PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS, Vec, Vec);

758: PETSC_EXTERN PetscErrorCode DMTSSetBoundaryLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *);
759: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM, TSRHSFunctionFn *, void *);
760: PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM, TSRHSFunctionFn **, void **);
761: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionContextDestroy(DM, PetscErrorCode (*)(void *));
762: PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM, TSRHSJacobianFn *, void *);
763: PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM, TSRHSJacobianFn **, void **);
764: PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobianContextDestroy(DM, PetscErrorCode (*)(void *));
765: PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM, TSIFunctionFn *, void *);
766: PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM, TSIFunctionFn **, void **);
767: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionContextDestroy(DM, PetscErrorCode (*)(void *));
768: PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM, TSIJacobianFn *, void *);
769: PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM, TSIJacobianFn **, void **);
770: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianContextDestroy(DM, PetscErrorCode (*)(void *));
771: PETSC_EXTERN PetscErrorCode DMTSSetI2Function(DM, TSI2FunctionFn *, void *);
772: PETSC_EXTERN PetscErrorCode DMTSGetI2Function(DM, TSI2FunctionFn **, void **);
773: PETSC_EXTERN PetscErrorCode DMTSSetI2FunctionContextDestroy(DM, PetscErrorCode (*)(void *));
774: PETSC_EXTERN PetscErrorCode DMTSSetI2Jacobian(DM, TSI2JacobianFn *, void *);
775: PETSC_EXTERN PetscErrorCode DMTSGetI2Jacobian(DM, TSI2JacobianFn **, void **);
776: PETSC_EXTERN PetscErrorCode DMTSSetI2JacobianContextDestroy(DM, PetscErrorCode (*)(void *));

778: /*S
779:   TSTransientVariableFn - A prototype of a function to transform from state to transient variables that would be passed to `TSSetTransientVariable()`

781:   Calling Sequence:
782: + ts  - timestep context
783: . p   - input vector (primitive form)
784: . c   - output vector, transient variables (conservative form)
785: - ctx - [optional] user-defined function context

787:   Level: advanced

789:   Note:
790:   The deprecated `TSTransientVariable` still works as a replacement for `TSTransientVariableFn` *.

792: .seealso: [](ch_ts), `TS`, `TSSetTransientVariable()`, `DMTSSetTransientVariable()`
793: S*/
794: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSTransientVariableFn)(TS ts, Vec p, Vec c, void *ctx);

796: PETSC_EXTERN_TYPEDEF typedef TSTransientVariableFn *TSTransientVariable;

798: PETSC_EXTERN PetscErrorCode TSSetTransientVariable(TS, TSTransientVariableFn *, void *);
799: PETSC_EXTERN PetscErrorCode DMTSSetTransientVariable(DM, TSTransientVariableFn *, void *);
800: PETSC_EXTERN PetscErrorCode DMTSGetTransientVariable(DM, TSTransientVariableFn **, void *);
801: PETSC_EXTERN PetscErrorCode TSComputeTransientVariable(TS, Vec, Vec);
802: PETSC_EXTERN PetscErrorCode TSHasTransientVariable(TS, PetscBool *);

804: PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM, TSSolutionFn *, void *);
805: PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM, TSSolutionFn **, void **);
806: PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM, TSForcingFn *, void *);
807: PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM, TSForcingFn **, void **);
808: PETSC_EXTERN PetscErrorCode DMTSCheckResidual(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscReal *);
809: PETSC_EXTERN PetscErrorCode DMTSCheckJacobian(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscBool *, PetscReal *);
810: PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec);

812: PETSC_EXTERN PetscErrorCode DMTSGetIFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, Vec, void *), void **);
813: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, Vec, void *), void *);
814: PETSC_EXTERN PetscErrorCode DMTSGetIJacobianLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void **);
815: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void *);
816: PETSC_EXTERN PetscErrorCode DMTSGetRHSFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, void *), void **);
817: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *);
818: PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrix(DM);
819: PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrixLumped(DM);
820: PETSC_EXTERN PetscErrorCode DMTSDestroyRHSMassMatrix(DM);

822: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM, PetscErrorCode (*)(void *, PetscViewer), PetscErrorCode (*)(void **, PetscViewer));
823: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM, PetscErrorCode (*)(void *, PetscViewer), PetscErrorCode (*)(void **, PetscViewer));

825: /*S
826:   DMDATSRHSFunctionLocalFn - A prototype of a local `TS` right-hand side residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetRHSFunctionLocal()`

828:   Calling Sequence:
829: + info - defines the subdomain to evaluate the residual on
830: . t    - time at which to evaluate residual
831: . x    - array of local state information
832: . f    - output array of local residual information
833: - ctx  - optional user context

835:   Level: beginner

837:   Note:
838:   The deprecated `DMDATSRHSFunctionLocal` still works as a replacement for `DMDATSRHSFunctionLocalFn` *.

840: .seealso: `DMDA`, `DMDATSSetRHSFunctionLocal()`, `TSRHSFunctionFn`, `DMDATSRHSJacobianLocalFn`, `DMDATSIJacobianLocalFn`, `DMDATSIFunctionLocalFn`
841: S*/
842: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDATSRHSFunctionLocalFn)(DMDALocalInfo *info, PetscReal t, void *x, void *f, void *ctx);

844: PETSC_EXTERN_TYPEDEF typedef DMDATSRHSFunctionLocalFn *DMDATSRHSFunctionLocal;

846: /*S
847:   DMDATSRHSJacobianLocalFn - A prototype of a local residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetRHSJacobianLocal()`

849:   Calling Sequence:
850: + info - defines the subdomain to evaluate the residual on
851: . t    - time at which to evaluate residual
852: . x    - array of local state information
853: . J    - Jacobian matrix
854: . B    - matrix from which to construct the preconditioner; often same as `J`
855: - ctx  - optional context

857:   Level: beginner

859:   Note:
860:   The deprecated `DMDATSRHSJacobianLocal` still works as a replacement for `DMDATSRHSJacobianLocalFn` *.

862: .seealso: `DMDA`, `DMDATSSetRHSJacobianLocal()`, `TSRHSJacobianFn`, `DMDATSRHSFunctionLocalFn`, `DMDATSIJacobianLocalFn`, `DMDATSIFunctionLocalFn`
863: S*/
864: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDATSRHSJacobianLocalFn)(DMDALocalInfo *info, PetscReal t, void *x, Mat J, Mat B, void *ctx);

866: PETSC_EXTERN_TYPEDEF typedef DMDATSRHSJacobianLocalFn *DMDATSRHSJacobianLocal;

868: /*S
869:   DMDATSIFunctionLocalFn - A prototype of a local residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetIFunctionLocal()`

871:   Calling Sequence:
872: + info  - defines the subdomain to evaluate the residual on
873: . t     - time at which to evaluate residual
874: . x     - array of local state information
875: . xdot  - array of local time derivative information
876: . imode - output array of local function evaluation information
877: - ctx   - optional context

879:   Level: beginner

881:   Note:
882:   The deprecated `DMDATSIFunctionLocal` still works as a replacement for `DMDATSIFunctionLocalFn` *.

884: .seealso: `DMDA`, `DMDATSSetIFunctionLocal()`, `DMDATSIJacobianLocalFn`, `TSIFunctionFn`
885: S*/
886: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDATSIFunctionLocalFn)(DMDALocalInfo *info, PetscReal t, void *x, void *xdot, void *imode, void *ctx);

888: PETSC_EXTERN_TYPEDEF typedef DMDATSIFunctionLocalFn *DMDATSIFunctionLocal;

890: /*S
891:   DMDATSIJacobianLocalFn - A prototype of a local residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetIJacobianLocal()`

893:   Calling Sequence:
894: + info  - defines the subdomain to evaluate the residual on
895: . t     - time at which to evaluate the jacobian
896: . x     - array of local state information
897: . xdot  - time derivative at this state
898: . shift - see `TSSetIJacobian()` for the meaning of this parameter
899: . J     - Jacobian matrix
900: . B     - matrix from which to construct the preconditioner; often same as `J`
901: - ctx   - optional context

903:   Level: beginner

905:   Note:
906:   The deprecated `DMDATSIJacobianLocal` still works as a replacement for `DMDATSIJacobianLocalFn` *.

908: .seealso: `DMDA` `DMDATSSetIJacobianLocal()`, `TSIJacobianFn`, `DMDATSIFunctionLocalFn`, `DMDATSRHSFunctionLocalFn`,  `DMDATSRHSJacobianlocal()`
909: S*/
910: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDATSIJacobianLocalFn)(DMDALocalInfo *info, PetscReal t, void *x, void *xdot, PetscReal shift, Mat J, Mat B, void *ctx);

912: PETSC_EXTERN_TYPEDEF typedef DMDATSIJacobianLocalFn *DMDATSIJacobianLocal;

914: PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM, InsertMode, DMDATSRHSFunctionLocalFn *, void *);
915: PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM, DMDATSRHSJacobianLocalFn *, void *);
916: PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM, InsertMode, DMDATSIFunctionLocalFn *, void *);
917: PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM, DMDATSIJacobianLocalFn *, void *);

919: typedef struct _n_TSMonitorLGCtx *TSMonitorLGCtx;
920: typedef struct {
921:   Vec            ray;
922:   VecScatter     scatter;
923:   PetscViewer    viewer;
924:   TSMonitorLGCtx lgctx;
925: } TSMonitorDMDARayCtx;
926: PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void **);
927: PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS, PetscInt, PetscReal, Vec, void *);
928: PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS, PetscInt, PetscReal, Vec, void *);

930: /* Dynamic creation and loading functions */
931: PETSC_EXTERN PetscFunctionList TSList;
932: PETSC_EXTERN PetscErrorCode    TSGetType(TS, TSType *);
933: PETSC_EXTERN PetscErrorCode    TSSetType(TS, TSType);
934: PETSC_EXTERN PetscErrorCode    TSRegister(const char[], PetscErrorCode (*)(TS));

936: PETSC_EXTERN PetscErrorCode TSGetSNES(TS, SNES *);
937: PETSC_EXTERN PetscErrorCode TSSetSNES(TS, SNES);
938: PETSC_EXTERN PetscErrorCode TSGetKSP(TS, KSP *);

940: PETSC_EXTERN PetscErrorCode TSView(TS, PetscViewer);
941: PETSC_EXTERN PetscErrorCode TSLoad(TS, PetscViewer);
942: PETSC_EXTERN PetscErrorCode TSViewFromOptions(TS, PetscObject, const char[]);
943: PETSC_EXTERN PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory, PetscObject, const char[]);

945: #define TS_FILE_CLASSID 1211225

947: PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS, void *);
948: PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS, void *);

950: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtx *);
951: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *);
952: PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS, PetscInt, PetscReal, Vec, void *);
953: PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS, PetscInt, PetscReal, Vec, void *);
954: PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS, const char *const *);
955: PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS, const char *const **);
956: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx, const char *const *);
957: PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS, const char *const *);
958: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx, const char *const *);
959: PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS, PetscErrorCode (*)(void *, Vec, Vec *), PetscErrorCode (*)(void *), void *);
960: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx, PetscErrorCode (*)(void *, Vec, Vec *), PetscErrorCode (*)(void *), void *);
961: PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS, PetscInt, PetscReal, Vec, void *);
962: PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS, PetscInt, PetscReal, Vec, void *);
963: PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS, PetscInt, PetscReal, Vec, void *);
964: PETSC_EXTERN PetscErrorCode TSMonitorError(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
965: PETSC_EXTERN PetscErrorCode TSDMSwarmMonitorMoments(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);

967: struct _n_TSMonitorLGCtxNetwork {
968:   PetscInt     nlg;
969:   PetscDrawLG *lg;
970:   PetscBool    semilogy;
971:   PetscInt     howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
972: };
973: typedef struct _n_TSMonitorLGCtxNetwork *TSMonitorLGCtxNetwork;
974: PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkDestroy(TSMonitorLGCtxNetwork *);
975: PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkCreate(TS, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtxNetwork *);
976: PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkSolution(TS, PetscInt, PetscReal, Vec, void *);

978: typedef struct _n_TSMonitorEnvelopeCtx *TSMonitorEnvelopeCtx;
979: PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeCtxCreate(TS, TSMonitorEnvelopeCtx *);
980: PETSC_EXTERN PetscErrorCode             TSMonitorEnvelope(TS, PetscInt, PetscReal, Vec, void *);
981: PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeGetBounds(TS, Vec *, Vec *);
982: PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *);

984: typedef struct _n_TSMonitorSPEigCtx *TSMonitorSPEigCtx;
985: PETSC_EXTERN PetscErrorCode          TSMonitorSPEigCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorSPEigCtx *);
986: PETSC_EXTERN PetscErrorCode          TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *);
987: PETSC_EXTERN PetscErrorCode          TSMonitorSPEig(TS, PetscInt, PetscReal, Vec, void *);

989: typedef struct _n_TSMonitorSPCtx *TSMonitorSPCtx;
990: PETSC_EXTERN PetscErrorCode       TSMonitorSPCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscBool, PetscBool, TSMonitorSPCtx *);
991: PETSC_EXTERN PetscErrorCode       TSMonitorSPCtxDestroy(TSMonitorSPCtx *);
992: PETSC_EXTERN PetscErrorCode       TSMonitorSPSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);

994: typedef struct _n_TSMonitorHGCtx *TSMonitorHGCtx;
995: PETSC_EXTERN PetscErrorCode       TSMonitorHGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscInt, PetscBool, TSMonitorHGCtx *);
996: PETSC_EXTERN PetscErrorCode       TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
997: PETSC_EXTERN PetscErrorCode       TSMonitorHGCtxDestroy(TSMonitorHGCtx *);
998: PETSC_EXTERN PetscErrorCode       TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);

1000: PETSC_EXTERN PetscErrorCode TSSetEventHandler(TS, PetscInt, PetscInt[], PetscBool[], PetscErrorCode (*)(TS, PetscReal, Vec, PetscReal[], void *), PetscErrorCode (*)(TS, PetscInt, PetscInt[], PetscReal, Vec, PetscBool, void *), void *);
1001: PETSC_EXTERN PetscErrorCode TSSetPostEventStep(TS, PetscReal);
1002: PETSC_EXTERN PetscErrorCode TSSetPostEventSecondStep(TS, PetscReal);
1003: PETSC_DEPRECATED_FUNCTION(3, 21, 0, "TSSetPostEventSecondStep()", ) static inline PetscErrorCode TSSetPostEventIntervalStep(TS ts, PetscReal dt)
1004: {
1005:   return TSSetPostEventSecondStep(ts, dt);
1006: }
1007: PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS, PetscReal, PetscReal[]);
1008: PETSC_EXTERN PetscErrorCode TSGetNumEvents(TS, PetscInt *);

1010: /*J
1011:    TSSSPType - string with the name of a `TSSSP` scheme.

1013:    Level: beginner

1015: .seealso: [](ch_ts), `TSSSPSetType()`, `TS`, `TSSSP`
1016: J*/
1017: typedef const char *TSSSPType;
1018: #define TSSSPRKS2  "rks2"
1019: #define TSSSPRKS3  "rks3"
1020: #define TSSSPRK104 "rk104"

1022: PETSC_EXTERN PetscErrorCode    TSSSPSetType(TS, TSSSPType);
1023: PETSC_EXTERN PetscErrorCode    TSSSPGetType(TS, TSSSPType *);
1024: PETSC_EXTERN PetscErrorCode    TSSSPSetNumStages(TS, PetscInt);
1025: PETSC_EXTERN PetscErrorCode    TSSSPGetNumStages(TS, PetscInt *);
1026: PETSC_EXTERN PetscErrorCode    TSSSPInitializePackage(void);
1027: PETSC_EXTERN PetscErrorCode    TSSSPFinalizePackage(void);
1028: PETSC_EXTERN PetscFunctionList TSSSPList;

1030: /*S
1031:    TSAdapt - Abstract object that manages time-step adaptivity

1033:    Level: beginner

1035: .seealso: [](ch_ts), `TS`, `TSGetAdapt()`, `TSAdaptCreate()`, `TSAdaptType`
1036: S*/
1037: typedef struct _p_TSAdapt *TSAdapt;

1039: /*J
1040:    TSAdaptType - String with the name of `TSAdapt` scheme.

1042:    Level: beginner

1044: .seealso: [](ch_ts), `TSGetAdapt()`, `TSAdaptSetType()`, `TS`, `TSAdapt`
1045: J*/
1046: typedef const char *TSAdaptType;
1047: #define TSADAPTNONE    "none"
1048: #define TSADAPTBASIC   "basic"
1049: #define TSADAPTDSP     "dsp"
1050: #define TSADAPTCFL     "cfl"
1051: #define TSADAPTGLEE    "glee"
1052: #define TSADAPTHISTORY "history"

1054: PETSC_EXTERN PetscErrorCode TSGetAdapt(TS, TSAdapt *);
1055: PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[], PetscErrorCode (*)(TSAdapt));
1056: PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void);
1057: PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
1058: PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm, TSAdapt *);
1059: PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt, TSAdaptType);
1060: PETSC_EXTERN PetscErrorCode TSAdaptGetType(TSAdapt, TSAdaptType *);
1061: PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt, const char[]);
1062: PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
1063: PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt, const char[], PetscInt, PetscInt, PetscReal, PetscReal, PetscBool);
1064: PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt, PetscInt *, const PetscInt **, const PetscInt **, const PetscReal **, const PetscReal **);
1065: PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt, TS, PetscReal, PetscInt *, PetscReal *, PetscBool *);
1066: PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt, TS, PetscReal, Vec, PetscBool *);
1067: PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt, PetscViewer);
1068: PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt, PetscViewer);
1069: PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(TSAdapt, PetscOptionItems *);
1070: PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt);
1071: PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt *);
1072: PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt, PetscBool);
1073: PETSC_EXTERN PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt, PetscBool);
1074: PETSC_EXTERN PetscErrorCode TSAdaptSetSafety(TSAdapt, PetscReal, PetscReal);
1075: PETSC_EXTERN PetscErrorCode TSAdaptGetSafety(TSAdapt, PetscReal *, PetscReal *);
1076: PETSC_EXTERN PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt, PetscReal);
1077: PETSC_EXTERN PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt, PetscReal *);
1078: PETSC_EXTERN PetscErrorCode TSAdaptSetClip(TSAdapt, PetscReal, PetscReal);
1079: PETSC_EXTERN PetscErrorCode TSAdaptGetClip(TSAdapt, PetscReal *, PetscReal *);
1080: PETSC_EXTERN PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt, PetscReal);
1081: PETSC_EXTERN PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt, PetscReal *);
1082: PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt, PetscReal, PetscReal);
1083: PETSC_EXTERN PetscErrorCode TSAdaptGetStepLimits(TSAdapt, PetscReal *, PetscReal *);
1084: PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt, PetscErrorCode (*)(TSAdapt, TS, PetscReal, Vec, PetscBool *));
1085: PETSC_EXTERN PetscErrorCode TSAdaptHistorySetHistory(TSAdapt, PetscInt n, PetscReal hist[], PetscBool);
1086: PETSC_EXTERN PetscErrorCode TSAdaptHistorySetTrajectory(TSAdapt, TSTrajectory, PetscBool);
1087: PETSC_EXTERN PetscErrorCode TSAdaptHistoryGetStep(TSAdapt, PetscInt, PetscReal *, PetscReal *);
1088: PETSC_EXTERN PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt, PetscInt);
1089: PETSC_EXTERN PetscErrorCode TSAdaptDSPSetFilter(TSAdapt, const char *);
1090: PETSC_EXTERN PetscErrorCode TSAdaptDSPSetPID(TSAdapt, PetscReal, PetscReal, PetscReal);

1092: /*S
1093:    TSGLLEAdapt - Abstract object that manages time-step adaptivity for `TSGLLE`

1095:    Level: beginner

1097:    Developer Note:
1098:    This functionality should be replaced by the `TSAdapt`.

1100: .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLEAdaptCreate()`, `TSGLLEAdaptType`
1101: S*/
1102: typedef struct _p_TSGLLEAdapt *TSGLLEAdapt;

1104: /*J
1105:    TSGLLEAdaptType - String with the name of `TSGLLEAdapt` scheme

1107:    Level: beginner

1109:    Developer Note:
1110:    This functionality should be replaced by the `TSAdaptType`.

1112: .seealso: [](ch_ts), `TSGLLEAdaptSetType()`, `TS`
1113: J*/
1114: typedef const char *TSGLLEAdaptType;
1115: #define TSGLLEADAPT_NONE "none"
1116: #define TSGLLEADAPT_SIZE "size"
1117: #define TSGLLEADAPT_BOTH "both"

1119: PETSC_EXTERN PetscErrorCode TSGLLEAdaptRegister(const char[], PetscErrorCode (*)(TSGLLEAdapt));
1120: PETSC_EXTERN PetscErrorCode TSGLLEAdaptInitializePackage(void);
1121: PETSC_EXTERN PetscErrorCode TSGLLEAdaptFinalizePackage(void);
1122: PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate(MPI_Comm, TSGLLEAdapt *);
1123: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt, TSGLLEAdaptType);
1124: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt, const char[]);
1125: PETSC_EXTERN PetscErrorCode TSGLLEAdaptChoose(TSGLLEAdapt, PetscInt, const PetscInt[], const PetscReal[], const PetscReal[], PetscInt, PetscReal, PetscReal, PetscInt *, PetscReal *, PetscBool *);
1126: PETSC_EXTERN PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt, PetscViewer);
1127: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetFromOptions(TSGLLEAdapt, PetscOptionItems *);
1128: PETSC_EXTERN PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt *);

1130: /*J
1131:    TSGLLEAcceptType - String with the name of `TSGLLEAccept` scheme

1133:    Level: beginner

1135: .seealso: [](ch_ts), `TSGLLESetAcceptType()`, `TS`, `TSGLLEAccept`
1136: J*/
1137: typedef const char *TSGLLEAcceptType;
1138: #define TSGLLEACCEPT_ALWAYS "always"

1140: /*S
1141:   TSGLLEAcceptFn - A prototype of a `TS` accept function that would be passed to `TSGLLEAcceptRegister()`

1143:   Calling Sequence:
1144: + ts  - timestep context
1145: . nt - time to end of solution time
1146: . h - the proposed step-size
1147: . enorm - unknown
1148: - accept - output, if the proposal is accepted

1150:   Level: beginner

1152:   Note:
1153:   The deprecated `TSGLLEAcceptFunction` still works as a replacement for `TSGLLEAcceptFn` *

1155: .seealso: [](ch_ts), `TS`, `TSSetRHSFunction()`, `DMTSSetRHSFunction()`, `TSIFunctionFn`,
1156: `TSIJacobianFn`, `TSRHSJacobianFn`, `TSGLLEAcceptRegister()`
1157: S*/
1158: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSGLLEAcceptFn)(TS ts, PetscReal nt, PetscReal h, const PetscReal enorm[], PetscBool *accept);

1160: PETSC_EXTERN_TYPEDEF typedef TSGLLEAcceptFn *TSGLLEAcceptFunction;

1162: PETSC_EXTERN PetscErrorCode TSGLLEAcceptRegister(const char[], TSGLLEAcceptFn *);

1164: /*J
1165:   TSGLLEType - string with the name of a General Linear `TSGLLE` type

1167:   Level: beginner

1169: .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLESetType()`, `TSGLLERegister()`, `TSGLLEAccept`
1170: J*/
1171: typedef const char *TSGLLEType;
1172: #define TSGLLE_IRKS "irks"

1174: PETSC_EXTERN PetscErrorCode TSGLLERegister(const char[], PetscErrorCode (*)(TS));
1175: PETSC_EXTERN PetscErrorCode TSGLLEInitializePackage(void);
1176: PETSC_EXTERN PetscErrorCode TSGLLEFinalizePackage(void);
1177: PETSC_EXTERN PetscErrorCode TSGLLESetType(TS, TSGLLEType);
1178: PETSC_EXTERN PetscErrorCode TSGLLEGetAdapt(TS, TSGLLEAdapt *);
1179: PETSC_EXTERN PetscErrorCode TSGLLESetAcceptType(TS, TSGLLEAcceptType);

1181: /*J
1182:    TSEIMEXType - String with the name of an Extrapolated IMEX `TSEIMEX` type

1184:    Level: beginner

1186: .seealso: [](ch_ts), `TSEIMEXSetType()`, `TS`, `TSEIMEX`, `TSEIMEXRegister()`
1187: J*/
1188: #define TSEIMEXType char *

1190: PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts, PetscInt);
1191: PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts, PetscInt, PetscInt);
1192: PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS, PetscBool);

1194: /*J
1195:    TSRKType - String with the name of a Runge-Kutta `TSRK` type

1197:    Level: beginner

1199: .seealso: [](ch_ts), `TS`, `TSRKSetType()`, `TS`, `TSRK`, `TSRKRegister()`
1200: J*/
1201: typedef const char *TSRKType;
1202: #define TSRK1FE "1fe"
1203: #define TSRK2A  "2a"
1204: #define TSRK2B  "2b"
1205: #define TSRK3   "3"
1206: #define TSRK3BS "3bs"
1207: #define TSRK4   "4"
1208: #define TSRK5F  "5f"
1209: #define TSRK5DP "5dp"
1210: #define TSRK5BS "5bs"
1211: #define TSRK6VR "6vr"
1212: #define TSRK7VR "7vr"
1213: #define TSRK8VR "8vr"

1215: PETSC_EXTERN PetscErrorCode TSRKGetOrder(TS, PetscInt *);
1216: PETSC_EXTERN PetscErrorCode TSRKGetType(TS, TSRKType *);
1217: PETSC_EXTERN PetscErrorCode TSRKSetType(TS, TSRKType);
1218: PETSC_EXTERN PetscErrorCode TSRKGetTableau(TS, PetscInt *, const PetscReal **, const PetscReal **, const PetscReal **, const PetscReal **, PetscInt *, const PetscReal **, PetscBool *);
1219: PETSC_EXTERN PetscErrorCode TSRKSetMultirate(TS, PetscBool);
1220: PETSC_EXTERN PetscErrorCode TSRKGetMultirate(TS, PetscBool *);
1221: PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
1222: PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void);
1223: PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void);
1224: PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void);

1226: /*J
1227:    TSMPRKType - String with the name of a partitioned Runge-Kutta `TSMPRK` type

1229:    Level: beginner

1231: .seealso: [](ch_ts), `TSMPRKSetType()`, `TS`, `TSMPRK`, `TSMPRKRegister()`
1232: J*/
1233: typedef const char *TSMPRKType;
1234: #define TSMPRK2A22 "2a22"
1235: #define TSMPRK2A23 "2a23"
1236: #define TSMPRK2A32 "2a32"
1237: #define TSMPRK2A33 "2a33"
1238: #define TSMPRKP2   "p2"
1239: #define TSMPRKP3   "p3"

1241: PETSC_EXTERN PetscErrorCode TSMPRKGetType(TS ts, TSMPRKType *);
1242: PETSC_EXTERN PetscErrorCode TSMPRKSetType(TS ts, TSMPRKType);
1243: 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[]);
1244: PETSC_EXTERN PetscErrorCode TSMPRKInitializePackage(void);
1245: PETSC_EXTERN PetscErrorCode TSMPRKFinalizePackage(void);
1246: PETSC_EXTERN PetscErrorCode TSMPRKRegisterDestroy(void);

1248: /*J
1249:    TSIRKType - String with the name of an implicit Runge-Kutta `TSIRK` type

1251:    Level: beginner

1253: .seealso: [](ch_ts), `TSIRKSetType()`, `TS`, `TSIRK`, `TSIRKRegister()`
1254: J*/
1255: typedef const char *TSIRKType;
1256: #define TSIRKGAUSS "gauss"

1258: PETSC_EXTERN PetscErrorCode TSIRKGetType(TS, TSIRKType *);
1259: PETSC_EXTERN PetscErrorCode TSIRKSetType(TS, TSIRKType);
1260: PETSC_EXTERN PetscErrorCode TSIRKGetNumStages(TS, PetscInt *);
1261: PETSC_EXTERN PetscErrorCode TSIRKSetNumStages(TS, PetscInt);
1262: PETSC_EXTERN PetscErrorCode TSIRKRegister(const char[], PetscErrorCode (*function)(TS));
1263: PETSC_EXTERN PetscErrorCode TSIRKTableauCreate(TS, PetscInt, const PetscReal *, const PetscReal *, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, const PetscScalar *);
1264: PETSC_EXTERN PetscErrorCode TSIRKInitializePackage(void);
1265: PETSC_EXTERN PetscErrorCode TSIRKFinalizePackage(void);
1266: PETSC_EXTERN PetscErrorCode TSIRKRegisterDestroy(void);

1268: /*J
1269:    TSGLEEType - String with the name of a General Linear with Error Estimation `TSGLEE` type

1271:    Level: beginner

1273: .seealso: [](ch_ts), `TSGLEESetType()`, `TS`, `TSGLEE`, `TSGLEERegister()`
1274: J*/
1275: typedef const char *TSGLEEType;
1276: #define TSGLEEi1      "BE1"
1277: #define TSGLEE23      "23"
1278: #define TSGLEE24      "24"
1279: #define TSGLEE25I     "25i"
1280: #define TSGLEE35      "35"
1281: #define TSGLEEEXRK2A  "exrk2a"
1282: #define TSGLEERK32G1  "rk32g1"
1283: #define TSGLEERK285EX "rk285ex"

1285: /*J
1286:    TSGLEEMode - String with the mode of error estimation for a General Linear with Error Estimation `TSGLEE` type

1288:    Level: beginner

1290: .seealso: [](ch_ts), `TSGLEESetMode()`, `TS`, `TSGLEE`, `TSGLEERegister()`
1291: J*/
1292: PETSC_EXTERN PetscErrorCode TSGLEEGetType(TS ts, TSGLEEType *);
1293: PETSC_EXTERN PetscErrorCode TSGLEESetType(TS ts, TSGLEEType);
1294: 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[]);
1295: PETSC_EXTERN PetscErrorCode TSGLEERegisterAll(void);
1296: PETSC_EXTERN PetscErrorCode TSGLEEFinalizePackage(void);
1297: PETSC_EXTERN PetscErrorCode TSGLEEInitializePackage(void);
1298: PETSC_EXTERN PetscErrorCode TSGLEERegisterDestroy(void);

1300: /*J
1301:    TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX `TSARKIMEX` type

1303:    Level: beginner

1305: .seealso: [](ch_ts), `TSARKIMEXSetType()`, `TS`, `TSARKIMEX`, `TSARKIMEXRegister()`
1306: J*/
1307: typedef const char *TSARKIMEXType;
1308: #define TSARKIMEX1BEE   "1bee"
1309: #define TSARKIMEXA2     "a2"
1310: #define TSARKIMEXL2     "l2"
1311: #define TSARKIMEXARS122 "ars122"
1312: #define TSARKIMEX2C     "2c"
1313: #define TSARKIMEX2D     "2d"
1314: #define TSARKIMEX2E     "2e"
1315: #define TSARKIMEXPRSSP2 "prssp2"
1316: #define TSARKIMEX3      "3"
1317: #define TSARKIMEXBPR3   "bpr3"
1318: #define TSARKIMEXARS443 "ars443"
1319: #define TSARKIMEX4      "4"
1320: #define TSARKIMEX5      "5"
1321: PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts, TSARKIMEXType *);
1322: PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts, TSARKIMEXType);
1323: PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS, PetscBool);
1324: PETSC_EXTERN PetscErrorCode TSARKIMEXGetFullyImplicit(TS, PetscBool *);
1325: 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[]);
1326: PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void);
1327: PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
1328: PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);

1330: /*J
1331:    TSDIRKType - String with the name of a Diagonally Implicit Runge-Kutta `TSDIRK` type

1333:    Level: beginner

1335: .seealso: [](ch_ts), `TSDIRKSetType()`, `TS`, `TSDIRK`, `TSDIRKRegister()`
1336: J*/
1337: typedef const char *TSDIRKType;
1338: #define TSDIRKS212      "s212"
1339: #define TSDIRKES122SAL  "es122sal"
1340: #define TSDIRKES213SAL  "es213sal"
1341: #define TSDIRKES324SAL  "es324sal"
1342: #define TSDIRKES325SAL  "es325sal"
1343: #define TSDIRK657A      "657a"
1344: #define TSDIRKES648SA   "es648sa"
1345: #define TSDIRK658A      "658a"
1346: #define TSDIRKS659A     "s659a"
1347: #define TSDIRK7510SAL   "7510sal"
1348: #define TSDIRKES7510SA  "es7510sa"
1349: #define TSDIRK759A      "759a"
1350: #define TSDIRKS7511SAL  "s7511sal"
1351: #define TSDIRK8614A     "8614a"
1352: #define TSDIRK8616SAL   "8616sal"
1353: #define TSDIRKES8516SAL "es8516sal"
1354: PETSC_EXTERN PetscErrorCode TSDIRKGetType(TS ts, TSDIRKType *);
1355: PETSC_EXTERN PetscErrorCode TSDIRKSetType(TS ts, TSDIRKType);
1356: PETSC_EXTERN PetscErrorCode TSDIRKRegister(TSDIRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);

1358: /*J
1359:    TSRosWType - String with the name of a Rosenbrock-W `TSROSW` type

1361:    Level: beginner

1363: .seealso: [](ch_ts), `TSRosWSetType()`, `TS`, `TSROSW`, `TSRosWRegister()`
1364: J*/
1365: typedef const char *TSRosWType;
1366: #define TSROSW2M          "2m"
1367: #define TSROSW2P          "2p"
1368: #define TSROSWRA3PW       "ra3pw"
1369: #define TSROSWRA34PW2     "ra34pw2"
1370: #define TSROSWR34PRW      "r34prw"
1371: #define TSROSWR3PRL2      "r3prl2"
1372: #define TSROSWRODAS3      "rodas3"
1373: #define TSROSWRODASPR     "rodaspr"
1374: #define TSROSWRODASPR2    "rodaspr2"
1375: #define TSROSWSANDU3      "sandu3"
1376: #define TSROSWASSP3P3S1C  "assp3p3s1c"
1377: #define TSROSWLASSP3P4S2C "lassp3p4s2c"
1378: #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
1379: #define TSROSWARK3        "ark3"
1380: #define TSROSWTHETA1      "theta1"
1381: #define TSROSWTHETA2      "theta2"
1382: #define TSROSWGRK4T       "grk4t"
1383: #define TSROSWSHAMP4      "shamp4"
1384: #define TSROSWVELDD4      "veldd4"
1385: #define TSROSW4L          "4l"

1387: PETSC_EXTERN PetscErrorCode TSRosWGetType(TS, TSRosWType *);
1388: PETSC_EXTERN PetscErrorCode TSRosWSetType(TS, TSRosWType);
1389: PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS, PetscBool);
1390: PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
1391: PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal);
1392: PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void);
1393: PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
1394: PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);

1396: PETSC_EXTERN PetscErrorCode TSBDFSetOrder(TS, PetscInt);
1397: PETSC_EXTERN PetscErrorCode TSBDFGetOrder(TS, PetscInt *);

1399: /*J
1400:   TSBasicSymplecticType - String with the name of a basic symplectic integration `TSBASICSYMPLECTIC` type

1402:   Level: beginner

1404: .seealso: [](ch_ts), `TSBasicSymplecticSetType()`, `TS`, `TSBASICSYMPLECTIC`, `TSBasicSymplecticRegister()`
1405: J*/
1406: typedef const char *TSBasicSymplecticType;
1407: #define TSBASICSYMPLECTICSIEULER   "1"
1408: #define TSBASICSYMPLECTICVELVERLET "2"
1409: #define TSBASICSYMPLECTIC3         "3"
1410: #define TSBASICSYMPLECTIC4         "4"
1411: PETSC_EXTERN PetscErrorCode TSBasicSymplecticSetType(TS, TSBasicSymplecticType);
1412: PETSC_EXTERN PetscErrorCode TSBasicSymplecticGetType(TS, TSBasicSymplecticType *);
1413: PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegister(TSBasicSymplecticType, PetscInt, PetscInt, PetscReal[], PetscReal[]);
1414: PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterAll(void);
1415: PETSC_EXTERN PetscErrorCode TSBasicSymplecticInitializePackage(void);
1416: PETSC_EXTERN PetscErrorCode TSBasicSymplecticFinalizePackage(void);
1417: PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterDestroy(void);

1419: /*J
1420:   TSDISCGRAD - The Discrete Gradient integrator is a timestepper for Hamiltonian systems designed to conserve the first integral (energy),
1421:   but also has the property for some systems of monotonicity in a functional.

1423:   Level: beginner

1425: .seealso: [](ch_ts), `TS`, TSDiscGradSetFormulation()`, `TSDiscGradGetFormulation()`
1426: J*/
1427: PETSC_EXTERN PetscErrorCode TSDiscGradSetFormulation(TS, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (*)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (*)(TS, PetscReal, Vec, Vec, void *), void *);
1428: PETSC_EXTERN PetscErrorCode TSDiscGradGetFormulation(TS, PetscErrorCode (**)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (**)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (**)(TS, PetscReal, Vec, Vec, void *), void *);
1429: PETSC_EXTERN PetscErrorCode TSDiscGradIsGonzalez(TS, PetscBool *);
1430: PETSC_EXTERN PetscErrorCode TSDiscGradUseGonzalez(TS, PetscBool);

1432: /*
1433:        PETSc interface to Sundials
1434: */
1435: #ifdef PETSC_HAVE_SUNDIALS2
1436: typedef enum {
1437:   SUNDIALS_ADAMS = 1,
1438:   SUNDIALS_BDF   = 2
1439: } TSSundialsLmmType;
1440: PETSC_EXTERN const char *const TSSundialsLmmTypes[];
1441: typedef enum {
1442:   SUNDIALS_MODIFIED_GS  = 1,
1443:   SUNDIALS_CLASSICAL_GS = 2
1444: } TSSundialsGramSchmidtType;
1445: PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];

1447: PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS, TSSundialsLmmType);
1448: PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS, PC *);
1449: PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS, PetscReal, PetscReal);
1450: PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS, PetscReal);
1451: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS, PetscReal);
1452: PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS, PetscInt *, PetscInt *);
1453: PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS, TSSundialsGramSchmidtType);
1454: PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS, PetscInt);
1455: PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS, PetscReal);
1456: PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS, PetscBool);
1457: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS, PetscInt);
1458: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxord(TS, PetscInt);
1459: PETSC_EXTERN PetscErrorCode TSSundialsSetUseDense(TS, PetscBool);
1460: #endif

1462: PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS, PetscReal);
1463: PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS, PetscReal *);
1464: PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS, PetscBool *);
1465: PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS, PetscBool);

1467: PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS, PetscReal);
1468: PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS, PetscReal, PetscReal, PetscReal);
1469: PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS, PetscReal *, PetscReal *, PetscReal *);

1471: PETSC_EXTERN PetscErrorCode TSAlpha2SetRadius(TS, PetscReal);
1472: PETSC_EXTERN PetscErrorCode TSAlpha2SetParams(TS, PetscReal, PetscReal, PetscReal, PetscReal);
1473: PETSC_EXTERN PetscErrorCode TSAlpha2GetParams(TS, PetscReal *, PetscReal *, PetscReal *, PetscReal *);

1475: /*S
1476:   TSAlpha2PredictorFn - A callback to set the predictor (i.e., the initial guess for the nonlinear solver) in
1477:   a second-order generalized-alpha time integrator.

1479:   Calling Sequence:
1480: + ts   - the `TS` context obtained from `TSCreate()`
1481: . X0   - the previous time step's state vector
1482: . V0   - the previous time step's first derivative of the state vector
1483: . A0   - the previous time step's second derivative of the state vector
1484: . X1   - the vector into which the initial guess for the current time step will be written
1485: - ctx  - [optional] user-defined context for the predictor evaluation routine (may be `NULL`)

1487:   Level: intermediate

1489:   Note:
1490:   The deprecated `TSAlpha2Predictor` still works as a replacement for `TSAlpha2PredictorFn` *.

1492: .seealso: [](ch_ts), `TS`, `TSAlpha2SetPredictor()`
1493: S*/
1494: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSAlpha2PredictorFn)(TS ts, Vec X0, Vec V0, Vec A0, Vec X1, void *ctx);

1496: PETSC_EXTERN_TYPEDEF typedef TSAlpha2PredictorFn *TSAlpha2Predictor;

1498: PETSC_EXTERN PetscErrorCode TSAlpha2SetPredictor(TS, TSAlpha2PredictorFn *, void *ctx);

1500: PETSC_EXTERN PetscErrorCode TSSetDM(TS, DM);
1501: PETSC_EXTERN PetscErrorCode TSGetDM(TS, DM *);

1503: PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES, Vec, Vec, void *);
1504: PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES, Vec, Mat, Mat, void *);

1506: PETSC_EXTERN PetscErrorCode TSRHSJacobianTest(TS, PetscBool *);
1507: PETSC_EXTERN PetscErrorCode TSRHSJacobianTestTranspose(TS, PetscBool *);

1509: PETSC_EXTERN PetscErrorCode TSGetComputeInitialCondition(TS, PetscErrorCode (**)(TS, Vec));
1510: PETSC_EXTERN PetscErrorCode TSSetComputeInitialCondition(TS, PetscErrorCode (*)(TS, Vec));
1511: PETSC_EXTERN PetscErrorCode TSComputeInitialCondition(TS, Vec);
1512: PETSC_EXTERN PetscErrorCode TSGetComputeExactError(TS, PetscErrorCode (**)(TS, Vec, Vec));
1513: PETSC_EXTERN PetscErrorCode TSSetComputeExactError(TS, PetscErrorCode (*)(TS, Vec, Vec));
1514: PETSC_EXTERN PetscErrorCode TSComputeExactError(TS, Vec, Vec);
1515: PETSC_EXTERN PetscErrorCode PetscConvEstUseTS(PetscConvEst, PetscBool);

1517: PETSC_EXTERN PetscErrorCode TSSetMatStructure(TS, MatStructure);