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 integrating an ODE.

 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. These are all the time/ODE integrators that PETSc provides.

 26:    Level: beginner

 28:    Note:
 29:    Use `TSSetType()` or the options database key `-ts_type` to set the ODE integrator method to use with a given `TS` object

 31: .seealso: [](integrator_table), [](ch_ts), `TSSetType()`, `TS`, `TSRegister()`
 32: J*/
 33: typedef const char *TSType;
 34: #define TSEULER           "euler"
 35: #define TSBEULER          "beuler"
 36: #define TSBASICSYMPLECTIC "basicsymplectic"
 37: #define TSPSEUDO          "pseudo"
 38: #define TSCN              "cn"
 39: #define TSSUNDIALS        "sundials"
 40: #define TSRK              "rk"
 41: #define TSPYTHON          "python"
 42: #define TSTHETA           "theta"
 43: #define TSALPHA           "alpha"
 44: #define TSALPHA2          "alpha2"
 45: #define TSGLLE            "glle"
 46: #define TSGLEE            "glee"
 47: #define TSSSP             "ssp"
 48: #define TSARKIMEX         "arkimex"
 49: #define TSROSW            "rosw"
 50: #define TSEIMEX           "eimex"
 51: #define TSMIMEX           "mimex"
 52: #define TSBDF             "bdf"
 53: #define TSRADAU5          "radau5"
 54: #define TSMPRK            "mprk"
 55: #define TSDISCGRAD        "discgrad"
 56: #define TSIRK             "irk"
 57: #define TSDIRK            "dirk"

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

 62:    Values:
 63:  + `TS_LINEAR`    - a linear ODE or DAE
 64:  - `TS_NONLINEAR` - a nonlinear ODE or DAE

 66:    Level: beginner

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

 75: /*E
 76:    TSEquationType - type of `TS` problem that is solved

 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

 83:    Level: beginner

 85: .seealso: [](ch_ts), `TS`, `TSGetEquationType()`, `TSSetEquationType()`
 86: E*/
 87: typedef enum {
 88:   TS_EQ_UNSPECIFIED               = -1,
 89:   TS_EQ_EXPLICIT                  = 0,
 90:   TS_EQ_ODE_EXPLICIT              = 1,
 91:   TS_EQ_DAE_SEMI_EXPLICIT_INDEX1  = 100,
 92:   TS_EQ_DAE_SEMI_EXPLICIT_INDEX2  = 200,
 93:   TS_EQ_DAE_SEMI_EXPLICIT_INDEX3  = 300,
 94:   TS_EQ_DAE_SEMI_EXPLICIT_INDEXHI = 500,
 95:   TS_EQ_IMPLICIT                  = 1000,
 96:   TS_EQ_ODE_IMPLICIT              = 1001,
 97:   TS_EQ_DAE_IMPLICIT_INDEX1       = 1100,
 98:   TS_EQ_DAE_IMPLICIT_INDEX2       = 1200,
 99:   TS_EQ_DAE_IMPLICIT_INDEX3       = 1300,
100:   TS_EQ_DAE_IMPLICIT_INDEXHI      = 1500
101: } TSEquationType;
102: PETSC_EXTERN const char *const *TSEquationTypes;

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

107:    Values:
108: +  `TS_CONVERGED_ITERATING`          - this only occurs if `TSGetConvergedReason()` is called during the `TSSolve()`
109: .  `TS_CONVERGED_TIME`               - the final time was reached
110: .  `TS_CONVERGED_ITS`                - the maximum number of iterations (time-steps) was reached prior to the final time
111: .  `TS_CONVERGED_USER`               - user requested termination
112: .  `TS_CONVERGED_EVENT`              - user requested termination on event detection
113: .  `TS_CONVERGED_PSEUDO_FATOL`       - stops when function norm decreased by a set amount, used only for `TSPSEUDO`
114: .  `TS_CONVERGED_PSEUDO_FRTOL`       - stops when function norm decreases below a set amount, used only for `TSPSEUDO`
115: .  `TS_DIVERGED_NONLINEAR_SOLVE`     - too many nonlinear solve failures have occurred
116: .  `TS_DIVERGED_STEP_REJECTED`       - too many steps were rejected
117: .  `TSFORWARD_DIVERGED_LINEAR_SOLVE` - tangent linear solve failed
118: -  `TSADJOINT_DIVERGED_LINEAR_SOLVE` - transposed linear solve failed

120:    Level: beginner

122: .seealso: [](ch_ts), `TS`, `TSGetConvergedReason()`
123: E*/
124: typedef enum {
125:   TS_CONVERGED_ITERATING          = 0,
126:   TS_CONVERGED_TIME               = 1,
127:   TS_CONVERGED_ITS                = 2,
128:   TS_CONVERGED_USER               = 3,
129:   TS_CONVERGED_EVENT              = 4,
130:   TS_CONVERGED_PSEUDO_FATOL       = 5,
131:   TS_CONVERGED_PSEUDO_FRTOL       = 6,
132:   TS_DIVERGED_NONLINEAR_SOLVE     = -1,
133:   TS_DIVERGED_STEP_REJECTED       = -2,
134:   TSFORWARD_DIVERGED_LINEAR_SOLVE = -3,
135:   TSADJOINT_DIVERGED_LINEAR_SOLVE = -4
136: } TSConvergedReason;
137: PETSC_EXTERN const char *const *TSConvergedReasons;

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

142:    Level: beginner

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

147: /*MC
148:    TS_CONVERGED_TIME - the final time was reached

150:    Level: beginner

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

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

158:    Level: beginner

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

163: /*MC
164:    TS_CONVERGED_USER - user requested termination

166:    Level: beginner

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

171: /*MC
172:    TS_CONVERGED_EVENT - user requested termination on event detection

174:    Level: beginner

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

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

182:    Options Database Key:
183: .   -ts_pseudo_frtol rtol - use specified rtol

185:    Level: beginner

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

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

193:    Options Database Key:
194: .   -ts_pseudo_fatol atol - use specified atol

196:    Level: beginner

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

201: /*MC
202:    TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed

204:    Level: beginner

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

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

212: /*MC
213:    TS_DIVERGED_STEP_REJECTED - too many steps were rejected

215:    Level: beginner

217:    Notes:
218:    See `TSSetMaxStepRejections()` for how to allow more step rejections.

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

223: /*E
224:    TSExactFinalTimeOption - option for handling of final time step

226:    Values:
227: +  `TS_EXACTFINALTIME_STEPOVER`    - Don't do anything if requested final time is exceeded
228: .  `TS_EXACTFINALTIME_INTERPOLATE` - Interpolate back to final time
229: -  `TS_EXACTFINALTIME_MATCHSTEP`   - Adapt final time step to match the final time requested

231:    Level: beginner

233: .seealso: [](ch_ts), `TS`, `TSGetConvergedReason()`, `TSSetExactFinalTime()`, `TSGetExactFinalTime()`
234: E*/
235: typedef enum {
236:   TS_EXACTFINALTIME_UNSPECIFIED = 0,
237:   TS_EXACTFINALTIME_STEPOVER    = 1,
238:   TS_EXACTFINALTIME_INTERPOLATE = 2,
239:   TS_EXACTFINALTIME_MATCHSTEP   = 3
240: } TSExactFinalTimeOption;
241: PETSC_EXTERN const char *const TSExactFinalTimeOptions[];

243: /* Logging support */
244: PETSC_EXTERN PetscClassId TS_CLASSID;
245: PETSC_EXTERN PetscClassId DMTS_CLASSID;
246: PETSC_EXTERN PetscClassId TSADAPT_CLASSID;

248: PETSC_EXTERN PetscErrorCode TSInitializePackage(void);
249: PETSC_EXTERN PetscErrorCode TSFinalizePackage(void);

251: PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm, TS *);
252: PETSC_EXTERN PetscErrorCode TSClone(TS, TS *);
253: PETSC_EXTERN PetscErrorCode TSDestroy(TS *);

255: PETSC_EXTERN PetscErrorCode TSSetProblemType(TS, TSProblemType);
256: PETSC_EXTERN PetscErrorCode TSGetProblemType(TS, TSProblemType *);
257: PETSC_EXTERN PetscErrorCode TSMonitor(TS, PetscInt, PetscReal, Vec);
258: PETSC_EXTERN PetscErrorCode TSMonitorSet(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscCtx), PetscCtx, PetscCtxDestroyFn *);
259: PETSC_EXTERN PetscErrorCode TSMonitorSetFromOptions(TS, const char[], const char[], const char[], PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *), PetscErrorCode (*)(TS, PetscViewerAndFormat *));
260: PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS);

262: PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS, const char[]);
263: PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS, const char[]);
264: PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS, const char *[]);
265: PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS);
266: PETSC_EXTERN PetscErrorCode TSSetUp(TS);
267: PETSC_EXTERN PetscErrorCode TSReset(TS);

269: PETSC_EXTERN PetscErrorCode TSSetSolution(TS, Vec);
270: PETSC_EXTERN PetscErrorCode TSGetSolution(TS, Vec *);

272: PETSC_EXTERN PetscErrorCode TS2SetSolution(TS, Vec, Vec);
273: PETSC_EXTERN PetscErrorCode TS2GetSolution(TS, Vec *, Vec *);

275: PETSC_EXTERN PetscErrorCode TSGetSolutionComponents(TS, PetscInt *, Vec *);
276: PETSC_EXTERN PetscErrorCode TSGetAuxSolution(TS, Vec *);
277: PETSC_EXTERN PetscErrorCode TSGetTimeError(TS, PetscInt, Vec *);
278: PETSC_EXTERN PetscErrorCode TSSetTimeError(TS, Vec);

280: PETSC_EXTERN PetscErrorCode TSSetRHSJacobianP(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, PetscCtx), PetscCtx);
281: PETSC_EXTERN PetscErrorCode TSGetRHSJacobianP(TS, Mat *, PetscErrorCode (**)(TS, PetscReal, Vec, Mat, PetscCtx), PetscCtxRt);
282: PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianP(TS, PetscReal, Vec, Mat);
283: PETSC_EXTERN PetscErrorCode TSSetIJacobianP(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Vec, PetscReal, Mat, PetscCtx), PetscCtx);
284: PETSC_EXTERN PetscErrorCode TSGetIJacobianP(TS, Mat *, PetscErrorCode (**)(TS, PetscReal, Vec, Vec, PetscReal, Mat, PetscCtx), PetscCtxRt);
285: PETSC_EXTERN PetscErrorCode TSComputeIJacobianP(TS, PetscReal, Vec, Vec, PetscReal, Mat, PetscBool);
286: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS() then TSComputeRHSJacobianP()", ) PetscErrorCode TSComputeDRDPFunction(TS, PetscReal, Vec, Vec *);
287: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS() then TSComputeRHSJacobian()", ) PetscErrorCode TSComputeDRDUFunction(TS, PetscReal, Vec, Vec *);
288: PETSC_EXTERN PetscErrorCode TSSetIHessianProduct(TS, Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, PetscCtx), Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, PetscCtx), Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, PetscCtx), Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, PetscCtx), PetscCtx);
289: PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUU(TS, PetscReal, Vec, Vec[], Vec, Vec[]);
290: PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUP(TS, PetscReal, Vec, Vec[], Vec, Vec[]);
291: PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPU(TS, PetscReal, Vec, Vec[], Vec, Vec[]);
292: PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPP(TS, PetscReal, Vec, Vec[], Vec, Vec[]);
293: PETSC_EXTERN PetscErrorCode TSSetRHSHessianProduct(TS, Vec[], PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, PetscCtx), Vec[], PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, PetscCtx), Vec[], PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, PetscCtx), Vec[], PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, PetscCtx), PetscCtx);
294: PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS, PetscReal, Vec, Vec[], Vec, Vec[]);
295: PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS, PetscReal, Vec, Vec[], Vec, Vec[]);
296: PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS, PetscReal, Vec, Vec[], Vec, Vec[]);
297: PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS, PetscReal, Vec, Vec[], Vec, Vec[]);
298: PETSC_EXTERN PetscErrorCode TSSetCostHessianProducts(TS, PetscInt, Vec[], Vec[], Vec);
299: PETSC_EXTERN PetscErrorCode TSGetCostHessianProducts(TS, PetscInt *, Vec *[], Vec *[], Vec *);
300: PETSC_EXTERN PetscErrorCode TSComputeSNESJacobian(TS, Vec, Mat, Mat);

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

305:    Level: advanced

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

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

314:    Level: intermediate

316: .seealso: [](ch_ts), `TS`, `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectoryDestroy()`
317: J*/
318: typedef const char *TSTrajectoryType;
319: #define TSTRAJECTORYBASIC         "basic"
320: #define TSTRAJECTORYSINGLEFILE    "singlefile"
321: #define TSTRAJECTORYMEMORY        "memory"
322: #define TSTRAJECTORYVISUALIZATION "visualization"

324: PETSC_EXTERN PetscFunctionList TSTrajectoryList;
325: PETSC_EXTERN PetscClassId      TSTRAJECTORY_CLASSID;
326: PETSC_EXTERN PetscBool         TSTrajectoryRegisterAllCalled;

328: PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS);
329: PETSC_EXTERN PetscErrorCode TSResetTrajectory(TS);
330: PETSC_EXTERN PetscErrorCode TSRemoveTrajectory(TS);

332: PETSC_EXTERN PetscErrorCode TSTrajectoryCreate(MPI_Comm, TSTrajectory *);
333: PETSC_EXTERN PetscErrorCode TSTrajectoryReset(TSTrajectory);
334: PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy(TSTrajectory *);
335: PETSC_EXTERN PetscErrorCode TSTrajectoryView(TSTrajectory, PetscViewer);
336: PETSC_EXTERN PetscErrorCode TSTrajectorySetType(TSTrajectory, TS, TSTrajectoryType);
337: PETSC_EXTERN PetscErrorCode TSTrajectoryGetType(TSTrajectory, TS, TSTrajectoryType *);
338: PETSC_EXTERN PetscErrorCode TSTrajectorySet(TSTrajectory, TS, PetscInt, PetscReal, Vec);
339: PETSC_EXTERN PetscErrorCode TSTrajectoryGet(TSTrajectory, TS, PetscInt, PetscReal *);
340: PETSC_EXTERN PetscErrorCode TSTrajectoryGetVecs(TSTrajectory, TS, PetscInt, PetscReal *, Vec, Vec);
341: PETSC_EXTERN PetscErrorCode TSTrajectoryGetUpdatedHistoryVecs(TSTrajectory, TS, PetscReal, Vec *, Vec *);
342: PETSC_EXTERN PetscErrorCode TSTrajectoryGetNumSteps(TSTrajectory, PetscInt *);
343: PETSC_EXTERN PetscErrorCode TSTrajectoryRestoreUpdatedHistoryVecs(TSTrajectory, Vec *, Vec *);
344: PETSC_EXTERN PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory, TS);
345: PETSC_EXTERN PetscErrorCode TSTrajectoryRegister(const char[], PetscErrorCode (*)(TSTrajectory, TS));
346: PETSC_EXTERN PetscErrorCode TSTrajectoryRegisterAll(void);
347: PETSC_EXTERN PetscErrorCode TSTrajectorySetUp(TSTrajectory, TS);
348: PETSC_EXTERN PetscErrorCode TSTrajectorySetUseHistory(TSTrajectory, PetscBool);
349: PETSC_EXTERN PetscErrorCode TSTrajectorySetMonitor(TSTrajectory, PetscBool);
350: PETSC_EXTERN PetscErrorCode TSTrajectorySetVariableNames(TSTrajectory, const char *const *);
351: PETSC_EXTERN PetscErrorCode TSTrajectorySetTransform(TSTrajectory, PetscErrorCode (*)(PetscCtx, Vec, Vec *), PetscCtxDestroyFn *, PetscCtx);
352: PETSC_EXTERN PetscErrorCode TSTrajectorySetSolutionOnly(TSTrajectory, PetscBool);
353: PETSC_EXTERN PetscErrorCode TSTrajectoryGetSolutionOnly(TSTrajectory, PetscBool *);
354: PETSC_EXTERN PetscErrorCode TSTrajectorySetKeepFiles(TSTrajectory, PetscBool);
355: PETSC_EXTERN PetscErrorCode TSTrajectorySetDirname(TSTrajectory, const char[]);
356: PETSC_EXTERN PetscErrorCode TSTrajectorySetFiletemplate(TSTrajectory, const char[]);
357: PETSC_EXTERN PetscErrorCode TSGetTrajectory(TS, TSTrajectory *);

359: typedef enum {
360:   TJ_REVOLVE,
361:   TJ_CAMS,
362:   TJ_PETSC
363: } TSTrajectoryMemoryType;
364: PETSC_EXTERN const char *const TSTrajectoryMemoryTypes[];

366: PETSC_EXTERN PetscErrorCode TSTrajectoryMemorySetType(TSTrajectory, TSTrajectoryMemoryType);
367: PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxCpsRAM(TSTrajectory, PetscInt);
368: PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxCpsDisk(TSTrajectory, PetscInt);
369: PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxUnitsRAM(TSTrajectory, PetscInt);
370: PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxUnitsDisk(TSTrajectory, PetscInt);

372: PETSC_EXTERN PetscErrorCode TSSetCostGradients(TS, PetscInt, Vec[], Vec[]);
373: PETSC_EXTERN PetscErrorCode TSGetCostGradients(TS, PetscInt *, Vec *[], Vec *[]);
374: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 12, 0, "TSCreateQuadratureTS() and TSForwardSetSensitivities()", ) PetscErrorCode TSSetCostIntegrand(TS, PetscInt, Vec, PetscErrorCode (*)(TS, PetscReal, Vec, Vec, PetscCtx), PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, PetscCtx), PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, PetscCtx), PetscBool, PetscCtx);
375: PETSC_EXTERN PetscErrorCode TSGetCostIntegral(TS, Vec *);
376: PETSC_EXTERN PetscErrorCode TSComputeCostIntegrand(TS, PetscReal, Vec, Vec);
377: PETSC_EXTERN PetscErrorCode TSCreateQuadratureTS(TS, PetscBool, TS *);
378: PETSC_EXTERN PetscErrorCode TSGetQuadratureTS(TS, PetscBool *, TS *);

380: PETSC_EXTERN PetscErrorCode TSAdjointSetFromOptions(TS, PetscOptionItems);
381: PETSC_EXTERN PetscErrorCode TSAdjointMonitor(TS, PetscInt, PetscReal, Vec, PetscInt, Vec[], Vec[]);
382: PETSC_EXTERN PetscErrorCode TSAdjointMonitorSet(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, PetscCtx), PetscCtx, PetscCtxDestroyFn *);
383: PETSC_EXTERN PetscErrorCode TSAdjointMonitorCancel(TS);
384: PETSC_EXTERN PetscErrorCode TSAdjointMonitorSetFromOptions(TS, const char[], const char[], const char[], PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, PetscViewerAndFormat *), PetscErrorCode (*)(TS, PetscViewerAndFormat *));

386: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSSetRHSJacobianP()", ) PetscErrorCode TSAdjointSetRHSJacobian(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, PetscCtx), PetscCtx);
387: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSComputeRHSJacobianP()", ) PetscErrorCode TSAdjointComputeRHSJacobian(TS, PetscReal, Vec, Mat);
388: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS()", ) PetscErrorCode TSAdjointComputeDRDPFunction(TS, PetscReal, Vec, Vec *);
389: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS()", ) PetscErrorCode TSAdjointComputeDRDYFunction(TS, PetscReal, Vec, Vec *);
390: PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS);
391: PETSC_EXTERN PetscErrorCode TSAdjointSetSteps(TS, PetscInt);

393: PETSC_EXTERN PetscErrorCode TSAdjointStep(TS);
394: PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS);
395: PETSC_EXTERN PetscErrorCode TSAdjointReset(TS);
396: PETSC_EXTERN PetscErrorCode TSAdjointCostIntegral(TS);
397: PETSC_EXTERN PetscErrorCode TSAdjointSetForward(TS, Mat);
398: PETSC_EXTERN PetscErrorCode TSAdjointResetForward(TS);

400: PETSC_EXTERN PetscErrorCode TSForwardSetSensitivities(TS, PetscInt, Mat);
401: PETSC_EXTERN PetscErrorCode TSForwardGetSensitivities(TS, PetscInt *, Mat *);
402: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 12, 0, "TSCreateQuadratureTS()", ) PetscErrorCode TSForwardSetIntegralGradients(TS, PetscInt, Vec[]);
403: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 12, 0, "TSForwardGetSensitivities()", ) PetscErrorCode TSForwardGetIntegralGradients(TS, PetscInt *, Vec *[]);
404: PETSC_EXTERN PetscErrorCode TSForwardSetUp(TS);
405: PETSC_EXTERN PetscErrorCode TSForwardReset(TS);
406: PETSC_EXTERN PetscErrorCode TSForwardCostIntegral(TS);
407: PETSC_EXTERN PetscErrorCode TSForwardStep(TS);
408: PETSC_EXTERN PetscErrorCode TSForwardSetInitialSensitivities(TS, Mat);
409: PETSC_EXTERN PetscErrorCode TSForwardGetStages(TS, PetscInt *, Mat *[]);

411: PETSC_EXTERN PetscErrorCode TSSetMaxSteps(TS, PetscInt);
412: PETSC_EXTERN PetscErrorCode TSGetMaxSteps(TS, PetscInt *);
413: PETSC_EXTERN PetscErrorCode TSSetRunSteps(TS, PetscInt);
414: PETSC_EXTERN PetscErrorCode TSGetRunSteps(TS, PetscInt *);
415: PETSC_EXTERN PetscErrorCode TSSetMaxTime(TS, PetscReal);
416: PETSC_EXTERN PetscErrorCode TSGetMaxTime(TS, PetscReal *);
417: PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS, TSExactFinalTimeOption);
418: PETSC_EXTERN PetscErrorCode TSGetExactFinalTime(TS, TSExactFinalTimeOption *);
419: PETSC_EXTERN PetscErrorCode TSSetEvaluationTimes(TS, PetscInt, PetscReal[]);
420: PETSC_EXTERN PetscErrorCode TSGetEvaluationTimes(TS, PetscInt *, const PetscReal *[]);
421: PETSC_EXTERN PetscErrorCode TSGetEvaluationSolutions(TS, PetscInt *, const PetscReal *[], Vec *[]);
422: PETSC_EXTERN PetscErrorCode TSSetTimeSpan(TS, PetscInt, PetscReal[]);

424: /*@C
425:   TSGetTimeSpan - gets the time span set with `TSSetTimeSpan()`

427:   Not Collective

429:   Input Parameter:
430: . ts - the time-stepper

432:   Output Parameters:
433: + n          - number of the time points (>=2)
434: - span_times - array of the time points. The first element and the last element are the initial time and the final time respectively.

436:   Level: deprecated

438:   Note:
439:   Deprecated, use `TSGetEvaluationTimes()`.

441:   The values obtained are valid until the `TS` object is destroyed.

443:   Both `n` and `span_times` can be `NULL`.

445: .seealso: [](ch_ts), `TS`, `TSGetEvaluationTimes()`, `TSSetTimeSpan()`, `TSSetEvaluationTimes()`, `TSGetEvaluationSolutions()`
446:  @*/
447: PETSC_DEPRECATED_FUNCTION(3, 23, 0, "TSGetEvaluationTimes()", ) static inline PetscErrorCode TSGetTimeSpan(TS ts, PetscInt *n, const PetscReal *span_times[])
448: {
449:   return TSGetEvaluationTimes(ts, n, span_times);
450: }

452: /*@C
453:   TSGetTimeSpanSolutions - Get the number of solutions and the solutions at the time points specified by the time span.

455:   Input Parameter:
456: . ts - the `TS` context obtained from `TSCreate()`

458:   Output Parameters:
459: + nsol - the number of solutions
460: - Sols - the solution vectors

462:   Level: deprecated

464:   Notes:
465:   Deprecated, use `TSGetEvaluationSolutions()`.

467:   Both `nsol` and `Sols` can be `NULL`.

469:   Some time points in the time span may be skipped by `TS` so that `nsol` is less than the number of points specified by `TSSetTimeSpan()`.
470:   For example, manipulating the step size, especially with a reduced precision, may cause `TS` to step over certain points in the span.
471:   This issue is alleviated in `TSGetEvaluationSolutions()` by returning the solution times that `Sols` were recorded at.

473: .seealso: [](ch_ts), `TS`, `TSGetEvaluationSolutions()`, `TSSetTimeSpan()`, `TSGetEvaluationTimes()`, `TSSetEvaluationTimes()`
474:  @*/
475: PETSC_DEPRECATED_FUNCTION(3, 23, 0, "TSGetEvaluationSolutions()", ) static inline PetscErrorCode TSGetTimeSpanSolutions(TS ts, PetscInt *nsol, Vec **Sols)
476: {
477:   return TSGetEvaluationSolutions(ts, nsol, NULL, Sols);
478: }

480: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSSetTime()", ) PetscErrorCode TSSetInitialTimeStep(TS, PetscReal, PetscReal);
481: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSSetMax()", ) PetscErrorCode TSSetDuration(TS, PetscInt, PetscReal);
482: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetMax()", ) PetscErrorCode TSGetDuration(TS, PetscInt *, PetscReal *);
483: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetStepNumber()", ) PetscErrorCode TSGetTimeStepNumber(TS, PetscInt *);
484: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetStepNumber()", ) PetscErrorCode TSGetTotalSteps(TS, PetscInt *);

486: PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
487: PETSC_EXTERN PetscErrorCode TSMonitorWallClockTime(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
488: PETSC_EXTERN PetscErrorCode TSMonitorWallClockTimeSetUp(TS, PetscViewerAndFormat *);
489: PETSC_EXTERN PetscErrorCode TSMonitorExtreme(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);

491: typedef struct _n_TSMonitorDrawCtx *TSMonitorDrawCtx;
492: PETSC_EXTERN PetscErrorCode         TSMonitorDrawCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorDrawCtx *);
493: PETSC_EXTERN PetscErrorCode         TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *);
494: PETSC_EXTERN PetscErrorCode         TSMonitorDrawSolution(TS, PetscInt, PetscReal, Vec, PetscCtx);
495: PETSC_EXTERN PetscErrorCode         TSMonitorDrawSolutionPhase(TS, PetscInt, PetscReal, Vec, PetscCtx);
496: PETSC_EXTERN PetscErrorCode         TSMonitorDrawError(TS, PetscInt, PetscReal, Vec, PetscCtx);
497: PETSC_EXTERN PetscErrorCode         TSMonitorDrawSolutionFunction(TS, PetscInt, PetscReal, Vec, PetscCtx);

499: PETSC_EXTERN PetscErrorCode TSAdjointMonitorDefault(TS, PetscInt, PetscReal, Vec, PetscInt, Vec[], Vec[], PetscViewerAndFormat *);
500: PETSC_EXTERN PetscErrorCode TSAdjointMonitorDrawSensi(TS, PetscInt, PetscReal, Vec, PetscInt, Vec[], Vec[], PetscCtx);

502: typedef struct _n_TSMonitorSolutionCtx *TSMonitorSolutionCtx;
503: PETSC_EXTERN PetscErrorCode             TSMonitorSolution(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
504: PETSC_EXTERN PetscErrorCode             TSMonitorSolutionSetup(TS, PetscViewerAndFormat *);

506: typedef struct _n_TSMonitorVTKCtx *TSMonitorVTKCtx;
507: PETSC_EXTERN PetscErrorCode        TSMonitorSolutionVTK(TS, PetscInt, PetscReal, Vec, TSMonitorVTKCtx);
508: PETSC_EXTERN PetscErrorCode        TSMonitorSolutionVTKDestroy(TSMonitorVTKCtx *);
509: PETSC_EXTERN PetscErrorCode        TSMonitorSolutionVTKCtxCreate(const char *, TSMonitorVTKCtx *);

511: PETSC_EXTERN PetscErrorCode TSStep(TS);
512: PETSC_EXTERN PetscErrorCode TSEvaluateWLTE(TS, NormType, PetscInt *, PetscReal *);
513: PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS, PetscInt, Vec, PetscBool *);
514: PETSC_EXTERN PetscErrorCode TSSolve(TS, Vec);
515: PETSC_EXTERN PetscErrorCode TSGetEquationType(TS, TSEquationType *);
516: PETSC_EXTERN PetscErrorCode TSSetEquationType(TS, TSEquationType);
517: PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS, TSConvergedReason *);
518: PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS, TSConvergedReason);
519: PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS, PetscReal *);
520: PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS, PetscInt *);
521: PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS, PetscInt *);
522: PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS, PetscInt *);
523: PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS, PetscInt);
524: PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS, PetscInt *);
525: PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS, PetscInt);
526: PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS, PetscBool);
527: PETSC_EXTERN PetscErrorCode TSRestartStep(TS);
528: PETSC_EXTERN PetscErrorCode TSRollBack(TS);
529: PETSC_EXTERN PetscErrorCode TSGetStepRollBack(TS, PetscBool *);

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

533: PETSC_EXTERN PetscErrorCode TSGetTime(TS, PetscReal *);
534: PETSC_EXTERN PetscErrorCode TSSetTime(TS, PetscReal);
535: PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS, PetscReal *);
536: PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS, PetscReal *);
537: PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS, PetscReal);
538: PETSC_EXTERN PetscErrorCode TSGetStepNumber(TS, PetscInt *);
539: PETSC_EXTERN PetscErrorCode TSSetStepNumber(TS, PetscInt);

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

544:   Calling Sequence:
545: + ts  - timestep context
546: . t   - current time
547: . u   - input vector
548: . F   - function vector
549: - ctx - [optional] user-defined function context

551:   Level: beginner

553:   Note:
554:   The deprecated `TSRHSFunction` still works as a replacement for `TSRHSFunctionFn` *.

556: .seealso: [](ch_ts), `TS`, `TSSetRHSFunction()`, `DMTSSetRHSFunction()`, `TSIFunctionFn`,
557: `TSIJacobianFn`, `TSRHSJacobianFn`
558: S*/
559: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSRHSFunctionFn(TS ts, PetscReal t, Vec u, Vec F, PetscCtx ctx);

561: PETSC_EXTERN_TYPEDEF typedef TSRHSFunctionFn *TSRHSFunction;

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

566:   Calling Sequence:
567: + ts   - the `TS` context obtained from `TSCreate()`
568: . t    - current time
569: . u    - input vector
570: . Amat - (approximate) Jacobian matrix
571: . Pmat - matrix from which preconditioner is to be constructed (usually the same as `Amat`)
572: - ctx  - [optional] user-defined context for matrix evaluation routine

574:   Level: beginner

576:   Note:
577:   The deprecated `TSRHSJacobian` still works as a replacement for `TSRHSJacobianFn` *.

579: .seealso: [](ch_ts), `TS`, `TSSetRHSJacobian()`, `DMTSSetRHSJacobian()`, `TSRHSFunctionFn`,
580: `TSIFunctionFn`, `TSIJacobianFn`
581: S*/
582: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSRHSJacobianFn(TS ts, PetscReal t, Vec u, Mat Amat, Mat Pmat, PetscCtx ctx);

584: PETSC_EXTERN_TYPEDEF typedef TSRHSJacobianFn *TSRHSJacobian;

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

590:   Calling Sequence:
591: + ts  - the `TS` context
592: . t   - current timestep
593: . U   - input vector (current ODE solution)
594: . A   - output matrix
595: - ctx - [optional] user-defined function context

597:   Level: beginner

599:   Note:
600:   The deprecated `TSRHSJacobianP` still works as a replacement for `TSRHSJacobianPFn` *.

602: .seealso: [](ch_ts), `TS`, `TSSetRHSJacobianP()`, `TSGetRHSJacobianP()`
603: S*/
604: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSRHSJacobianPFn(TS ts, PetscReal t, Vec U, Mat A, PetscCtx ctx);

606: PETSC_EXTERN_TYPEDEF typedef TSRHSJacobianPFn *TSRHSJacobianP;

608: PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS, Vec, TSRHSFunctionFn *, PetscCtx);
609: PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS, Vec *, TSRHSFunctionFn **, PetscCtxRt);
610: PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS, Mat, Mat, TSRHSJacobianFn *, PetscCtx);
611: PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS, Mat *, Mat *, TSRHSJacobianFn **, PetscCtxRt);
612: PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS, PetscBool);

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

617:   Calling Sequence:
618: + ts  - timestep context
619: . t   - current time
620: . u   - output vector
621: - ctx - [optional] user-defined function context

623:   Level: advanced

625:   Note:
626:   The deprecated `TSSolutionFunction` still works as a replacement for `TSSolutionFn` *.

628: .seealso: [](ch_ts), `TS`, `TSSetSolutionFunction()`, `DMTSSetSolutionFunction()`
629: S*/
630: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSSolutionFn(TS ts, PetscReal t, Vec u, PetscCtx ctx);

632: PETSC_EXTERN_TYPEDEF typedef TSSolutionFn *TSSolutionFunction;

634: PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS, TSSolutionFn *, PetscCtx);

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

639:   Calling Sequence:
640: + ts  - timestep context
641: . t   - current time
642: . f   - output vector
643: - ctx - [optional] user-defined function context

645:   Level: advanced

647:   Note:
648:   The deprecated `TSForcingFunction` still works as a replacement for `TSForcingFn` *.

650: .seealso: [](ch_ts), `TS`, `TSSetForcingFunction()`, `DMTSSetForcingFunction()`
651: S*/
652: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSForcingFn(TS ts, PetscReal t, Vec f, PetscCtx ctx);

654: PETSC_EXTERN_TYPEDEF typedef TSForcingFn *TSForcingFunction;

656: PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS, TSForcingFn *, PetscCtx);

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

661:   Calling Sequence:
662: + ts  - the `TS` context obtained from `TSCreate()`
663: . t   - time at step/stage being solved
664: . U   - state vector
665: . U_t - time derivative of state vector
666: . F   - function vector
667: - ctx - [optional] user-defined context for function

669:   Level: beginner

671:   Note:
672:   The deprecated `TSIFunction` still works as a replacement for `TSIFunctionFn` *.

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

678: PETSC_EXTERN_TYPEDEF typedef TSIFunctionFn *TSIFunction;

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

683:   Calling Sequence:
684: + ts   - the `TS` context obtained from `TSCreate()`
685: . t    - time at step/stage being solved
686: . U    - state vector
687: . U_t  - time derivative of state vector
688: . a    - shift
689: . Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
690: . Pmat - matrix used for constructing preconditioner, usually the same as `Amat`
691: - ctx  - [optional] user-defined context for Jacobian evaluation routine

693:   Level: beginner

695:   Note:
696:   The deprecated `TSIJacobian` still works as a replacement for `TSIJacobianFn` *.

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

702: PETSC_EXTERN_TYPEDEF typedef TSIJacobianFn *TSIJacobian;

704: PETSC_EXTERN PetscErrorCode TSSetIFunction(TS, Vec, TSIFunctionFn *, PetscCtx);
705: PETSC_EXTERN PetscErrorCode TSGetIFunction(TS, Vec *, TSIFunctionFn **, PetscCtxRt);
706: PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS, Mat, Mat, TSIJacobianFn *, PetscCtx);
707: PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS, Mat *, Mat *, TSIJacobianFn **, PetscCtxRt);

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

712:   Calling Sequence:
713: + ts   - the `TS` context obtained from `TSCreate()`
714: . t    - time at step/stage being solved
715: . U    - state vector
716: . U_t  - time derivative of state vector
717: . U_tt - second time derivative of state vector
718: . F    - function vector
719: - ctx  - [optional] user-defined context for matrix evaluation routine (may be `NULL`)

721:   Level: advanced

723:   Note:
724:   The deprecated `TSI2Function` still works as a replacement for `TSI2FunctionFn` *.

726: .seealso: [](ch_ts), `TS`, `TSSetI2Function()`, `DMTSSetI2Function()`, `TSIFunctionFn`
727: S*/
728: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSI2FunctionFn(TS ts, PetscReal t, Vec U, Vec U_t, Vec U_tt, Vec F, PetscCtx ctx);

730: PETSC_EXTERN_TYPEDEF typedef TSI2FunctionFn *TSI2Function;

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

735:   Calling Sequence:
736: + ts   - the `TS` context obtained from `TSCreate()`
737: . t    - time at step/stage being solved
738: . U    - state vector
739: . U_t  - time derivative of state vector
740: . U_tt - second time derivative of state vector
741: . v    - shift for U_t
742: . a    - shift for U_tt
743: . 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
744: . jac  - matrix from which to construct the preconditioner, may be same as `J`
745: - ctx  - [optional] user-defined context for matrix evaluation routine

747:   Level: advanced

749:   Note:
750:   The deprecated `TSI2Jacobian` still works as a replacement for `TSI2JacobianFn` *.

752: .seealso: [](ch_ts), `TS`, `TSSetI2Jacobian()`, `DMTSSetI2Jacobian()`, `TSIFunctionFn`, `TSIJacobianFn`, `TSRHSFunctionFn`, `TSRHSJacobianFn`
753: S*/
754: 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, PetscCtx ctx);

756: PETSC_EXTERN_TYPEDEF typedef TSI2JacobianFn *TSI2Jacobian;

758: PETSC_EXTERN PetscErrorCode TSSetI2Function(TS, Vec, TSI2FunctionFn *, PetscCtx);
759: PETSC_EXTERN PetscErrorCode TSGetI2Function(TS, Vec *, TSI2FunctionFn **, PetscCtxRt);
760: PETSC_EXTERN PetscErrorCode TSSetI2Jacobian(TS, Mat, Mat, TSI2JacobianFn *, PetscCtx);
761: PETSC_EXTERN PetscErrorCode TSGetI2Jacobian(TS, Mat *, Mat *, TSI2JacobianFn **, PetscCtxRt);

763: PETSC_EXTERN PetscErrorCode TSRHSSplitSetIS(TS, const char[], IS);
764: PETSC_EXTERN PetscErrorCode TSRHSSplitGetIS(TS, const char[], IS *);
765: PETSC_EXTERN PetscErrorCode TSRHSSplitSetRHSFunction(TS, const char[], Vec, TSRHSFunctionFn *, PetscCtx);
766: PETSC_EXTERN PetscErrorCode TSRHSSplitSetIFunction(TS, const char[], Vec, TSIFunctionFn *, PetscCtx);
767: PETSC_EXTERN PetscErrorCode TSRHSSplitSetIJacobian(TS, const char[], Mat, Mat, TSIJacobianFn *, PetscCtx);
768: PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTS(TS, const char[], TS *);
769: PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTSs(TS, PetscInt *, TS *[]);
770: PETSC_EXTERN PetscErrorCode TSSetUseSplitRHSFunction(TS, PetscBool);
771: PETSC_EXTERN PetscErrorCode TSGetUseSplitRHSFunction(TS, PetscBool *);
772: PETSC_EXTERN PetscErrorCode TSRHSSplitGetSNES(TS, SNES *);
773: PETSC_EXTERN PetscErrorCode TSRHSSplitSetSNES(TS, SNES);

775: PETSC_EXTERN TSRHSFunctionFn TSComputeRHSFunctionLinear;
776: PETSC_EXTERN TSRHSJacobianFn TSComputeRHSJacobianConstant;
777: PETSC_EXTERN PetscErrorCode  TSComputeIFunctionLinear(TS, PetscReal, Vec, Vec, Vec, PetscCtx);
778: PETSC_EXTERN PetscErrorCode  TSComputeIJacobianConstant(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscCtx);
779: PETSC_EXTERN PetscErrorCode  TSComputeSolutionFunction(TS, PetscReal, Vec);
780: PETSC_EXTERN PetscErrorCode  TSComputeForcingFunction(TS, PetscReal, Vec);
781: PETSC_EXTERN PetscErrorCode  TSComputeIJacobianDefaultColor(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscCtx);
782: PETSC_EXTERN PetscErrorCode  TSPruneIJacobianColor(TS, Mat, Mat);

784: PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
785: PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS, PetscReal));
786: PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS, PetscReal, PetscInt, Vec *));
787: PETSC_EXTERN PetscErrorCode TSSetPostEvaluate(TS, PetscErrorCode (*)(TS));
788: PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
789: PETSC_EXTERN PetscErrorCode TSSetResize(TS, PetscBool, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscBool *, PetscCtx), PetscErrorCode (*)(TS, PetscInt, Vec[], Vec[], PetscCtx), PetscCtx);
790: PETSC_EXTERN PetscErrorCode TSPreStep(TS);
791: PETSC_EXTERN PetscErrorCode TSPreStage(TS, PetscReal);
792: PETSC_EXTERN PetscErrorCode TSPostStage(TS, PetscReal, PetscInt, Vec[]);
793: PETSC_EXTERN PetscErrorCode TSPostEvaluate(TS);
794: PETSC_EXTERN PetscErrorCode TSPostStep(TS);
795: PETSC_EXTERN PetscErrorCode TSResize(TS);
796: PETSC_EXTERN PetscErrorCode TSResizeRetrieveVec(TS, const char *, Vec *);
797: PETSC_EXTERN PetscErrorCode TSResizeRegisterVec(TS, const char *, Vec);
798: PETSC_EXTERN PetscErrorCode TSGetStepResize(TS, PetscBool *);

800: PETSC_EXTERN PetscErrorCode TSInterpolate(TS, PetscReal, Vec);
801: PETSC_EXTERN PetscErrorCode TSSetTolerances(TS, PetscReal, Vec, PetscReal, Vec);
802: PETSC_EXTERN PetscErrorCode TSGetTolerances(TS, PetscReal *, Vec *, PetscReal *, Vec *);
803: PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *);
804: PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm(TS, Vec, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *);
805: PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS, PetscReal);
806: PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS, PetscReal *);
807: PETSC_EXTERN PetscErrorCode TSSetFunctionDomainError(TS, PetscErrorCode (*)(TS, PetscReal, Vec, PetscBool *));
808: PETSC_EXTERN PetscErrorCode TSFunctionDomainError(TS, PetscReal, Vec, PetscBool *);

810: PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS, PetscErrorCode (*)(TS, PetscReal *, PetscCtx), PetscCtx);
811: PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS, PetscReal *, PetscCtx);
812: PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS, PetscReal);
813: PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS, PetscErrorCode (*)(TS, Vec, PetscCtx, PetscReal *, PetscBool *), PetscCtx);
814: PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS, PetscReal);
815: PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);
816: PETSC_EXTERN PetscErrorCode TSPseudoComputeFunction(TS, Vec, Vec *, PetscReal *);

818: PETSC_EXTERN PetscErrorCode TSPythonSetType(TS, const char[]);
819: PETSC_EXTERN PetscErrorCode TSPythonGetType(TS, const char *[]);

821: PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS, PetscReal, Vec, Vec);
822: PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS, PetscReal, Vec, Mat, Mat);
823: PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS, PetscReal, Vec, Vec, Vec, PetscBool);
824: PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscBool);
825: PETSC_EXTERN PetscErrorCode TSComputeI2Function(TS, PetscReal, Vec, Vec, Vec, Vec);
826: PETSC_EXTERN PetscErrorCode TSComputeI2Jacobian(TS, PetscReal, Vec, Vec, Vec, PetscReal, PetscReal, Mat, Mat);
827: PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS, PetscReal, PetscReal, PetscReal *, PetscReal *);

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

831: PETSC_EXTERN PetscErrorCode DMTSSetBoundaryLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, PetscCtx), PetscCtx);
832: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM, TSRHSFunctionFn *, PetscCtx);
833: PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM, TSRHSFunctionFn **, PetscCtxRt);
834: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionContextDestroy(DM, PetscCtxDestroyFn *);
835: PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM, TSRHSJacobianFn *, PetscCtx);
836: PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM, TSRHSJacobianFn **, PetscCtxRt);
837: PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobianContextDestroy(DM, PetscCtxDestroyFn *);
838: PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM, TSIFunctionFn *, PetscCtx);
839: PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM, TSIFunctionFn **, PetscCtxRt);
840: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionContextDestroy(DM, PetscCtxDestroyFn *);
841: PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM, TSIJacobianFn *, PetscCtx);
842: PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM, TSIJacobianFn **, PetscCtxRt);
843: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianContextDestroy(DM, PetscCtxDestroyFn *);
844: PETSC_EXTERN PetscErrorCode DMTSSetI2Function(DM, TSI2FunctionFn *, PetscCtx);
845: PETSC_EXTERN PetscErrorCode DMTSGetI2Function(DM, TSI2FunctionFn **, PetscCtxRt);
846: PETSC_EXTERN PetscErrorCode DMTSSetI2FunctionContextDestroy(DM, PetscCtxDestroyFn *);
847: PETSC_EXTERN PetscErrorCode DMTSSetI2Jacobian(DM, TSI2JacobianFn *, PetscCtx);
848: PETSC_EXTERN PetscErrorCode DMTSGetI2Jacobian(DM, TSI2JacobianFn **, PetscCtxRt);
849: PETSC_EXTERN PetscErrorCode DMTSSetI2JacobianContextDestroy(DM, PetscCtxDestroyFn *);

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

854:   Calling Sequence:
855: + ts  - timestep context
856: . p   - input vector (primitive form)
857: . c   - output vector, transient variables (conservative form)
858: - ctx - [optional] user-defined function context

860:   Level: advanced

862:   Note:
863:   The deprecated `TSTransientVariable` still works as a replacement for `TSTransientVariableFn` *.

865: .seealso: [](ch_ts), `TS`, `TSSetTransientVariable()`, `DMTSSetTransientVariable()`
866: S*/
867: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSTransientVariableFn(TS ts, Vec p, Vec c, PetscCtx ctx);

869: PETSC_EXTERN_TYPEDEF typedef TSTransientVariableFn *TSTransientVariable;

871: PETSC_EXTERN PetscErrorCode TSSetTransientVariable(TS, TSTransientVariableFn *, PetscCtx);
872: PETSC_EXTERN PetscErrorCode DMTSSetTransientVariable(DM, TSTransientVariableFn *, PetscCtx);
873: PETSC_EXTERN PetscErrorCode DMTSGetTransientVariable(DM, TSTransientVariableFn **, PetscCtx);
874: PETSC_EXTERN PetscErrorCode TSComputeTransientVariable(TS, Vec, Vec);
875: PETSC_EXTERN PetscErrorCode TSHasTransientVariable(TS, PetscBool *);

877: PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM, TSSolutionFn *, PetscCtx);
878: PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM, TSSolutionFn **, PetscCtxRt);
879: PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM, TSForcingFn *, PetscCtx);
880: PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM, TSForcingFn **, PetscCtxRt);
881: PETSC_EXTERN PetscErrorCode DMTSCheckResidual(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscReal *);
882: PETSC_EXTERN PetscErrorCode DMTSCheckJacobian(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscBool *, PetscReal *);
883: PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec);

885: PETSC_EXTERN PetscErrorCode DMTSGetIFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, Vec, PetscCtx), PetscCtxRt);
886: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, Vec, PetscCtx), PetscCtx);
887: PETSC_EXTERN PetscErrorCode DMTSGetIJacobianLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscCtx), PetscCtxRt);
888: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscCtx), PetscCtx);
889: PETSC_EXTERN PetscErrorCode DMTSGetRHSFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, PetscCtx), PetscCtxRt);
890: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, PetscCtx), PetscCtx);
891: PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrix(DM);
892: PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrixLumped(DM);
893: PETSC_EXTERN PetscErrorCode DMTSDestroyRHSMassMatrix(DM);

895: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM, PetscErrorCode (*)(PetscCtx, PetscViewer), PetscErrorCode (*)(PetscCtx *, PetscViewer));
896: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM, PetscErrorCode (*)(PetscCtx, PetscViewer), PetscErrorCode (*)(PetscCtx *, PetscViewer));

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

901:   Calling Sequence:
902: + info - defines the subdomain to evaluate the residual on
903: . t    - time at which to evaluate residual
904: . x    - array of local state information
905: . f    - output array of local residual information
906: - ctx  - optional user context

908:   Level: beginner

910:   Note:
911:   The deprecated `DMDATSRHSFunctionLocal` still works as a replacement for `DMDATSRHSFunctionLocalFn` *.

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

917: PETSC_EXTERN_TYPEDEF typedef DMDATSRHSFunctionLocalFn *DMDATSRHSFunctionLocal;

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

922:   Calling Sequence:
923: + info - defines the subdomain to evaluate the residual on
924: . t    - time at which to evaluate residual
925: . x    - array of local state information
926: . J    - Jacobian matrix
927: . B    - matrix from which to construct the preconditioner; often same as `J`
928: - ctx  - optional context

930:   Level: beginner

932:   Note:
933:   The deprecated `DMDATSRHSJacobianLocal` still works as a replacement for `DMDATSRHSJacobianLocalFn` *.

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

939: PETSC_EXTERN_TYPEDEF typedef DMDATSRHSJacobianLocalFn *DMDATSRHSJacobianLocal;

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

944:   Calling Sequence:
945: + info  - defines the subdomain to evaluate the residual on
946: . t     - time at which to evaluate residual
947: . x     - array of local state information
948: . xdot  - array of local time derivative information
949: . imode - output array of local function evaluation information
950: - ctx   - optional context

952:   Level: beginner

954:   Note:
955:   The deprecated `DMDATSIFunctionLocal` still works as a replacement for `DMDATSIFunctionLocalFn` *.

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

961: PETSC_EXTERN_TYPEDEF typedef DMDATSIFunctionLocalFn *DMDATSIFunctionLocal;

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

966:   Calling Sequence:
967: + info  - defines the subdomain to evaluate the residual on
968: . t     - time at which to evaluate the jacobian
969: . x     - array of local state information
970: . xdot  - time derivative at this state
971: . shift - see `TSSetIJacobian()` for the meaning of this parameter
972: . J     - Jacobian matrix
973: . B     - matrix from which to construct the preconditioner; often same as `J`
974: - ctx   - optional context

976:   Level: beginner

978:   Note:
979:   The deprecated `DMDATSIJacobianLocal` still works as a replacement for `DMDATSIJacobianLocalFn` *.

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

985: PETSC_EXTERN_TYPEDEF typedef DMDATSIJacobianLocalFn *DMDATSIJacobianLocal;

987: PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM, InsertMode, DMDATSRHSFunctionLocalFn *, void *);
988: PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM, DMDATSRHSJacobianLocalFn *, void *);
989: PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM, InsertMode, DMDATSIFunctionLocalFn *, void *);
990: PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM, DMDATSIJacobianLocalFn *, void *);

992: typedef struct _n_TSMonitorLGCtx *TSMonitorLGCtx;
993: typedef struct {
994:   Vec            ray;
995:   VecScatter     scatter;
996:   PetscViewer    viewer;
997:   TSMonitorLGCtx lgctx;
998: } TSMonitorDMDARayCtx;
999: PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(PetscCtxRt);
1000: PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS, PetscInt, PetscReal, Vec, void *);
1001: PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS, PetscInt, PetscReal, Vec, void *);

1003: /* Dynamic creation and loading functions */
1004: PETSC_EXTERN PetscFunctionList TSList;
1005: PETSC_EXTERN PetscErrorCode    TSGetType(TS, TSType *);
1006: PETSC_EXTERN PetscErrorCode    TSSetType(TS, TSType);
1007: PETSC_EXTERN PetscErrorCode    TSRegister(const char[], PetscErrorCode (*)(TS));

1009: PETSC_EXTERN PetscErrorCode TSGetSNES(TS, SNES *);
1010: PETSC_EXTERN PetscErrorCode TSSetSNES(TS, SNES);
1011: PETSC_EXTERN PetscErrorCode TSGetKSP(TS, KSP *);

1013: PETSC_EXTERN PetscErrorCode TSView(TS, PetscViewer);
1014: PETSC_EXTERN PetscErrorCode TSLoad(TS, PetscViewer);
1015: PETSC_EXTERN PetscErrorCode TSViewFromOptions(TS, PetscObject, const char[]);
1016: PETSC_EXTERN PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory, PetscObject, const char[]);

1018: #define TS_FILE_CLASSID 1211225

1020: PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS, PetscCtx);
1021: PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS, PetscCtxRt);

1023: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtx *);
1024: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *);
1025: PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS, PetscInt, PetscReal, Vec, void *);
1026: PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS, PetscInt, PetscReal, Vec, void *);
1027: PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS, const char *const *);
1028: PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS, const char *const **);
1029: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx, const char *const *);
1030: PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS, const char *const *);
1031: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx, const char *const *);
1032: PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS, PetscErrorCode (*)(void *, Vec, Vec *), PetscCtxDestroyFn *, void *);
1033: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx, PetscErrorCode (*)(void *, Vec, Vec *), PetscCtxDestroyFn *, void *);
1034: PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS, PetscInt, PetscReal, Vec, void *);
1035: PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS, PetscInt, PetscReal, Vec, void *);
1036: PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS, PetscInt, PetscReal, Vec, void *);
1037: PETSC_EXTERN PetscErrorCode TSMonitorError(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
1038: PETSC_EXTERN PetscErrorCode TSDMSwarmMonitorMoments(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);

1040: struct _n_TSMonitorLGCtxNetwork {
1041:   PetscInt     nlg;
1042:   PetscDrawLG *lg;
1043:   PetscBool    semilogy;
1044:   PetscInt     howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
1045: };
1046: typedef struct _n_TSMonitorLGCtxNetwork *TSMonitorLGCtxNetwork;
1047: PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkDestroy(TSMonitorLGCtxNetwork *);
1048: PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkCreate(TS, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtxNetwork *);
1049: PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkSolution(TS, PetscInt, PetscReal, Vec, void *);

1051: typedef struct _n_TSMonitorEnvelopeCtx *TSMonitorEnvelopeCtx;
1052: PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeCtxCreate(TS, TSMonitorEnvelopeCtx *);
1053: PETSC_EXTERN PetscErrorCode             TSMonitorEnvelope(TS, PetscInt, PetscReal, Vec, void *);
1054: PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeGetBounds(TS, Vec *, Vec *);
1055: PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *);

1057: typedef struct _n_TSMonitorSPEigCtx *TSMonitorSPEigCtx;
1058: PETSC_EXTERN PetscErrorCode          TSMonitorSPEigCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorSPEigCtx *);
1059: PETSC_EXTERN PetscErrorCode          TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *);
1060: PETSC_EXTERN PetscErrorCode          TSMonitorSPEig(TS, PetscInt, PetscReal, Vec, void *);

1062: typedef struct _n_TSMonitorSPCtx *TSMonitorSPCtx;
1063: PETSC_EXTERN PetscErrorCode       TSMonitorSPCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscBool, PetscBool, TSMonitorSPCtx *);
1064: PETSC_EXTERN PetscErrorCode       TSMonitorSPCtxDestroy(TSMonitorSPCtx *);
1065: PETSC_EXTERN PetscErrorCode       TSMonitorSPSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);

1067: typedef struct _n_TSMonitorHGCtx *TSMonitorHGCtx;
1068: PETSC_EXTERN PetscErrorCode       TSMonitorHGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscInt, PetscBool, TSMonitorHGCtx *);
1069: PETSC_EXTERN PetscErrorCode       TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
1070: PETSC_EXTERN PetscErrorCode       TSMonitorHGCtxDestroy(TSMonitorHGCtx *);
1071: PETSC_EXTERN PetscErrorCode       TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);

1073: /*M
1074:    TSEvent - Abstract object to handle event detection in `TS` time integrator

1076:    Level: intermediate

1078:    Note:
1079:    See `TSSetEventHandler()` for the management of events.

1081: .seealso: [](sec_ts_event), `TS`, `TSSetEventHandler()`, `TSSetPostEventStep()`, `TSSetPostEventSecondStep()`, `TSSetEventTolerances()`, `TSGetNumEvents()`
1082: M*/
1083: typedef struct _n_TSEvent *TSEvent;

1085: PETSC_EXTERN PetscErrorCode TSSetEventHandler(TS, PetscInt, PetscInt[], PetscBool[], PetscErrorCode (*)(TS, PetscReal, Vec, PetscReal[], void *), PetscErrorCode (*)(TS, PetscInt, PetscInt[], PetscReal, Vec, PetscBool, void *), void *);
1086: PETSC_EXTERN PetscErrorCode TSSetPostEventStep(TS, PetscReal);
1087: PETSC_EXTERN PetscErrorCode TSSetPostEventSecondStep(TS, PetscReal);
1088: PETSC_DEPRECATED_FUNCTION(3, 21, 0, "TSSetPostEventSecondStep()", ) static inline PetscErrorCode TSSetPostEventIntervalStep(TS ts, PetscReal dt)
1089: {
1090:   return TSSetPostEventSecondStep(ts, dt);
1091: }
1092: PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS, PetscReal, PetscReal[]);
1093: PETSC_EXTERN PetscErrorCode TSGetNumEvents(TS, PetscInt *);

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

1098:    Level: beginner

1100: .seealso: [](ch_ts), `TSSSPSetType()`, `TS`, `TSSSP`
1101: J*/
1102: typedef const char *TSSSPType;
1103: #define TSSSPRKS2  "rks2"
1104: #define TSSSPRKS3  "rks3"
1105: #define TSSSPRK104 "rk104"

1107: PETSC_EXTERN PetscErrorCode    TSSSPSetType(TS, TSSSPType);
1108: PETSC_EXTERN PetscErrorCode    TSSSPGetType(TS, TSSSPType *);
1109: PETSC_EXTERN PetscErrorCode    TSSSPSetNumStages(TS, PetscInt);
1110: PETSC_EXTERN PetscErrorCode    TSSSPGetNumStages(TS, PetscInt *);
1111: PETSC_EXTERN PetscErrorCode    TSSSPInitializePackage(void);
1112: PETSC_EXTERN PetscErrorCode    TSSSPFinalizePackage(void);
1113: PETSC_EXTERN PetscFunctionList TSSSPList;

1115: /*S
1116:    TSAdapt - Abstract object that manages time-step adaptivity

1118:    Level: beginner

1120: .seealso: [](ch_ts), [](sec_ts_error_control), `TS`, `TSGetAdapt()`, `TSAdaptCreate()`, `TSAdaptType`
1121: S*/
1122: typedef struct _p_TSAdapt *TSAdapt;

1124: /*J
1125:    TSAdaptType - String with the name of `TSAdapt` scheme.

1127:    Level: beginner

1129: .seealso: [](ch_ts), [](sec_ts_error_control), `TSGetAdapt()`, `TSAdaptSetType()`, `TS`, `TSAdapt`
1130: J*/
1131: typedef const char *TSAdaptType;
1132: #define TSADAPTNONE    "none"
1133: #define TSADAPTBASIC   "basic"
1134: #define TSADAPTDSP     "dsp"
1135: #define TSADAPTCFL     "cfl"
1136: #define TSADAPTGLEE    "glee"
1137: #define TSADAPTHISTORY "history"

1139: PETSC_EXTERN PetscErrorCode TSGetAdapt(TS, TSAdapt *);
1140: PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[], PetscErrorCode (*)(TSAdapt));
1141: PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void);
1142: PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
1143: PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm, TSAdapt *);
1144: PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt, TSAdaptType);
1145: PETSC_EXTERN PetscErrorCode TSAdaptGetType(TSAdapt, TSAdaptType *);
1146: PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt, const char[]);
1147: PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
1148: PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt, const char[], PetscInt, PetscInt, PetscReal, PetscReal, PetscBool);
1149: PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt, PetscInt *, const PetscInt **, const PetscInt **, const PetscReal **, const PetscReal **);
1150: PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt, TS, PetscReal, PetscInt *, PetscReal *, PetscBool *);
1151: PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt, TS, PetscReal, Vec, PetscBool *);
1152: PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt, PetscViewer);
1153: PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt, PetscViewer);
1154: PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(TSAdapt, PetscOptionItems);
1155: PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt);
1156: PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt *);
1157: PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt, PetscBool);
1158: PETSC_EXTERN PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt, PetscBool);
1159: PETSC_EXTERN PetscErrorCode TSAdaptSetSafety(TSAdapt, PetscReal, PetscReal);
1160: PETSC_EXTERN PetscErrorCode TSAdaptGetSafety(TSAdapt, PetscReal *, PetscReal *);
1161: PETSC_EXTERN PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt, PetscReal);
1162: PETSC_EXTERN PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt, PetscReal *);
1163: PETSC_EXTERN PetscErrorCode TSAdaptSetClip(TSAdapt, PetscReal, PetscReal);
1164: PETSC_EXTERN PetscErrorCode TSAdaptGetClip(TSAdapt, PetscReal *, PetscReal *);
1165: PETSC_EXTERN PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt, PetscReal);
1166: PETSC_EXTERN PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt, PetscReal *);
1167: PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt, PetscReal, PetscReal);
1168: PETSC_EXTERN PetscErrorCode TSAdaptGetStepLimits(TSAdapt, PetscReal *, PetscReal *);
1169: PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt, PetscErrorCode (*)(TSAdapt, TS, PetscReal, Vec, PetscBool *));
1170: PETSC_EXTERN PetscErrorCode TSAdaptHistorySetHistory(TSAdapt, PetscInt, PetscReal[], PetscBool);
1171: PETSC_EXTERN PetscErrorCode TSAdaptHistorySetTrajectory(TSAdapt, TSTrajectory, PetscBool);
1172: PETSC_EXTERN PetscErrorCode TSAdaptHistoryGetStep(TSAdapt, PetscInt, PetscReal *, PetscReal *);
1173: PETSC_EXTERN PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt, PetscInt);
1174: PETSC_EXTERN PetscErrorCode TSAdaptDSPSetFilter(TSAdapt, const char *);
1175: PETSC_EXTERN PetscErrorCode TSAdaptDSPSetPID(TSAdapt, PetscReal, PetscReal, PetscReal);

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

1180:    Level: beginner

1182:    Developer Note:
1183:    This functionality should be replaced by the `TSAdapt`.

1185: .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLEAdaptCreate()`, `TSGLLEAdaptType`
1186: S*/
1187: typedef struct _p_TSGLLEAdapt *TSGLLEAdapt;

1189: /*J
1190:    TSGLLEAdaptType - String with the name of `TSGLLEAdapt` scheme

1192:    Level: beginner

1194:    Developer Note:
1195:    This functionality should be replaced by the `TSAdaptType`.

1197: .seealso: [](ch_ts), `TSGLLEAdaptSetType()`, `TS`
1198: J*/
1199: typedef const char *TSGLLEAdaptType;
1200: #define TSGLLEADAPT_NONE "none"
1201: #define TSGLLEADAPT_SIZE "size"
1202: #define TSGLLEADAPT_BOTH "both"

1204: PETSC_EXTERN PetscErrorCode TSGLLEAdaptRegister(const char[], PetscErrorCode (*)(TSGLLEAdapt));
1205: PETSC_EXTERN PetscErrorCode TSGLLEAdaptInitializePackage(void);
1206: PETSC_EXTERN PetscErrorCode TSGLLEAdaptFinalizePackage(void);
1207: PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate(MPI_Comm, TSGLLEAdapt *);
1208: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt, TSGLLEAdaptType);
1209: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt, const char[]);
1210: PETSC_EXTERN PetscErrorCode TSGLLEAdaptChoose(TSGLLEAdapt, PetscInt, const PetscInt[], const PetscReal[], const PetscReal[], PetscInt, PetscReal, PetscReal, PetscInt *, PetscReal *, PetscBool *);
1211: PETSC_EXTERN PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt, PetscViewer);
1212: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetFromOptions(TSGLLEAdapt, PetscOptionItems);
1213: PETSC_EXTERN PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt *);

1215: /*J
1216:    TSGLLEAcceptType - String with the name of `TSGLLEAccept` scheme

1218:    Level: beginner

1220: .seealso: [](ch_ts), `TSGLLESetAcceptType()`, `TS`, `TSGLLEAccept`
1221: J*/
1222: typedef const char *TSGLLEAcceptType;
1223: #define TSGLLEACCEPT_ALWAYS "always"

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

1228:   Calling Sequence:
1229: + ts     - timestep context
1230: . nt     - time to end of solution time
1231: . h      - the proposed step-size
1232: . enorm  - unknown
1233: - accept - output, if the proposal is accepted

1235:   Level: beginner

1237:   Note:
1238:   The deprecated `TSGLLEAcceptFunction` still works as a replacement for `TSGLLEAcceptFn` *

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

1245: PETSC_EXTERN_TYPEDEF typedef TSGLLEAcceptFn *TSGLLEAcceptFunction;

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

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

1252:   Level: beginner

1254: .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLESetType()`, `TSGLLERegister()`, `TSGLLEAccept`
1255: J*/
1256: typedef const char *TSGLLEType;
1257: #define TSGLLE_IRKS "irks"

1259: PETSC_EXTERN PetscErrorCode TSGLLERegister(const char[], PetscErrorCode (*)(TS));
1260: PETSC_EXTERN PetscErrorCode TSGLLEInitializePackage(void);
1261: PETSC_EXTERN PetscErrorCode TSGLLEFinalizePackage(void);
1262: PETSC_EXTERN PetscErrorCode TSGLLESetType(TS, TSGLLEType);
1263: PETSC_EXTERN PetscErrorCode TSGLLEGetAdapt(TS, TSGLLEAdapt *);
1264: PETSC_EXTERN PetscErrorCode TSGLLESetAcceptType(TS, TSGLLEAcceptType);

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

1269:    Level: beginner

1271: .seealso: [](ch_ts), `TSEIMEXSetType()`, `TS`, `TSEIMEX`, `TSEIMEXRegister()`
1272: J*/
1273: #define TSEIMEXType char *

1275: PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS, PetscInt);
1276: PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS, PetscInt, PetscInt);
1277: PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS, PetscBool);

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

1282:    Level: beginner

1284: .seealso: [](ch_ts), `TS`, `TSRKSetType()`, `TSRK`, `TSRKRegister()`
1285: J*/
1286: typedef const char *TSRKType;
1287: #define TSRK1FE "1fe"
1288: #define TSRK2A  "2a"
1289: #define TSRK2B  "2b"
1290: #define TSRK3   "3"
1291: #define TSRK3BS "3bs"
1292: #define TSRK4   "4"
1293: #define TSRK5F  "5f"
1294: #define TSRK5DP "5dp"
1295: #define TSRK5BS "5bs"
1296: #define TSRK6VR "6vr"
1297: #define TSRK7VR "7vr"
1298: #define TSRK8VR "8vr"

1300: PETSC_EXTERN PetscErrorCode TSRKGetOrder(TS, PetscInt *);
1301: PETSC_EXTERN PetscErrorCode TSRKGetType(TS, TSRKType *);
1302: PETSC_EXTERN PetscErrorCode TSRKSetType(TS, TSRKType);
1303: PETSC_EXTERN PetscErrorCode TSRKGetTableau(TS, PetscInt *, const PetscReal *[], const PetscReal *[], const PetscReal *[], const PetscReal *[], PetscInt *, const PetscReal *[], PetscBool *);
1304: PETSC_EXTERN PetscErrorCode TSRKSetMultirate(TS, PetscBool);
1305: PETSC_EXTERN PetscErrorCode TSRKGetMultirate(TS, PetscBool *);
1306: PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
1307: PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void);
1308: PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void);
1309: PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void);

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

1314:    Level: beginner

1316: .seealso: [](ch_ts), `TSMPRKSetType()`, `TS`, `TSMPRK`, `TSMPRKRegister()`
1317: J*/
1318: typedef const char *TSMPRKType;
1319: #define TSMPRK2A22 "2a22"
1320: #define TSMPRK2A23 "2a23"
1321: #define TSMPRK2A32 "2a32"
1322: #define TSMPRK2A33 "2a33"
1323: #define TSMPRKP2   "p2"
1324: #define TSMPRKP3   "p3"

1326: PETSC_EXTERN PetscErrorCode TSMPRKGetType(TS, TSMPRKType *);
1327: PETSC_EXTERN PetscErrorCode TSMPRKSetType(TS, TSMPRKType);
1328: 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[]);
1329: PETSC_EXTERN PetscErrorCode TSMPRKInitializePackage(void);
1330: PETSC_EXTERN PetscErrorCode TSMPRKFinalizePackage(void);
1331: PETSC_EXTERN PetscErrorCode TSMPRKRegisterDestroy(void);

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

1336:    Level: beginner

1338: .seealso: [](ch_ts), `TSIRKSetType()`, `TS`, `TSIRK`, `TSIRKRegister()`
1339: J*/
1340: typedef const char *TSIRKType;
1341: #define TSIRKGAUSS "gauss"

1343: PETSC_EXTERN PetscErrorCode TSIRKGetType(TS, TSIRKType *);
1344: PETSC_EXTERN PetscErrorCode TSIRKSetType(TS, TSIRKType);
1345: PETSC_EXTERN PetscErrorCode TSIRKGetNumStages(TS, PetscInt *);
1346: PETSC_EXTERN PetscErrorCode TSIRKSetNumStages(TS, PetscInt);
1347: PETSC_EXTERN PetscErrorCode TSIRKRegister(const char[], PetscErrorCode (*function)(TS));
1348: PETSC_EXTERN PetscErrorCode TSIRKTableauCreate(TS, PetscInt, const PetscReal *, const PetscReal *, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, const PetscScalar *);
1349: PETSC_EXTERN PetscErrorCode TSIRKInitializePackage(void);
1350: PETSC_EXTERN PetscErrorCode TSIRKFinalizePackage(void);
1351: PETSC_EXTERN PetscErrorCode TSIRKRegisterDestroy(void);

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

1356:    Level: beginner

1358: .seealso: [](ch_ts), `TSGLEESetType()`, `TS`, `TSGLEE`, `TSGLEERegister()`
1359: J*/
1360: typedef const char *TSGLEEType;
1361: #define TSGLEEi1      "BE1"
1362: #define TSGLEE23      "23"
1363: #define TSGLEE24      "24"
1364: #define TSGLEE25I     "25i"
1365: #define TSGLEE35      "35"
1366: #define TSGLEEEXRK2A  "exrk2a"
1367: #define TSGLEERK32G1  "rk32g1"
1368: #define TSGLEERK285EX "rk285ex"

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

1373:    Level: beginner

1375: .seealso: [](ch_ts), `TSGLEESetMode()`, `TS`, `TSGLEE`, `TSGLEERegister()`
1376: J*/
1377: PETSC_EXTERN PetscErrorCode TSGLEEGetType(TS, TSGLEEType *);
1378: PETSC_EXTERN PetscErrorCode TSGLEESetType(TS, TSGLEEType);
1379: 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[]);
1380: PETSC_EXTERN PetscErrorCode TSGLEERegisterAll(void);
1381: PETSC_EXTERN PetscErrorCode TSGLEEFinalizePackage(void);
1382: PETSC_EXTERN PetscErrorCode TSGLEEInitializePackage(void);
1383: PETSC_EXTERN PetscErrorCode TSGLEERegisterDestroy(void);

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

1388:    Level: beginner

1390: .seealso: [](ch_ts), `TSARKIMEXSetType()`, `TS`, `TSARKIMEX`, `TSARKIMEXRegister()`
1391: J*/
1392: typedef const char *TSARKIMEXType;
1393: #define TSARKIMEX1BEE   "1bee"
1394: #define TSARKIMEXA2     "a2"
1395: #define TSARKIMEXL2     "l2"
1396: #define TSARKIMEXARS122 "ars122"
1397: #define TSARKIMEX2C     "2c"
1398: #define TSARKIMEX2D     "2d"
1399: #define TSARKIMEX2E     "2e"
1400: #define TSARKIMEXPRSSP2 "prssp2"
1401: #define TSARKIMEX3      "3"
1402: #define TSARKIMEXBPR3   "bpr3"
1403: #define TSARKIMEXARS443 "ars443"
1404: #define TSARKIMEX4      "4"
1405: #define TSARKIMEX5      "5"

1407: PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS, TSARKIMEXType *);
1408: PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS, TSARKIMEXType);
1409: PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS, PetscBool);
1410: PETSC_EXTERN PetscErrorCode TSARKIMEXGetFullyImplicit(TS, PetscBool *);
1411: PETSC_EXTERN PetscErrorCode TSARKIMEXSetFastSlowSplit(TS, PetscBool);
1412: PETSC_EXTERN PetscErrorCode TSARKIMEXGetFastSlowSplit(TS, PetscBool *);
1413: 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[]);
1414: PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void);
1415: PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
1416: PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);

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

1421:    Level: beginner

1423: .seealso: [](ch_ts), `TSDIRKSetType()`, `TS`, `TSDIRK`, `TSDIRKRegister()`
1424: J*/
1425: typedef const char *TSDIRKType;
1426: #define TSDIRKS212      "s212"
1427: #define TSDIRKES122SAL  "es122sal"
1428: #define TSDIRKES213SAL  "es213sal"
1429: #define TSDIRKES324SAL  "es324sal"
1430: #define TSDIRKES325SAL  "es325sal"
1431: #define TSDIRK657A      "657a"
1432: #define TSDIRKES648SA   "es648sa"
1433: #define TSDIRK658A      "658a"
1434: #define TSDIRKS659A     "s659a"
1435: #define TSDIRK7510SAL   "7510sal"
1436: #define TSDIRKES7510SA  "es7510sa"
1437: #define TSDIRK759A      "759a"
1438: #define TSDIRKS7511SAL  "s7511sal"
1439: #define TSDIRK8614A     "8614a"
1440: #define TSDIRK8616SAL   "8616sal"
1441: #define TSDIRKES8516SAL "es8516sal"

1443: PETSC_EXTERN PetscErrorCode TSDIRKGetType(TS, TSDIRKType *);
1444: PETSC_EXTERN PetscErrorCode TSDIRKSetType(TS, TSDIRKType);
1445: PETSC_EXTERN PetscErrorCode TSDIRKRegister(TSDIRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);

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

1450:    Level: beginner

1452: .seealso: [](ch_ts), `TSRosWSetType()`, `TS`, `TSROSW`, `TSRosWRegister()`
1453: J*/
1454: typedef const char *TSRosWType;
1455: #define TSROSW2M          "2m"
1456: #define TSROSW2P          "2p"
1457: #define TSROSWRA3PW       "ra3pw"
1458: #define TSROSWRA34PW2     "ra34pw2"
1459: #define TSROSWR34PRW      "r34prw"
1460: #define TSROSWR3PRL2      "r3prl2"
1461: #define TSROSWRODAS3      "rodas3"
1462: #define TSROSWRODASPR     "rodaspr"
1463: #define TSROSWRODASPR2    "rodaspr2"
1464: #define TSROSWSANDU3      "sandu3"
1465: #define TSROSWASSP3P3S1C  "assp3p3s1c"
1466: #define TSROSWLASSP3P4S2C "lassp3p4s2c"
1467: #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
1468: #define TSROSWARK3        "ark3"
1469: #define TSROSWTHETA1      "theta1"
1470: #define TSROSWTHETA2      "theta2"
1471: #define TSROSWGRK4T       "grk4t"
1472: #define TSROSWSHAMP4      "shamp4"
1473: #define TSROSWVELDD4      "veldd4"
1474: #define TSROSW4L          "4l"

1476: PETSC_EXTERN PetscErrorCode TSRosWGetType(TS, TSRosWType *);
1477: PETSC_EXTERN PetscErrorCode TSRosWSetType(TS, TSRosWType);
1478: PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS, PetscBool);
1479: PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
1480: PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal);
1481: PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void);
1482: PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
1483: PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);

1485: PETSC_EXTERN PetscErrorCode TSBDFSetOrder(TS, PetscInt);
1486: PETSC_EXTERN PetscErrorCode TSBDFGetOrder(TS, PetscInt *);

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

1491:   Level: beginner

1493: .seealso: [](ch_ts), `TSBasicSymplecticSetType()`, `TS`, `TSBASICSYMPLECTIC`, `TSBasicSymplecticRegister()`
1494: J*/
1495: typedef const char *TSBasicSymplecticType;
1496: #define TSBASICSYMPLECTICSIEULER   "1"
1497: #define TSBASICSYMPLECTICVELVERLET "2"
1498: #define TSBASICSYMPLECTIC3         "3"
1499: #define TSBASICSYMPLECTIC4         "4"

1501: PETSC_EXTERN PetscErrorCode TSBasicSymplecticSetType(TS, TSBasicSymplecticType);
1502: PETSC_EXTERN PetscErrorCode TSBasicSymplecticGetType(TS, TSBasicSymplecticType *);
1503: PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegister(TSBasicSymplecticType, PetscInt, PetscInt, PetscReal[], PetscReal[]);
1504: PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterAll(void);
1505: PETSC_EXTERN PetscErrorCode TSBasicSymplecticInitializePackage(void);
1506: PETSC_EXTERN PetscErrorCode TSBasicSymplecticFinalizePackage(void);
1507: PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterDestroy(void);

1509: typedef enum {
1510:   TS_DG_GONZALEZ,
1511:   TS_DG_AVERAGE,
1512:   TS_DG_NONE
1513: } TSDGType;
1514: PETSC_EXTERN PetscErrorCode TSDiscGradSetFormulation(TS, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (*)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (*)(TS, PetscReal, Vec, Vec, void *), void *);
1515: PETSC_EXTERN PetscErrorCode TSDiscGradGetFormulation(TS, PetscErrorCode (**)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (**)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (**)(TS, PetscReal, Vec, Vec, void *), void *);
1516: PETSC_EXTERN PetscErrorCode TSDiscGradSetType(TS, TSDGType);
1517: PETSC_EXTERN PetscErrorCode TSDiscGradGetType(TS, TSDGType *);

1519: /*
1520:        PETSc interface to Sundials
1521: */
1522: #ifdef PETSC_HAVE_SUNDIALS2
1523: typedef enum {
1524:   SUNDIALS_ADAMS = 1,
1525:   SUNDIALS_BDF   = 2
1526: } TSSundialsLmmType;
1527: PETSC_EXTERN const char *const TSSundialsLmmTypes[];
1528: typedef enum {
1529:   SUNDIALS_MODIFIED_GS  = 1,
1530:   SUNDIALS_CLASSICAL_GS = 2
1531: } TSSundialsGramSchmidtType;
1532: PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];

1534: PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS, TSSundialsLmmType);
1535: PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS, PC *);
1536: PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS, PetscReal, PetscReal);
1537: PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS, PetscReal);
1538: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS, PetscReal);
1539: PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS, PetscInt *, PetscInt *);
1540: PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS, TSSundialsGramSchmidtType);
1541: PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS, PetscInt);
1542: PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS, PetscReal);
1543: PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS, PetscBool);
1544: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS, PetscInt);
1545: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxord(TS, PetscInt);
1546: PETSC_EXTERN PetscErrorCode TSSundialsSetUseDense(TS, PetscBool);
1547: #endif

1549: PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS, PetscReal);
1550: PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS, PetscReal *);
1551: PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS, PetscBool *);
1552: PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS, PetscBool);

1554: PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS, PetscReal);
1555: PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS, PetscReal, PetscReal, PetscReal);
1556: PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS, PetscReal *, PetscReal *, PetscReal *);

1558: PETSC_EXTERN PetscErrorCode TSAlpha2SetRadius(TS, PetscReal);
1559: PETSC_EXTERN PetscErrorCode TSAlpha2SetParams(TS, PetscReal, PetscReal, PetscReal, PetscReal);
1560: PETSC_EXTERN PetscErrorCode TSAlpha2GetParams(TS, PetscReal *, PetscReal *, PetscReal *, PetscReal *);

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

1566:   Calling Sequence:
1567: + ts   - the `TS` context obtained from `TSCreate()`
1568: . X0   - the previous time step's state vector
1569: . V0   - the previous time step's first derivative of the state vector
1570: . A0   - the previous time step's second derivative of the state vector
1571: . X1   - the vector into which the initial guess for the current time step will be written
1572: - ctx  - [optional] user-defined context for the predictor evaluation routine (may be `NULL`)

1574:   Level: intermediate

1576:   Note:
1577:   The deprecated `TSAlpha2Predictor` still works as a replacement for `TSAlpha2PredictorFn` *.

1579: .seealso: [](ch_ts), `TS`, `TSAlpha2SetPredictor()`
1580: S*/
1581: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSAlpha2PredictorFn(TS ts, Vec X0, Vec V0, Vec A0, Vec X1, PetscCtx ctx);

1583: PETSC_EXTERN_TYPEDEF typedef TSAlpha2PredictorFn *TSAlpha2Predictor;

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

1587: PETSC_EXTERN PetscErrorCode TSSetDM(TS, DM);
1588: PETSC_EXTERN PetscErrorCode TSGetDM(TS, DM *);

1590: PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES, Vec, Vec, void *);
1591: PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES, Vec, Mat, Mat, void *);

1593: PETSC_EXTERN PetscErrorCode TSRHSJacobianTest(TS, PetscBool *);
1594: PETSC_EXTERN PetscErrorCode TSRHSJacobianTestTranspose(TS, PetscBool *);

1596: PETSC_EXTERN PetscErrorCode TSGetComputeInitialCondition(TS, PetscErrorCode (**)(TS, Vec));
1597: PETSC_EXTERN PetscErrorCode TSSetComputeInitialCondition(TS, PetscErrorCode (*)(TS, Vec));
1598: PETSC_EXTERN PetscErrorCode TSComputeInitialCondition(TS, Vec);
1599: PETSC_EXTERN PetscErrorCode TSGetComputeExactError(TS, PetscErrorCode (**)(TS, Vec, Vec));
1600: PETSC_EXTERN PetscErrorCode TSSetComputeExactError(TS, PetscErrorCode (*)(TS, Vec, Vec));
1601: PETSC_EXTERN PetscErrorCode TSComputeExactError(TS, Vec, Vec);
1602: PETSC_EXTERN PetscErrorCode PetscConvEstUseTS(PetscConvEst, PetscBool);

1604: PETSC_EXTERN PetscErrorCode TSSetMatStructure(TS, MatStructure);