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: /*S
492: TSMonitorDrawCtx - Context object for the `TS` graphical monitor routines that draw the solution, phase plot or error using a `PetscDraw`
494: Level: developer
496: .seealso: `TS`, `TSMonitorDrawCtxCreate()`, `TSMonitorDrawCtxDestroy()`, `TSMonitorDrawSolution()`, `TSMonitorDrawSolutionPhase()`, `TSMonitorDrawError()`, `TSMonitorDrawSolutionFunction()`
497: S*/
498: typedef struct _n_TSMonitorDrawCtx *TSMonitorDrawCtx;
499: PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorDrawCtx *);
500: PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *);
501: PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS, PetscInt, PetscReal, Vec, PetscCtx);
502: PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionPhase(TS, PetscInt, PetscReal, Vec, PetscCtx);
503: PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS, PetscInt, PetscReal, Vec, PetscCtx);
504: PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionFunction(TS, PetscInt, PetscReal, Vec, PetscCtx);
506: PETSC_EXTERN PetscErrorCode TSAdjointMonitorDefault(TS, PetscInt, PetscReal, Vec, PetscInt, Vec[], Vec[], PetscViewerAndFormat *);
507: PETSC_EXTERN PetscErrorCode TSAdjointMonitorDrawSensi(TS, PetscInt, PetscReal, Vec, PetscInt, Vec[], Vec[], PetscCtx);
509: /*S
510: TSMonitorSolutionCtx - Context object for the `TS` `TSMonitorSolution()` monitor that views the solution at each time step using a `PetscViewer`
512: Level: developer
514: .seealso: `TS`, `TSMonitorSet()`, `TSMonitorSolution()`, `TSMonitorSolutionSetup()`
515: S*/
516: typedef struct _n_TSMonitorSolutionCtx *TSMonitorSolutionCtx;
517: PETSC_EXTERN PetscErrorCode TSMonitorSolution(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
518: PETSC_EXTERN PetscErrorCode TSMonitorSolutionSetup(TS, PetscViewerAndFormat *);
520: /*S
521: TSMonitorVTKCtx - Context object for the `TS` `TSMonitorSolutionVTK()` monitor that dumps the solution to VTK files at each time step
523: Level: developer
525: .seealso: `TS`, `TSMonitorSet()`, `TSMonitorSolutionVTK()`, `TSMonitorSolutionVTKCtxCreate()`, `TSMonitorSolutionVTKDestroy()`
526: S*/
527: typedef struct _n_TSMonitorVTKCtx *TSMonitorVTKCtx;
528: PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS, PetscInt, PetscReal, Vec, TSMonitorVTKCtx);
529: PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(TSMonitorVTKCtx *);
530: PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKCtxCreate(const char *, TSMonitorVTKCtx *);
532: PETSC_EXTERN PetscErrorCode TSStep(TS);
533: PETSC_EXTERN PetscErrorCode TSEvaluateWLTE(TS, NormType, PetscInt *, PetscReal *);
534: PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS, PetscInt, Vec, PetscBool *);
535: PETSC_EXTERN PetscErrorCode TSSolve(TS, Vec);
536: PETSC_EXTERN PetscErrorCode TSGetEquationType(TS, TSEquationType *);
537: PETSC_EXTERN PetscErrorCode TSSetEquationType(TS, TSEquationType);
538: PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS, TSConvergedReason *);
539: PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS, TSConvergedReason);
540: PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS, PetscReal *);
541: PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS, PetscInt *);
542: PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS, PetscInt *);
543: PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS, PetscInt *);
544: PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS, PetscInt);
545: PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS, PetscInt *);
546: PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS, PetscInt);
547: PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS, PetscBool);
548: PETSC_EXTERN PetscErrorCode TSRestartStep(TS);
549: PETSC_EXTERN PetscErrorCode TSRollBack(TS);
550: PETSC_EXTERN PetscErrorCode TSGetStepRollBack(TS, PetscBool *);
552: PETSC_EXTERN PetscErrorCode TSGetStages(TS, PetscInt *, Vec *[]);
554: PETSC_EXTERN PetscErrorCode TSGetTime(TS, PetscReal *);
555: PETSC_EXTERN PetscErrorCode TSSetTime(TS, PetscReal);
556: PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS, PetscReal *);
557: PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS, PetscReal *);
558: PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS, PetscReal);
559: PETSC_EXTERN PetscErrorCode TSGetStepNumber(TS, PetscInt *);
560: PETSC_EXTERN PetscErrorCode TSSetStepNumber(TS, PetscInt);
562: /*S
563: TSRHSFunctionFn - A prototype of a `TS` right-hand-side evaluation function that would be passed to `TSSetRHSFunction()`
565: Calling Sequence:
566: + ts - timestep context
567: . t - current time
568: . u - input vector
569: . F - function vector
570: - ctx - [optional] user-defined function context
572: Level: beginner
574: Note:
575: The deprecated `TSRHSFunction` still works as a replacement for `TSRHSFunctionFn` *.
577: .seealso: [](ch_ts), `TS`, `TSSetRHSFunction()`, `DMTSSetRHSFunction()`, `TSIFunctionFn`,
578: `TSIJacobianFn`, `TSRHSJacobianFn`
579: S*/
580: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSRHSFunctionFn(TS ts, PetscReal t, Vec u, Vec F, PetscCtx ctx);
582: PETSC_EXTERN_TYPEDEF typedef TSRHSFunctionFn *TSRHSFunction;
584: /*S
585: TSRHSJacobianFn - A prototype of a `TS` right-hand-side Jacobian evaluation function that would be passed to `TSSetRHSJacobian()`
587: Calling Sequence:
588: + ts - the `TS` context obtained from `TSCreate()`
589: . t - current time
590: . u - input vector
591: . Amat - (approximate) Jacobian matrix
592: . Pmat - matrix from which preconditioner is to be constructed (usually the same as `Amat`)
593: - ctx - [optional] user-defined context for matrix evaluation routine
595: Level: beginner
597: Note:
598: The deprecated `TSRHSJacobian` still works as a replacement for `TSRHSJacobianFn` *.
600: .seealso: [](ch_ts), `TS`, `TSSetRHSJacobian()`, `DMTSSetRHSJacobian()`, `TSRHSFunctionFn`,
601: `TSIFunctionFn`, `TSIJacobianFn`
602: S*/
603: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSRHSJacobianFn(TS ts, PetscReal t, Vec u, Mat Amat, Mat Pmat, PetscCtx ctx);
605: PETSC_EXTERN_TYPEDEF typedef TSRHSJacobianFn *TSRHSJacobian;
607: /*S
608: TSRHSJacobianPFn - A prototype of a function that computes the Jacobian of G w.r.t. the parameters P where
609: U_t = G(U,P,t), as well as the location to store the matrix that would be passed to `TSSetRHSJacobianP()`
611: Calling Sequence:
612: + ts - the `TS` context
613: . t - current timestep
614: . U - input vector (current ODE solution)
615: . A - output matrix
616: - ctx - [optional] user-defined function context
618: Level: beginner
620: Note:
621: The deprecated `TSRHSJacobianP` still works as a replacement for `TSRHSJacobianPFn` *.
623: .seealso: [](ch_ts), `TS`, `TSSetRHSJacobianP()`, `TSGetRHSJacobianP()`
624: S*/
625: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSRHSJacobianPFn(TS ts, PetscReal t, Vec U, Mat A, PetscCtx ctx);
627: PETSC_EXTERN_TYPEDEF typedef TSRHSJacobianPFn *TSRHSJacobianP;
629: PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS, Vec, TSRHSFunctionFn *, PetscCtx);
630: PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS, Vec *, TSRHSFunctionFn **, PetscCtxRt);
631: PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS, Mat, Mat, TSRHSJacobianFn *, PetscCtx);
632: PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS, Mat *, Mat *, TSRHSJacobianFn **, PetscCtxRt);
633: PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS, PetscBool);
635: /*S
636: TSSolutionFn - A prototype of a `TS` solution evaluation function that would be passed to `TSSetSolutionFunction()`
638: Calling Sequence:
639: + ts - timestep context
640: . t - current time
641: . u - output vector
642: - ctx - [optional] user-defined function context
644: Level: advanced
646: Note:
647: The deprecated `TSSolutionFunction` still works as a replacement for `TSSolutionFn` *.
649: .seealso: [](ch_ts), `TS`, `TSSetSolutionFunction()`, `DMTSSetSolutionFunction()`
650: S*/
651: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSSolutionFn(TS ts, PetscReal t, Vec u, PetscCtx ctx);
653: PETSC_EXTERN_TYPEDEF typedef TSSolutionFn *TSSolutionFunction;
655: PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS, TSSolutionFn *, PetscCtx);
657: /*S
658: TSForcingFn - A prototype of a `TS` forcing function evaluation function that would be passed to `TSSetForcingFunction()`
660: Calling Sequence:
661: + ts - timestep context
662: . t - current time
663: . f - output vector
664: - ctx - [optional] user-defined function context
666: Level: advanced
668: Note:
669: The deprecated `TSForcingFunction` still works as a replacement for `TSForcingFn` *.
671: .seealso: [](ch_ts), `TS`, `TSSetForcingFunction()`, `DMTSSetForcingFunction()`
672: S*/
673: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSForcingFn(TS ts, PetscReal t, Vec f, PetscCtx ctx);
675: PETSC_EXTERN_TYPEDEF typedef TSForcingFn *TSForcingFunction;
677: PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS, TSForcingFn *, PetscCtx);
679: /*S
680: TSIFunctionFn - A prototype of a `TS` implicit function evaluation function that would be passed to `TSSetIFunction()
682: Calling Sequence:
683: + ts - the `TS` context obtained from `TSCreate()`
684: . t - time at step/stage being solved
685: . U - state vector
686: . U_t - time derivative of state vector
687: . F - function vector
688: - ctx - [optional] user-defined context for function
690: Level: beginner
692: Note:
693: The deprecated `TSIFunction` still works as a replacement for `TSIFunctionFn` *.
695: .seealso: [](ch_ts), `TS`, `TSSetIFunction()`, `DMTSSetIFunction()`, `TSIJacobianFn`, `TSRHSFunctionFn`, `TSRHSJacobianFn`
696: S*/
697: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSIFunctionFn(TS ts, PetscReal t, Vec U, Vec U_t, Vec F, PetscCtx ctx);
699: PETSC_EXTERN_TYPEDEF typedef TSIFunctionFn *TSIFunction;
701: /*S
702: TSIJacobianFn - A prototype of a `TS` Jacobian evaluation function that would be passed to `TSSetIJacobian()`
704: Calling Sequence:
705: + ts - the `TS` context obtained from `TSCreate()`
706: . t - time at step/stage being solved
707: . U - state vector
708: . U_t - time derivative of state vector
709: . a - shift
710: . Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
711: . Pmat - matrix used for constructing preconditioner, usually the same as `Amat`
712: - ctx - [optional] user-defined context for Jacobian evaluation routine
714: Level: beginner
716: Note:
717: The deprecated `TSIJacobian` still works as a replacement for `TSIJacobianFn` *.
719: .seealso: [](ch_ts), `TSSetIJacobian()`, `DMTSSetIJacobian()`, `TSIFunctionFn`, `TSRHSFunctionFn`, `TSRHSJacobianFn`
720: S*/
721: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSIJacobianFn(TS ts, PetscReal t, Vec U, Vec U_t, PetscReal a, Mat Amat, Mat Pmat, PetscCtx ctx);
723: PETSC_EXTERN_TYPEDEF typedef TSIJacobianFn *TSIJacobian;
725: PETSC_EXTERN PetscErrorCode TSSetIFunction(TS, Vec, TSIFunctionFn *, PetscCtx);
726: PETSC_EXTERN PetscErrorCode TSGetIFunction(TS, Vec *, TSIFunctionFn **, PetscCtxRt);
727: PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS, Mat, Mat, TSIJacobianFn *, PetscCtx);
728: PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS, Mat *, Mat *, TSIJacobianFn **, PetscCtxRt);
730: /*S
731: TSI2FunctionFn - A prototype of a `TS` implicit function evaluation function for 2nd order systems that would be passed to `TSSetI2Function()`
733: Calling Sequence:
734: + ts - the `TS` context obtained from `TSCreate()`
735: . t - time at step/stage being solved
736: . U - state vector
737: . U_t - time derivative of state vector
738: . U_tt - second time derivative of state vector
739: . F - function vector
740: - ctx - [optional] user-defined context for matrix evaluation routine (may be `NULL`)
742: Level: advanced
744: Note:
745: The deprecated `TSI2Function` still works as a replacement for `TSI2FunctionFn` *.
747: .seealso: [](ch_ts), `TS`, `TSSetI2Function()`, `DMTSSetI2Function()`, `TSIFunctionFn`
748: S*/
749: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSI2FunctionFn(TS ts, PetscReal t, Vec U, Vec U_t, Vec U_tt, Vec F, PetscCtx ctx);
751: PETSC_EXTERN_TYPEDEF typedef TSI2FunctionFn *TSI2Function;
753: /*S
754: TSI2JacobianFn - A prototype of a `TS` implicit Jacobian evaluation function for 2nd order systems that would be passed to `TSSetI2Jacobian()`
756: Calling Sequence:
757: + ts - the `TS` context obtained from `TSCreate()`
758: . t - time at step/stage being solved
759: . U - state vector
760: . U_t - time derivative of state vector
761: . U_tt - second time derivative of state vector
762: . v - shift for U_t
763: . a - shift for U_tt
764: . 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
765: . jac - matrix from which to construct the preconditioner, may be same as `J`
766: - ctx - [optional] user-defined context for matrix evaluation routine
768: Level: advanced
770: Note:
771: The deprecated `TSI2Jacobian` still works as a replacement for `TSI2JacobianFn` *.
773: .seealso: [](ch_ts), `TS`, `TSSetI2Jacobian()`, `DMTSSetI2Jacobian()`, `TSIFunctionFn`, `TSIJacobianFn`, `TSRHSFunctionFn`, `TSRHSJacobianFn`
774: S*/
775: 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);
777: PETSC_EXTERN_TYPEDEF typedef TSI2JacobianFn *TSI2Jacobian;
779: PETSC_EXTERN PetscErrorCode TSSetI2Function(TS, Vec, TSI2FunctionFn *, PetscCtx);
780: PETSC_EXTERN PetscErrorCode TSGetI2Function(TS, Vec *, TSI2FunctionFn **, PetscCtxRt);
781: PETSC_EXTERN PetscErrorCode TSSetI2Jacobian(TS, Mat, Mat, TSI2JacobianFn *, PetscCtx);
782: PETSC_EXTERN PetscErrorCode TSGetI2Jacobian(TS, Mat *, Mat *, TSI2JacobianFn **, PetscCtxRt);
784: PETSC_EXTERN PetscErrorCode TSRHSSplitSetIS(TS, const char[], IS);
785: PETSC_EXTERN PetscErrorCode TSRHSSplitGetIS(TS, const char[], IS *);
786: PETSC_EXTERN PetscErrorCode TSRHSSplitSetRHSFunction(TS, const char[], Vec, TSRHSFunctionFn *, PetscCtx);
787: PETSC_EXTERN PetscErrorCode TSRHSSplitSetIFunction(TS, const char[], Vec, TSIFunctionFn *, PetscCtx);
788: PETSC_EXTERN PetscErrorCode TSRHSSplitSetIJacobian(TS, const char[], Mat, Mat, TSIJacobianFn *, PetscCtx);
789: PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTS(TS, const char[], TS *);
790: PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTSs(TS, PetscInt *, TS *[]);
791: PETSC_EXTERN PetscErrorCode TSSetUseSplitRHSFunction(TS, PetscBool);
792: PETSC_EXTERN PetscErrorCode TSGetUseSplitRHSFunction(TS, PetscBool *);
793: PETSC_EXTERN PetscErrorCode TSRHSSplitGetSNES(TS, SNES *);
794: PETSC_EXTERN PetscErrorCode TSRHSSplitSetSNES(TS, SNES);
796: PETSC_EXTERN TSRHSFunctionFn TSComputeRHSFunctionLinear;
797: PETSC_EXTERN TSRHSJacobianFn TSComputeRHSJacobianConstant;
798: PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS, PetscReal, Vec, Vec, Vec, PetscCtx);
799: PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscCtx);
800: PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS, PetscReal, Vec);
801: PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS, PetscReal, Vec);
802: PETSC_EXTERN PetscErrorCode TSComputeIJacobianDefaultColor(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscCtx);
803: PETSC_EXTERN PetscErrorCode TSPruneIJacobianColor(TS, Mat, Mat);
805: PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
806: PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS, PetscReal));
807: PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS, PetscReal, PetscInt, Vec *));
808: PETSC_EXTERN PetscErrorCode TSSetPostEvaluate(TS, PetscErrorCode (*)(TS));
809: PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
810: PETSC_EXTERN PetscErrorCode TSSetResize(TS, PetscBool, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscBool *, PetscCtx), PetscErrorCode (*)(TS, PetscInt, Vec[], Vec[], PetscCtx), PetscCtx);
811: PETSC_EXTERN PetscErrorCode TSPreStep(TS);
812: PETSC_EXTERN PetscErrorCode TSPreStage(TS, PetscReal);
813: PETSC_EXTERN PetscErrorCode TSPostStage(TS, PetscReal, PetscInt, Vec[]);
814: PETSC_EXTERN PetscErrorCode TSPostEvaluate(TS);
815: PETSC_EXTERN PetscErrorCode TSPostStep(TS);
816: PETSC_EXTERN PetscErrorCode TSResize(TS);
817: PETSC_EXTERN PetscErrorCode TSResizeRetrieveVec(TS, const char *, Vec *);
818: PETSC_EXTERN PetscErrorCode TSResizeRegisterVec(TS, const char *, Vec);
819: PETSC_EXTERN PetscErrorCode TSGetStepResize(TS, PetscBool *);
821: PETSC_EXTERN PetscErrorCode TSInterpolate(TS, PetscReal, Vec);
822: PETSC_EXTERN PetscErrorCode TSSetTolerances(TS, PetscReal, Vec, PetscReal, Vec);
823: PETSC_EXTERN PetscErrorCode TSGetTolerances(TS, PetscReal *, Vec *, PetscReal *, Vec *);
824: PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *);
825: PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm(TS, Vec, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *);
826: PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS, PetscReal);
827: PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS, PetscReal *);
828: PETSC_EXTERN PetscErrorCode TSSetFunctionDomainError(TS, PetscErrorCode (*)(TS, PetscReal, Vec, PetscBool *));
829: PETSC_EXTERN PetscErrorCode TSFunctionDomainError(TS, PetscReal, Vec, PetscBool *);
831: PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS, PetscErrorCode (*)(TS, PetscReal *, PetscCtx), PetscCtx);
832: PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS, PetscReal *, PetscCtx);
833: PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS, PetscReal);
834: PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS, PetscErrorCode (*)(TS, Vec, PetscCtx, PetscReal *, PetscBool *), PetscCtx);
835: PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS, PetscReal);
836: PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);
837: PETSC_EXTERN PetscErrorCode TSPseudoComputeFunction(TS, Vec, Vec *, PetscReal *);
839: PETSC_EXTERN PetscErrorCode TSPythonSetType(TS, const char[]);
840: PETSC_EXTERN PetscErrorCode TSPythonGetType(TS, const char *[]);
842: PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS, PetscReal, Vec, Vec);
843: PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS, PetscReal, Vec, Mat, Mat);
844: PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS, PetscReal, Vec, Vec, Vec, PetscBool);
845: PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscBool);
846: PETSC_EXTERN PetscErrorCode TSComputeI2Function(TS, PetscReal, Vec, Vec, Vec, Vec);
847: PETSC_EXTERN PetscErrorCode TSComputeI2Jacobian(TS, PetscReal, Vec, Vec, Vec, PetscReal, PetscReal, Mat, Mat);
848: PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS, PetscReal, PetscReal, PetscReal *, PetscReal *);
850: PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS, Vec, Vec);
852: PETSC_EXTERN PetscErrorCode DMTSSetBoundaryLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, PetscCtx), PetscCtx);
853: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionPre(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, PetscCtx), PetscCtx);
854: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM, TSRHSFunctionFn *, PetscCtx);
855: PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM, TSRHSFunctionFn **, PetscCtxRt);
856: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionContextDestroy(DM, PetscCtxDestroyFn *);
857: PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM, TSRHSJacobianFn *, PetscCtx);
858: PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM, TSRHSJacobianFn **, PetscCtxRt);
859: PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobianContextDestroy(DM, PetscCtxDestroyFn *);
860: PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM, TSIFunctionFn *, PetscCtx);
861: PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM, TSIFunctionFn **, PetscCtxRt);
862: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionContextDestroy(DM, PetscCtxDestroyFn *);
863: PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM, TSIJacobianFn *, PetscCtx);
864: PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM, TSIJacobianFn **, PetscCtxRt);
865: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianContextDestroy(DM, PetscCtxDestroyFn *);
866: PETSC_EXTERN PetscErrorCode DMTSSetI2Function(DM, TSI2FunctionFn *, PetscCtx);
867: PETSC_EXTERN PetscErrorCode DMTSGetI2Function(DM, TSI2FunctionFn **, PetscCtxRt);
868: PETSC_EXTERN PetscErrorCode DMTSSetI2FunctionContextDestroy(DM, PetscCtxDestroyFn *);
869: PETSC_EXTERN PetscErrorCode DMTSSetI2Jacobian(DM, TSI2JacobianFn *, PetscCtx);
870: PETSC_EXTERN PetscErrorCode DMTSGetI2Jacobian(DM, TSI2JacobianFn **, PetscCtxRt);
871: PETSC_EXTERN PetscErrorCode DMTSSetI2JacobianContextDestroy(DM, PetscCtxDestroyFn *);
873: /*S
874: TSTransientVariableFn - A prototype of a function to transform from state to transient variables that would be passed to `TSSetTransientVariable()`
876: Calling Sequence:
877: + ts - timestep context
878: . p - input vector (primitive form)
879: . c - output vector, transient variables (conservative form)
880: - ctx - [optional] user-defined function context
882: Level: advanced
884: Note:
885: The deprecated `TSTransientVariable` still works as a replacement for `TSTransientVariableFn` *.
887: .seealso: [](ch_ts), `TS`, `TSSetTransientVariable()`, `DMTSSetTransientVariable()`
888: S*/
889: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSTransientVariableFn(TS ts, Vec p, Vec c, PetscCtx ctx);
891: PETSC_EXTERN_TYPEDEF typedef TSTransientVariableFn *TSTransientVariable;
893: PETSC_EXTERN PetscErrorCode TSSetTransientVariable(TS, TSTransientVariableFn *, PetscCtx);
894: PETSC_EXTERN PetscErrorCode DMTSSetTransientVariable(DM, TSTransientVariableFn *, PetscCtx);
895: PETSC_EXTERN PetscErrorCode DMTSGetTransientVariable(DM, TSTransientVariableFn **, PetscCtx);
896: PETSC_EXTERN PetscErrorCode TSComputeTransientVariable(TS, Vec, Vec);
897: PETSC_EXTERN PetscErrorCode TSHasTransientVariable(TS, PetscBool *);
899: PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM, TSSolutionFn *, PetscCtx);
900: PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM, TSSolutionFn **, PetscCtxRt);
901: PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM, TSForcingFn *, PetscCtx);
902: PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM, TSForcingFn **, PetscCtxRt);
903: PETSC_EXTERN PetscErrorCode DMTSCheckResidual(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscReal *);
904: PETSC_EXTERN PetscErrorCode DMTSCheckJacobian(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscBool *, PetscReal *);
905: PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec);
907: PETSC_EXTERN PetscErrorCode DMTSGetIFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, Vec, PetscCtx), PetscCtxRt);
908: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, Vec, PetscCtx), PetscCtx);
909: PETSC_EXTERN PetscErrorCode DMTSGetIJacobianLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscCtx), PetscCtxRt);
910: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscCtx), PetscCtx);
911: PETSC_EXTERN PetscErrorCode DMTSGetRHSFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, PetscCtx), PetscCtxRt);
912: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, PetscCtx), PetscCtx);
913: PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrix(DM);
914: PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrixLumped(DM);
915: PETSC_EXTERN PetscErrorCode DMTSDestroyRHSMassMatrix(DM);
917: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM, PetscErrorCode (*)(PetscCtx, PetscViewer), PetscErrorCode (*)(PetscCtx *, PetscViewer));
918: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM, PetscErrorCode (*)(PetscCtx, PetscViewer), PetscErrorCode (*)(PetscCtx *, PetscViewer));
920: /*S
921: DMDATSRHSFunctionLocalFn - A prototype of a local `TS` right-hand side residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetRHSFunctionLocal()`
923: Calling Sequence:
924: + info - defines the subdomain to evaluate the residual on
925: . t - time at which to evaluate residual
926: . x - array of local state information
927: . f - output array of local residual information
928: - ctx - optional user context
930: Level: beginner
932: Note:
933: The deprecated `DMDATSRHSFunctionLocal` still works as a replacement for `DMDATSRHSFunctionLocalFn` *.
935: .seealso: `DMDA`, `DMDATSSetRHSFunctionLocal()`, `TSRHSFunctionFn`, `DMDATSRHSJacobianLocalFn`, `DMDATSIJacobianLocalFn`, `DMDATSIFunctionLocalFn`
936: S*/
937: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDATSRHSFunctionLocalFn(DMDALocalInfo *info, PetscReal t, void *x, void *f, PetscCtx ctx);
939: PETSC_EXTERN_TYPEDEF typedef DMDATSRHSFunctionLocalFn *DMDATSRHSFunctionLocal;
941: /*S
942: DMDATSRHSJacobianLocalFn - A prototype of a local residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetRHSJacobianLocal()`
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: . J - Jacobian matrix
949: . B - matrix from which to construct the preconditioner; often same as `J`
950: - ctx - optional context
952: Level: beginner
954: Note:
955: The deprecated `DMDATSRHSJacobianLocal` still works as a replacement for `DMDATSRHSJacobianLocalFn` *.
957: .seealso: `DMDA`, `DMDATSSetRHSJacobianLocal()`, `TSRHSJacobianFn`, `DMDATSRHSFunctionLocalFn`, `DMDATSIJacobianLocalFn`, `DMDATSIFunctionLocalFn`
958: S*/
959: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDATSRHSJacobianLocalFn(DMDALocalInfo *info, PetscReal t, void *x, Mat J, Mat B, PetscCtx ctx);
961: PETSC_EXTERN_TYPEDEF typedef DMDATSRHSJacobianLocalFn *DMDATSRHSJacobianLocal;
963: /*S
964: DMDATSIFunctionLocalFn - A prototype of a local residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetIFunctionLocal()`
966: Calling Sequence:
967: + info - defines the subdomain to evaluate the residual on
968: . t - time at which to evaluate residual
969: . x - array of local state information
970: . xdot - array of local time derivative information
971: . imode - output array of local function evaluation information
972: - ctx - optional context
974: Level: beginner
976: Note:
977: The deprecated `DMDATSIFunctionLocal` still works as a replacement for `DMDATSIFunctionLocalFn` *.
979: .seealso: `DMDA`, `DMDATSSetIFunctionLocal()`, `DMDATSIJacobianLocalFn`, `TSIFunctionFn`
980: S*/
981: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDATSIFunctionLocalFn(DMDALocalInfo *info, PetscReal t, void *x, void *xdot, void *imode, PetscCtx ctx);
983: PETSC_EXTERN_TYPEDEF typedef DMDATSIFunctionLocalFn *DMDATSIFunctionLocal;
985: /*S
986: DMDATSIJacobianLocalFn - A prototype of a local residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetIJacobianLocal()`
988: Calling Sequence:
989: + info - defines the subdomain to evaluate the residual on
990: . t - time at which to evaluate the jacobian
991: . x - array of local state information
992: . xdot - time derivative at this state
993: . shift - see `TSSetIJacobian()` for the meaning of this parameter
994: . J - Jacobian matrix
995: . B - matrix from which to construct the preconditioner; often same as `J`
996: - ctx - optional context
998: Level: beginner
1000: Note:
1001: The deprecated `DMDATSIJacobianLocal` still works as a replacement for `DMDATSIJacobianLocalFn` *.
1003: .seealso: `DMDA`, `DMDATSSetIJacobianLocal()`, `TSIJacobianFn`, `DMDATSIFunctionLocalFn`, `DMDATSRHSFunctionLocalFn`, `DMDATSRHSJacobianlocal()`
1004: S*/
1005: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDATSIJacobianLocalFn(DMDALocalInfo *info, PetscReal t, void *x, void *xdot, PetscReal shift, Mat J, Mat B, PetscCtx ctx);
1007: PETSC_EXTERN_TYPEDEF typedef DMDATSIJacobianLocalFn *DMDATSIJacobianLocal;
1009: PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM, InsertMode, DMDATSRHSFunctionLocalFn *, void *);
1010: PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM, DMDATSRHSJacobianLocalFn *, void *);
1011: PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM, InsertMode, DMDATSIFunctionLocalFn *, void *);
1012: PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM, DMDATSIJacobianLocalFn *, void *);
1014: /*S
1015: TSMonitorLGCtx - Context object for `TS` line-graph monitor routines that plot residuals, iteration counts or solution components at each time step on a `PetscDrawLG`
1017: Level: developer
1019: .seealso: `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorLGCtxDestroy()`, `TSMonitorLGSolution()`, `TSMonitorLGTimeStep()`, `TSMonitorLGError()`
1020: S*/
1021: typedef struct _n_TSMonitorLGCtx *TSMonitorLGCtx;
1023: /*S
1024: TSMonitorDMDARayCtx - Context object for `TSMonitorDMDARay()`
1026: Level: developer
1028: .seealso: `TS`, `TSSetMonitor()`, `TSMonitorDMDARay()`, `TSMonitorDMDARayDestroy()`
1029: S*/
1030: typedef struct {
1031: Vec ray;
1032: VecScatter scatter;
1033: PetscViewer viewer;
1034: TSMonitorLGCtx lgctx;
1035: } TSMonitorDMDARayCtx;
1036: PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(PetscCtxRt);
1037: PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS, PetscInt, PetscReal, Vec, void *);
1038: PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS, PetscInt, PetscReal, Vec, void *);
1040: /* Dynamic creation and loading functions */
1041: PETSC_EXTERN PetscFunctionList TSList;
1042: PETSC_EXTERN PetscErrorCode TSGetType(TS, TSType *);
1043: PETSC_EXTERN PetscErrorCode TSSetType(TS, TSType);
1044: PETSC_EXTERN PetscErrorCode TSRegister(const char[], PetscErrorCode (*)(TS));
1046: PETSC_EXTERN PetscErrorCode TSGetSNES(TS, SNES *);
1047: PETSC_EXTERN PetscErrorCode TSSetSNES(TS, SNES);
1048: PETSC_EXTERN PetscErrorCode TSGetKSP(TS, KSP *);
1049: PETSC_EXTERN PetscErrorCode TSIsImplicit(TS, PetscBool *);
1051: PETSC_EXTERN PetscErrorCode TSView(TS, PetscViewer);
1052: PETSC_EXTERN PetscErrorCode TSLoad(TS, PetscViewer);
1053: PETSC_EXTERN PetscErrorCode TSViewFromOptions(TS, PetscObject, const char[]);
1054: PETSC_EXTERN PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory, PetscObject, const char[]);
1056: #define TS_FILE_CLASSID 1211225
1058: PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS, PetscCtx);
1059: PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS, PetscCtxRt);
1061: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtx *);
1062: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *);
1063: PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS, PetscInt, PetscReal, Vec, void *);
1064: PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS, PetscInt, PetscReal, Vec, void *);
1065: PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS, const char *const *);
1066: PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS, const char *const **);
1067: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx, const char *const *);
1068: PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS, const char *const *);
1069: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx, const char *const *);
1070: PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS, PetscErrorCode (*)(void *, Vec, Vec *), PetscCtxDestroyFn *, void *);
1071: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx, PetscErrorCode (*)(void *, Vec, Vec *), PetscCtxDestroyFn *, void *);
1072: PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS, PetscInt, PetscReal, Vec, void *);
1073: PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS, PetscInt, PetscReal, Vec, void *);
1074: PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS, PetscInt, PetscReal, Vec, void *);
1075: PETSC_EXTERN PetscErrorCode TSMonitorError(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
1076: PETSC_EXTERN PetscErrorCode TSDMSwarmMonitorMoments(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
1078: struct _n_TSMonitorLGCtxNetwork {
1079: PetscInt nlg;
1080: PetscDrawLG *lg;
1081: PetscBool semilogy;
1082: PetscInt howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
1083: };
1084: /*S
1085: TSMonitorLGCtxNetwork - Context object for the `TSMonitorLGCtxNetworkSolution()` line-graph monitor that plots solution components on each subnetwork of a `DMNETWORK`
1087: Level: developer
1089: .seealso: `TS`, `DMNETWORK`, `TSMonitorLGCtxNetworkCreate()`, `TSMonitorLGCtxNetworkDestroy()`, `TSMonitorLGCtxNetworkSolution()`
1090: S*/
1091: typedef struct _n_TSMonitorLGCtxNetwork *TSMonitorLGCtxNetwork;
1092: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxNetworkDestroy(TSMonitorLGCtxNetwork *);
1093: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxNetworkCreate(TS, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtxNetwork *);
1094: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxNetworkSolution(TS, PetscInt, PetscReal, Vec, void *);
1096: /*S
1097: TSMonitorEnvelopeCtx - Context object for the `TSMonitorEnvelope()` monitor that tracks the per-component min/max envelope of the solution over a time integration
1099: Level: developer
1101: .seealso: `TS`, `TSMonitorEnvelopeCtxCreate()`, `TSMonitorEnvelopeCtxDestroy()`, `TSMonitorEnvelope()`, `TSMonitorEnvelopeGetBounds()`
1102: S*/
1103: typedef struct _n_TSMonitorEnvelopeCtx *TSMonitorEnvelopeCtx;
1104: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxCreate(TS, TSMonitorEnvelopeCtx *);
1105: PETSC_EXTERN PetscErrorCode TSMonitorEnvelope(TS, PetscInt, PetscReal, Vec, void *);
1106: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeGetBounds(TS, Vec *, Vec *);
1107: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *);
1109: /*S
1110: TSMonitorSPEigCtx - Context object for the `TSMonitorSPEig()` monitor that displays an estimate of the spectrum of the operator using a `PetscDrawSP` scatter plot
1112: Level: developer
1114: .seealso: `TS`, `TSMonitorSPEigCtxCreate()`, `TSMonitorSPEigCtxDestroy()`, `TSMonitorSPEig()`
1115: S*/
1116: typedef struct _n_TSMonitorSPEigCtx *TSMonitorSPEigCtx;
1117: PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorSPEigCtx *);
1118: PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *);
1119: PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS, PetscInt, PetscReal, Vec, void *);
1121: /*S
1122: TSMonitorSPCtx - Context object for the `TSMonitorSPSwarmSolution()` scatter-plot monitor that draws the swarm particle positions at each time step on a `PetscDrawSP`
1124: Level: developer
1126: .seealso: `TS`, `DMSWARM`, `TSMonitorSPCtxCreate()`, `TSMonitorSPCtxDestroy()`, `TSMonitorSPSwarmSolution()`
1127: S*/
1128: typedef struct _n_TSMonitorSPCtx *TSMonitorSPCtx;
1129: PETSC_EXTERN PetscErrorCode TSMonitorSPCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscBool, PetscBool, TSMonitorSPCtx *);
1130: PETSC_EXTERN PetscErrorCode TSMonitorSPCtxDestroy(TSMonitorSPCtx *);
1131: PETSC_EXTERN PetscErrorCode TSMonitorSPSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
1133: /*S
1134: TSMonitorHGCtx - Context object for the `TSMonitorHGSwarmSolution()` histogram monitor that displays a histogram of `DMSWARM` particle quantities at each time step
1136: Level: developer
1138: .seealso: `TS`, `DMSWARM`, `TSMonitorHGCtxCreate()`, `TSMonitorHGCtxDestroy()`, `TSMonitorHGSwarmSolution()`
1139: S*/
1140: typedef struct _n_TSMonitorHGCtx *TSMonitorHGCtx;
1141: PETSC_EXTERN PetscErrorCode TSMonitorHGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscInt, PetscBool, TSMonitorHGCtx *);
1142: PETSC_EXTERN PetscErrorCode TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
1143: PETSC_EXTERN PetscErrorCode TSMonitorHGCtxDestroy(TSMonitorHGCtx *);
1144: PETSC_EXTERN PetscErrorCode TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
1146: /*M
1147: TSEvent - Abstract object to handle event detection in `TS` time integrator
1149: Level: intermediate
1151: Note:
1152: See `TSSetEventHandler()` for the management of events.
1154: .seealso: [](sec_ts_event), `TS`, `TSSetEventHandler()`, `TSSetPostEventStep()`, `TSSetPostEventSecondStep()`, `TSSetEventTolerances()`, `TSGetNumEvents()`
1155: M*/
1156: typedef struct _n_TSEvent *TSEvent;
1158: PETSC_EXTERN PetscErrorCode TSSetEventHandler(TS, PetscInt, PetscInt[], PetscBool[], PetscErrorCode (*)(TS, PetscReal, Vec, PetscReal[], void *), PetscErrorCode (*)(TS, PetscInt, PetscInt[], PetscReal, Vec, PetscBool, void *), void *);
1159: PETSC_EXTERN PetscErrorCode TSSetPostEventStep(TS, PetscReal);
1160: PETSC_EXTERN PetscErrorCode TSSetPostEventSecondStep(TS, PetscReal);
1161: PETSC_DEPRECATED_FUNCTION(3, 21, 0, "TSSetPostEventSecondStep()", ) static inline PetscErrorCode TSSetPostEventIntervalStep(TS ts, PetscReal dt)
1162: {
1163: return TSSetPostEventSecondStep(ts, dt);
1164: }
1165: PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS, PetscReal, PetscReal[]);
1166: PETSC_EXTERN PetscErrorCode TSGetNumEvents(TS, PetscInt *);
1168: /*J
1169: TSSSPType - string with the name of a `TSSSP` scheme.
1171: Level: beginner
1173: .seealso: [](ch_ts), `TSSSPSetType()`, `TS`, `TSSSP`
1174: J*/
1175: typedef const char *TSSSPType;
1176: #define TSSSPRKS2 "rks2"
1177: #define TSSSPRKS3 "rks3"
1178: #define TSSSPRK104 "rk104"
1180: PETSC_EXTERN PetscErrorCode TSSSPSetType(TS, TSSSPType);
1181: PETSC_EXTERN PetscErrorCode TSSSPGetType(TS, TSSSPType *);
1182: PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS, PetscInt);
1183: PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS, PetscInt *);
1184: PETSC_EXTERN PetscErrorCode TSSSPInitializePackage(void);
1185: PETSC_EXTERN PetscErrorCode TSSSPFinalizePackage(void);
1186: PETSC_EXTERN PetscFunctionList TSSSPList;
1188: /*S
1189: TSAdapt - Abstract object that manages time-step adaptivity
1191: Level: beginner
1193: .seealso: [](ch_ts), [](sec_ts_error_control), `TS`, `TSGetAdapt()`, `TSAdaptCreate()`, `TSAdaptType`
1194: S*/
1195: typedef struct _p_TSAdapt *TSAdapt;
1197: /*J
1198: TSAdaptType - String with the name of `TSAdapt` scheme.
1200: Level: beginner
1202: .seealso: [](ch_ts), [](sec_ts_error_control), `TSGetAdapt()`, `TSAdaptSetType()`, `TS`, `TSAdapt`
1203: J*/
1204: typedef const char *TSAdaptType;
1205: #define TSADAPTNONE "none"
1206: #define TSADAPTBASIC "basic"
1207: #define TSADAPTDSP "dsp"
1208: #define TSADAPTCFL "cfl"
1209: #define TSADAPTGLEE "glee"
1210: #define TSADAPTHISTORY "history"
1212: PETSC_EXTERN PetscErrorCode TSGetAdapt(TS, TSAdapt *);
1213: PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[], PetscErrorCode (*)(TSAdapt));
1214: PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void);
1215: PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
1216: PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm, TSAdapt *);
1217: PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt, TSAdaptType);
1218: PETSC_EXTERN PetscErrorCode TSAdaptGetType(TSAdapt, TSAdaptType *);
1219: PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt, const char[]);
1220: PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
1221: PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt, const char[], PetscInt, PetscInt, PetscReal, PetscReal, PetscBool);
1222: PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt, PetscInt *, const PetscInt **, const PetscInt **, const PetscReal **, const PetscReal **);
1223: PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt, TS, PetscReal, PetscInt *, PetscReal *, PetscBool *);
1224: PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt, TS, PetscReal, Vec, PetscBool *);
1225: PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt, PetscViewer);
1226: PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt, PetscViewer);
1227: PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(TSAdapt, PetscOptionItems);
1228: PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt);
1229: PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt *);
1230: PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt, PetscBool);
1231: PETSC_EXTERN PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt, PetscBool);
1232: PETSC_EXTERN PetscErrorCode TSAdaptSetSafety(TSAdapt, PetscReal, PetscReal);
1233: PETSC_EXTERN PetscErrorCode TSAdaptGetSafety(TSAdapt, PetscReal *, PetscReal *);
1234: PETSC_EXTERN PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt, PetscReal);
1235: PETSC_EXTERN PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt, PetscReal *);
1236: PETSC_EXTERN PetscErrorCode TSAdaptSetClip(TSAdapt, PetscReal, PetscReal);
1237: PETSC_EXTERN PetscErrorCode TSAdaptGetClip(TSAdapt, PetscReal *, PetscReal *);
1238: PETSC_EXTERN PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt, PetscReal);
1239: PETSC_EXTERN PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt, PetscReal *);
1240: PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt, PetscReal, PetscReal);
1241: PETSC_EXTERN PetscErrorCode TSAdaptGetStepLimits(TSAdapt, PetscReal *, PetscReal *);
1242: PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt, PetscErrorCode (*)(TSAdapt, TS, PetscReal, Vec, PetscBool *));
1243: PETSC_EXTERN PetscErrorCode TSAdaptHistorySetHistory(TSAdapt, PetscInt, PetscReal[], PetscBool);
1244: PETSC_EXTERN PetscErrorCode TSAdaptHistorySetTrajectory(TSAdapt, TSTrajectory, PetscBool);
1245: PETSC_EXTERN PetscErrorCode TSAdaptHistoryGetStep(TSAdapt, PetscInt, PetscReal *, PetscReal *);
1246: PETSC_EXTERN PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt, PetscInt);
1247: PETSC_EXTERN PetscErrorCode TSAdaptDSPSetFilter(TSAdapt, const char *);
1248: PETSC_EXTERN PetscErrorCode TSAdaptDSPSetPID(TSAdapt, PetscReal, PetscReal, PetscReal);
1250: /*S
1251: TSGLLEAdapt - Abstract object that manages time-step adaptivity for `TSGLLE`
1253: Level: beginner
1255: Developer Note:
1256: This functionality should be replaced by the `TSAdapt`.
1258: .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLEAdaptCreate()`, `TSGLLEAdaptType`
1259: S*/
1260: typedef struct _p_TSGLLEAdapt *TSGLLEAdapt;
1262: /*J
1263: TSGLLEAdaptType - String with the name of `TSGLLEAdapt` scheme
1265: Level: beginner
1267: Developer Note:
1268: This functionality should be replaced by the `TSAdaptType`.
1270: .seealso: [](ch_ts), `TSGLLEAdaptSetType()`, `TS`
1271: J*/
1272: typedef const char *TSGLLEAdaptType;
1273: #define TSGLLEADAPT_NONE "none"
1274: #define TSGLLEADAPT_SIZE "size"
1275: #define TSGLLEADAPT_BOTH "both"
1277: PETSC_EXTERN PetscErrorCode TSGLLEAdaptRegister(const char[], PetscErrorCode (*)(TSGLLEAdapt));
1278: PETSC_EXTERN PetscErrorCode TSGLLEAdaptInitializePackage(void);
1279: PETSC_EXTERN PetscErrorCode TSGLLEAdaptFinalizePackage(void);
1280: PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate(MPI_Comm, TSGLLEAdapt *);
1281: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt, TSGLLEAdaptType);
1282: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt, const char[]);
1283: PETSC_EXTERN PetscErrorCode TSGLLEAdaptChoose(TSGLLEAdapt, PetscInt, const PetscInt[], const PetscReal[], const PetscReal[], PetscInt, PetscReal, PetscReal, PetscInt *, PetscReal *, PetscBool *);
1284: PETSC_EXTERN PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt, PetscViewer);
1285: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetFromOptions(TSGLLEAdapt, PetscOptionItems);
1286: PETSC_EXTERN PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt *);
1288: /*J
1289: TSGLLEAcceptType - String with the name of `TSGLLEAccept` scheme
1291: Level: beginner
1293: .seealso: [](ch_ts), `TSGLLESetAcceptType()`, `TS`, `TSGLLEAccept`
1294: J*/
1295: typedef const char *TSGLLEAcceptType;
1296: #define TSGLLEACCEPT_ALWAYS "always"
1298: /*S
1299: TSGLLEAcceptFn - A prototype of a `TS` accept function that would be passed to `TSGLLEAcceptRegister()`
1301: Calling Sequence:
1302: + ts - timestep context
1303: . nt - time to end of solution time
1304: . h - the proposed step-size
1305: . enorm - unknown
1306: - accept - output, if the proposal is accepted
1308: Level: beginner
1310: Note:
1311: The deprecated `TSGLLEAcceptFunction` still works as a replacement for `TSGLLEAcceptFn` *
1313: .seealso: [](ch_ts), `TS`, `TSSetRHSFunction()`, `DMTSSetRHSFunction()`, `TSIFunctionFn`,
1314: `TSIJacobianFn`, `TSRHSJacobianFn`, `TSGLLEAcceptRegister()`
1315: S*/
1316: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSGLLEAcceptFn(TS ts, PetscReal nt, PetscReal h, const PetscReal enorm[], PetscBool *accept);
1318: PETSC_EXTERN_TYPEDEF typedef TSGLLEAcceptFn *TSGLLEAcceptFunction;
1320: PETSC_EXTERN PetscErrorCode TSGLLEAcceptRegister(const char[], TSGLLEAcceptFn *);
1322: /*J
1323: TSGLLEType - string with the name of a General Linear `TSGLLE` type
1325: Level: beginner
1327: .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLESetType()`, `TSGLLERegister()`, `TSGLLEAccept`
1328: J*/
1329: typedef const char *TSGLLEType;
1330: #define TSGLLE_IRKS "irks"
1332: PETSC_EXTERN PetscErrorCode TSGLLERegister(const char[], PetscErrorCode (*)(TS));
1333: PETSC_EXTERN PetscErrorCode TSGLLEInitializePackage(void);
1334: PETSC_EXTERN PetscErrorCode TSGLLEFinalizePackage(void);
1335: PETSC_EXTERN PetscErrorCode TSGLLESetType(TS, TSGLLEType);
1336: PETSC_EXTERN PetscErrorCode TSGLLEGetAdapt(TS, TSGLLEAdapt *);
1337: PETSC_EXTERN PetscErrorCode TSGLLESetAcceptType(TS, TSGLLEAcceptType);
1339: /*J
1340: TSEIMEXType - String with the name of an Extrapolated IMEX `TSEIMEX` type
1342: Level: beginner
1344: .seealso: [](ch_ts), `TSEIMEXSetType()`, `TS`, `TSEIMEX`, `TSEIMEXRegister()`
1345: J*/
1346: #define TSEIMEXType char *
1348: PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS, PetscInt);
1349: PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS, PetscInt, PetscInt);
1350: PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS, PetscBool);
1352: /*J
1353: TSRKType - String with the name of a Runge-Kutta `TSRK` type
1355: Level: beginner
1357: .seealso: [](ch_ts), `TS`, `TSRKSetType()`, `TSRK`, `TSRKRegister()`
1358: J*/
1359: typedef const char *TSRKType;
1360: #define TSRK1FE "1fe"
1361: #define TSRK2A "2a"
1362: #define TSRK2B "2b"
1363: #define TSRK3 "3"
1364: #define TSRK3BS "3bs"
1365: #define TSRK4 "4"
1366: #define TSRK5F "5f"
1367: #define TSRK5DP "5dp"
1368: #define TSRK5BS "5bs"
1369: #define TSRK6VR "6vr"
1370: #define TSRK7VR "7vr"
1371: #define TSRK8VR "8vr"
1373: PETSC_EXTERN PetscErrorCode TSRKGetOrder(TS, PetscInt *);
1374: PETSC_EXTERN PetscErrorCode TSRKGetType(TS, TSRKType *);
1375: PETSC_EXTERN PetscErrorCode TSRKSetType(TS, TSRKType);
1376: PETSC_EXTERN PetscErrorCode TSRKGetTableau(TS, PetscInt *, const PetscReal *[], const PetscReal *[], const PetscReal *[], const PetscReal *[], PetscInt *, const PetscReal *[], PetscBool *);
1377: PETSC_EXTERN PetscErrorCode TSRKSetMultirate(TS, PetscBool);
1378: PETSC_EXTERN PetscErrorCode TSRKGetMultirate(TS, PetscBool *);
1379: PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
1380: PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void);
1381: PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void);
1382: PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void);
1384: /*J
1385: TSMPRKType - String with the name of a partitioned Runge-Kutta `TSMPRK` type
1387: Level: beginner
1389: .seealso: [](ch_ts), `TSMPRKSetType()`, `TS`, `TSMPRK`, `TSMPRKRegister()`
1390: J*/
1391: typedef const char *TSMPRKType;
1392: #define TSMPRK2A22 "2a22"
1393: #define TSMPRK2A23 "2a23"
1394: #define TSMPRK2A32 "2a32"
1395: #define TSMPRK2A33 "2a33"
1396: #define TSMPRKP2 "p2"
1397: #define TSMPRKP3 "p3"
1399: PETSC_EXTERN PetscErrorCode TSMPRKGetType(TS, TSMPRKType *);
1400: PETSC_EXTERN PetscErrorCode TSMPRKSetType(TS, TSMPRKType);
1401: 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[]);
1402: PETSC_EXTERN PetscErrorCode TSMPRKInitializePackage(void);
1403: PETSC_EXTERN PetscErrorCode TSMPRKFinalizePackage(void);
1404: PETSC_EXTERN PetscErrorCode TSMPRKRegisterDestroy(void);
1406: /*J
1407: TSIRKType - String with the name of an implicit Runge-Kutta `TSIRK` type
1409: Level: beginner
1411: .seealso: [](ch_ts), `TSIRKSetType()`, `TS`, `TSIRK`, `TSIRKRegister()`
1412: J*/
1413: typedef const char *TSIRKType;
1414: #define TSIRKGAUSS "gauss"
1416: PETSC_EXTERN PetscErrorCode TSIRKGetType(TS, TSIRKType *);
1417: PETSC_EXTERN PetscErrorCode TSIRKSetType(TS, TSIRKType);
1418: PETSC_EXTERN PetscErrorCode TSIRKGetNumStages(TS, PetscInt *);
1419: PETSC_EXTERN PetscErrorCode TSIRKSetNumStages(TS, PetscInt);
1420: PETSC_EXTERN PetscErrorCode TSIRKRegister(const char[], PetscErrorCode (*function)(TS));
1421: PETSC_EXTERN PetscErrorCode TSIRKTableauCreate(TS, PetscInt, const PetscReal *, const PetscReal *, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, const PetscScalar *);
1422: PETSC_EXTERN PetscErrorCode TSIRKInitializePackage(void);
1423: PETSC_EXTERN PetscErrorCode TSIRKFinalizePackage(void);
1424: PETSC_EXTERN PetscErrorCode TSIRKRegisterDestroy(void);
1426: /*J
1427: TSGLEEType - String with the name of a General Linear with Error Estimation `TSGLEE` type
1429: Level: beginner
1431: .seealso: [](ch_ts), `TSGLEESetType()`, `TS`, `TSGLEE`, `TSGLEERegister()`
1432: J*/
1433: typedef const char *TSGLEEType;
1434: #define TSGLEEi1 "BE1"
1435: #define TSGLEE23 "23"
1436: #define TSGLEE24 "24"
1437: #define TSGLEE25I "25i"
1438: #define TSGLEE35 "35"
1439: #define TSGLEEEXRK2A "exrk2a"
1440: #define TSGLEERK32G1 "rk32g1"
1441: #define TSGLEERK285EX "rk285ex"
1443: /*J
1444: TSGLEEMode - String with the mode of error estimation for a General Linear with Error Estimation `TSGLEE` type
1446: Level: beginner
1448: .seealso: [](ch_ts), `TSGLEESetMode()`, `TS`, `TSGLEE`, `TSGLEERegister()`
1449: J*/
1450: PETSC_EXTERN PetscErrorCode TSGLEEGetType(TS, TSGLEEType *);
1451: PETSC_EXTERN PetscErrorCode TSGLEESetType(TS, TSGLEEType);
1452: 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[]);
1453: PETSC_EXTERN PetscErrorCode TSGLEERegisterAll(void);
1454: PETSC_EXTERN PetscErrorCode TSGLEEFinalizePackage(void);
1455: PETSC_EXTERN PetscErrorCode TSGLEEInitializePackage(void);
1456: PETSC_EXTERN PetscErrorCode TSGLEERegisterDestroy(void);
1458: /*J
1459: TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX `TSARKIMEX` type
1461: Options Database Key:
1462: . -ts_arkimex_type (1bee|a2|l2|ars122|2c|2d|2e|prssp2|3|bpr3|ars443|4|5) - set `TSARKIMEX` scheme type, see `TSARKIMEXType`
1464: Level: beginner
1466: .seealso: [](ch_ts), `TSARKIMEXSetType()`, `TS`, `TSARKIMEX`, `TSARKIMEXRegister()`
1467: J*/
1468: typedef const char *TSARKIMEXType;
1469: #define TSARKIMEX1BEE "1bee"
1470: #define TSARKIMEXA2 "a2"
1471: #define TSARKIMEXL2 "l2"
1472: #define TSARKIMEXARS122 "ars122"
1473: #define TSARKIMEX2C "2c"
1474: #define TSARKIMEX2D "2d"
1475: #define TSARKIMEX2E "2e"
1476: #define TSARKIMEXPRSSP2 "prssp2"
1477: #define TSARKIMEX3 "3"
1478: #define TSARKIMEXBPR3 "bpr3"
1479: #define TSARKIMEXARS443 "ars443"
1480: #define TSARKIMEX4 "4"
1481: #define TSARKIMEX5 "5"
1483: PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS, TSARKIMEXType *);
1484: PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS, TSARKIMEXType);
1485: PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS, PetscBool);
1486: PETSC_EXTERN PetscErrorCode TSARKIMEXGetFullyImplicit(TS, PetscBool *);
1487: PETSC_EXTERN PetscErrorCode TSARKIMEXSetFastSlowSplit(TS, PetscBool);
1488: PETSC_EXTERN PetscErrorCode TSARKIMEXGetFastSlowSplit(TS, PetscBool *);
1489: 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[]);
1490: PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void);
1491: PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
1492: PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);
1494: /*J
1495: TSDIRKType - String with the name of a Diagonally Implicit Runge-Kutta `TSDIRK` type
1497: Level: beginner
1499: .seealso: [](ch_ts), `TSDIRKSetType()`, `TS`, `TSDIRK`, `TSDIRKRegister()`
1500: J*/
1501: typedef const char *TSDIRKType;
1502: #define TSDIRKS212 "s212"
1503: #define TSDIRKES122SAL "es122sal"
1504: #define TSDIRKES213SAL "es213sal"
1505: #define TSDIRKES324SAL "es324sal"
1506: #define TSDIRKES325SAL "es325sal"
1507: #define TSDIRK657A "657a"
1508: #define TSDIRKES648SA "es648sa"
1509: #define TSDIRK658A "658a"
1510: #define TSDIRKS659A "s659a"
1511: #define TSDIRK7510SAL "7510sal"
1512: #define TSDIRKES7510SA "es7510sa"
1513: #define TSDIRK759A "759a"
1514: #define TSDIRKS7511SAL "s7511sal"
1515: #define TSDIRK8614A "8614a"
1516: #define TSDIRK8616SAL "8616sal"
1517: #define TSDIRKES8516SAL "es8516sal"
1519: PETSC_EXTERN PetscErrorCode TSDIRKGetType(TS, TSDIRKType *);
1520: PETSC_EXTERN PetscErrorCode TSDIRKSetType(TS, TSDIRKType);
1521: PETSC_EXTERN PetscErrorCode TSDIRKRegister(TSDIRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
1523: /*J
1524: TSRosWType - String with the name of a Rosenbrock-W `TSROSW` type
1526: Level: beginner
1528: .seealso: [](ch_ts), `TSRosWSetType()`, `TS`, `TSROSW`, `TSRosWRegister()`
1529: J*/
1530: typedef const char *TSRosWType;
1531: #define TSROSW2M "2m"
1532: #define TSROSW2P "2p"
1533: #define TSROSWRA3PW "ra3pw"
1534: #define TSROSWRA34PW2 "ra34pw2"
1535: #define TSROSWR34PRW "r34prw"
1536: #define TSROSWR3PRL2 "r3prl2"
1537: #define TSROSWRODAS3 "rodas3"
1538: #define TSROSWRODASPR "rodaspr"
1539: #define TSROSWRODASPR2 "rodaspr2"
1540: #define TSROSWSANDU3 "sandu3"
1541: #define TSROSWASSP3P3S1C "assp3p3s1c"
1542: #define TSROSWLASSP3P4S2C "lassp3p4s2c"
1543: #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
1544: #define TSROSWARK3 "ark3"
1545: #define TSROSWTHETA1 "theta1"
1546: #define TSROSWTHETA2 "theta2"
1547: #define TSROSWGRK4T "grk4t"
1548: #define TSROSWSHAMP4 "shamp4"
1549: #define TSROSWVELDD4 "veldd4"
1550: #define TSROSW4L "4l"
1552: PETSC_EXTERN PetscErrorCode TSRosWGetType(TS, TSRosWType *);
1553: PETSC_EXTERN PetscErrorCode TSRosWSetType(TS, TSRosWType);
1554: PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS, PetscBool);
1555: PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
1556: PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal);
1557: PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void);
1558: PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
1559: PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);
1561: PETSC_EXTERN PetscErrorCode TSBDFSetOrder(TS, PetscInt);
1562: PETSC_EXTERN PetscErrorCode TSBDFGetOrder(TS, PetscInt *);
1564: /*J
1565: TSBasicSymplecticType - String with the name of a basic symplectic integration `TSBASICSYMPLECTIC` type
1567: Level: beginner
1569: .seealso: [](ch_ts), `TSBasicSymplecticSetType()`, `TS`, `TSBASICSYMPLECTIC`, `TSBasicSymplecticRegister()`
1570: J*/
1571: typedef const char *TSBasicSymplecticType;
1572: #define TSBASICSYMPLECTICSIEULER "1"
1573: #define TSBASICSYMPLECTICVELVERLET "2"
1574: #define TSBASICSYMPLECTIC3 "3"
1575: #define TSBASICSYMPLECTIC4 "4"
1577: PETSC_EXTERN PetscErrorCode TSBasicSymplecticSetType(TS, TSBasicSymplecticType);
1578: PETSC_EXTERN PetscErrorCode TSBasicSymplecticGetType(TS, TSBasicSymplecticType *);
1579: PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegister(TSBasicSymplecticType, PetscInt, PetscInt, PetscReal[], PetscReal[]);
1580: PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterAll(void);
1581: PETSC_EXTERN PetscErrorCode TSBasicSymplecticInitializePackage(void);
1582: PETSC_EXTERN PetscErrorCode TSBasicSymplecticFinalizePackage(void);
1583: PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterDestroy(void);
1585: typedef enum {
1586: TS_DG_GONZALEZ,
1587: TS_DG_AVERAGE,
1588: TS_DG_NONE
1589: } TSDGType;
1590: PETSC_EXTERN PetscErrorCode TSDiscGradSetFormulation(TS, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (*)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (*)(TS, PetscReal, Vec, Vec, void *), void *);
1591: PETSC_EXTERN PetscErrorCode TSDiscGradGetFormulation(TS, PetscErrorCode (**)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (**)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (**)(TS, PetscReal, Vec, Vec, void *), void *);
1592: PETSC_EXTERN PetscErrorCode TSDiscGradSetImplicitFormulation(TS, PetscErrorCode (*)(TS, PetscReal, Vec, Vec, Vec, void *), PetscErrorCode (*)(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *));
1593: PETSC_EXTERN PetscErrorCode TSDiscGradSetType(TS, TSDGType);
1594: PETSC_EXTERN PetscErrorCode TSDiscGradGetType(TS, TSDGType *);
1595: PETSC_EXTERN PetscErrorCode TSDiscGradGetX0AndXdot(TS, DM, Vec *, Vec *);
1596: PETSC_EXTERN PetscErrorCode TSDiscGradRestoreX0AndXdot(TS, DM, Vec *, Vec *);
1598: /*
1599: PETSc interface to Sundials
1600: */
1601: #ifdef PETSC_HAVE_SUNDIALS2
1602: typedef enum {
1603: SUNDIALS_ADAMS = 1,
1604: SUNDIALS_BDF = 2
1605: } TSSundialsLmmType;
1606: PETSC_EXTERN const char *const TSSundialsLmmTypes[];
1607: typedef enum {
1608: SUNDIALS_MODIFIED_GS = 1,
1609: SUNDIALS_CLASSICAL_GS = 2
1610: } TSSundialsGramSchmidtType;
1611: PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
1613: PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS, TSSundialsLmmType);
1614: PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS, PC *);
1615: PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS, PetscReal, PetscReal);
1616: PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS, PetscReal);
1617: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS, PetscReal);
1618: PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS, PetscInt *, PetscInt *);
1619: PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS, TSSundialsGramSchmidtType);
1620: PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS, PetscInt);
1621: PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS, PetscReal);
1622: PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS, PetscBool);
1623: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS, PetscInt);
1624: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxord(TS, PetscInt);
1625: PETSC_EXTERN PetscErrorCode TSSundialsSetUseDense(TS, PetscBool);
1626: #endif
1628: PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS, PetscReal);
1629: PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS, PetscReal *);
1630: PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS, PetscBool *);
1631: PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS, PetscBool);
1633: PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS, PetscReal);
1634: PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS, PetscReal, PetscReal, PetscReal);
1635: PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS, PetscReal *, PetscReal *, PetscReal *);
1637: PETSC_EXTERN PetscErrorCode TSAlpha2SetRadius(TS, PetscReal);
1638: PETSC_EXTERN PetscErrorCode TSAlpha2SetParams(TS, PetscReal, PetscReal, PetscReal, PetscReal);
1639: PETSC_EXTERN PetscErrorCode TSAlpha2GetParams(TS, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
1641: /*S
1642: TSAlpha2PredictorFn - A callback to set the predictor (i.e., the initial guess for the nonlinear solver) in
1643: a second-order generalized-alpha time integrator.
1645: Calling Sequence:
1646: + ts - the `TS` context obtained from `TSCreate()`
1647: . X0 - the previous time step's state vector
1648: . V0 - the previous time step's first derivative of the state vector
1649: . A0 - the previous time step's second derivative of the state vector
1650: . X1 - the vector into which the initial guess for the current time step will be written
1651: - ctx - [optional] user-defined context for the predictor evaluation routine (may be `NULL`)
1653: Level: intermediate
1655: Note:
1656: The deprecated `TSAlpha2Predictor` still works as a replacement for `TSAlpha2PredictorFn` *.
1658: .seealso: [](ch_ts), `TS`, `TSAlpha2SetPredictor()`
1659: S*/
1660: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSAlpha2PredictorFn(TS ts, Vec X0, Vec V0, Vec A0, Vec X1, PetscCtx ctx);
1662: PETSC_EXTERN_TYPEDEF typedef TSAlpha2PredictorFn *TSAlpha2Predictor;
1664: PETSC_EXTERN PetscErrorCode TSAlpha2SetPredictor(TS, TSAlpha2PredictorFn *, void *);
1666: PETSC_EXTERN PetscErrorCode TSSetDM(TS, DM);
1667: PETSC_EXTERN PetscErrorCode TSGetDM(TS, DM *);
1669: PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES, Vec, Vec, void *);
1670: PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES, Vec, Mat, Mat, void *);
1672: PETSC_EXTERN PetscErrorCode TSRHSJacobianTest(TS, PetscBool *);
1673: PETSC_EXTERN PetscErrorCode TSRHSJacobianTestTranspose(TS, PetscBool *);
1675: PETSC_EXTERN PetscErrorCode TSGetComputeInitialCondition(TS, PetscErrorCode (**)(TS, Vec));
1676: PETSC_EXTERN PetscErrorCode TSSetComputeInitialCondition(TS, PetscErrorCode (*)(TS, Vec));
1677: PETSC_EXTERN PetscErrorCode TSComputeInitialCondition(TS, Vec);
1678: PETSC_EXTERN PetscErrorCode TSGetComputeExactError(TS, PetscErrorCode (**)(TS, Vec, Vec));
1679: PETSC_EXTERN PetscErrorCode TSSetComputeExactError(TS, PetscErrorCode (*)(TS, Vec, Vec));
1680: PETSC_EXTERN PetscErrorCode TSComputeExactError(TS, Vec, Vec);
1681: PETSC_EXTERN PetscErrorCode PetscConvEstUseTS(PetscConvEst, PetscBool);
1683: PETSC_EXTERN PetscErrorCode TSSetMatStructure(TS, MatStructure);