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: #ifndef PETSCTS_H
6: #define PETSCTS_H
8: #include <petscsnes.h>
9: #include <petscconvest.h>
11: /* SUBMANSEC = TS */
13: /*S
14: TS - Abstract PETSc object that manages all time-steppers (ODE integrators)
16: Level: beginner
18: .seealso: [](integrator_table), [](ch_ts), `TSCreate()`, `TSSetType()`, `TSType`, `SNES`, `KSP`, `PC`, `TSDestroy()`
19: S*/
20: typedef struct _p_TS *TS;
22: /*J
23: TSType - String with the name of a PETSc `TS` method.
25: Level: beginner
27: .seealso: [](integrator_table), [](ch_ts), `TSSetType()`, `TS`, `TSRegister()`
28: J*/
29: typedef const char *TSType;
30: #define TSEULER "euler"
31: #define TSBEULER "beuler"
32: #define TSBASICSYMPLECTIC "basicsymplectic"
33: #define TSPSEUDO "pseudo"
34: #define TSCN "cn"
35: #define TSSUNDIALS "sundials"
36: #define TSRK "rk"
37: #define TSPYTHON "python"
38: #define TSTHETA "theta"
39: #define TSALPHA "alpha"
40: #define TSALPHA2 "alpha2"
41: #define TSGLLE "glle"
42: #define TSGLEE "glee"
43: #define TSSSP "ssp"
44: #define TSARKIMEX "arkimex"
45: #define TSROSW "rosw"
46: #define TSEIMEX "eimex"
47: #define TSMIMEX "mimex"
48: #define TSBDF "bdf"
49: #define TSRADAU5 "radau5"
50: #define TSMPRK "mprk"
51: #define TSDISCGRAD "discgrad"
52: #define TSIRK "irk"
54: /*E
55: TSProblemType - Determines the type of problem this `TS` object is to be used to solve
57: Values:
58: + `TS_LINEAR` - a linear ODE or DAE
59: - `TS_NONLINEAR` - a nonlinear ODE or DAE
61: Level: beginner
63: .seealso: [](ch_ts), `TS`, `TSCreate()`
64: E*/
65: typedef enum {
66: TS_LINEAR,
67: TS_NONLINEAR
68: } TSProblemType;
70: /*E
71: TSEquationType - type of `TS` problem that is solved
73: Level: beginner
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: .seealso: [](ch_ts), `TS`, `TSGetEquationType()`, `TSSetEquationType()`
81: E*/
82: typedef enum {
83: TS_EQ_UNSPECIFIED = -1,
84: TS_EQ_EXPLICIT = 0,
85: TS_EQ_ODE_EXPLICIT = 1,
86: TS_EQ_DAE_SEMI_EXPLICIT_INDEX1 = 100,
87: TS_EQ_DAE_SEMI_EXPLICIT_INDEX2 = 200,
88: TS_EQ_DAE_SEMI_EXPLICIT_INDEX3 = 300,
89: TS_EQ_DAE_SEMI_EXPLICIT_INDEXHI = 500,
90: TS_EQ_IMPLICIT = 1000,
91: TS_EQ_ODE_IMPLICIT = 1001,
92: TS_EQ_DAE_IMPLICIT_INDEX1 = 1100,
93: TS_EQ_DAE_IMPLICIT_INDEX2 = 1200,
94: TS_EQ_DAE_IMPLICIT_INDEX3 = 1300,
95: TS_EQ_DAE_IMPLICIT_INDEXHI = 1500
96: } TSEquationType;
97: PETSC_EXTERN const char *const *TSEquationTypes;
99: /*E
100: TSConvergedReason - reason a `TS` method has converged (integrated to the requested time) or not
102: Values:
103: + `TS_CONVERGED_ITERATING` - this only occurs if `TSGetConvergedReason()` is called during the `TSSolve()`
104: . `TS_CONVERGED_TIME` - the final time was reached
105: . `TS_CONVERGED_ITS` - the maximum number of iterations (time-steps) was reached prior to the final time
106: . `TS_CONVERGED_USER` - user requested termination
107: . `TS_CONVERGED_EVENT` - user requested termination on event detection
108: . `TS_CONVERGED_PSEUDO_FATOL` - stops when function norm decreased by a set amount, used only for `TSPSEUDO`
109: . `TS_CONVERGED_PSEUDO_FRTOL` - stops when function norm decreases below a set amount, used only for `TSPSEUDO`
110: . `TS_DIVERGED_NONLINEAR_SOLVE` - too many nonlinear solve failures have occurred
111: . `TS_DIVERGED_STEP_REJECTED` - too many steps were rejected
112: . `TSFORWARD_DIVERGED_LINEAR_SOLVE` - tangent linear solve failed
113: - `TSADJOINT_DIVERGED_LINEAR_SOLVE` - transposed linear solve failed
115: Level: beginner
117: .seealso: [](ch_ts), `TS`, `TSGetConvergedReason()`
118: E*/
119: typedef enum {
120: TS_CONVERGED_ITERATING = 0,
121: TS_CONVERGED_TIME = 1,
122: TS_CONVERGED_ITS = 2,
123: TS_CONVERGED_USER = 3,
124: TS_CONVERGED_EVENT = 4,
125: TS_CONVERGED_PSEUDO_FATOL = 5,
126: TS_CONVERGED_PSEUDO_FRTOL = 6,
127: TS_DIVERGED_NONLINEAR_SOLVE = -1,
128: TS_DIVERGED_STEP_REJECTED = -2,
129: TSFORWARD_DIVERGED_LINEAR_SOLVE = -3,
130: TSADJOINT_DIVERGED_LINEAR_SOLVE = -4
131: } TSConvergedReason;
132: PETSC_EXTERN const char *const *TSConvergedReasons;
134: /*MC
135: TS_CONVERGED_ITERATING - this only occurs if `TSGetConvergedReason()` is called during the `TSSolve()`
137: Level: beginner
139: .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`
140: M*/
142: /*MC
143: TS_CONVERGED_TIME - the final time was reached
145: Level: beginner
147: .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxTime()`, `TSGetMaxTime()`, `TSGetSolveTime()`
148: M*/
150: /*MC
151: TS_CONVERGED_ITS - the maximum number of iterations (time-steps) was reached prior to the final time
153: Level: beginner
155: .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxSteps()`, `TSGetMaxSteps()`
156: M*/
158: /*MC
159: TS_CONVERGED_USER - user requested termination
161: Level: beginner
163: .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`
164: M*/
166: /*MC
167: TS_CONVERGED_EVENT - user requested termination on event detection
169: Level: beginner
171: .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`
172: M*/
174: /*MC
175: TS_CONVERGED_PSEUDO_FRTOL - stops when function norm decreased by a set amount, used only for `TSPSEUDO`
177: Options Database Key:
178: . -ts_pseudo_frtol <rtol> - use specified rtol
180: Level: beginner
182: .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`, `TS_CONVERGED_PSEUDO_FATOL`
183: M*/
185: /*MC
186: TS_CONVERGED_PSEUDO_FATOL - stops when function norm decreases below a set amount, used only for `TSPSEUDO`
188: Options Database Key:
189: . -ts_pseudo_fatol <atol> - use specified atol
191: Level: beginner
193: .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`, `TS_CONVERGED_PSEUDO_FRTOL`
194: M*/
196: /*MC
197: TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed
199: Level: beginner
201: Note:
202: See `TSSetMaxSNESFailures()` for how to allow more nonlinear solver failures.
204: .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSGetSNES()`, `SNESGetConvergedReason()`, `TSSetMaxSNESFailures()`
205: M*/
207: /*MC
208: TS_DIVERGED_STEP_REJECTED - too many steps were rejected
210: Level: beginner
212: Notes:
213: See `TSSetMaxStepRejections()` for how to allow more step rejections.
215: .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxStepRejections()`
216: M*/
218: /*E
219: TSExactFinalTimeOption - option for handling of final time step
221: Values:
222: + `TS_EXACTFINALTIME_STEPOVER` - Don't do anything if requested final time is exceeded
223: . `TS_EXACTFINALTIME_INTERPOLATE` - Interpolate back to final time
224: - `TS_EXACTFINALTIME_MATCHSTEP` - Adapt final time step to match the final time requested
226: Level: beginner
228: .seealso: [](ch_ts), `TS`, `TSGetConvergedReason()`, `TSSetExactFinalTime()`, `TSGetExactFinalTime()`
229: E*/
230: typedef enum {
231: TS_EXACTFINALTIME_UNSPECIFIED = 0,
232: TS_EXACTFINALTIME_STEPOVER = 1,
233: TS_EXACTFINALTIME_INTERPOLATE = 2,
234: TS_EXACTFINALTIME_MATCHSTEP = 3
235: } TSExactFinalTimeOption;
236: PETSC_EXTERN const char *const TSExactFinalTimeOptions[];
238: /* Logging support */
239: PETSC_EXTERN PetscClassId TS_CLASSID;
240: PETSC_EXTERN PetscClassId DMTS_CLASSID;
241: PETSC_EXTERN PetscClassId TSADAPT_CLASSID;
243: PETSC_EXTERN PetscErrorCode TSInitializePackage(void);
244: PETSC_EXTERN PetscErrorCode TSFinalizePackage(void);
246: PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm, TS *);
247: PETSC_EXTERN PetscErrorCode TSClone(TS, TS *);
248: PETSC_EXTERN PetscErrorCode TSDestroy(TS *);
250: PETSC_EXTERN PetscErrorCode TSSetProblemType(TS, TSProblemType);
251: PETSC_EXTERN PetscErrorCode TSGetProblemType(TS, TSProblemType *);
252: PETSC_EXTERN PetscErrorCode TSMonitor(TS, PetscInt, PetscReal, Vec);
253: PETSC_EXTERN PetscErrorCode TSMonitorSet(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, void *), void *, PetscErrorCode (*)(void **));
254: PETSC_EXTERN PetscErrorCode TSMonitorSetFromOptions(TS, const char[], const char[], const char[], PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *), PetscErrorCode (*)(TS, PetscViewerAndFormat *));
255: PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS);
257: PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS, const char[]);
258: PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS, const char[]);
259: PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS, const char *[]);
260: PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS);
261: PETSC_EXTERN PetscErrorCode TSSetUp(TS);
262: PETSC_EXTERN PetscErrorCode TSReset(TS);
264: PETSC_EXTERN PetscErrorCode TSSetSolution(TS, Vec);
265: PETSC_EXTERN PetscErrorCode TSGetSolution(TS, Vec *);
267: PETSC_EXTERN PetscErrorCode TS2SetSolution(TS, Vec, Vec);
268: PETSC_EXTERN PetscErrorCode TS2GetSolution(TS, Vec *, Vec *);
270: PETSC_EXTERN PetscErrorCode TSGetSolutionComponents(TS, PetscInt *, Vec *);
271: PETSC_EXTERN PetscErrorCode TSGetAuxSolution(TS, Vec *);
272: PETSC_EXTERN PetscErrorCode TSGetTimeError(TS, PetscInt, Vec *);
273: PETSC_EXTERN PetscErrorCode TSSetTimeError(TS, Vec);
275: PETSC_EXTERN PetscErrorCode TSSetRHSJacobianP(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), void *);
276: PETSC_EXTERN PetscErrorCode TSGetRHSJacobianP(TS, Mat *, PetscErrorCode (**)(TS, PetscReal, Vec, Mat, void *), void **);
277: PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianP(TS, PetscReal, Vec, Mat);
278: PETSC_EXTERN PetscErrorCode TSSetIJacobianP(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Vec, PetscReal, Mat, void *), void *);
279: PETSC_EXTERN PetscErrorCode TSComputeIJacobianP(TS, PetscReal, Vec, Vec, PetscReal, Mat, PetscBool);
280: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetQuadratureTS then TSComputeRHSJacobianP") PetscErrorCode TSComputeDRDPFunction(TS, PetscReal, Vec, Vec *);
281: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetQuadratureTS then TSComputeRHSJacobian") PetscErrorCode TSComputeDRDUFunction(TS, PetscReal, Vec, Vec *);
282: 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 *);
283: PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUU(TS, PetscReal, Vec, Vec *, Vec, Vec *);
284: PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUP(TS, PetscReal, Vec, Vec *, Vec, Vec *);
285: PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPU(TS, PetscReal, Vec, Vec *, Vec, Vec *);
286: PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPP(TS, PetscReal, Vec, Vec *, Vec, Vec *);
287: 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 *);
288: PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS, PetscReal, Vec, Vec *, Vec, Vec *);
289: PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS, PetscReal, Vec, Vec *, Vec, Vec *);
290: PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS, PetscReal, Vec, Vec *, Vec, Vec *);
291: PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS, PetscReal, Vec, Vec *, Vec, Vec *);
292: PETSC_EXTERN PetscErrorCode TSSetCostHessianProducts(TS, PetscInt, Vec *, Vec *, Vec);
293: PETSC_EXTERN PetscErrorCode TSGetCostHessianProducts(TS, PetscInt *, Vec **, Vec **, Vec *);
294: PETSC_EXTERN PetscErrorCode TSComputeSNESJacobian(TS, Vec, Mat, Mat);
296: /*S
297: TSTrajectory - Abstract PETSc object that stores the trajectory (solution of ODE/DAE at each time step)
299: Level: advanced
301: .seealso: [](ch_ts), `TS`, `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectorySetType()`, `TSTrajectoryDestroy()`, `TSTrajectoryReset()`
302: S*/
303: typedef struct _p_TSTrajectory *TSTrajectory;
305: /*J
306: TSTrajectoryType - String with the name of a PETSc `TS` trajectory storage method
308: Level: intermediate
310: .seealso: [](ch_ts), `TS`, `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectoryDestroy()`
311: J*/
312: typedef const char *TSTrajectoryType;
313: #define TSTRAJECTORYBASIC "basic"
314: #define TSTRAJECTORYSINGLEFILE "singlefile"
315: #define TSTRAJECTORYMEMORY "memory"
316: #define TSTRAJECTORYVISUALIZATION "visualization"
318: PETSC_EXTERN PetscFunctionList TSTrajectoryList;
319: PETSC_EXTERN PetscClassId TSTRAJECTORY_CLASSID;
320: PETSC_EXTERN PetscBool TSTrajectoryRegisterAllCalled;
322: PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS);
323: PETSC_EXTERN PetscErrorCode TSResetTrajectory(TS);
324: PETSC_EXTERN PetscErrorCode TSRemoveTrajectory(TS);
326: PETSC_EXTERN PetscErrorCode TSTrajectoryCreate(MPI_Comm, TSTrajectory *);
327: PETSC_EXTERN PetscErrorCode TSTrajectoryReset(TSTrajectory);
328: PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy(TSTrajectory *);
329: PETSC_EXTERN PetscErrorCode TSTrajectoryView(TSTrajectory, PetscViewer);
330: PETSC_EXTERN PetscErrorCode TSTrajectorySetType(TSTrajectory, TS, TSTrajectoryType);
331: PETSC_EXTERN PetscErrorCode TSTrajectoryGetType(TSTrajectory, TS, TSTrajectoryType *);
332: PETSC_EXTERN PetscErrorCode TSTrajectorySet(TSTrajectory, TS, PetscInt, PetscReal, Vec);
333: PETSC_EXTERN PetscErrorCode TSTrajectoryGet(TSTrajectory, TS, PetscInt, PetscReal *);
334: PETSC_EXTERN PetscErrorCode TSTrajectoryGetVecs(TSTrajectory, TS, PetscInt, PetscReal *, Vec, Vec);
335: PETSC_EXTERN PetscErrorCode TSTrajectoryGetUpdatedHistoryVecs(TSTrajectory, TS, PetscReal, Vec *, Vec *);
336: PETSC_EXTERN PetscErrorCode TSTrajectoryGetNumSteps(TSTrajectory, PetscInt *);
337: PETSC_EXTERN PetscErrorCode TSTrajectoryRestoreUpdatedHistoryVecs(TSTrajectory, Vec *, Vec *);
338: PETSC_EXTERN PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory, TS);
339: PETSC_EXTERN PetscErrorCode TSTrajectoryRegister(const char[], PetscErrorCode (*)(TSTrajectory, TS));
340: PETSC_EXTERN PetscErrorCode TSTrajectoryRegisterAll(void);
341: PETSC_EXTERN PetscErrorCode TSTrajectorySetUp(TSTrajectory, TS);
342: PETSC_EXTERN PetscErrorCode TSTrajectorySetUseHistory(TSTrajectory, PetscBool);
343: PETSC_EXTERN PetscErrorCode TSTrajectorySetMonitor(TSTrajectory, PetscBool);
344: PETSC_EXTERN PetscErrorCode TSTrajectorySetVariableNames(TSTrajectory, const char *const *);
345: PETSC_EXTERN PetscErrorCode TSTrajectorySetTransform(TSTrajectory, PetscErrorCode (*)(void *, Vec, Vec *), PetscErrorCode (*)(void *), void *);
346: PETSC_EXTERN PetscErrorCode TSTrajectorySetSolutionOnly(TSTrajectory, PetscBool);
347: PETSC_EXTERN PetscErrorCode TSTrajectoryGetSolutionOnly(TSTrajectory, PetscBool *);
348: PETSC_EXTERN PetscErrorCode TSTrajectorySetKeepFiles(TSTrajectory, PetscBool);
349: PETSC_EXTERN PetscErrorCode TSTrajectorySetDirname(TSTrajectory, const char[]);
350: PETSC_EXTERN PetscErrorCode TSTrajectorySetFiletemplate(TSTrajectory, const char[]);
351: PETSC_EXTERN PetscErrorCode TSGetTrajectory(TS, TSTrajectory *);
353: PETSC_EXTERN PetscErrorCode TSSetCostGradients(TS, PetscInt, Vec *, Vec *);
354: PETSC_EXTERN PetscErrorCode TSGetCostGradients(TS, PetscInt *, Vec **, Vec **);
355: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSCreateQuadratureTS() then set up the sub-TS (since version 3.12)") 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 *);
356: PETSC_EXTERN PetscErrorCode TSGetCostIntegral(TS, Vec *);
357: PETSC_EXTERN PetscErrorCode TSComputeCostIntegrand(TS, PetscReal, Vec, Vec);
358: PETSC_EXTERN PetscErrorCode TSCreateQuadratureTS(TS, PetscBool, TS *);
359: PETSC_EXTERN PetscErrorCode TSGetQuadratureTS(TS, PetscBool *, TS *);
361: PETSC_EXTERN PetscErrorCode TSAdjointSetFromOptions(TS, PetscOptionItems *);
362: PETSC_EXTERN PetscErrorCode TSAdjointMonitor(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *);
363: PETSC_EXTERN PetscErrorCode TSAdjointMonitorSet(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *), void *, PetscErrorCode (*)(void **));
364: PETSC_EXTERN PetscErrorCode TSAdjointMonitorCancel(TS);
365: PETSC_EXTERN PetscErrorCode TSAdjointMonitorSetFromOptions(TS, const char[], const char[], const char[], PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, PetscViewerAndFormat *), PetscErrorCode (*)(TS, PetscViewerAndFormat *));
367: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSSetRHSJacobianP()") PetscErrorCode TSAdjointSetRHSJacobian(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), void *);
368: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSComputeRHSJacobianP()") PetscErrorCode TSAdjointComputeRHSJacobian(TS, PetscReal, Vec, Mat);
369: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetQuadratureTS then TSComputeRHSJacobianP") PetscErrorCode TSAdjointComputeDRDPFunction(TS, PetscReal, Vec, Vec *);
370: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetQuadratureTS then TSComputeRHSJacobian") PetscErrorCode TSAdjointComputeDRDYFunction(TS, PetscReal, Vec, Vec *);
371: PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS);
372: PETSC_EXTERN PetscErrorCode TSAdjointSetSteps(TS, PetscInt);
374: PETSC_EXTERN PetscErrorCode TSAdjointStep(TS);
375: PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS);
376: PETSC_EXTERN PetscErrorCode TSAdjointReset(TS);
377: PETSC_EXTERN PetscErrorCode TSAdjointCostIntegral(TS);
378: PETSC_EXTERN PetscErrorCode TSAdjointSetForward(TS, Mat);
379: PETSC_EXTERN PetscErrorCode TSAdjointResetForward(TS);
381: PETSC_EXTERN PetscErrorCode TSForwardSetSensitivities(TS, PetscInt, Mat);
382: PETSC_EXTERN PetscErrorCode TSForwardGetSensitivities(TS, PetscInt *, Mat *);
383: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSCreateQuadratureTS() and TSForwardSetSensitivities() (since version 3.12)") PetscErrorCode TSForwardSetIntegralGradients(TS, PetscInt, Vec *);
384: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSForwardGetSensitivities()") PetscErrorCode TSForwardGetIntegralGradients(TS, PetscInt *, Vec **);
385: PETSC_EXTERN PetscErrorCode TSForwardSetUp(TS);
386: PETSC_EXTERN PetscErrorCode TSForwardReset(TS);
387: PETSC_EXTERN PetscErrorCode TSForwardCostIntegral(TS);
388: PETSC_EXTERN PetscErrorCode TSForwardStep(TS);
389: PETSC_EXTERN PetscErrorCode TSForwardSetInitialSensitivities(TS, Mat);
390: PETSC_EXTERN PetscErrorCode TSForwardGetStages(TS, PetscInt *, Mat *[]);
392: PETSC_EXTERN PetscErrorCode TSSetMaxSteps(TS, PetscInt);
393: PETSC_EXTERN PetscErrorCode TSGetMaxSteps(TS, PetscInt *);
394: PETSC_EXTERN PetscErrorCode TSSetMaxTime(TS, PetscReal);
395: PETSC_EXTERN PetscErrorCode TSGetMaxTime(TS, PetscReal *);
396: PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS, TSExactFinalTimeOption);
397: PETSC_EXTERN PetscErrorCode TSGetExactFinalTime(TS, TSExactFinalTimeOption *);
398: PETSC_EXTERN PetscErrorCode TSSetTimeSpan(TS, PetscInt, PetscReal *);
399: PETSC_EXTERN PetscErrorCode TSGetTimeSpan(TS, PetscInt *, const PetscReal **);
400: PETSC_EXTERN PetscErrorCode TSGetTimeSpanSolutions(TS, PetscInt *, Vec **);
402: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSSetTime[Step]() (since version 3.8)") PetscErrorCode TSSetInitialTimeStep(TS, PetscReal, PetscReal);
403: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSSetMax{Steps|Time}() (since version 3.8)") PetscErrorCode TSSetDuration(TS, PetscInt, PetscReal);
404: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetMax{Steps|Time}() (since version 3.8)") PetscErrorCode TSGetDuration(TS, PetscInt *, PetscReal *);
405: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetStepNumber() (since version 3.8)") PetscErrorCode TSGetTimeStepNumber(TS, PetscInt *);
406: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetStepNumber() (since version 3.8)") PetscErrorCode TSGetTotalSteps(TS, PetscInt *);
408: PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
409: PETSC_EXTERN PetscErrorCode TSMonitorExtreme(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
411: typedef struct _n_TSMonitorDrawCtx *TSMonitorDrawCtx;
412: PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorDrawCtx *);
413: PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *);
414: PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS, PetscInt, PetscReal, Vec, void *);
415: PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionPhase(TS, PetscInt, PetscReal, Vec, void *);
416: PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS, PetscInt, PetscReal, Vec, void *);
417: PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionFunction(TS, PetscInt, PetscReal, Vec, void *);
419: PETSC_EXTERN PetscErrorCode TSAdjointMonitorDefault(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, PetscViewerAndFormat *);
420: PETSC_EXTERN PetscErrorCode TSAdjointMonitorDrawSensi(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *);
422: PETSC_EXTERN PetscErrorCode TSMonitorSolution(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
423: PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS, PetscInt, PetscReal, Vec, void *);
424: PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void *);
426: PETSC_EXTERN PetscErrorCode TSStep(TS);
427: PETSC_EXTERN PetscErrorCode TSEvaluateWLTE(TS, NormType, PetscInt *, PetscReal *);
428: PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS, PetscInt, Vec, PetscBool *);
429: PETSC_EXTERN PetscErrorCode TSSolve(TS, Vec);
430: PETSC_EXTERN PetscErrorCode TSGetEquationType(TS, TSEquationType *);
431: PETSC_EXTERN PetscErrorCode TSSetEquationType(TS, TSEquationType);
432: PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS, TSConvergedReason *);
433: PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS, TSConvergedReason);
434: PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS, PetscReal *);
435: PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS, PetscInt *);
436: PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS, PetscInt *);
437: PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS, PetscInt *);
438: PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS, PetscInt);
439: PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS, PetscInt *);
440: PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS, PetscInt);
441: PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS, PetscBool);
442: PETSC_EXTERN PetscErrorCode TSRestartStep(TS);
443: PETSC_EXTERN PetscErrorCode TSRollBack(TS);
445: PETSC_EXTERN PetscErrorCode TSGetStages(TS, PetscInt *, Vec *[]);
447: PETSC_EXTERN PetscErrorCode TSGetTime(TS, PetscReal *);
448: PETSC_EXTERN PetscErrorCode TSSetTime(TS, PetscReal);
449: PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS, PetscReal *);
450: PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS, PetscReal *);
451: PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS, PetscReal);
452: PETSC_EXTERN PetscErrorCode TSGetStepNumber(TS, PetscInt *);
453: PETSC_EXTERN PetscErrorCode TSSetStepNumber(TS, PetscInt);
455: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS, PetscReal, Vec, Vec, void *);
456: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS, PetscReal, Vec, Mat, Mat, void *);
457: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobianP)(TS, PetscReal, Vec, Mat, void *);
458: PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS, Vec, TSRHSFunction, void *);
459: PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS, Vec *, TSRHSFunction *, void **);
460: PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS, Mat, Mat, TSRHSJacobian, void *);
461: PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS, Mat *, Mat *, TSRHSJacobian *, void **);
462: PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS, PetscBool);
464: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS, PetscReal, Vec, void *);
465: PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS, TSSolutionFunction, void *);
466: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSForcingFunction)(TS, PetscReal, Vec, void *);
467: PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS, TSForcingFunction, void *);
469: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS, PetscReal, Vec, Vec, Vec, void *);
470: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
471: PETSC_EXTERN PetscErrorCode TSSetIFunction(TS, Vec, TSIFunction, void *);
472: PETSC_EXTERN PetscErrorCode TSGetIFunction(TS, Vec *, TSIFunction *, void **);
473: PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS, Mat, Mat, TSIJacobian, void *);
474: PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS, Mat *, Mat *, TSIJacobian *, void **);
476: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Function)(TS, PetscReal, Vec, Vec, Vec, Vec, void *);
477: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Jacobian)(TS, PetscReal, Vec, Vec, Vec, PetscReal, PetscReal, Mat, Mat, void *);
478: PETSC_EXTERN PetscErrorCode TSSetI2Function(TS, Vec, TSI2Function, void *);
479: PETSC_EXTERN PetscErrorCode TSGetI2Function(TS, Vec *, TSI2Function *, void **);
480: PETSC_EXTERN PetscErrorCode TSSetI2Jacobian(TS, Mat, Mat, TSI2Jacobian, void *);
481: PETSC_EXTERN PetscErrorCode TSGetI2Jacobian(TS, Mat *, Mat *, TSI2Jacobian *, void **);
483: PETSC_EXTERN PetscErrorCode TSRHSSplitSetIS(TS, const char[], IS);
484: PETSC_EXTERN PetscErrorCode TSRHSSplitGetIS(TS, const char[], IS *);
485: PETSC_EXTERN PetscErrorCode TSRHSSplitSetRHSFunction(TS, const char[], Vec, TSRHSFunction, void *);
486: PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTS(TS, const char[], TS *);
487: PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTSs(TS, PetscInt *, TS *[]);
488: PETSC_EXTERN PetscErrorCode TSSetUseSplitRHSFunction(TS, PetscBool);
489: PETSC_EXTERN PetscErrorCode TSGetUseSplitRHSFunction(TS, PetscBool *);
491: PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS, PetscReal, Vec, Vec, void *);
492: PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS, PetscReal, Vec, Mat, Mat, void *);
493: PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS, PetscReal, Vec, Vec, Vec, void *);
494: PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
495: PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS, PetscReal, Vec);
496: PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS, PetscReal, Vec);
497: PETSC_EXTERN PetscErrorCode TSComputeIJacobianDefaultColor(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
498: PETSC_EXTERN PetscErrorCode TSPruneIJacobianColor(TS, Mat, Mat);
500: PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
501: PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS, PetscReal));
502: PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS, PetscReal, PetscInt, Vec *));
503: PETSC_EXTERN PetscErrorCode TSSetPostEvaluate(TS, PetscErrorCode (*)(TS));
504: PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
505: PETSC_EXTERN PetscErrorCode TSPreStep(TS);
506: PETSC_EXTERN PetscErrorCode TSPreStage(TS, PetscReal);
507: PETSC_EXTERN PetscErrorCode TSPostStage(TS, PetscReal, PetscInt, Vec *);
508: PETSC_EXTERN PetscErrorCode TSPostEvaluate(TS);
509: PETSC_EXTERN PetscErrorCode TSPostStep(TS);
510: PETSC_EXTERN PetscErrorCode TSInterpolate(TS, PetscReal, Vec);
511: PETSC_EXTERN PetscErrorCode TSSetTolerances(TS, PetscReal, Vec, PetscReal, Vec);
512: PETSC_EXTERN PetscErrorCode TSGetTolerances(TS, PetscReal *, Vec *, PetscReal *, Vec *);
513: PETSC_EXTERN PetscErrorCode TSErrorWeightedNormInfinity(TS, Vec, Vec, PetscReal *, PetscReal *, PetscReal *);
514: PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm2(TS, Vec, Vec, PetscReal *, PetscReal *, PetscReal *);
515: PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *);
516: PETSC_EXTERN PetscErrorCode TSErrorWeightedENormInfinity(TS, Vec, Vec, Vec, PetscReal *, PetscReal *, PetscReal *);
517: PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm2(TS, Vec, Vec, Vec, PetscReal *, PetscReal *, PetscReal *);
518: PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm(TS, Vec, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *);
519: PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS, PetscReal);
520: PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS, PetscReal *);
521: PETSC_EXTERN PetscErrorCode TSSetFunctionDomainError(TS, PetscErrorCode (*)(TS, PetscReal, Vec, PetscBool *));
522: PETSC_EXTERN PetscErrorCode TSFunctionDomainError(TS, PetscReal, Vec, PetscBool *);
524: PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS, PetscErrorCode (*)(TS, PetscReal *, void *), void *);
525: PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS, PetscReal *, void *);
526: PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS, PetscReal *);
527: PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS, PetscReal);
528: PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS, PetscErrorCode (*)(TS, Vec, void *, PetscReal *, PetscBool *), void *);
529: PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS, Vec, void *, PetscReal *, PetscBool *);
530: PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS, Vec, PetscReal *, PetscBool *);
531: PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS, PetscReal);
532: PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);
534: PETSC_EXTERN PetscErrorCode TSPythonSetType(TS, const char[]);
535: PETSC_EXTERN PetscErrorCode TSPythonGetType(TS, const char *[]);
537: PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS, PetscReal, Vec, Vec);
538: PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS, PetscReal, Vec, Mat, Mat);
539: PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS, PetscReal, Vec, Vec, Vec, PetscBool);
540: PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscBool);
541: PETSC_EXTERN PetscErrorCode TSComputeI2Function(TS, PetscReal, Vec, Vec, Vec, Vec);
542: PETSC_EXTERN PetscErrorCode TSComputeI2Jacobian(TS, PetscReal, Vec, Vec, Vec, PetscReal, PetscReal, Mat, Mat);
543: PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS, PetscReal, PetscReal, PetscReal *, PetscReal *);
545: PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS, Vec, Vec);
547: PETSC_EXTERN PetscErrorCode DMTSSetBoundaryLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *);
548: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM, TSRHSFunction, void *);
549: PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM, TSRHSFunction *, void **);
550: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionContextDestroy(DM, PetscErrorCode (*)(void *));
551: PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM, TSRHSJacobian, void *);
552: PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM, TSRHSJacobian *, void **);
553: PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobianContextDestroy(DM, PetscErrorCode (*)(void *));
554: PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM, TSIFunction, void *);
555: PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM, TSIFunction *, void **);
556: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionContextDestroy(DM, PetscErrorCode (*)(void *));
557: PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM, TSIJacobian, void *);
558: PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM, TSIJacobian *, void **);
559: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianContextDestroy(DM, PetscErrorCode (*)(void *));
560: PETSC_EXTERN PetscErrorCode DMTSSetI2Function(DM, TSI2Function, void *);
561: PETSC_EXTERN PetscErrorCode DMTSGetI2Function(DM, TSI2Function *, void **);
562: PETSC_EXTERN PetscErrorCode DMTSSetI2FunctionContextDestroy(DM, PetscErrorCode (*)(void *));
563: PETSC_EXTERN PetscErrorCode DMTSSetI2Jacobian(DM, TSI2Jacobian, void *);
564: PETSC_EXTERN PetscErrorCode DMTSGetI2Jacobian(DM, TSI2Jacobian *, void **);
565: PETSC_EXTERN PetscErrorCode DMTSSetI2JacobianContextDestroy(DM, PetscErrorCode (*)(void *));
567: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSTransientVariable)(TS, Vec, Vec, void *);
568: PETSC_EXTERN PetscErrorCode TSSetTransientVariable(TS, TSTransientVariable, void *);
569: PETSC_EXTERN PetscErrorCode DMTSSetTransientVariable(DM, TSTransientVariable, void *);
570: PETSC_EXTERN PetscErrorCode DMTSGetTransientVariable(DM, TSTransientVariable *, void *);
571: PETSC_EXTERN PetscErrorCode TSComputeTransientVariable(TS, Vec, Vec);
572: PETSC_EXTERN PetscErrorCode TSHasTransientVariable(TS, PetscBool *);
574: PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM, TSSolutionFunction, void *);
575: PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM, TSSolutionFunction *, void **);
576: PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM, TSForcingFunction, void *);
577: PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM, TSForcingFunction *, void **);
578: PETSC_EXTERN PetscErrorCode DMTSCheckResidual(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscReal *);
579: PETSC_EXTERN PetscErrorCode DMTSCheckJacobian(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscBool *, PetscReal *);
580: PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec);
582: PETSC_EXTERN PetscErrorCode DMTSGetIFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, Vec, void *), void **);
583: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, Vec, void *), void *);
584: PETSC_EXTERN PetscErrorCode DMTSGetIJacobianLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void **);
585: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void *);
586: PETSC_EXTERN PetscErrorCode DMTSGetRHSFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, void *), void **);
587: PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *);
588: PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrix(DM);
589: PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrixLumped(DM);
590: PETSC_EXTERN PetscErrorCode DMTSDestroyRHSMassMatrix(DM);
592: PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM, PetscErrorCode (*)(void *, PetscViewer), PetscErrorCode (*)(void **, PetscViewer));
593: PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM, PetscErrorCode (*)(void *, PetscViewer), PetscErrorCode (*)(void **, PetscViewer));
595: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo *, PetscReal, void *, void *, void *);
596: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo *, PetscReal, void *, Mat, Mat, void *);
597: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo *, PetscReal, void *, void *, void *, void *);
598: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo *, PetscReal, void *, void *, PetscReal, Mat, Mat, void *);
600: PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM, InsertMode, PetscErrorCode (*)(DMDALocalInfo *, PetscReal, void *, void *, void *), void *);
601: PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM, PetscErrorCode (*)(DMDALocalInfo *, PetscReal, void *, Mat, Mat, void *), void *);
602: PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM, InsertMode, PetscErrorCode (*)(DMDALocalInfo *, PetscReal, void *, void *, void *, void *), void *);
603: PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM, PetscErrorCode (*)(DMDALocalInfo *, PetscReal, void *, void *, PetscReal, Mat, Mat, void *), void *);
605: PETSC_EXTERN PetscErrorCode DMPlexTSGetGeometryFVM(DM, Vec *, Vec *, PetscReal *);
607: typedef struct _n_TSMonitorLGCtx *TSMonitorLGCtx;
608: typedef struct {
609: Vec ray;
610: VecScatter scatter;
611: PetscViewer viewer;
612: TSMonitorLGCtx lgctx;
613: } TSMonitorDMDARayCtx;
614: PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void **);
615: PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS, PetscInt, PetscReal, Vec, void *);
616: PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS, PetscInt, PetscReal, Vec, void *);
618: /* Dynamic creation and loading functions */
619: PETSC_EXTERN PetscFunctionList TSList;
620: PETSC_EXTERN PetscErrorCode TSGetType(TS, TSType *);
621: PETSC_EXTERN PetscErrorCode TSSetType(TS, TSType);
622: PETSC_EXTERN PetscErrorCode TSRegister(const char[], PetscErrorCode (*)(TS));
624: PETSC_EXTERN PetscErrorCode TSGetSNES(TS, SNES *);
625: PETSC_EXTERN PetscErrorCode TSSetSNES(TS, SNES);
626: PETSC_EXTERN PetscErrorCode TSGetKSP(TS, KSP *);
628: PETSC_EXTERN PetscErrorCode TSView(TS, PetscViewer);
629: PETSC_EXTERN PetscErrorCode TSLoad(TS, PetscViewer);
630: PETSC_EXTERN PetscErrorCode TSViewFromOptions(TS, PetscObject, const char[]);
631: PETSC_EXTERN PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory, PetscObject, const char[]);
633: #define TS_FILE_CLASSID 1211225
635: PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS, void *);
636: PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS, void *);
638: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtx *);
639: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *);
640: PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS, PetscInt, PetscReal, Vec, void *);
641: PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS, PetscInt, PetscReal, Vec, void *);
642: PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS, const char *const *);
643: PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS, const char *const **);
644: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx, const char *const *);
645: PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS, const char *const *);
646: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx, const char *const *);
647: PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS, PetscErrorCode (*)(void *, Vec, Vec *), PetscErrorCode (*)(void *), void *);
648: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx, PetscErrorCode (*)(void *, Vec, Vec *), PetscErrorCode (*)(void *), void *);
649: PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS, PetscInt, PetscReal, Vec, void *);
650: PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS, PetscInt, PetscReal, Vec, void *);
651: PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS, PetscInt, PetscReal, Vec, void *);
652: PETSC_EXTERN PetscErrorCode TSMonitorError(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
653: PETSC_EXTERN PetscErrorCode TSDMSwarmMonitorMoments(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
655: struct _n_TSMonitorLGCtxNetwork {
656: PetscInt nlg;
657: PetscDrawLG *lg;
658: PetscBool semilogy;
659: PetscInt howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
660: };
661: typedef struct _n_TSMonitorLGCtxNetwork *TSMonitorLGCtxNetwork;
662: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxNetworkDestroy(TSMonitorLGCtxNetwork *);
663: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxNetworkCreate(TS, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtxNetwork *);
664: PETSC_EXTERN PetscErrorCode TSMonitorLGCtxNetworkSolution(TS, PetscInt, PetscReal, Vec, void *);
666: typedef struct _n_TSMonitorEnvelopeCtx *TSMonitorEnvelopeCtx;
667: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxCreate(TS, TSMonitorEnvelopeCtx *);
668: PETSC_EXTERN PetscErrorCode TSMonitorEnvelope(TS, PetscInt, PetscReal, Vec, void *);
669: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeGetBounds(TS, Vec *, Vec *);
670: PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *);
672: typedef struct _n_TSMonitorSPEigCtx *TSMonitorSPEigCtx;
673: PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorSPEigCtx *);
674: PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *);
675: PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS, PetscInt, PetscReal, Vec, void *);
677: typedef struct _n_TSMonitorSPCtx *TSMonitorSPCtx;
678: PETSC_EXTERN PetscErrorCode TSMonitorSPCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscBool, PetscBool, TSMonitorSPCtx *);
679: PETSC_EXTERN PetscErrorCode TSMonitorSPCtxDestroy(TSMonitorSPCtx *);
680: PETSC_EXTERN PetscErrorCode TSMonitorSPSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
682: typedef struct _n_TSMonitorHGCtx *TSMonitorHGCtx;
683: PETSC_EXTERN PetscErrorCode TSMonitorHGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscInt, PetscBool, TSMonitorHGCtx *);
684: PETSC_EXTERN PetscErrorCode TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
685: PETSC_EXTERN PetscErrorCode TSMonitorHGCtxDestroy(TSMonitorHGCtx *);
686: PETSC_EXTERN PetscErrorCode TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
688: PETSC_EXTERN PetscErrorCode TSSetEventHandler(TS, PetscInt, PetscInt[], PetscBool[], PetscErrorCode (*)(TS, PetscReal, Vec, PetscScalar[], void *), PetscErrorCode (*)(TS, PetscInt, PetscInt[], PetscReal, Vec, PetscBool, void *), void *);
689: PETSC_EXTERN PetscErrorCode TSSetPostEventIntervalStep(TS, PetscReal);
690: PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS, PetscReal, PetscReal[]);
691: PETSC_EXTERN PetscErrorCode TSGetNumEvents(TS, PetscInt *);
693: /*J
694: TSSSPType - string with the name of a `TSSSP` scheme.
696: Level: beginner
698: .seealso: [](ch_ts), `TSSSPSetType()`, `TS`, `TSSSP`
699: J*/
700: typedef const char *TSSSPType;
701: #define TSSSPRKS2 "rks2"
702: #define TSSSPRKS3 "rks3"
703: #define TSSSPRK104 "rk104"
705: PETSC_EXTERN PetscErrorCode TSSSPSetType(TS, TSSSPType);
706: PETSC_EXTERN PetscErrorCode TSSSPGetType(TS, TSSSPType *);
707: PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS, PetscInt);
708: PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS, PetscInt *);
709: PETSC_EXTERN PetscErrorCode TSSSPInitializePackage(void);
710: PETSC_EXTERN PetscErrorCode TSSSPFinalizePackage(void);
711: PETSC_EXTERN PetscFunctionList TSSSPList;
713: /*S
714: TSAdapt - Abstract object that manages time-step adaptivity
716: Level: beginner
718: .seealso: [](ch_ts), `TS`, `TSAdaptCreate()`, `TSAdaptType`
719: S*/
720: typedef struct _p_TSAdapt *TSAdapt;
722: /*J
723: TSAdaptType - String with the name of `TSAdapt` scheme.
725: Level: beginner
727: .seealso: [](ch_ts), `TSAdaptSetType()`, `TS`, `TSAdapt`
728: J*/
729: typedef const char *TSAdaptType;
730: #define TSADAPTNONE "none"
731: #define TSADAPTBASIC "basic"
732: #define TSADAPTDSP "dsp"
733: #define TSADAPTCFL "cfl"
734: #define TSADAPTGLEE "glee"
735: #define TSADAPTHISTORY "history"
737: PETSC_EXTERN PetscErrorCode TSGetAdapt(TS, TSAdapt *);
738: PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[], PetscErrorCode (*)(TSAdapt));
739: PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void);
740: PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
741: PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm, TSAdapt *);
742: PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt, TSAdaptType);
743: PETSC_EXTERN PetscErrorCode TSAdaptGetType(TSAdapt, TSAdaptType *);
744: PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt, const char[]);
745: PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
746: PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt, const char[], PetscInt, PetscInt, PetscReal, PetscReal, PetscBool);
747: PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt, PetscInt *, const PetscInt **, const PetscInt **, const PetscReal **, const PetscReal **);
748: PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt, TS, PetscReal, PetscInt *, PetscReal *, PetscBool *);
749: PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt, TS, PetscReal, Vec, PetscBool *);
750: PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt, PetscViewer);
751: PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt, PetscViewer);
752: PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(TSAdapt, PetscOptionItems *);
753: PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt);
754: PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt *);
755: PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt, PetscBool);
756: PETSC_EXTERN PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt, PetscBool);
757: PETSC_EXTERN PetscErrorCode TSAdaptSetSafety(TSAdapt, PetscReal, PetscReal);
758: PETSC_EXTERN PetscErrorCode TSAdaptGetSafety(TSAdapt, PetscReal *, PetscReal *);
759: PETSC_EXTERN PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt, PetscReal);
760: PETSC_EXTERN PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt, PetscReal *);
761: PETSC_EXTERN PetscErrorCode TSAdaptSetClip(TSAdapt, PetscReal, PetscReal);
762: PETSC_EXTERN PetscErrorCode TSAdaptGetClip(TSAdapt, PetscReal *, PetscReal *);
763: PETSC_EXTERN PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt, PetscReal);
764: PETSC_EXTERN PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt, PetscReal *);
765: PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt, PetscReal, PetscReal);
766: PETSC_EXTERN PetscErrorCode TSAdaptGetStepLimits(TSAdapt, PetscReal *, PetscReal *);
767: PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt, PetscErrorCode (*)(TSAdapt, TS, PetscReal, Vec, PetscBool *));
768: PETSC_EXTERN PetscErrorCode TSAdaptHistorySetHistory(TSAdapt, PetscInt n, PetscReal hist[], PetscBool);
769: PETSC_EXTERN PetscErrorCode TSAdaptHistorySetTrajectory(TSAdapt, TSTrajectory, PetscBool);
770: PETSC_EXTERN PetscErrorCode TSAdaptHistoryGetStep(TSAdapt, PetscInt, PetscReal *, PetscReal *);
771: PETSC_EXTERN PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt, PetscInt);
772: PETSC_EXTERN PetscErrorCode TSAdaptDSPSetFilter(TSAdapt, const char *);
773: PETSC_EXTERN PetscErrorCode TSAdaptDSPSetPID(TSAdapt, PetscReal, PetscReal, PetscReal);
775: /*S
776: TSGLLEAdapt - Abstract object that manages time-step adaptivity for `TSGLLE`
778: Level: beginner
780: Developer Note:
781: This functionality should be replaced by the `TSAdapt`.
783: .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLEAdaptCreate()`, `TSGLLEAdaptType`
784: S*/
785: typedef struct _p_TSGLLEAdapt *TSGLLEAdapt;
787: /*J
788: TSGLLEAdaptType - String with the name of `TSGLLEAdapt` scheme
790: Level: beginner
792: Developer Note:
793: This functionality should be replaced by the `TSAdaptType`.
795: .seealso: [](ch_ts), `TSGLLEAdaptSetType()`, `TS`
796: J*/
797: typedef const char *TSGLLEAdaptType;
798: #define TSGLLEADAPT_NONE "none"
799: #define TSGLLEADAPT_SIZE "size"
800: #define TSGLLEADAPT_BOTH "both"
802: PETSC_EXTERN PetscErrorCode TSGLLEAdaptRegister(const char[], PetscErrorCode (*)(TSGLLEAdapt));
803: PETSC_EXTERN PetscErrorCode TSGLLEAdaptInitializePackage(void);
804: PETSC_EXTERN PetscErrorCode TSGLLEAdaptFinalizePackage(void);
805: PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate(MPI_Comm, TSGLLEAdapt *);
806: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt, TSGLLEAdaptType);
807: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt, const char[]);
808: PETSC_EXTERN PetscErrorCode TSGLLEAdaptChoose(TSGLLEAdapt, PetscInt, const PetscInt[], const PetscReal[], const PetscReal[], PetscInt, PetscReal, PetscReal, PetscInt *, PetscReal *, PetscBool *);
809: PETSC_EXTERN PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt, PetscViewer);
810: PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetFromOptions(TSGLLEAdapt, PetscOptionItems *);
811: PETSC_EXTERN PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt *);
813: /*J
814: TSGLLEAcceptType - String with the name of `TSGLLEAccept` scheme
816: Level: beginner
818: .seealso: [](ch_ts), `TSGLLESetAcceptType()`, `TS`, `TSGLLEAccept`
819: J*/
820: typedef const char *TSGLLEAcceptType;
821: #define TSGLLEACCEPT_ALWAYS "always"
823: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLLEAcceptFunction)(TS, PetscReal, PetscReal, const PetscReal[], PetscBool *);
824: PETSC_EXTERN PetscErrorCode TSGLLEAcceptRegister(const char[], TSGLLEAcceptFunction);
826: /*J
827: TSGLLEType - string with the name of a General Linear `TSGLLE` type
829: Level: beginner
831: .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLESetType()`, `TSGLLERegister()`, `TSGLLEAccept`
832: J*/
833: typedef const char *TSGLLEType;
834: #define TSGLLE_IRKS "irks"
836: PETSC_EXTERN PetscErrorCode TSGLLERegister(const char[], PetscErrorCode (*)(TS));
837: PETSC_EXTERN PetscErrorCode TSGLLEInitializePackage(void);
838: PETSC_EXTERN PetscErrorCode TSGLLEFinalizePackage(void);
839: PETSC_EXTERN PetscErrorCode TSGLLESetType(TS, TSGLLEType);
840: PETSC_EXTERN PetscErrorCode TSGLLEGetAdapt(TS, TSGLLEAdapt *);
841: PETSC_EXTERN PetscErrorCode TSGLLESetAcceptType(TS, TSGLLEAcceptType);
843: /*J
844: TSEIMEXType - String with the name of an Extrapolated IMEX `TSEIMEX` type
846: Level: beginner
848: .seealso: [](ch_ts), `TSEIMEXSetType()`, `TS`, `TSEIMEX`, `TSEIMEXRegister()`
849: J*/
850: #define TSEIMEXType char *
852: PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts, PetscInt);
853: PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts, PetscInt, PetscInt);
854: PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS, PetscBool);
856: /*J
857: TSRKType - String with the name of a Runge-Kutta `TSRK` type
859: Level: beginner
861: .seealso: [](ch_ts), `TS`, `TSRKSetType()`, `TS`, `TSRK`, `TSRKRegister()`
862: J*/
863: typedef const char *TSRKType;
864: #define TSRK1FE "1fe"
865: #define TSRK2A "2a"
866: #define TSRK2B "2b"
867: #define TSRK3 "3"
868: #define TSRK3BS "3bs"
869: #define TSRK4 "4"
870: #define TSRK5F "5f"
871: #define TSRK5DP "5dp"
872: #define TSRK5BS "5bs"
873: #define TSRK6VR "6vr"
874: #define TSRK7VR "7vr"
875: #define TSRK8VR "8vr"
877: PETSC_EXTERN PetscErrorCode TSRKGetOrder(TS, PetscInt *);
878: PETSC_EXTERN PetscErrorCode TSRKGetType(TS, TSRKType *);
879: PETSC_EXTERN PetscErrorCode TSRKSetType(TS, TSRKType);
880: PETSC_EXTERN PetscErrorCode TSRKGetTableau(TS, PetscInt *, const PetscReal **, const PetscReal **, const PetscReal **, const PetscReal **, PetscInt *, const PetscReal **, PetscBool *);
881: PETSC_EXTERN PetscErrorCode TSRKSetMultirate(TS, PetscBool);
882: PETSC_EXTERN PetscErrorCode TSRKGetMultirate(TS, PetscBool *);
883: PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
884: PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void);
885: PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void);
886: PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void);
888: /*J
889: TSMPRKType - String with the name of a Partitioned Runge-Kutta `TSMPRK` type
891: Level: beginner
893: .seealso: [](ch_ts), `TSMPRKSetType()`, `TS`, `TSMPRK`, `TSMPRKRegister()`
894: J*/
895: typedef const char *TSMPRKType;
896: #define TSMPRK2A22 "2a22"
897: #define TSMPRK2A23 "2a23"
898: #define TSMPRK2A32 "2a32"
899: #define TSMPRK2A33 "2a33"
900: #define TSMPRKP2 "p2"
901: #define TSMPRKP3 "p3"
903: PETSC_EXTERN PetscErrorCode TSMPRKGetType(TS ts, TSMPRKType *);
904: PETSC_EXTERN PetscErrorCode TSMPRKSetType(TS ts, TSMPRKType);
905: 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[]);
906: PETSC_EXTERN PetscErrorCode TSMPRKInitializePackage(void);
907: PETSC_EXTERN PetscErrorCode TSMPRKFinalizePackage(void);
908: PETSC_EXTERN PetscErrorCode TSMPRKRegisterDestroy(void);
910: /*J
911: TSIRKType - String with the name of an implicit Runge-Kutta `TSIRK` type
913: Level: beginner
915: .seealso: [](ch_ts), `TSIRKSetType()`, `TS`, `TSIRK`, `TSIRKRegister()`
916: J*/
917: typedef const char *TSIRKType;
918: #define TSIRKGAUSS "gauss"
920: PETSC_EXTERN PetscErrorCode TSIRKGetOrder(TS, PetscInt *);
921: PETSC_EXTERN PetscErrorCode TSIRKGetType(TS, TSIRKType *);
922: PETSC_EXTERN PetscErrorCode TSIRKSetType(TS, TSIRKType);
923: PETSC_EXTERN PetscErrorCode TSIRKGetNumStages(TS, PetscInt *);
924: PETSC_EXTERN PetscErrorCode TSIRKSetNumStages(TS, PetscInt);
925: PETSC_EXTERN PetscErrorCode TSIRKRegister(const char[], PetscErrorCode (*function)(TS));
926: PETSC_EXTERN PetscErrorCode TSIRKTableauCreate(TS, PetscInt, const PetscReal *, const PetscReal *, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, const PetscScalar *);
927: PETSC_EXTERN PetscErrorCode TSIRKInitializePackage(void);
928: PETSC_EXTERN PetscErrorCode TSIRKFinalizePackage(void);
929: PETSC_EXTERN PetscErrorCode TSIRKRegisterDestroy(void);
931: /*J
932: TSGLEEType - String with the name of a General Linear with Error Estimation `TSGLEE` type
934: Level: beginner
936: .seealso: [](ch_ts), `TSGLEESetType()`, `TS`, `TSGLEE`, `TSGLEERegister()`
937: J*/
938: typedef const char *TSGLEEType;
939: #define TSGLEEi1 "BE1"
940: #define TSGLEE23 "23"
941: #define TSGLEE24 "24"
942: #define TSGLEE25I "25i"
943: #define TSGLEE35 "35"
944: #define TSGLEEEXRK2A "exrk2a"
945: #define TSGLEERK32G1 "rk32g1"
946: #define TSGLEERK285EX "rk285ex"
948: /*J
949: TSGLEEMode - String with the mode of error estimation for a General Linear with Error Estimation `TSGLEE` type
951: Level: beginner
953: .seealso: [](ch_ts), `TSGLEESetMode()`, `TS`, `TSGLEE`, `TSGLEERegister()`
954: J*/
955: PETSC_EXTERN PetscErrorCode TSGLEEGetType(TS ts, TSGLEEType *);
956: PETSC_EXTERN PetscErrorCode TSGLEESetType(TS ts, TSGLEEType);
957: 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[]);
958: PETSC_EXTERN PetscErrorCode TSGLEEFinalizePackage(void);
959: PETSC_EXTERN PetscErrorCode TSGLEEInitializePackage(void);
960: PETSC_EXTERN PetscErrorCode TSGLEERegisterDestroy(void);
962: /*J
963: TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX `TSARKIMEX` type
965: Level: beginner
967: .seealso: [](ch_ts), `TSARKIMEXSetType()`, `TS`, `TSARKIMEX`, `TSARKIMEXRegister()`
968: J*/
969: typedef const char *TSARKIMEXType;
970: #define TSARKIMEX1BEE "1bee"
971: #define TSARKIMEXA2 "a2"
972: #define TSARKIMEXL2 "l2"
973: #define TSARKIMEXARS122 "ars122"
974: #define TSARKIMEX2C "2c"
975: #define TSARKIMEX2D "2d"
976: #define TSARKIMEX2E "2e"
977: #define TSARKIMEXPRSSP2 "prssp2"
978: #define TSARKIMEX3 "3"
979: #define TSARKIMEXBPR3 "bpr3"
980: #define TSARKIMEXARS443 "ars443"
981: #define TSARKIMEX4 "4"
982: #define TSARKIMEX5 "5"
983: PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts, TSARKIMEXType *);
984: PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts, TSARKIMEXType);
985: PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS, PetscBool);
986: PETSC_EXTERN PetscErrorCode TSARKIMEXGetFullyImplicit(TS, PetscBool *);
987: 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[]);
988: PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void);
989: PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
990: PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);
992: /*J
993: TSRosWType - String with the name of a Rosenbrock-W `TSROSW` type
995: Level: beginner
997: .seealso: [](ch_ts), `TSRosWSetType()`, `TS`, `TSROSW`, `TSRosWRegister()`
998: J*/
999: typedef const char *TSRosWType;
1000: #define TSROSW2M "2m"
1001: #define TSROSW2P "2p"
1002: #define TSROSWRA3PW "ra3pw"
1003: #define TSROSWRA34PW2 "ra34pw2"
1004: #define TSROSWRODAS3 "rodas3"
1005: #define TSROSWSANDU3 "sandu3"
1006: #define TSROSWASSP3P3S1C "assp3p3s1c"
1007: #define TSROSWLASSP3P4S2C "lassp3p4s2c"
1008: #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
1009: #define TSROSWARK3 "ark3"
1010: #define TSROSWTHETA1 "theta1"
1011: #define TSROSWTHETA2 "theta2"
1012: #define TSROSWGRK4T "grk4t"
1013: #define TSROSWSHAMP4 "shamp4"
1014: #define TSROSWVELDD4 "veldd4"
1015: #define TSROSW4L "4l"
1017: PETSC_EXTERN PetscErrorCode TSRosWGetType(TS, TSRosWType *);
1018: PETSC_EXTERN PetscErrorCode TSRosWSetType(TS, TSRosWType);
1019: PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS, PetscBool);
1020: PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
1021: PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal);
1022: PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void);
1023: PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
1024: PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);
1026: PETSC_EXTERN PetscErrorCode TSBDFSetOrder(TS, PetscInt);
1027: PETSC_EXTERN PetscErrorCode TSBDFGetOrder(TS, PetscInt *);
1029: /*J
1030: TSBasicSymplecticType - String with the name of a basic symplectic integration `TSBASICSYMPLECTIC` type
1032: Level: beginner
1034: .seealso: [](ch_ts), `TSBasicSymplecticSetType()`, `TS`, `TSBASICSYMPLECTIC`, `TSBasicSymplecticRegister()`
1035: J*/
1036: typedef const char *TSBasicSymplecticType;
1037: #define TSBASICSYMPLECTICSIEULER "1"
1038: #define TSBASICSYMPLECTICVELVERLET "2"
1039: #define TSBASICSYMPLECTIC3 "3"
1040: #define TSBASICSYMPLECTIC4 "4"
1041: PETSC_EXTERN PetscErrorCode TSBasicSymplecticSetType(TS, TSBasicSymplecticType);
1042: PETSC_EXTERN PetscErrorCode TSBasicSymplecticGetType(TS, TSBasicSymplecticType *);
1043: PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegister(TSBasicSymplecticType, PetscInt, PetscInt, PetscReal[], PetscReal[]);
1044: PETSC_EXTERN PetscErrorCode TSBasicSymplecticInitializePackage(void);
1045: PETSC_EXTERN PetscErrorCode TSBasicSymplecticFinalizePackage(void);
1046: PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterDestroy(void);
1048: /*J
1049: TSDISCGRAD - The Discrete Gradient integrator is a timestepper for Hamiltonian systems designed to conserve the first integral (energy),
1050: but also has the property for some systems of monotonicity in a functional.
1052: Level: beginner
1054: .seealso: [](ch_ts), `TS`, TSDiscGradSetFormulation()`, `TSDiscGradGetFormulation()`
1055: J*/
1056: PETSC_EXTERN PetscErrorCode TSDiscGradSetFormulation(TS, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (*)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (*)(TS, PetscReal, Vec, Vec, void *), void *);
1057: PETSC_EXTERN PetscErrorCode TSDiscGradIsGonzalez(TS, PetscBool *);
1058: PETSC_EXTERN PetscErrorCode TSDiscGradUseGonzalez(TS, PetscBool);
1060: /*
1061: PETSc interface to Sundials
1062: */
1063: #ifdef PETSC_HAVE_SUNDIALS2
1064: typedef enum {
1065: SUNDIALS_ADAMS = 1,
1066: SUNDIALS_BDF = 2
1067: } TSSundialsLmmType;
1068: PETSC_EXTERN const char *const TSSundialsLmmTypes[];
1069: typedef enum {
1070: SUNDIALS_MODIFIED_GS = 1,
1071: SUNDIALS_CLASSICAL_GS = 2
1072: } TSSundialsGramSchmidtType;
1073: PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
1074: PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS, TSSundialsLmmType);
1075: PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS, PC *);
1076: PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS, PetscReal, PetscReal);
1077: PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS, PetscReal);
1078: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS, PetscReal);
1079: PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS, PetscInt *, PetscInt *);
1080: PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS, TSSundialsGramSchmidtType);
1081: PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS, PetscInt);
1082: PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS, PetscReal);
1083: PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS, PetscBool);
1084: PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS, PetscInt *, long *[], double *[]);
1085: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS, PetscInt);
1086: PETSC_EXTERN PetscErrorCode TSSundialsSetMaxord(TS, PetscInt);
1087: PETSC_EXTERN PetscErrorCode TSSundialsSetUseDense(TS, PetscBool);
1088: #endif
1090: PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS, PetscReal);
1091: PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS, PetscReal *);
1092: PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS, PetscBool *);
1093: PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS, PetscBool);
1095: PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS, PetscReal);
1096: PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS, PetscReal, PetscReal, PetscReal);
1097: PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS, PetscReal *, PetscReal *, PetscReal *);
1099: PETSC_EXTERN PetscErrorCode TSAlpha2SetRadius(TS, PetscReal);
1100: PETSC_EXTERN PetscErrorCode TSAlpha2SetParams(TS, PetscReal, PetscReal, PetscReal, PetscReal);
1101: PETSC_EXTERN PetscErrorCode TSAlpha2GetParams(TS, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
1103: PETSC_EXTERN PetscErrorCode TSSetDM(TS, DM);
1104: PETSC_EXTERN PetscErrorCode TSGetDM(TS, DM *);
1106: PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES, Vec, Vec, void *);
1107: PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES, Vec, Mat, Mat, void *);
1109: PETSC_EXTERN PetscErrorCode TSRHSJacobianTest(TS, PetscBool *);
1110: PETSC_EXTERN PetscErrorCode TSRHSJacobianTestTranspose(TS, PetscBool *);
1112: PETSC_EXTERN PetscErrorCode TSGetComputeInitialCondition(TS, PetscErrorCode (**)(TS, Vec));
1113: PETSC_EXTERN PetscErrorCode TSSetComputeInitialCondition(TS, PetscErrorCode (*)(TS, Vec));
1114: PETSC_EXTERN PetscErrorCode TSComputeInitialCondition(TS, Vec);
1115: PETSC_EXTERN PetscErrorCode TSGetComputeExactError(TS, PetscErrorCode (**)(TS, Vec, Vec));
1116: PETSC_EXTERN PetscErrorCode TSSetComputeExactError(TS, PetscErrorCode (*)(TS, Vec, Vec));
1117: PETSC_EXTERN PetscErrorCode TSComputeExactError(TS, Vec, Vec);
1118: PETSC_EXTERN PetscErrorCode PetscConvEstUseTS(PetscConvEst, PetscBool);
1120: PETSC_EXTERN PetscErrorCode TSSetMatStructure(TS, MatStructure);
1121: #endif