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 *, PetscCtxDestroyFn *);
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 *, PetscCtxDestroyFn *);
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 *);
439: typedef struct _n_TSMonitorVTKCtx *TSMonitorVTKCtx;
440: PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS, PetscInt, PetscReal, Vec, TSMonitorVTKCtx);
441: PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(TSMonitorVTKCtx *);
442: PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKCtxCreate(const char *, TSMonitorVTKCtx *);
444: PETSC_EXTERN PetscErrorCode TSStep(TS);
445: PETSC_EXTERN PetscErrorCode TSEvaluateWLTE(TS, NormType, PetscInt *, PetscReal *);
446: PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS, PetscInt, Vec, PetscBool *);
447: PETSC_EXTERN PetscErrorCode TSSolve(TS, Vec);
448: PETSC_EXTERN PetscErrorCode TSGetEquationType(TS, TSEquationType *);
449: PETSC_EXTERN PetscErrorCode TSSetEquationType(TS, TSEquationType);
450: PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS, TSConvergedReason *);
451: PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS, TSConvergedReason);
452: PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS, PetscReal *);
453: PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS, PetscInt *);
454: PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS, PetscInt *);
455: PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS, PetscInt *);
456: PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS, PetscInt);
457: PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS, PetscInt *);
458: PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS, PetscInt);
459: PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS, PetscBool);
460: PETSC_EXTERN PetscErrorCode TSRestartStep(TS);
461: PETSC_EXTERN PetscErrorCode TSRollBack(TS);
462: PETSC_EXTERN PetscErrorCode TSGetStepRollBack(TS, PetscBool *);
464: PETSC_EXTERN PetscErrorCode TSGetStages(TS, PetscInt *, Vec *[]);
466: PETSC_EXTERN PetscErrorCode TSGetTime(TS, PetscReal *);
467: PETSC_EXTERN PetscErrorCode TSSetTime(TS, PetscReal);
468: PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS, PetscReal *);
469: PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS, PetscReal *);
470: PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS, PetscReal);
471: PETSC_EXTERN PetscErrorCode TSGetStepNumber(TS, PetscInt *);
472: PETSC_EXTERN PetscErrorCode TSSetStepNumber(TS, PetscInt);
474: /*S
475: TSRHSFunctionFn - A prototype of a `TS` right-hand-side evaluation function that would be passed to `TSSetRHSFunction()`
477: Calling Sequence:
478: + ts - timestep context
479: . t - current time
480: . u - input vector
481: . F - function vector
482: - ctx - [optional] user-defined function context
484: Level: beginner
486: Note:
487: The deprecated `TSRHSFunction` still works as a replacement for `TSRHSFunctionFn` *.
489: .seealso: [](ch_ts), `TS`, `TSSetRHSFunction()`, `DMTSSetRHSFunction()`, `TSIFunctionFn`,
490: `TSIJacobianFn`, `TSRHSJacobianFn`
491: S*/
492: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSRHSFunctionFn)(TS ts, PetscReal t, Vec u, Vec F, void *ctx);
494: PETSC_EXTERN_TYPEDEF typedef TSRHSFunctionFn *TSRHSFunction;
496: /*S
497: TSRHSJacobianFn - A prototype of a `TS` right-hand-side Jacobian evaluation function that would be passed to `TSSetRHSJacobian()`
499: Calling Sequence:
500: + ts - the `TS` context obtained from `TSCreate()`
501: . t - current time
502: . u - input vector
503: . Amat - (approximate) Jacobian matrix
504: . Pmat - matrix from which preconditioner is to be constructed (usually the same as `Amat`)
505: - ctx - [optional] user-defined context for matrix evaluation routine
507: Level: beginner
509: Note:
510: The deprecated `TSRHSJacobian` still works as a replacement for `TSRHSJacobianFn` *.
512: .seealso: [](ch_ts), `TS`, `TSSetRHSJacobian()`, `DMTSSetRHSJacobian()`, `TSRHSFunctionFn`,
513: `TSIFunctionFn`, `TSIJacobianFn`
514: S*/
515: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSRHSJacobianFn)(TS ts, PetscReal t, Vec u, Mat Amat, Mat Pmat, void *ctx);
517: PETSC_EXTERN_TYPEDEF typedef TSRHSJacobianFn *TSRHSJacobian;
519: /*S
520: TSRHSJacobianPFn - A prototype of a function that computes the Jacobian of G w.r.t. the parameters P where
521: U_t = G(U,P,t), as well as the location to store the matrix that would be passed to `TSSetRHSJacobianP()`
523: Calling Sequence:
524: + ts - the `TS` context
525: . t - current timestep
526: . U - input vector (current ODE solution)
527: . A - output matrix
528: - ctx - [optional] user-defined function context
530: Level: beginner
532: Note:
533: The deprecated `TSRHSJacobianP` still works as a replacement for `TSRHSJacobianPFn` *.
535: .seealso: [](ch_ts), `TS`, `TSSetRHSJacobianP()`, `TSGetRHSJacobianP()`
536: S*/
537: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSRHSJacobianPFn)(TS ts, PetscReal t, Vec U, Mat A, void *ctx);
539: PETSC_EXTERN_TYPEDEF typedef TSRHSJacobianPFn *TSRHSJacobianP;
541: PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS, Vec, TSRHSFunctionFn *, void *);
542: PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS, Vec *, TSRHSFunctionFn **, void **);
543: PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS, Mat, Mat, TSRHSJacobianFn *, void *);
544: PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS, Mat *, Mat *, TSRHSJacobianFn **, void **);
545: PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS, PetscBool);
547: /*S
548: TSSolutionFn - A prototype of a `TS` solution evaluation function that would be passed to `TSSetSolutionFunction()`
550: Calling Sequence:
551: + ts - timestep context
552: . t - current time
553: . u - output vector
554: - ctx - [optional] user-defined function context
556: Level: advanced
558: Note:
559: The deprecated `TSSolutionFunction` still works as a replacement for `TSSolutionFn` *.
561: .seealso: [](ch_ts), `TS`, `TSSetSolutionFunction()`, `DMTSSetSolutionFunction()`
562: S*/
563: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSSolutionFn)(TS ts, PetscReal t, Vec u, void *ctx);
565: PETSC_EXTERN_TYPEDEF typedef TSSolutionFn *TSSolutionFunction;
567: PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS, TSSolutionFn *, void *);
569: /*S
570: TSForcingFn - A prototype of a `TS` forcing function evaluation function that would be passed to `TSSetForcingFunction()`
572: Calling Sequence:
573: + ts - timestep context
574: . t - current time
575: . f - output vector
576: - ctx - [optional] user-defined function context
578: Level: advanced
580: Note:
581: The deprecated `TSForcingFunction` still works as a replacement for `TSForcingFn` *.
583: .seealso: [](ch_ts), `TS`, `TSSetForcingFunction()`, `DMTSSetForcingFunction()`
584: S*/
585: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSForcingFn)(TS ts, PetscReal t, Vec f, void *ctx);
587: PETSC_EXTERN_TYPEDEF typedef TSForcingFn *TSForcingFunction;
589: PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS, TSForcingFn *, void *);
591: /*S
592: TSIFunctionFn - A prototype of a `TS` implicit function evaluation function that would be passed to `TSSetIFunction()
594: Calling Sequence:
595: + ts - the `TS` context obtained from `TSCreate()`
596: . t - time at step/stage being solved
597: . U - state vector
598: . U_t - time derivative of state vector
599: . F - function vector
600: - ctx - [optional] user-defined context for function
602: Level: beginner
604: Note:
605: The deprecated `TSIFunction` still works as a replacement for `TSIFunctionFn` *.
607: .seealso: [](ch_ts), `TS`, `TSSetIFunction()`, `DMTSSetIFunction()`, `TSIJacobianFn`, `TSRHSFunctionFn`, `TSRHSJacobianFn`
608: S*/
609: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSIFunctionFn)(TS ts, PetscReal t, Vec U, Vec U_t, Vec F, void *ctx);
611: PETSC_EXTERN_TYPEDEF typedef TSIFunctionFn *TSIFunction;
613: /*S
614: TSIJacobianFn - A prototype of a `TS` Jacobian evaluation function that would be passed to `TSSetIJacobian()`
616: Calling Sequence:
617: + ts - the `TS` context obtained from `TSCreate()`
618: . t - time at step/stage being solved
619: . U - state vector
620: . U_t - time derivative of state vector
621: . a - shift
622: . Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
623: . Pmat - matrix used for constructing preconditioner, usually the same as `Amat`
624: - ctx - [optional] user-defined context for Jacobian evaluation routine
626: Level: beginner
628: Note:
629: The deprecated `TSIJacobian` still works as a replacement for `TSIJacobianFn` *.
631: .seealso: [](ch_ts), `TSSetIJacobian()`, `DMTSSetIJacobian()`, `TSIFunctionFn`, `TSRHSFunctionFn`, `TSRHSJacobianFn`
632: S*/
633: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSIJacobianFn)(TS ts, PetscReal t, Vec U, Vec U_t, PetscReal a, Mat Amat, Mat Pmat, void *ctx);
635: PETSC_EXTERN_TYPEDEF typedef TSIJacobianFn *TSIJacobian;
637: PETSC_EXTERN PetscErrorCode TSSetIFunction(TS, Vec, TSIFunctionFn *, void *);
638: PETSC_EXTERN PetscErrorCode TSGetIFunction(TS, Vec *, TSIFunctionFn **, void **);
639: PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS, Mat, Mat, TSIJacobianFn *, void *);
640: PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS, Mat *, Mat *, TSIJacobianFn **, void **);
642: /*S
643: TSI2FunctionFn - A prototype of a `TS` implicit function evaluation function for 2nd order systems that would be passed to `TSSetI2Function()`
645: Calling Sequence:
646: + ts - the `TS` context obtained from `TSCreate()`
647: . t - time at step/stage being solved
648: . U - state vector
649: . U_t - time derivative of state vector
650: . U_tt - second time derivative of state vector
651: . F - function vector
652: - ctx - [optional] user-defined context for matrix evaluation routine (may be `NULL`)
654: Level: advanced
656: Note:
657: The deprecated `TSI2Function` still works as a replacement for `TSI2FunctionFn` *.
659: .seealso: [](ch_ts), `TS`, `TSSetI2Function()`, `DMTSSetI2Function()`, `TSIFunctionFn`
660: S*/
661: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSI2FunctionFn)(TS ts, PetscReal t, Vec U, Vec U_t, Vec U_tt, Vec F, void *ctx);
663: PETSC_EXTERN_TYPEDEF typedef TSI2FunctionFn *TSI2Function;
665: /*S
666: TSI2JacobianFn - A prototype of a `TS` implicit Jacobian evaluation function for 2nd order systems that would be passed to `TSSetI2Jacobian()`
668: Calling Sequence:
669: + ts - the `TS` context obtained from `TSCreate()`
670: . t - time at step/stage being solved
671: . U - state vector
672: . U_t - time derivative of state vector
673: . U_tt - second time derivative of state vector
674: . v - shift for U_t
675: . a - shift for U_tt
676: . 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
677: . jac - matrix from which to construct the preconditioner, may be same as `J`
678: - ctx - [optional] user-defined context for matrix evaluation routine
680: Level: advanced
682: Note:
683: The deprecated `TSI2Jacobian` still works as a replacement for `TSI2JacobianFn` *.
685: .seealso: [](ch_ts), `TS`, `TSSetI2Jacobian()`, `DMTSSetI2Jacobian()`, `TSIFunctionFn`, `TSIJacobianFn`, `TSRHSFunctionFn`, `TSRHSJacobianFn`
686: S*/
687: 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);
689: PETSC_EXTERN_TYPEDEF typedef TSI2JacobianFn *TSI2Jacobian;
691: PETSC_EXTERN PetscErrorCode TSSetI2Function(TS, Vec, TSI2FunctionFn *, void *);
692: PETSC_EXTERN PetscErrorCode TSGetI2Function(TS, Vec *, TSI2FunctionFn **, void **);
693: PETSC_EXTERN PetscErrorCode TSSetI2Jacobian(TS, Mat, Mat, TSI2JacobianFn *, void *);
694: PETSC_EXTERN PetscErrorCode TSGetI2Jacobian(TS, Mat *, Mat *, TSI2JacobianFn **, void **);
696: PETSC_EXTERN PetscErrorCode TSRHSSplitSetIS(TS, const char[], IS);
697: PETSC_EXTERN PetscErrorCode TSRHSSplitGetIS(TS, const char[], IS *);
698: PETSC_EXTERN PetscErrorCode TSRHSSplitSetRHSFunction(TS, const char[], Vec, TSRHSFunctionFn *, void *);
699: PETSC_EXTERN PetscErrorCode TSRHSSplitSetIFunction(TS, const char[], Vec, TSIFunctionFn *, void *);
700: PETSC_EXTERN PetscErrorCode TSRHSSplitSetIJacobian(TS, const char[], Mat, Mat, TSIJacobianFn *, void *);
701: PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTS(TS, const char[], TS *);
702: PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTSs(TS, PetscInt *, TS *[]);
703: PETSC_EXTERN PetscErrorCode TSSetUseSplitRHSFunction(TS, PetscBool);
704: PETSC_EXTERN PetscErrorCode TSGetUseSplitRHSFunction(TS, PetscBool *);
705: PETSC_EXTERN PetscErrorCode TSRHSSplitGetSNES(TS, SNES *);
706: PETSC_EXTERN PetscErrorCode TSRHSSplitSetSNES(TS, SNES);
708: PETSC_EXTERN TSRHSFunctionFn TSComputeRHSFunctionLinear;
709: PETSC_EXTERN TSRHSJacobianFn TSComputeRHSJacobianConstant;
710: PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS, PetscReal, Vec, Vec, Vec, void *);
711: PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
712: PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS, PetscReal, Vec);
713: PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS, PetscReal, Vec);
714: PETSC_EXTERN PetscErrorCode TSComputeIJacobianDefaultColor(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
715: PETSC_EXTERN PetscErrorCode TSPruneIJacobianColor(TS, Mat, Mat);
717: PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
718: PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS, PetscReal));
719: PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS, PetscReal, PetscInt, Vec *));
720: PETSC_EXTERN PetscErrorCode TSSetPostEvaluate(TS, PetscErrorCode (*)(TS));
721: PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
722: PETSC_EXTERN PetscErrorCode TSSetResize(TS, PetscBool, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscBool *, void *), PetscErrorCode (*)(TS, PetscInt, Vec[], Vec[], void *), void *);
723: PETSC_EXTERN PetscErrorCode TSPreStep(TS);
724: PETSC_EXTERN PetscErrorCode TSPreStage(TS, PetscReal);
725: PETSC_EXTERN PetscErrorCode TSPostStage(TS, PetscReal, PetscInt, Vec *);
726: PETSC_EXTERN PetscErrorCode TSPostEvaluate(TS);
727: PETSC_EXTERN PetscErrorCode TSPostStep(TS);
728: PETSC_EXTERN PetscErrorCode TSResize(TS);
729: PETSC_EXTERN PetscErrorCode TSResizeRetrieveVec(TS, const char *, Vec *);
730: PETSC_EXTERN PetscErrorCode TSResizeRegisterVec(TS, const char *, Vec);
731: PETSC_EXTERN PetscErrorCode TSGetStepResize(TS, PetscBool *);
733: PETSC_EXTERN PetscErrorCode TSInterpolate(TS, PetscReal, Vec);
734: PETSC_EXTERN PetscErrorCode TSSetTolerances(TS, PetscReal, Vec, PetscReal, Vec);
735: PETSC_EXTERN PetscErrorCode TSGetTolerances(TS, PetscReal *, Vec *, PetscReal *, Vec *);
736: PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *);
737: PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm(TS, Vec, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *);
738: PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS, PetscReal);
739: PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS, PetscReal *);
740: PETSC_EXTERN PetscErrorCode TSSetFunctionDomainError(TS, PetscErrorCode (*)(TS, PetscReal, Vec, PetscBool *));
741: PETSC_EXTERN PetscErrorCode TSFunctionDomainError(TS, PetscReal, Vec, PetscBool *);
743: PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS, PetscErrorCode (*)(TS, PetscReal *, void *), void *);
744: PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS, PetscReal *, void *);
745: PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS, PetscReal *);
746: PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS, PetscReal);
747: PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS, PetscErrorCode (*)(TS, Vec, void *, PetscReal *, PetscBool *), void *);
748: PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS, Vec, void *, PetscReal *, PetscBool *);
749: PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS, Vec, PetscReal *, PetscBool *);
750: PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS, PetscReal);
751: PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);
753: PETSC_EXTERN PetscErrorCode TSPythonSetType(TS, const char[]);
754: PETSC_EXTERN PetscErrorCode TSPythonGetType(TS, const char *[]);
756: PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS, PetscReal, Vec, Vec);
757: PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS, PetscReal, Vec, Mat, Mat);
758: PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS, PetscReal, Vec, Vec, Vec, PetscBool);
759: PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscBool);
760: PETSC_EXTERN PetscErrorCode TSComputeI2Function(TS, PetscReal, Vec, Vec, Vec, Vec);
761: PETSC_EXTERN PetscErrorCode TSComputeI2Jacobian(TS, PetscReal, Vec, Vec, Vec, PetscReal, PetscReal, Mat, Mat);
762: PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS, PetscReal, PetscReal, PetscReal *, PetscReal *);
764: PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS, Vec, Vec);
766: PETSC_EXTERN PetscErrorCode DMTSSetBoundaryLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *);
767: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM, TSRHSFunctionFn *, void *);
768: PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM, TSRHSFunctionFn **, void **);
769: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionContextDestroy(DM, PetscCtxDestroyFn *);
770: PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM, TSRHSJacobianFn *, void *);
771: PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM, TSRHSJacobianFn **, void **);
772: PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobianContextDestroy(DM, PetscCtxDestroyFn *);
773: PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM, TSIFunctionFn *, void *);
774: PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM, TSIFunctionFn **, void **);
775: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionContextDestroy(DM, PetscCtxDestroyFn *);
776: PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM, TSIJacobianFn *, void *);
777: PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM, TSIJacobianFn **, void **);
778: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianContextDestroy(DM, PetscCtxDestroyFn *);
779: PETSC_EXTERN PetscErrorCode DMTSSetI2Function(DM, TSI2FunctionFn *, void *);
780: PETSC_EXTERN PetscErrorCode DMTSGetI2Function(DM, TSI2FunctionFn **, void **);
781: PETSC_EXTERN PetscErrorCode DMTSSetI2FunctionContextDestroy(DM, PetscCtxDestroyFn *);
782: PETSC_EXTERN PetscErrorCode DMTSSetI2Jacobian(DM, TSI2JacobianFn *, void *);
783: PETSC_EXTERN PetscErrorCode DMTSGetI2Jacobian(DM, TSI2JacobianFn **, void **);
784: PETSC_EXTERN PetscErrorCode DMTSSetI2JacobianContextDestroy(DM, PetscCtxDestroyFn *);
786: /*S
787: TSTransientVariableFn - A prototype of a function to transform from state to transient variables that would be passed to `TSSetTransientVariable()`
789: Calling Sequence:
790: + ts - timestep context
791: . p - input vector (primitive form)
792: . c - output vector, transient variables (conservative form)
793: - ctx - [optional] user-defined function context
795: Level: advanced
797: Note:
798: The deprecated `TSTransientVariable` still works as a replacement for `TSTransientVariableFn` *.
800: .seealso: [](ch_ts), `TS`, `TSSetTransientVariable()`, `DMTSSetTransientVariable()`
801: S*/
802: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSTransientVariableFn)(TS ts, Vec p, Vec c, void *ctx);
804: PETSC_EXTERN_TYPEDEF typedef TSTransientVariableFn *TSTransientVariable;
806: PETSC_EXTERN PetscErrorCode TSSetTransientVariable(TS, TSTransientVariableFn *, void *);
807: PETSC_EXTERN PetscErrorCode DMTSSetTransientVariable(DM, TSTransientVariableFn *, void *);
808: PETSC_EXTERN PetscErrorCode DMTSGetTransientVariable(DM, TSTransientVariableFn **, void *);
809: PETSC_EXTERN PetscErrorCode TSComputeTransientVariable(TS, Vec, Vec);
810: PETSC_EXTERN PetscErrorCode TSHasTransientVariable(TS, PetscBool *);
812: PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM, TSSolutionFn *, void *);
813: PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM, TSSolutionFn **, void **);
814: PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM, TSForcingFn *, void *);
815: PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM, TSForcingFn **, void **);
816: PETSC_EXTERN PetscErrorCode DMTSCheckResidual(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscReal *);
817: PETSC_EXTERN PetscErrorCode DMTSCheckJacobian(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscBool *, PetscReal *);
818: PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec);
820: PETSC_EXTERN PetscErrorCode DMTSGetIFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, Vec, void *), void **);
821: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, Vec, void *), void *);
822: PETSC_EXTERN PetscErrorCode DMTSGetIJacobianLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void **);
823: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void *);
824: PETSC_EXTERN PetscErrorCode DMTSGetRHSFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, void *), void **);
825: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *);
826: PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrix(DM);
827: PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrixLumped(DM);
828: PETSC_EXTERN PetscErrorCode DMTSDestroyRHSMassMatrix(DM);
830: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM, PetscErrorCode (*)(void *, PetscViewer), PetscErrorCode (*)(void **, PetscViewer));
831: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM, PetscErrorCode (*)(void *, PetscViewer), PetscErrorCode (*)(void **, PetscViewer));
833: /*S
834: DMDATSRHSFunctionLocalFn - A prototype of a local `TS` right-hand side residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetRHSFunctionLocal()`
836: Calling Sequence:
837: + info - defines the subdomain to evaluate the residual on
838: . t - time at which to evaluate residual
839: . x - array of local state information
840: . f - output array of local residual information
841: - ctx - optional user context
843: Level: beginner
845: Note:
846: The deprecated `DMDATSRHSFunctionLocal` still works as a replacement for `DMDATSRHSFunctionLocalFn` *.
848: .seealso: `DMDA`, `DMDATSSetRHSFunctionLocal()`, `TSRHSFunctionFn`, `DMDATSRHSJacobianLocalFn`, `DMDATSIJacobianLocalFn`, `DMDATSIFunctionLocalFn`
849: S*/
850: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDATSRHSFunctionLocalFn)(DMDALocalInfo *info, PetscReal t, void *x, void *f, void *ctx);
852: PETSC_EXTERN_TYPEDEF typedef DMDATSRHSFunctionLocalFn *DMDATSRHSFunctionLocal;
854: /*S
855: DMDATSRHSJacobianLocalFn - A prototype of a local residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetRHSJacobianLocal()`
857: Calling Sequence:
858: + info - defines the subdomain to evaluate the residual on
859: . t - time at which to evaluate residual
860: . x - array of local state information
861: . J - Jacobian matrix
862: . B - matrix from which to construct the preconditioner; often same as `J`
863: - ctx - optional context
865: Level: beginner
867: Note:
868: The deprecated `DMDATSRHSJacobianLocal` still works as a replacement for `DMDATSRHSJacobianLocalFn` *.
870: .seealso: `DMDA`, `DMDATSSetRHSJacobianLocal()`, `TSRHSJacobianFn`, `DMDATSRHSFunctionLocalFn`, `DMDATSIJacobianLocalFn`, `DMDATSIFunctionLocalFn`
871: S*/
872: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDATSRHSJacobianLocalFn)(DMDALocalInfo *info, PetscReal t, void *x, Mat J, Mat B, void *ctx);
874: PETSC_EXTERN_TYPEDEF typedef DMDATSRHSJacobianLocalFn *DMDATSRHSJacobianLocal;
876: /*S
877: DMDATSIFunctionLocalFn - A prototype of a local residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetIFunctionLocal()`
879: Calling Sequence:
880: + info - defines the subdomain to evaluate the residual on
881: . t - time at which to evaluate residual
882: . x - array of local state information
883: . xdot - array of local time derivative information
884: . imode - output array of local function evaluation information
885: - ctx - optional context
887: Level: beginner
889: Note:
890: The deprecated `DMDATSIFunctionLocal` still works as a replacement for `DMDATSIFunctionLocalFn` *.
892: .seealso: `DMDA`, `DMDATSSetIFunctionLocal()`, `DMDATSIJacobianLocalFn`, `TSIFunctionFn`
893: S*/
894: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDATSIFunctionLocalFn)(DMDALocalInfo *info, PetscReal t, void *x, void *xdot, void *imode, void *ctx);
896: PETSC_EXTERN_TYPEDEF typedef DMDATSIFunctionLocalFn *DMDATSIFunctionLocal;
898: /*S
899: DMDATSIJacobianLocalFn - A prototype of a local residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetIJacobianLocal()`
901: Calling Sequence:
902: + info - defines the subdomain to evaluate the residual on
903: . t - time at which to evaluate the jacobian
904: . x - array of local state information
905: . xdot - time derivative at this state
906: . shift - see `TSSetIJacobian()` for the meaning of this parameter
907: . J - Jacobian matrix
908: . B - matrix from which to construct the preconditioner; often same as `J`
909: - ctx - optional context
911: Level: beginner
913: Note:
914: The deprecated `DMDATSIJacobianLocal` still works as a replacement for `DMDATSIJacobianLocalFn` *.
916: .seealso: `DMDA` `DMDATSSetIJacobianLocal()`, `TSIJacobianFn`, `DMDATSIFunctionLocalFn`, `DMDATSRHSFunctionLocalFn`, `DMDATSRHSJacobianlocal()`
917: S*/
918: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDATSIJacobianLocalFn)(DMDALocalInfo *info, PetscReal t, void *x, void *xdot, PetscReal shift, Mat J, Mat B, void *ctx);
920: PETSC_EXTERN_TYPEDEF typedef DMDATSIJacobianLocalFn *DMDATSIJacobianLocal;
922: PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM, InsertMode, DMDATSRHSFunctionLocalFn *, void *);
923: PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM, DMDATSRHSJacobianLocalFn *, void *);
924: PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM, InsertMode, DMDATSIFunctionLocalFn *, void *);
925: PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM, DMDATSIJacobianLocalFn *, void *);
927: typedef struct _n_TSMonitorLGCtx *TSMonitorLGCtx;
928: typedef struct {
929: Vec ray;
930: VecScatter scatter;
931: PetscViewer viewer;
932: TSMonitorLGCtx lgctx;
933: } TSMonitorDMDARayCtx;
934: PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void **);
935: PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS, PetscInt, PetscReal, Vec, void *);
936: PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS, PetscInt, PetscReal, Vec, void *);
938: /* Dynamic creation and loading functions */
939: PETSC_EXTERN PetscFunctionList TSList;
940: PETSC_EXTERN PetscErrorCode TSGetType(TS, TSType *);
941: PETSC_EXTERN PetscErrorCode TSSetType(TS, TSType);
942: PETSC_EXTERN PetscErrorCode TSRegister(const char[], PetscErrorCode (*)(TS));
944: PETSC_EXTERN PetscErrorCode TSGetSNES(TS, SNES *);
945: PETSC_EXTERN PetscErrorCode TSSetSNES(TS, SNES);
946: PETSC_EXTERN PetscErrorCode TSGetKSP(TS, KSP *);
948: PETSC_EXTERN PetscErrorCode TSView(TS, PetscViewer);
949: PETSC_EXTERN PetscErrorCode TSLoad(TS, PetscViewer);
950: PETSC_EXTERN PetscErrorCode TSViewFromOptions(TS, PetscObject, const char[]);
951: PETSC_EXTERN PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory, PetscObject, const char[]);
953: #define TS_FILE_CLASSID 1211225
955: PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS, void *);
956: PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS, void *);
958: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtx *);
959: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *);
960: PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS, PetscInt, PetscReal, Vec, void *);
961: PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS, PetscInt, PetscReal, Vec, void *);
962: PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS, const char *const *);
963: PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS, const char *const **);
964: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx, const char *const *);
965: PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS, const char *const *);
966: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx, const char *const *);
967: PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS, PetscErrorCode (*)(void *, Vec, Vec *), PetscCtxDestroyFn *, void *);
968: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx, PetscErrorCode (*)(void *, Vec, Vec *), PetscCtxDestroyFn *, void *);
969: PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS, PetscInt, PetscReal, Vec, void *);
970: PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS, PetscInt, PetscReal, Vec, void *);
971: PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS, PetscInt, PetscReal, Vec, void *);
972: PETSC_EXTERN PetscErrorCode TSMonitorError(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
973: PETSC_EXTERN PetscErrorCode TSDMSwarmMonitorMoments(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
975: struct _n_TSMonitorLGCtxNetwork {
976: PetscInt nlg;
977: PetscDrawLG *lg;
978: PetscBool semilogy;
979: PetscInt howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
980: };
981: typedef struct _n_TSMonitorLGCtxNetwork *TSMonitorLGCtxNetwork;
982: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxNetworkDestroy(TSMonitorLGCtxNetwork *);
983: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxNetworkCreate(TS, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtxNetwork *);
984: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxNetworkSolution(TS, PetscInt, PetscReal, Vec, void *);
986: typedef struct _n_TSMonitorEnvelopeCtx *TSMonitorEnvelopeCtx;
987: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxCreate(TS, TSMonitorEnvelopeCtx *);
988: PETSC_EXTERN PetscErrorCode TSMonitorEnvelope(TS, PetscInt, PetscReal, Vec, void *);
989: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeGetBounds(TS, Vec *, Vec *);
990: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *);
992: typedef struct _n_TSMonitorSPEigCtx *TSMonitorSPEigCtx;
993: PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorSPEigCtx *);
994: PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *);
995: PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS, PetscInt, PetscReal, Vec, void *);
997: typedef struct _n_TSMonitorSPCtx *TSMonitorSPCtx;
998: PETSC_EXTERN PetscErrorCode TSMonitorSPCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscBool, PetscBool, TSMonitorSPCtx *);
999: PETSC_EXTERN PetscErrorCode TSMonitorSPCtxDestroy(TSMonitorSPCtx *);
1000: PETSC_EXTERN PetscErrorCode TSMonitorSPSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
1002: typedef struct _n_TSMonitorHGCtx *TSMonitorHGCtx;
1003: PETSC_EXTERN PetscErrorCode TSMonitorHGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscInt, PetscBool, TSMonitorHGCtx *);
1004: PETSC_EXTERN PetscErrorCode TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
1005: PETSC_EXTERN PetscErrorCode TSMonitorHGCtxDestroy(TSMonitorHGCtx *);
1006: PETSC_EXTERN PetscErrorCode TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
1008: PETSC_EXTERN PetscErrorCode TSSetEventHandler(TS, PetscInt, PetscInt[], PetscBool[], PetscErrorCode (*)(TS, PetscReal, Vec, PetscReal[], void *), PetscErrorCode (*)(TS, PetscInt, PetscInt[], PetscReal, Vec, PetscBool, void *), void *);
1009: PETSC_EXTERN PetscErrorCode TSSetPostEventStep(TS, PetscReal);
1010: PETSC_EXTERN PetscErrorCode TSSetPostEventSecondStep(TS, PetscReal);
1011: PETSC_DEPRECATED_FUNCTION(3, 21, 0, "TSSetPostEventSecondStep()", ) static inline PetscErrorCode TSSetPostEventIntervalStep(TS ts, PetscReal dt)
1012: {
1013: return TSSetPostEventSecondStep(ts, dt);
1014: }
1015: PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS, PetscReal, PetscReal[]);
1016: PETSC_EXTERN PetscErrorCode TSGetNumEvents(TS, PetscInt *);
1018: /*J
1019: TSSSPType - string with the name of a `TSSSP` scheme.
1021: Level: beginner
1023: .seealso: [](ch_ts), `TSSSPSetType()`, `TS`, `TSSSP`
1024: J*/
1025: typedef const char *TSSSPType;
1026: #define TSSSPRKS2 "rks2"
1027: #define TSSSPRKS3 "rks3"
1028: #define TSSSPRK104 "rk104"
1030: PETSC_EXTERN PetscErrorCode TSSSPSetType(TS, TSSSPType);
1031: PETSC_EXTERN PetscErrorCode TSSSPGetType(TS, TSSSPType *);
1032: PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS, PetscInt);
1033: PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS, PetscInt *);
1034: PETSC_EXTERN PetscErrorCode TSSSPInitializePackage(void);
1035: PETSC_EXTERN PetscErrorCode TSSSPFinalizePackage(void);
1036: PETSC_EXTERN PetscFunctionList TSSSPList;
1038: /*S
1039: TSAdapt - Abstract object that manages time-step adaptivity
1041: Level: beginner
1043: .seealso: [](ch_ts), `TS`, `TSGetAdapt()`, `TSAdaptCreate()`, `TSAdaptType`
1044: S*/
1045: typedef struct _p_TSAdapt *TSAdapt;
1047: /*J
1048: TSAdaptType - String with the name of `TSAdapt` scheme.
1050: Level: beginner
1052: .seealso: [](ch_ts), `TSGetAdapt()`, `TSAdaptSetType()`, `TS`, `TSAdapt`
1053: J*/
1054: typedef const char *TSAdaptType;
1055: #define TSADAPTNONE "none"
1056: #define TSADAPTBASIC "basic"
1057: #define TSADAPTDSP "dsp"
1058: #define TSADAPTCFL "cfl"
1059: #define TSADAPTGLEE "glee"
1060: #define TSADAPTHISTORY "history"
1062: PETSC_EXTERN PetscErrorCode TSGetAdapt(TS, TSAdapt *);
1063: PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[], PetscErrorCode (*)(TSAdapt));
1064: PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void);
1065: PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
1066: PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm, TSAdapt *);
1067: PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt, TSAdaptType);
1068: PETSC_EXTERN PetscErrorCode TSAdaptGetType(TSAdapt, TSAdaptType *);
1069: PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt, const char[]);
1070: PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
1071: PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt, const char[], PetscInt, PetscInt, PetscReal, PetscReal, PetscBool);
1072: PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt, PetscInt *, const PetscInt **, const PetscInt **, const PetscReal **, const PetscReal **);
1073: PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt, TS, PetscReal, PetscInt *, PetscReal *, PetscBool *);
1074: PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt, TS, PetscReal, Vec, PetscBool *);
1075: PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt, PetscViewer);
1076: PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt, PetscViewer);
1077: PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(TSAdapt, PetscOptionItems *);
1078: PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt);
1079: PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt *);
1080: PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt, PetscBool);
1081: PETSC_EXTERN PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt, PetscBool);
1082: PETSC_EXTERN PetscErrorCode TSAdaptSetSafety(TSAdapt, PetscReal, PetscReal);
1083: PETSC_EXTERN PetscErrorCode TSAdaptGetSafety(TSAdapt, PetscReal *, PetscReal *);
1084: PETSC_EXTERN PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt, PetscReal);
1085: PETSC_EXTERN PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt, PetscReal *);
1086: PETSC_EXTERN PetscErrorCode TSAdaptSetClip(TSAdapt, PetscReal, PetscReal);
1087: PETSC_EXTERN PetscErrorCode TSAdaptGetClip(TSAdapt, PetscReal *, PetscReal *);
1088: PETSC_EXTERN PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt, PetscReal);
1089: PETSC_EXTERN PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt, PetscReal *);
1090: PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt, PetscReal, PetscReal);
1091: PETSC_EXTERN PetscErrorCode TSAdaptGetStepLimits(TSAdapt, PetscReal *, PetscReal *);
1092: PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt, PetscErrorCode (*)(TSAdapt, TS, PetscReal, Vec, PetscBool *));
1093: PETSC_EXTERN PetscErrorCode TSAdaptHistorySetHistory(TSAdapt, PetscInt n, PetscReal hist[], PetscBool);
1094: PETSC_EXTERN PetscErrorCode TSAdaptHistorySetTrajectory(TSAdapt, TSTrajectory, PetscBool);
1095: PETSC_EXTERN PetscErrorCode TSAdaptHistoryGetStep(TSAdapt, PetscInt, PetscReal *, PetscReal *);
1096: PETSC_EXTERN PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt, PetscInt);
1097: PETSC_EXTERN PetscErrorCode TSAdaptDSPSetFilter(TSAdapt, const char *);
1098: PETSC_EXTERN PetscErrorCode TSAdaptDSPSetPID(TSAdapt, PetscReal, PetscReal, PetscReal);
1100: /*S
1101: TSGLLEAdapt - Abstract object that manages time-step adaptivity for `TSGLLE`
1103: Level: beginner
1105: Developer Note:
1106: This functionality should be replaced by the `TSAdapt`.
1108: .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLEAdaptCreate()`, `TSGLLEAdaptType`
1109: S*/
1110: typedef struct _p_TSGLLEAdapt *TSGLLEAdapt;
1112: /*J
1113: TSGLLEAdaptType - String with the name of `TSGLLEAdapt` scheme
1115: Level: beginner
1117: Developer Note:
1118: This functionality should be replaced by the `TSAdaptType`.
1120: .seealso: [](ch_ts), `TSGLLEAdaptSetType()`, `TS`
1121: J*/
1122: typedef const char *TSGLLEAdaptType;
1123: #define TSGLLEADAPT_NONE "none"
1124: #define TSGLLEADAPT_SIZE "size"
1125: #define TSGLLEADAPT_BOTH "both"
1127: PETSC_EXTERN PetscErrorCode TSGLLEAdaptRegister(const char[], PetscErrorCode (*)(TSGLLEAdapt));
1128: PETSC_EXTERN PetscErrorCode TSGLLEAdaptInitializePackage(void);
1129: PETSC_EXTERN PetscErrorCode TSGLLEAdaptFinalizePackage(void);
1130: PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate(MPI_Comm, TSGLLEAdapt *);
1131: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt, TSGLLEAdaptType);
1132: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt, const char[]);
1133: PETSC_EXTERN PetscErrorCode TSGLLEAdaptChoose(TSGLLEAdapt, PetscInt, const PetscInt[], const PetscReal[], const PetscReal[], PetscInt, PetscReal, PetscReal, PetscInt *, PetscReal *, PetscBool *);
1134: PETSC_EXTERN PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt, PetscViewer);
1135: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetFromOptions(TSGLLEAdapt, PetscOptionItems *);
1136: PETSC_EXTERN PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt *);
1138: /*J
1139: TSGLLEAcceptType - String with the name of `TSGLLEAccept` scheme
1141: Level: beginner
1143: .seealso: [](ch_ts), `TSGLLESetAcceptType()`, `TS`, `TSGLLEAccept`
1144: J*/
1145: typedef const char *TSGLLEAcceptType;
1146: #define TSGLLEACCEPT_ALWAYS "always"
1148: /*S
1149: TSGLLEAcceptFn - A prototype of a `TS` accept function that would be passed to `TSGLLEAcceptRegister()`
1151: Calling Sequence:
1152: + ts - timestep context
1153: . nt - time to end of solution time
1154: . h - the proposed step-size
1155: . enorm - unknown
1156: - accept - output, if the proposal is accepted
1158: Level: beginner
1160: Note:
1161: The deprecated `TSGLLEAcceptFunction` still works as a replacement for `TSGLLEAcceptFn` *
1163: .seealso: [](ch_ts), `TS`, `TSSetRHSFunction()`, `DMTSSetRHSFunction()`, `TSIFunctionFn`,
1164: `TSIJacobianFn`, `TSRHSJacobianFn`, `TSGLLEAcceptRegister()`
1165: S*/
1166: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSGLLEAcceptFn)(TS ts, PetscReal nt, PetscReal h, const PetscReal enorm[], PetscBool *accept);
1168: PETSC_EXTERN_TYPEDEF typedef TSGLLEAcceptFn *TSGLLEAcceptFunction;
1170: PETSC_EXTERN PetscErrorCode TSGLLEAcceptRegister(const char[], TSGLLEAcceptFn *);
1172: /*J
1173: TSGLLEType - string with the name of a General Linear `TSGLLE` type
1175: Level: beginner
1177: .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLESetType()`, `TSGLLERegister()`, `TSGLLEAccept`
1178: J*/
1179: typedef const char *TSGLLEType;
1180: #define TSGLLE_IRKS "irks"
1182: PETSC_EXTERN PetscErrorCode TSGLLERegister(const char[], PetscErrorCode (*)(TS));
1183: PETSC_EXTERN PetscErrorCode TSGLLEInitializePackage(void);
1184: PETSC_EXTERN PetscErrorCode TSGLLEFinalizePackage(void);
1185: PETSC_EXTERN PetscErrorCode TSGLLESetType(TS, TSGLLEType);
1186: PETSC_EXTERN PetscErrorCode TSGLLEGetAdapt(TS, TSGLLEAdapt *);
1187: PETSC_EXTERN PetscErrorCode TSGLLESetAcceptType(TS, TSGLLEAcceptType);
1189: /*J
1190: TSEIMEXType - String with the name of an Extrapolated IMEX `TSEIMEX` type
1192: Level: beginner
1194: .seealso: [](ch_ts), `TSEIMEXSetType()`, `TS`, `TSEIMEX`, `TSEIMEXRegister()`
1195: J*/
1196: #define TSEIMEXType char *
1198: PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts, PetscInt);
1199: PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts, PetscInt, PetscInt);
1200: PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS, PetscBool);
1202: /*J
1203: TSRKType - String with the name of a Runge-Kutta `TSRK` type
1205: Level: beginner
1207: .seealso: [](ch_ts), `TS`, `TSRKSetType()`, `TS`, `TSRK`, `TSRKRegister()`
1208: J*/
1209: typedef const char *TSRKType;
1210: #define TSRK1FE "1fe"
1211: #define TSRK2A "2a"
1212: #define TSRK2B "2b"
1213: #define TSRK3 "3"
1214: #define TSRK3BS "3bs"
1215: #define TSRK4 "4"
1216: #define TSRK5F "5f"
1217: #define TSRK5DP "5dp"
1218: #define TSRK5BS "5bs"
1219: #define TSRK6VR "6vr"
1220: #define TSRK7VR "7vr"
1221: #define TSRK8VR "8vr"
1223: PETSC_EXTERN PetscErrorCode TSRKGetOrder(TS, PetscInt *);
1224: PETSC_EXTERN PetscErrorCode TSRKGetType(TS, TSRKType *);
1225: PETSC_EXTERN PetscErrorCode TSRKSetType(TS, TSRKType);
1226: PETSC_EXTERN PetscErrorCode TSRKGetTableau(TS, PetscInt *, const PetscReal **, const PetscReal **, const PetscReal **, const PetscReal **, PetscInt *, const PetscReal **, PetscBool *);
1227: PETSC_EXTERN PetscErrorCode TSRKSetMultirate(TS, PetscBool);
1228: PETSC_EXTERN PetscErrorCode TSRKGetMultirate(TS, PetscBool *);
1229: PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
1230: PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void);
1231: PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void);
1232: PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void);
1234: /*J
1235: TSMPRKType - String with the name of a partitioned Runge-Kutta `TSMPRK` type
1237: Level: beginner
1239: .seealso: [](ch_ts), `TSMPRKSetType()`, `TS`, `TSMPRK`, `TSMPRKRegister()`
1240: J*/
1241: typedef const char *TSMPRKType;
1242: #define TSMPRK2A22 "2a22"
1243: #define TSMPRK2A23 "2a23"
1244: #define TSMPRK2A32 "2a32"
1245: #define TSMPRK2A33 "2a33"
1246: #define TSMPRKP2 "p2"
1247: #define TSMPRKP3 "p3"
1249: PETSC_EXTERN PetscErrorCode TSMPRKGetType(TS ts, TSMPRKType *);
1250: PETSC_EXTERN PetscErrorCode TSMPRKSetType(TS ts, TSMPRKType);
1251: 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[]);
1252: PETSC_EXTERN PetscErrorCode TSMPRKInitializePackage(void);
1253: PETSC_EXTERN PetscErrorCode TSMPRKFinalizePackage(void);
1254: PETSC_EXTERN PetscErrorCode TSMPRKRegisterDestroy(void);
1256: /*J
1257: TSIRKType - String with the name of an implicit Runge-Kutta `TSIRK` type
1259: Level: beginner
1261: .seealso: [](ch_ts), `TSIRKSetType()`, `TS`, `TSIRK`, `TSIRKRegister()`
1262: J*/
1263: typedef const char *TSIRKType;
1264: #define TSIRKGAUSS "gauss"
1266: PETSC_EXTERN PetscErrorCode TSIRKGetType(TS, TSIRKType *);
1267: PETSC_EXTERN PetscErrorCode TSIRKSetType(TS, TSIRKType);
1268: PETSC_EXTERN PetscErrorCode TSIRKGetNumStages(TS, PetscInt *);
1269: PETSC_EXTERN PetscErrorCode TSIRKSetNumStages(TS, PetscInt);
1270: PETSC_EXTERN PetscErrorCode TSIRKRegister(const char[], PetscErrorCode (*function)(TS));
1271: PETSC_EXTERN PetscErrorCode TSIRKTableauCreate(TS, PetscInt, const PetscReal *, const PetscReal *, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, const PetscScalar *);
1272: PETSC_EXTERN PetscErrorCode TSIRKInitializePackage(void);
1273: PETSC_EXTERN PetscErrorCode TSIRKFinalizePackage(void);
1274: PETSC_EXTERN PetscErrorCode TSIRKRegisterDestroy(void);
1276: /*J
1277: TSGLEEType - String with the name of a General Linear with Error Estimation `TSGLEE` type
1279: Level: beginner
1281: .seealso: [](ch_ts), `TSGLEESetType()`, `TS`, `TSGLEE`, `TSGLEERegister()`
1282: J*/
1283: typedef const char *TSGLEEType;
1284: #define TSGLEEi1 "BE1"
1285: #define TSGLEE23 "23"
1286: #define TSGLEE24 "24"
1287: #define TSGLEE25I "25i"
1288: #define TSGLEE35 "35"
1289: #define TSGLEEEXRK2A "exrk2a"
1290: #define TSGLEERK32G1 "rk32g1"
1291: #define TSGLEERK285EX "rk285ex"
1293: /*J
1294: TSGLEEMode - String with the mode of error estimation for a General Linear with Error Estimation `TSGLEE` type
1296: Level: beginner
1298: .seealso: [](ch_ts), `TSGLEESetMode()`, `TS`, `TSGLEE`, `TSGLEERegister()`
1299: J*/
1300: PETSC_EXTERN PetscErrorCode TSGLEEGetType(TS ts, TSGLEEType *);
1301: PETSC_EXTERN PetscErrorCode TSGLEESetType(TS ts, TSGLEEType);
1302: 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[]);
1303: PETSC_EXTERN PetscErrorCode TSGLEERegisterAll(void);
1304: PETSC_EXTERN PetscErrorCode TSGLEEFinalizePackage(void);
1305: PETSC_EXTERN PetscErrorCode TSGLEEInitializePackage(void);
1306: PETSC_EXTERN PetscErrorCode TSGLEERegisterDestroy(void);
1308: /*J
1309: TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX `TSARKIMEX` type
1311: Level: beginner
1313: .seealso: [](ch_ts), `TSARKIMEXSetType()`, `TS`, `TSARKIMEX`, `TSARKIMEXRegister()`
1314: J*/
1315: typedef const char *TSARKIMEXType;
1316: #define TSARKIMEX1BEE "1bee"
1317: #define TSARKIMEXA2 "a2"
1318: #define TSARKIMEXL2 "l2"
1319: #define TSARKIMEXARS122 "ars122"
1320: #define TSARKIMEX2C "2c"
1321: #define TSARKIMEX2D "2d"
1322: #define TSARKIMEX2E "2e"
1323: #define TSARKIMEXPRSSP2 "prssp2"
1324: #define TSARKIMEX3 "3"
1325: #define TSARKIMEXBPR3 "bpr3"
1326: #define TSARKIMEXARS443 "ars443"
1327: #define TSARKIMEX4 "4"
1328: #define TSARKIMEX5 "5"
1329: PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts, TSARKIMEXType *);
1330: PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts, TSARKIMEXType);
1331: PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS, PetscBool);
1332: PETSC_EXTERN PetscErrorCode TSARKIMEXGetFullyImplicit(TS, PetscBool *);
1333: PETSC_EXTERN PetscErrorCode TSARKIMEXSetFastSlowSplit(TS, PetscBool);
1334: PETSC_EXTERN PetscErrorCode TSARKIMEXGetFastSlowSplit(TS, PetscBool *);
1335: 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[]);
1336: PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void);
1337: PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
1338: PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);
1340: /*J
1341: TSDIRKType - String with the name of a Diagonally Implicit Runge-Kutta `TSDIRK` type
1343: Level: beginner
1345: .seealso: [](ch_ts), `TSDIRKSetType()`, `TS`, `TSDIRK`, `TSDIRKRegister()`
1346: J*/
1347: typedef const char *TSDIRKType;
1348: #define TSDIRKS212 "s212"
1349: #define TSDIRKES122SAL "es122sal"
1350: #define TSDIRKES213SAL "es213sal"
1351: #define TSDIRKES324SAL "es324sal"
1352: #define TSDIRKES325SAL "es325sal"
1353: #define TSDIRK657A "657a"
1354: #define TSDIRKES648SA "es648sa"
1355: #define TSDIRK658A "658a"
1356: #define TSDIRKS659A "s659a"
1357: #define TSDIRK7510SAL "7510sal"
1358: #define TSDIRKES7510SA "es7510sa"
1359: #define TSDIRK759A "759a"
1360: #define TSDIRKS7511SAL "s7511sal"
1361: #define TSDIRK8614A "8614a"
1362: #define TSDIRK8616SAL "8616sal"
1363: #define TSDIRKES8516SAL "es8516sal"
1364: PETSC_EXTERN PetscErrorCode TSDIRKGetType(TS ts, TSDIRKType *);
1365: PETSC_EXTERN PetscErrorCode TSDIRKSetType(TS ts, TSDIRKType);
1366: PETSC_EXTERN PetscErrorCode TSDIRKRegister(TSDIRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
1368: /*J
1369: TSRosWType - String with the name of a Rosenbrock-W `TSROSW` type
1371: Level: beginner
1373: .seealso: [](ch_ts), `TSRosWSetType()`, `TS`, `TSROSW`, `TSRosWRegister()`
1374: J*/
1375: typedef const char *TSRosWType;
1376: #define TSROSW2M "2m"
1377: #define TSROSW2P "2p"
1378: #define TSROSWRA3PW "ra3pw"
1379: #define TSROSWRA34PW2 "ra34pw2"
1380: #define TSROSWR34PRW "r34prw"
1381: #define TSROSWR3PRL2 "r3prl2"
1382: #define TSROSWRODAS3 "rodas3"
1383: #define TSROSWRODASPR "rodaspr"
1384: #define TSROSWRODASPR2 "rodaspr2"
1385: #define TSROSWSANDU3 "sandu3"
1386: #define TSROSWASSP3P3S1C "assp3p3s1c"
1387: #define TSROSWLASSP3P4S2C "lassp3p4s2c"
1388: #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
1389: #define TSROSWARK3 "ark3"
1390: #define TSROSWTHETA1 "theta1"
1391: #define TSROSWTHETA2 "theta2"
1392: #define TSROSWGRK4T "grk4t"
1393: #define TSROSWSHAMP4 "shamp4"
1394: #define TSROSWVELDD4 "veldd4"
1395: #define TSROSW4L "4l"
1397: PETSC_EXTERN PetscErrorCode TSRosWGetType(TS, TSRosWType *);
1398: PETSC_EXTERN PetscErrorCode TSRosWSetType(TS, TSRosWType);
1399: PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS, PetscBool);
1400: PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
1401: PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal);
1402: PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void);
1403: PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
1404: PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);
1406: PETSC_EXTERN PetscErrorCode TSBDFSetOrder(TS, PetscInt);
1407: PETSC_EXTERN PetscErrorCode TSBDFGetOrder(TS, PetscInt *);
1409: /*J
1410: TSBasicSymplecticType - String with the name of a basic symplectic integration `TSBASICSYMPLECTIC` type
1412: Level: beginner
1414: .seealso: [](ch_ts), `TSBasicSymplecticSetType()`, `TS`, `TSBASICSYMPLECTIC`, `TSBasicSymplecticRegister()`
1415: J*/
1416: typedef const char *TSBasicSymplecticType;
1417: #define TSBASICSYMPLECTICSIEULER "1"
1418: #define TSBASICSYMPLECTICVELVERLET "2"
1419: #define TSBASICSYMPLECTIC3 "3"
1420: #define TSBASICSYMPLECTIC4 "4"
1421: PETSC_EXTERN PetscErrorCode TSBasicSymplecticSetType(TS, TSBasicSymplecticType);
1422: PETSC_EXTERN PetscErrorCode TSBasicSymplecticGetType(TS, TSBasicSymplecticType *);
1423: PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegister(TSBasicSymplecticType, PetscInt, PetscInt, PetscReal[], PetscReal[]);
1424: PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterAll(void);
1425: PETSC_EXTERN PetscErrorCode TSBasicSymplecticInitializePackage(void);
1426: PETSC_EXTERN PetscErrorCode TSBasicSymplecticFinalizePackage(void);
1427: PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterDestroy(void);
1429: /*J
1430: TSDISCGRAD - The Discrete Gradient integrator is a timestepper for Hamiltonian systems designed to conserve the first integral (energy),
1431: but also has the property for some systems of monotonicity in a functional.
1433: Level: beginner
1435: .seealso: [](ch_ts), `TS`, TSDiscGradSetFormulation()`, `TSDiscGradGetFormulation()`
1436: J*/
1437: PETSC_EXTERN PetscErrorCode TSDiscGradSetFormulation(TS, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (*)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (*)(TS, PetscReal, Vec, Vec, void *), void *);
1438: PETSC_EXTERN PetscErrorCode TSDiscGradGetFormulation(TS, PetscErrorCode (**)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (**)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (**)(TS, PetscReal, Vec, Vec, void *), void *);
1439: PETSC_EXTERN PetscErrorCode TSDiscGradIsGonzalez(TS, PetscBool *);
1440: PETSC_EXTERN PetscErrorCode TSDiscGradUseGonzalez(TS, PetscBool);
1442: /*
1443: PETSc interface to Sundials
1444: */
1445: #ifdef PETSC_HAVE_SUNDIALS2
1446: typedef enum {
1447: SUNDIALS_ADAMS = 1,
1448: SUNDIALS_BDF = 2
1449: } TSSundialsLmmType;
1450: PETSC_EXTERN const char *const TSSundialsLmmTypes[];
1451: typedef enum {
1452: SUNDIALS_MODIFIED_GS = 1,
1453: SUNDIALS_CLASSICAL_GS = 2
1454: } TSSundialsGramSchmidtType;
1455: PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
1457: PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS, TSSundialsLmmType);
1458: PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS, PC *);
1459: PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS, PetscReal, PetscReal);
1460: PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS, PetscReal);
1461: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS, PetscReal);
1462: PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS, PetscInt *, PetscInt *);
1463: PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS, TSSundialsGramSchmidtType);
1464: PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS, PetscInt);
1465: PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS, PetscReal);
1466: PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS, PetscBool);
1467: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS, PetscInt);
1468: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxord(TS, PetscInt);
1469: PETSC_EXTERN PetscErrorCode TSSundialsSetUseDense(TS, PetscBool);
1470: #endif
1472: PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS, PetscReal);
1473: PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS, PetscReal *);
1474: PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS, PetscBool *);
1475: PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS, PetscBool);
1477: PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS, PetscReal);
1478: PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS, PetscReal, PetscReal, PetscReal);
1479: PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS, PetscReal *, PetscReal *, PetscReal *);
1481: PETSC_EXTERN PetscErrorCode TSAlpha2SetRadius(TS, PetscReal);
1482: PETSC_EXTERN PetscErrorCode TSAlpha2SetParams(TS, PetscReal, PetscReal, PetscReal, PetscReal);
1483: PETSC_EXTERN PetscErrorCode TSAlpha2GetParams(TS, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
1485: /*S
1486: TSAlpha2PredictorFn - A callback to set the predictor (i.e., the initial guess for the nonlinear solver) in
1487: a second-order generalized-alpha time integrator.
1489: Calling Sequence:
1490: + ts - the `TS` context obtained from `TSCreate()`
1491: . X0 - the previous time step's state vector
1492: . V0 - the previous time step's first derivative of the state vector
1493: . A0 - the previous time step's second derivative of the state vector
1494: . X1 - the vector into which the initial guess for the current time step will be written
1495: - ctx - [optional] user-defined context for the predictor evaluation routine (may be `NULL`)
1497: Level: intermediate
1499: Note:
1500: The deprecated `TSAlpha2Predictor` still works as a replacement for `TSAlpha2PredictorFn` *.
1502: .seealso: [](ch_ts), `TS`, `TSAlpha2SetPredictor()`
1503: S*/
1504: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSAlpha2PredictorFn)(TS ts, Vec X0, Vec V0, Vec A0, Vec X1, void *ctx);
1506: PETSC_EXTERN_TYPEDEF typedef TSAlpha2PredictorFn *TSAlpha2Predictor;
1508: PETSC_EXTERN PetscErrorCode TSAlpha2SetPredictor(TS, TSAlpha2PredictorFn *, void *ctx);
1510: PETSC_EXTERN PetscErrorCode TSSetDM(TS, DM);
1511: PETSC_EXTERN PetscErrorCode TSGetDM(TS, DM *);
1513: PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES, Vec, Vec, void *);
1514: PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES, Vec, Mat, Mat, void *);
1516: PETSC_EXTERN PetscErrorCode TSRHSJacobianTest(TS, PetscBool *);
1517: PETSC_EXTERN PetscErrorCode TSRHSJacobianTestTranspose(TS, PetscBool *);
1519: PETSC_EXTERN PetscErrorCode TSGetComputeInitialCondition(TS, PetscErrorCode (**)(TS, Vec));
1520: PETSC_EXTERN PetscErrorCode TSSetComputeInitialCondition(TS, PetscErrorCode (*)(TS, Vec));
1521: PETSC_EXTERN PetscErrorCode TSComputeInitialCondition(TS, Vec);
1522: PETSC_EXTERN PetscErrorCode TSGetComputeExactError(TS, PetscErrorCode (**)(TS, Vec, Vec));
1523: PETSC_EXTERN PetscErrorCode TSSetComputeExactError(TS, PetscErrorCode (*)(TS, Vec, Vec));
1524: PETSC_EXTERN PetscErrorCode TSComputeExactError(TS, Vec, Vec);
1525: PETSC_EXTERN PetscErrorCode PetscConvEstUseTS(PetscConvEst, PetscBool);
1527: PETSC_EXTERN PetscErrorCode TSSetMatStructure(TS, MatStructure);