xref: /petsc/include/petscts.h (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1818ad0c1SBarry Smith /*
2316643e7SJed Brown    User interface for the timestepping package. This package
3f64a0f93SLois Curfman McInnes    is for use in solving time-dependent PDEs.
4818ad0c1SBarry Smith */
5a4963045SJacob Faibussowitsch #pragma once
6ac09b921SBarry Smith 
72c8e378dSBarry Smith #include <petscsnes.h>
8900f6b5bSMatthew G. Knepley #include <petscconvest.h>
9818ad0c1SBarry Smith 
1014d0ab18SJacob Faibussowitsch /*I <petscts.h> I*/
1114d0ab18SJacob Faibussowitsch 
12ac09b921SBarry Smith /* SUBMANSEC = TS */
13ac09b921SBarry Smith 
14435da068SBarry Smith /*S
150b4b7b1cSBarry Smith    TS - Abstract PETSc object that manages integrating an ODE.
16435da068SBarry Smith 
17435da068SBarry Smith    Level: beginner
18435da068SBarry Smith 
191cc06b55SBarry Smith .seealso: [](integrator_table), [](ch_ts), `TSCreate()`, `TSSetType()`, `TSType`, `SNES`, `KSP`, `PC`, `TSDestroy()`
20435da068SBarry Smith S*/
21f09e8eb9SSatish Balay typedef struct _p_TS *TS;
22435da068SBarry Smith 
2376bdecfbSBarry Smith /*J
240b4b7b1cSBarry Smith    TSType - String with the name of a PETSc `TS` method. These are all the time/ODE integrators that PETSc provides.
25435da068SBarry Smith 
26435da068SBarry Smith    Level: beginner
27435da068SBarry Smith 
280b4b7b1cSBarry Smith    Note:
290b4b7b1cSBarry Smith    Use `TSSetType()` or the options database key `-ts_type` to set the ODE integrator method to use with a given `TS` object
300b4b7b1cSBarry Smith 
311cc06b55SBarry Smith .seealso: [](integrator_table), [](ch_ts), `TSSetType()`, `TS`, `TSRegister()`
3276bdecfbSBarry Smith J*/
3319fd82e9SBarry Smith typedef const char *TSType;
349596e0b4SJed Brown #define TSEULER           "euler"
359596e0b4SJed Brown #define TSBEULER          "beuler"
36e49d4f37SHong Zhang #define TSBASICSYMPLECTIC "basicsymplectic"
379596e0b4SJed Brown #define TSPSEUDO          "pseudo"
384d91e141SJed Brown #define TSCN              "cn"
399596e0b4SJed Brown #define TSSUNDIALS        "sundials"
404d91e141SJed Brown #define TSRK              "rk"
419596e0b4SJed Brown #define TSPYTHON          "python"
429596e0b4SJed Brown #define TSTHETA           "theta"
4388df8a41SLisandro Dalcin #define TSALPHA           "alpha"
44818efac9SLisandro Dalcin #define TSALPHA2          "alpha2"
4526d28e4eSEmil Constantinescu #define TSGLLE            "glle"
46b6a60446SDebojyoti Ghosh #define TSGLEE            "glee"
47ef7bb5aaSJed Brown #define TSSSP             "ssp"
488a381b04SJed Brown #define TSARKIMEX         "arkimex"
49e27a552bSJed Brown #define TSROSW            "rosw"
507d5697caSHong Zhang #define TSEIMEX           "eimex"
51abadfa0fSMatthew G. Knepley #define TSMIMEX           "mimex"
52211a84d6SLisandro Dalcin #define TSBDF             "bdf"
53d249a78fSBarry Smith #define TSRADAU5          "radau5"
544b84eec9SHong Zhang #define TSMPRK            "mprk"
5540be0ff1SMatthew G. Knepley #define TSDISCGRAD        "discgrad"
56d2567f34SHong Zhang #define TSIRK             "irk"
57d27334e2SStefano Zampini #define TSDIRK            "dirk"
5840be0ff1SMatthew G. Knepley 
59435da068SBarry Smith /*E
6087497f52SBarry Smith    TSProblemType - Determines the type of problem this `TS` object is to be used to solve
61435da068SBarry Smith 
62195e9b02SBarry Smith    Values:
63195e9b02SBarry Smith  + `TS_LINEAR`    - a linear ODE or DAE
64195e9b02SBarry Smith  - `TS_NONLINEAR` - a nonlinear ODE or DAE
65195e9b02SBarry Smith 
66435da068SBarry Smith    Level: beginner
67435da068SBarry Smith 
681cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSCreate()`
69435da068SBarry Smith E*/
709371c9d4SSatish Balay typedef enum {
719371c9d4SSatish Balay   TS_LINEAR,
729371c9d4SSatish Balay   TS_NONLINEAR
739371c9d4SSatish Balay } TSProblemType;
74818ad0c1SBarry Smith 
75b93fea0eSJed Brown /*E
7687497f52SBarry Smith    TSEquationType - type of `TS` problem that is solved
77e817cc15SEmil Constantinescu 
78195e9b02SBarry Smith    Values:
79195e9b02SBarry Smith +  `TS_EQ_UNSPECIFIED` - (default)
8016a05f60SBarry Smith .  `TS_EQ_EXPLICIT`    - {ODE and DAE index 1, 2, 3, HI} F(t,U,U_t) := M(t) U_t - G(U,t) = 0
8116a05f60SBarry Smith -  `TS_EQ_IMPLICIT`    - {ODE and DAE index 1, 2, 3, HI} F(t,U,U_t) = 0
82e817cc15SEmil Constantinescu 
8395bd0b28SBarry Smith    Level: beginner
8495bd0b28SBarry Smith 
851cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGetEquationType()`, `TSSetEquationType()`
86e817cc15SEmil Constantinescu E*/
87e817cc15SEmil Constantinescu typedef enum {
88e817cc15SEmil Constantinescu   TS_EQ_UNSPECIFIED               = -1,
89e817cc15SEmil Constantinescu   TS_EQ_EXPLICIT                  = 0,
90e817cc15SEmil Constantinescu   TS_EQ_ODE_EXPLICIT              = 1,
91e817cc15SEmil Constantinescu   TS_EQ_DAE_SEMI_EXPLICIT_INDEX1  = 100,
92e817cc15SEmil Constantinescu   TS_EQ_DAE_SEMI_EXPLICIT_INDEX2  = 200,
93e817cc15SEmil Constantinescu   TS_EQ_DAE_SEMI_EXPLICIT_INDEX3  = 300,
94e817cc15SEmil Constantinescu   TS_EQ_DAE_SEMI_EXPLICIT_INDEXHI = 500,
95e817cc15SEmil Constantinescu   TS_EQ_IMPLICIT                  = 1000,
96e817cc15SEmil Constantinescu   TS_EQ_ODE_IMPLICIT              = 1001,
97e817cc15SEmil Constantinescu   TS_EQ_DAE_IMPLICIT_INDEX1       = 1100,
98e817cc15SEmil Constantinescu   TS_EQ_DAE_IMPLICIT_INDEX2       = 1200,
99e817cc15SEmil Constantinescu   TS_EQ_DAE_IMPLICIT_INDEX3       = 1300,
10019436ca2SJed Brown   TS_EQ_DAE_IMPLICIT_INDEXHI      = 1500
101e817cc15SEmil Constantinescu } TSEquationType;
102e817cc15SEmil Constantinescu PETSC_EXTERN const char *const *TSEquationTypes;
103e817cc15SEmil Constantinescu 
104e817cc15SEmil Constantinescu /*E
10516a05f60SBarry Smith    TSConvergedReason - reason a `TS` method has converged (integrated to the requested time) or not
106b93fea0eSJed Brown 
107195e9b02SBarry Smith    Values:
108195e9b02SBarry Smith +  `TS_CONVERGED_ITERATING`          - this only occurs if `TSGetConvergedReason()` is called during the `TSSolve()`
109195e9b02SBarry Smith .  `TS_CONVERGED_TIME`               - the final time was reached
110195e9b02SBarry Smith .  `TS_CONVERGED_ITS`                - the maximum number of iterations (time-steps) was reached prior to the final time
111195e9b02SBarry Smith .  `TS_CONVERGED_USER`               - user requested termination
112195e9b02SBarry Smith .  `TS_CONVERGED_EVENT`              - user requested termination on event detection
113195e9b02SBarry Smith .  `TS_CONVERGED_PSEUDO_FATOL`       - stops when function norm decreased by a set amount, used only for `TSPSEUDO`
114195e9b02SBarry Smith .  `TS_CONVERGED_PSEUDO_FRTOL`       - stops when function norm decreases below a set amount, used only for `TSPSEUDO`
115195e9b02SBarry Smith .  `TS_DIVERGED_NONLINEAR_SOLVE`     - too many nonlinear solve failures have occurred
116195e9b02SBarry Smith .  `TS_DIVERGED_STEP_REJECTED`       - too many steps were rejected
117195e9b02SBarry Smith .  `TSFORWARD_DIVERGED_LINEAR_SOLVE` - tangent linear solve failed
118195e9b02SBarry Smith -  `TSADJOINT_DIVERGED_LINEAR_SOLVE` - transposed linear solve failed
119195e9b02SBarry Smith 
120b93fea0eSJed Brown    Level: beginner
121b93fea0eSJed Brown 
1221cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGetConvergedReason()`
123b93fea0eSJed Brown E*/
124193ac0bcSJed Brown typedef enum {
125193ac0bcSJed Brown   TS_CONVERGED_ITERATING          = 0,
126193ac0bcSJed Brown   TS_CONVERGED_TIME               = 1,
127193ac0bcSJed Brown   TS_CONVERGED_ITS                = 2,
128d6ad946cSShri Abhyankar   TS_CONVERGED_USER               = 3,
1292b7db910SShri Abhyankar   TS_CONVERGED_EVENT              = 4,
1303118ae5eSBarry Smith   TS_CONVERGED_PSEUDO_FATOL       = 5,
1313118ae5eSBarry Smith   TS_CONVERGED_PSEUDO_FRTOL       = 6,
132193ac0bcSJed Brown   TS_DIVERGED_NONLINEAR_SOLVE     = -1,
133b5b4017aSHong Zhang   TS_DIVERGED_STEP_REJECTED       = -2,
13405c9950eSHong Zhang   TSFORWARD_DIVERGED_LINEAR_SOLVE = -3,
13505c9950eSHong Zhang   TSADJOINT_DIVERGED_LINEAR_SOLVE = -4
136193ac0bcSJed Brown } TSConvergedReason;
137014dd563SJed Brown PETSC_EXTERN const char *const *TSConvergedReasons;
1381957e957SBarry Smith 
139b93fea0eSJed Brown /*MC
14087497f52SBarry Smith    TS_CONVERGED_ITERATING - this only occurs if `TSGetConvergedReason()` is called during the `TSSolve()`
141b93fea0eSJed Brown 
142b93fea0eSJed Brown    Level: beginner
143b93fea0eSJed Brown 
1441cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`
145b93fea0eSJed Brown M*/
146b93fea0eSJed Brown 
147b93fea0eSJed Brown /*MC
148b93fea0eSJed Brown    TS_CONVERGED_TIME - the final time was reached
149b93fea0eSJed Brown 
150b93fea0eSJed Brown    Level: beginner
151b93fea0eSJed Brown 
1521cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxTime()`, `TSGetMaxTime()`, `TSGetSolveTime()`
153b93fea0eSJed Brown M*/
154b93fea0eSJed Brown 
155b93fea0eSJed Brown /*MC
1561957e957SBarry Smith    TS_CONVERGED_ITS - the maximum number of iterations (time-steps) was reached prior to the final time
157b93fea0eSJed Brown 
158b93fea0eSJed Brown    Level: beginner
159b93fea0eSJed Brown 
1601cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxSteps()`, `TSGetMaxSteps()`
161b93fea0eSJed Brown M*/
1621957e957SBarry Smith 
163d6ad946cSShri Abhyankar /*MC
164d6ad946cSShri Abhyankar    TS_CONVERGED_USER - user requested termination
165d6ad946cSShri Abhyankar 
166d6ad946cSShri Abhyankar    Level: beginner
167d6ad946cSShri Abhyankar 
1681cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`
169b93fea0eSJed Brown M*/
170b93fea0eSJed Brown 
171b93fea0eSJed Brown /*MC
172d76fc68cSShri Abhyankar    TS_CONVERGED_EVENT - user requested termination on event detection
173d76fc68cSShri Abhyankar 
174d76fc68cSShri Abhyankar    Level: beginner
175d76fc68cSShri Abhyankar 
1761cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`
177d76fc68cSShri Abhyankar M*/
178d76fc68cSShri Abhyankar 
179d76fc68cSShri Abhyankar /*MC
18087497f52SBarry Smith    TS_CONVERGED_PSEUDO_FRTOL - stops when function norm decreased by a set amount, used only for `TSPSEUDO`
1813118ae5eSBarry Smith 
1823c7db156SBarry Smith    Options Database Key:
183d5b43468SJose E. Roman .   -ts_pseudo_frtol <rtol> - use specified rtol
1843118ae5eSBarry Smith 
185195e9b02SBarry Smith    Level: beginner
186195e9b02SBarry Smith 
1871cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`, `TS_CONVERGED_PSEUDO_FATOL`
1883118ae5eSBarry Smith M*/
1893118ae5eSBarry Smith 
1903118ae5eSBarry Smith /*MC
19187497f52SBarry Smith    TS_CONVERGED_PSEUDO_FATOL - stops when function norm decreases below a set amount, used only for `TSPSEUDO`
1923118ae5eSBarry Smith 
1933c7db156SBarry Smith    Options Database Key:
19467b8a455SSatish Balay .   -ts_pseudo_fatol <atol> - use specified atol
1953118ae5eSBarry Smith 
196195e9b02SBarry Smith    Level: beginner
197195e9b02SBarry Smith 
1981cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`, `TS_CONVERGED_PSEUDO_FRTOL`
1993118ae5eSBarry Smith M*/
2003118ae5eSBarry Smith 
2013118ae5eSBarry Smith /*MC
202b93fea0eSJed Brown    TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed
203b93fea0eSJed Brown 
204b93fea0eSJed Brown    Level: beginner
205b93fea0eSJed Brown 
206195e9b02SBarry Smith    Note:
207195e9b02SBarry Smith    See `TSSetMaxSNESFailures()` for how to allow more nonlinear solver failures.
2081957e957SBarry Smith 
2091cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSGetSNES()`, `SNESGetConvergedReason()`, `TSSetMaxSNESFailures()`
210b93fea0eSJed Brown M*/
211b93fea0eSJed Brown 
212b93fea0eSJed Brown /*MC
213b93fea0eSJed Brown    TS_DIVERGED_STEP_REJECTED - too many steps were rejected
214b93fea0eSJed Brown 
215b93fea0eSJed Brown    Level: beginner
216b93fea0eSJed Brown 
21795452b02SPatrick Sanan    Notes:
218195e9b02SBarry Smith    See `TSSetMaxStepRejections()` for how to allow more step rejections.
2191957e957SBarry Smith 
2201cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxStepRejections()`
221b93fea0eSJed Brown M*/
222b93fea0eSJed Brown 
223d6ad946cSShri Abhyankar /*E
224d6ad946cSShri Abhyankar    TSExactFinalTimeOption - option for handling of final time step
225d6ad946cSShri Abhyankar 
226195e9b02SBarry Smith    Values:
22716a05f60SBarry Smith +  `TS_EXACTFINALTIME_STEPOVER`    - Don't do anything if requested final time is exceeded
228195e9b02SBarry Smith .  `TS_EXACTFINALTIME_INTERPOLATE` - Interpolate back to final time
22916a05f60SBarry Smith -  `TS_EXACTFINALTIME_MATCHSTEP`   - Adapt final time step to match the final time requested
230195e9b02SBarry Smith 
231d6ad946cSShri Abhyankar    Level: beginner
232d6ad946cSShri Abhyankar 
2331cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGetConvergedReason()`, `TSSetExactFinalTime()`, `TSGetExactFinalTime()`
234d6ad946cSShri Abhyankar E*/
2359371c9d4SSatish Balay typedef enum {
2369371c9d4SSatish Balay   TS_EXACTFINALTIME_UNSPECIFIED = 0,
2379371c9d4SSatish Balay   TS_EXACTFINALTIME_STEPOVER    = 1,
2389371c9d4SSatish Balay   TS_EXACTFINALTIME_INTERPOLATE = 2,
2399371c9d4SSatish Balay   TS_EXACTFINALTIME_MATCHSTEP   = 3
2409371c9d4SSatish Balay } TSExactFinalTimeOption;
241d6ad946cSShri Abhyankar PETSC_EXTERN const char *const TSExactFinalTimeOptions[];
242d6ad946cSShri Abhyankar 
243000e7ae3SMatthew Knepley /* Logging support */
244014dd563SJed Brown PETSC_EXTERN PetscClassId TS_CLASSID;
245d74926cbSBarry Smith PETSC_EXTERN PetscClassId DMTS_CLASSID;
2461b9b13dfSLisandro Dalcin PETSC_EXTERN PetscClassId TSADAPT_CLASSID;
247000e7ae3SMatthew Knepley 
248607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSInitializePackage(void);
249560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSFinalizePackage(void);
2508ba1e511SMatthew Knepley 
251014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm, TS *);
252baa10174SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSClone(TS, TS *);
253014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSDestroy(TS *);
254818ad0c1SBarry Smith 
255014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetProblemType(TS, TSProblemType);
256014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetProblemType(TS, TSProblemType *);
257014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSMonitor(TS, PetscInt, PetscReal, Vec);
258*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorSet(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscCtx), PetscCtx, PetscCtxDestroyFn *);
259721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorSetFromOptions(TS, const char[], const char[], const char[], PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *), PetscErrorCode (*)(TS, PetscViewerAndFormat *));
260014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS);
261818ad0c1SBarry Smith 
262014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS, const char[]);
263014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS, const char[]);
264014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS, const char *[]);
265014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS);
266014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetUp(TS);
267014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSReset(TS);
268818ad0c1SBarry Smith 
269014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetSolution(TS, Vec);
270014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetSolution(TS, Vec *);
271818ad0c1SBarry Smith 
272efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TS2SetSolution(TS, Vec, Vec);
273efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TS2GetSolution(TS, Vec *, Vec *);
27457df6a1bSDebojyoti Ghosh 
275642f3702SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetSolutionComponents(TS, PetscInt *, Vec *);
276642f3702SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetAuxSolution(TS, Vec *);
2770a01e1b2SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetTimeError(TS, PetscInt, Vec *);
27857df6a1bSDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSSetTimeError(TS, Vec);
2796e2e232bSDebojyoti Ghosh 
280*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSSetRHSJacobianP(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, PetscCtx), PetscCtx);
281*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSGetRHSJacobianP(TS, Mat *, PetscErrorCode (**)(TS, PetscReal, Vec, Mat, PetscCtx), PetscCtxRt);
282a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianP(TS, PetscReal, Vec, Mat);
283*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSSetIJacobianP(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Vec, PetscReal, Mat, PetscCtx), PetscCtx);
284*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSGetIJacobianP(TS, Mat *, PetscErrorCode (**)(TS, PetscReal, Vec, Vec, PetscReal, Mat, PetscCtx), PetscCtxRt);
28533f72c5dSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIJacobianP(TS, PetscReal, Vec, Vec, PetscReal, Mat, PetscBool);
286edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS() then TSComputeRHSJacobianP()", ) PetscErrorCode TSComputeDRDPFunction(TS, PetscReal, Vec, Vec *);
287edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS() then TSComputeRHSJacobian()", ) PetscErrorCode TSComputeDRDUFunction(TS, PetscReal, Vec, Vec *);
288*2a8381b2SBarry Smith 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);
2890fc7ecceSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUU(TS, PetscReal, Vec, Vec[], Vec, Vec[]);
2900fc7ecceSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUP(TS, PetscReal, Vec, Vec[], Vec, Vec[]);
2910fc7ecceSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPU(TS, PetscReal, Vec, Vec[], Vec, Vec[]);
2920fc7ecceSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPP(TS, PetscReal, Vec, Vec[], Vec, Vec[]);
293*2a8381b2SBarry Smith 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);
2940fc7ecceSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS, PetscReal, Vec, Vec[], Vec, Vec[]);
2950fc7ecceSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS, PetscReal, Vec, Vec[], Vec, Vec[]);
2960fc7ecceSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS, PetscReal, Vec, Vec[], Vec, Vec[]);
2970fc7ecceSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS, PetscReal, Vec, Vec[], Vec, Vec[]);
2980fc7ecceSBarry Smith PETSC_EXTERN PetscErrorCode TSSetCostHessianProducts(TS, PetscInt, Vec[], Vec[], Vec);
2990fc7ecceSBarry Smith PETSC_EXTERN PetscErrorCode TSGetCostHessianProducts(TS, PetscInt *, Vec *[], Vec *[], Vec *);
300ba0675f6SHong Zhang PETSC_EXTERN PetscErrorCode TSComputeSNESJacobian(TS, Vec, Mat, Mat);
301a05bf03eSHong Zhang 
302bc952696SBarry Smith /*S
303e91eccc2SStefano Zampini    TSTrajectory - Abstract PETSc object that stores the trajectory (solution of ODE/DAE at each time step)
304bc952696SBarry Smith 
305bc952696SBarry Smith    Level: advanced
306bc952696SBarry Smith 
3071cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectorySetType()`, `TSTrajectoryDestroy()`, `TSTrajectoryReset()`
308bc952696SBarry Smith S*/
309bc952696SBarry Smith typedef struct _p_TSTrajectory *TSTrajectory;
310bc952696SBarry Smith 
311bc952696SBarry Smith /*J
312195e9b02SBarry Smith    TSTrajectoryType - String with the name of a PETSc `TS` trajectory storage method
313bc952696SBarry Smith 
314bc952696SBarry Smith    Level: intermediate
315bc952696SBarry Smith 
3161cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectoryDestroy()`
317bc952696SBarry Smith J*/
318bc952696SBarry Smith typedef const char *TSTrajectoryType;
319bc952696SBarry Smith #define TSTRAJECTORYBASIC         "basic"
3201c8c567eSBarry Smith #define TSTRAJECTORYSINGLEFILE    "singlefile"
3219a53571cSHong Zhang #define TSTRAJECTORYMEMORY        "memory"
3222b043167SHong Zhang #define TSTRAJECTORYVISUALIZATION "visualization"
323bc952696SBarry Smith 
324bc952696SBarry Smith PETSC_EXTERN PetscFunctionList TSTrajectoryList;
325bc952696SBarry Smith PETSC_EXTERN PetscClassId      TSTRAJECTORY_CLASSID;
326bc952696SBarry Smith PETSC_EXTERN PetscBool         TSTrajectoryRegisterAllCalled;
327bc952696SBarry Smith 
328bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS);
3292d29f1f2SStefano Zampini PETSC_EXTERN PetscErrorCode TSResetTrajectory(TS);
33067a3cfb0SHong Zhang PETSC_EXTERN PetscErrorCode TSRemoveTrajectory(TS);
331bc952696SBarry Smith 
332bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectoryCreate(MPI_Comm, TSTrajectory *);
3339a992471SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryReset(TSTrajectory);
334bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy(TSTrajectory *);
335560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSTrajectoryView(TSTrajectory, PetscViewer);
336fd9d3c67SJed Brown PETSC_EXTERN PetscErrorCode TSTrajectorySetType(TSTrajectory, TS, TSTrajectoryType);
337881c1a9bSHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryGetType(TSTrajectory, TS, TSTrajectoryType *);
338bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySet(TSTrajectory, TS, PetscInt, PetscReal, Vec);
339c679fc15SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryGet(TSTrajectory, TS, PetscInt, PetscReal *);
340fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetVecs(TSTrajectory, TS, PetscInt, PetscReal *, Vec, Vec);
341fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetUpdatedHistoryVecs(TSTrajectory, TS, PetscReal, Vec *, Vec *);
342fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetNumSteps(TSTrajectory, PetscInt *);
343fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryRestoreUpdatedHistoryVecs(TSTrajectory, Vec *, Vec *);
344972caf09SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory, TS);
345560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSTrajectoryRegister(const char[], PetscErrorCode (*)(TSTrajectory, TS));
346bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectoryRegisterAll(void);
34768bece0bSHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetUp(TSTrajectory, TS);
3489ffb3502SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetUseHistory(TSTrajectory, PetscBool);
3492bee684fSHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetMonitor(TSTrajectory, PetscBool);
35078fbdcc8SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySetVariableNames(TSTrajectory, const char *const *);
351*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySetTransform(TSTrajectory, PetscErrorCode (*)(PetscCtx, Vec, Vec *), PetscCtxDestroyFn *, PetscCtx);
352fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectorySetSolutionOnly(TSTrajectory, PetscBool);
353fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetSolutionOnly(TSTrajectory, PetscBool *);
35464fc91eeSBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySetKeepFiles(TSTrajectory, PetscBool);
35564e38db7SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetDirname(TSTrajectory, const char[]);
35664e38db7SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetFiletemplate(TSTrajectory, const char[]);
35778fbdcc8SBarry Smith PETSC_EXTERN PetscErrorCode TSGetTrajectory(TS, TSTrajectory *);
358bc952696SBarry Smith 
3594bf303faSJacob Faibussowitsch typedef enum {
3604bf303faSJacob Faibussowitsch   TJ_REVOLVE,
3614bf303faSJacob Faibussowitsch   TJ_CAMS,
3624bf303faSJacob Faibussowitsch   TJ_PETSC
3634bf303faSJacob Faibussowitsch } TSTrajectoryMemoryType;
3644bf303faSJacob Faibussowitsch PETSC_EXTERN const char *const TSTrajectoryMemoryTypes[];
3654bf303faSJacob Faibussowitsch 
3664bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectoryMemorySetType(TSTrajectory, TSTrajectoryMemoryType);
3674bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxCpsRAM(TSTrajectory, PetscInt);
3684bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxCpsDisk(TSTrajectory, PetscInt);
3694bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxUnitsRAM(TSTrajectory, PetscInt);
3704bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxUnitsDisk(TSTrajectory, PetscInt);
3714bf303faSJacob Faibussowitsch 
3720fc7ecceSBarry Smith PETSC_EXTERN PetscErrorCode TSSetCostGradients(TS, PetscInt, Vec[], Vec[]);
3730fc7ecceSBarry Smith PETSC_EXTERN PetscErrorCode TSGetCostGradients(TS, PetscInt *, Vec *[], Vec *[]);
374*2a8381b2SBarry Smith 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);
375dfb21088SHong Zhang PETSC_EXTERN PetscErrorCode TSGetCostIntegral(TS, Vec *);
376022c081aSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeCostIntegrand(TS, PetscReal, Vec, Vec);
377cd4cee2dSHong Zhang PETSC_EXTERN PetscErrorCode TSCreateQuadratureTS(TS, PetscBool, TS *);
378cd4cee2dSHong Zhang PETSC_EXTERN PetscErrorCode TSGetQuadratureTS(TS, PetscBool *, TS *);
379d4aa0a58SBarry Smith 
380ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSetFromOptions(TS, PetscOptionItems);
3810fc7ecceSBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointMonitor(TS, PetscInt, PetscReal, Vec, PetscInt, Vec[], Vec[]);
382*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointMonitorSet(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, PetscCtx), PetscCtx, PetscCtxDestroyFn *);
383a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitorCancel(TS);
384a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitorSetFromOptions(TS, const char[], const char[], const char[], PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, PetscViewerAndFormat *), PetscErrorCode (*)(TS, PetscViewerAndFormat *));
385a05bf03eSHong Zhang 
386*2a8381b2SBarry Smith PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSSetRHSJacobianP()", ) PetscErrorCode TSAdjointSetRHSJacobian(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, PetscCtx), PetscCtx);
387edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSComputeRHSJacobianP()", ) PetscErrorCode TSAdjointComputeRHSJacobian(TS, PetscReal, Vec, Mat);
388edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS()", ) PetscErrorCode TSAdjointComputeDRDPFunction(TS, PetscReal, Vec, Vec *);
389edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS()", ) PetscErrorCode TSAdjointComputeDRDYFunction(TS, PetscReal, Vec, Vec *);
3902c39e106SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS);
39173b18844SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSetSteps(TS, PetscInt);
392d4aa0a58SBarry Smith 
39373b18844SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointStep(TS);
394f6a906c0SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS);
395ece46509SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointReset(TS);
396999a3721SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointCostIntegral(TS);
397ecf68647SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointSetForward(TS, Mat);
398ecf68647SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointResetForward(TS);
399715f1b00SHong Zhang 
40013b1b0ffSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardSetSensitivities(TS, PetscInt, Mat);
40113b1b0ffSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardGetSensitivities(TS, PetscInt *, Mat *);
4020fc7ecceSBarry Smith PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 12, 0, "TSCreateQuadratureTS()", ) PetscErrorCode TSForwardSetIntegralGradients(TS, PetscInt, Vec[]);
4030fc7ecceSBarry Smith PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 12, 0, "TSForwardGetSensitivities()", ) PetscErrorCode TSForwardGetIntegralGradients(TS, PetscInt *, Vec *[]);
404715f1b00SHong Zhang PETSC_EXTERN PetscErrorCode TSForwardSetUp(TS);
4057adebddeSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardReset(TS);
406999a3721SHong Zhang PETSC_EXTERN PetscErrorCode TSForwardCostIntegral(TS);
407715f1b00SHong Zhang PETSC_EXTERN PetscErrorCode TSForwardStep(TS);
408b5b4017aSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardSetInitialSensitivities(TS, Mat);
409cc4f23bcSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardGetStages(TS, PetscInt *, Mat *[]);
410c235aa8dSHong Zhang 
411618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetMaxSteps(TS, PetscInt);
412618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetMaxSteps(TS, PetscInt *);
4138e562f8dSJames Wright PETSC_EXTERN PetscErrorCode TSSetRunSteps(TS, PetscInt);
4148e562f8dSJames Wright PETSC_EXTERN PetscErrorCode TSGetRunSteps(TS, PetscInt *);
415618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetMaxTime(TS, PetscReal);
416618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetMaxTime(TS, PetscReal *);
41749354f04SShri Abhyankar PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS, TSExactFinalTimeOption);
418f6953c82SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetExactFinalTime(TS, TSExactFinalTimeOption *);
4190fc7ecceSBarry Smith PETSC_EXTERN PetscErrorCode TSSetEvaluationTimes(TS, PetscInt, PetscReal[]);
4200fc7ecceSBarry Smith PETSC_EXTERN PetscErrorCode TSGetEvaluationTimes(TS, PetscInt *, const PetscReal *[]);
4210fc7ecceSBarry Smith PETSC_EXTERN PetscErrorCode TSGetEvaluationSolutions(TS, PetscInt *, const PetscReal *[], Vec *[]);
4220fc7ecceSBarry Smith PETSC_EXTERN PetscErrorCode TSSetTimeSpan(TS, PetscInt, PetscReal[]);
423c80d56d9SJames Wright 
424c80d56d9SJames Wright /*@C
425c80d56d9SJames Wright   TSGetTimeSpan - gets the time span set with `TSSetTimeSpan()`
426c80d56d9SJames Wright 
427c80d56d9SJames Wright   Not Collective
428c80d56d9SJames Wright 
429c80d56d9SJames Wright   Input Parameter:
430c80d56d9SJames Wright . ts - the time-stepper
431c80d56d9SJames Wright 
432c80d56d9SJames Wright   Output Parameters:
433c80d56d9SJames Wright + n          - number of the time points (>=2)
434c80d56d9SJames Wright - span_times - array of the time points. The first element and the last element are the initial time and the final time respectively.
435c80d56d9SJames Wright 
436c80d56d9SJames Wright   Level: deprecated
437c80d56d9SJames Wright 
438c80d56d9SJames Wright   Note:
439c80d56d9SJames Wright   Deprecated, use `TSGetEvaluationTimes()`.
440c80d56d9SJames Wright 
441c80d56d9SJames Wright   The values obtained are valid until the `TS` object is destroyed.
442c80d56d9SJames Wright 
443c80d56d9SJames Wright   Both `n` and `span_times` can be `NULL`.
444c80d56d9SJames Wright 
445c80d56d9SJames Wright .seealso: [](ch_ts), `TS`, `TSGetEvaluationTimes()`, `TSSetTimeSpan()`, `TSSetEvaluationTimes()`, `TSGetEvaluationSolutions()`
446c80d56d9SJames Wright  @*/
TSGetTimeSpan(TS ts,PetscInt * n,const PetscReal * span_times[])447c80d56d9SJames Wright PETSC_DEPRECATED_FUNCTION(3, 23, 0, "TSGetEvaluationTimes()", ) static inline PetscErrorCode TSGetTimeSpan(TS ts, PetscInt *n, const PetscReal *span_times[])
448c80d56d9SJames Wright {
449c80d56d9SJames Wright   return TSGetEvaluationTimes(ts, n, span_times);
450c80d56d9SJames Wright }
451c80d56d9SJames Wright 
452c80d56d9SJames Wright /*@C
453c80d56d9SJames Wright   TSGetTimeSpanSolutions - Get the number of solutions and the solutions at the time points specified by the time span.
454c80d56d9SJames Wright 
455c80d56d9SJames Wright   Input Parameter:
456c80d56d9SJames Wright . ts - the `TS` context obtained from `TSCreate()`
457c80d56d9SJames Wright 
458c80d56d9SJames Wright   Output Parameters:
459c80d56d9SJames Wright + nsol - the number of solutions
460c80d56d9SJames Wright - Sols - the solution vectors
461c80d56d9SJames Wright 
462c80d56d9SJames Wright   Level: deprecated
463c80d56d9SJames Wright 
464c80d56d9SJames Wright   Notes:
465c80d56d9SJames Wright   Deprecated, use `TSGetEvaluationSolutions()`.
466c80d56d9SJames Wright 
467c80d56d9SJames Wright   Both `nsol` and `Sols` can be `NULL`.
468c80d56d9SJames Wright 
469c80d56d9SJames Wright   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()`.
470c80d56d9SJames Wright   For example, manipulating the step size, especially with a reduced precision, may cause `TS` to step over certain points in the span.
471c80d56d9SJames Wright   This issue is alleviated in `TSGetEvaluationSolutions()` by returning the solution times that `Sols` were recorded at.
472c80d56d9SJames Wright 
473c80d56d9SJames Wright .seealso: [](ch_ts), `TS`, `TSGetEvaluationSolutions()`, `TSSetTimeSpan()`, `TSGetEvaluationTimes()`, `TSSetEvaluationTimes()`
474c80d56d9SJames Wright  @*/
TSGetTimeSpanSolutions(TS ts,PetscInt * nsol,Vec ** Sols)475c80d56d9SJames Wright PETSC_DEPRECATED_FUNCTION(3, 23, 0, "TSGetEvaluationSolutions()", ) static inline PetscErrorCode TSGetTimeSpanSolutions(TS ts, PetscInt *nsol, Vec **Sols)
476c80d56d9SJames Wright {
477c80d56d9SJames Wright   return TSGetEvaluationSolutions(ts, nsol, NULL, Sols);
478c80d56d9SJames Wright }
479818ad0c1SBarry Smith 
480edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSSetTime()", ) PetscErrorCode TSSetInitialTimeStep(TS, PetscReal, PetscReal);
481edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSSetMax()", ) PetscErrorCode TSSetDuration(TS, PetscInt, PetscReal);
482edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetMax()", ) PetscErrorCode TSGetDuration(TS, PetscInt *, PetscReal *);
483edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetStepNumber()", ) PetscErrorCode TSGetTimeStepNumber(TS, PetscInt *);
484edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetStepNumber()", ) PetscErrorCode TSGetTotalSteps(TS, PetscInt *);
48519eac22cSLisandro Dalcin 
486721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
487340b794cSJed Brown PETSC_EXTERN PetscErrorCode TSMonitorWallClockTime(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
488340b794cSJed Brown PETSC_EXTERN PetscErrorCode TSMonitorWallClockTimeSetUp(TS, PetscViewerAndFormat *);
489cc9c3a59SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorExtreme(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
49083a4ac43SBarry Smith 
49183a4ac43SBarry Smith typedef struct _n_TSMonitorDrawCtx *TSMonitorDrawCtx;
49283a4ac43SBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorDrawCtx *);
49383a4ac43SBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *);
494*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawSolution(TS, PetscInt, PetscReal, Vec, PetscCtx);
495*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawSolutionPhase(TS, PetscInt, PetscReal, Vec, PetscCtx);
496*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawError(TS, PetscInt, PetscReal, Vec, PetscCtx);
497*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawSolutionFunction(TS, PetscInt, PetscReal, Vec, PetscCtx);
49883a4ac43SBarry Smith 
4990fc7ecceSBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointMonitorDefault(TS, PetscInt, PetscReal, Vec, PetscInt, Vec[], Vec[], PetscViewerAndFormat *);
500*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointMonitorDrawSensi(TS, PetscInt, PetscReal, Vec, PetscInt, Vec[], Vec[], PetscCtx);
5019110b2e7SHong Zhang 
5028e562f8dSJames Wright typedef struct _n_TSMonitorSolutionCtx *TSMonitorSolutionCtx;
503721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode             TSMonitorSolution(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
5048e562f8dSJames Wright PETSC_EXTERN PetscErrorCode             TSMonitorSolutionSetup(TS, PetscViewerAndFormat *);
5057f27e910SStefano Zampini 
5067f27e910SStefano Zampini typedef struct _n_TSMonitorVTKCtx *TSMonitorVTKCtx;
5077f27e910SStefano Zampini PETSC_EXTERN PetscErrorCode        TSMonitorSolutionVTK(TS, PetscInt, PetscReal, Vec, TSMonitorVTKCtx);
5087f27e910SStefano Zampini PETSC_EXTERN PetscErrorCode        TSMonitorSolutionVTKDestroy(TSMonitorVTKCtx *);
5097f27e910SStefano Zampini PETSC_EXTERN PetscErrorCode        TSMonitorSolutionVTKCtxCreate(const char *, TSMonitorVTKCtx *);
510fb1732b5SBarry Smith 
511014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSStep(TS);
5127cbde773SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSEvaluateWLTE(TS, NormType, PetscInt *, PetscReal *);
513014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS, PetscInt, Vec, PetscBool *);
514cc708dedSBarry Smith PETSC_EXTERN PetscErrorCode TSSolve(TS, Vec);
515e817cc15SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetEquationType(TS, TSEquationType *);
516e817cc15SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSSetEquationType(TS, TSEquationType);
517014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS, TSConvergedReason *);
518d6ad946cSShri Abhyankar PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS, TSConvergedReason);
519487e0bb9SJed Brown PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS, PetscReal *);
5205ef26d82SJed Brown PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS, PetscInt *);
5215ef26d82SJed Brown PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS, PetscInt *);
522cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS, PetscInt *);
523cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS, PetscInt);
524cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS, PetscInt *);
525cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS, PetscInt);
526cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS, PetscBool);
527dcb233daSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRestartStep(TS);
52824655328SShri PETSC_EXTERN PetscErrorCode TSRollBack(TS);
5299dc50cb5SStefano Zampini PETSC_EXTERN PetscErrorCode TSGetStepRollBack(TS, PetscBool *);
530d2daff3dSHong Zhang 
531cc4f23bcSHong Zhang PETSC_EXTERN PetscErrorCode TSGetStages(TS, PetscInt *, Vec *[]);
532818ad0c1SBarry Smith 
533014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetTime(TS, PetscReal *);
534014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetTime(TS, PetscReal);
535e5e524a1SHong Zhang PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS, PetscReal *);
53680275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS, PetscReal *);
53780275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS, PetscReal);
53880275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetStepNumber(TS, PetscInt *);
53980275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetStepNumber(TS, PetscInt);
540818ad0c1SBarry Smith 
54114d0ab18SJacob Faibussowitsch /*S
5428434afd1SBarry Smith   TSRHSFunctionFn - A prototype of a `TS` right-hand-side evaluation function that would be passed to `TSSetRHSFunction()`
54314d0ab18SJacob Faibussowitsch 
54414d0ab18SJacob Faibussowitsch   Calling Sequence:
54514d0ab18SJacob Faibussowitsch + ts  - timestep context
54614d0ab18SJacob Faibussowitsch . t   - current time
54714d0ab18SJacob Faibussowitsch . u   - input vector
54814d0ab18SJacob Faibussowitsch . F   - function vector
54914d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context
55014d0ab18SJacob Faibussowitsch 
55114d0ab18SJacob Faibussowitsch   Level: beginner
55214d0ab18SJacob Faibussowitsch 
553d1c5d1fcSBarry Smith   Note:
5548434afd1SBarry Smith   The deprecated `TSRHSFunction` still works as a replacement for `TSRHSFunctionFn` *.
555d1c5d1fcSBarry Smith 
5568434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetRHSFunction()`, `DMTSSetRHSFunction()`, `TSIFunctionFn`,
5578434afd1SBarry Smith `TSIJacobianFn`, `TSRHSJacobianFn`
55814d0ab18SJacob Faibussowitsch S*/
559*2a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSRHSFunctionFn(TS ts, PetscReal t, Vec u, Vec F, PetscCtx ctx);
560d1c5d1fcSBarry Smith 
5618434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSRHSFunctionFn *TSRHSFunction;
56214d0ab18SJacob Faibussowitsch 
56314d0ab18SJacob Faibussowitsch /*S
5648434afd1SBarry Smith   TSRHSJacobianFn - A prototype of a `TS` right-hand-side Jacobian evaluation function that would be passed to `TSSetRHSJacobian()`
56514d0ab18SJacob Faibussowitsch 
56614d0ab18SJacob Faibussowitsch   Calling Sequence:
56714d0ab18SJacob Faibussowitsch + ts   - the `TS` context obtained from `TSCreate()`
56814d0ab18SJacob Faibussowitsch . t    - current time
56914d0ab18SJacob Faibussowitsch . u    - input vector
57014d0ab18SJacob Faibussowitsch . Amat - (approximate) Jacobian matrix
571af27ebaaSBarry Smith . Pmat - matrix from which preconditioner is to be constructed (usually the same as `Amat`)
57214d0ab18SJacob Faibussowitsch - ctx  - [optional] user-defined context for matrix evaluation routine
57314d0ab18SJacob Faibussowitsch 
57414d0ab18SJacob Faibussowitsch   Level: beginner
57514d0ab18SJacob Faibussowitsch 
576d1c5d1fcSBarry Smith   Note:
5778434afd1SBarry Smith   The deprecated `TSRHSJacobian` still works as a replacement for `TSRHSJacobianFn` *.
578d1c5d1fcSBarry Smith 
5798434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetRHSJacobian()`, `DMTSSetRHSJacobian()`, `TSRHSFunctionFn`,
5808434afd1SBarry Smith `TSIFunctionFn`, `TSIJacobianFn`
58114d0ab18SJacob Faibussowitsch S*/
582*2a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSRHSJacobianFn(TS ts, PetscReal t, Vec u, Mat Amat, Mat Pmat, PetscCtx ctx);
583d1c5d1fcSBarry Smith 
5848434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSRHSJacobianFn *TSRHSJacobian;
58514d0ab18SJacob Faibussowitsch 
58614d0ab18SJacob Faibussowitsch /*S
5878434afd1SBarry Smith   TSRHSJacobianPFn - A prototype of a function that computes the Jacobian of G w.r.t. the parameters P where
588af27ebaaSBarry Smith   U_t = G(U,P,t), as well as the location to store the matrix that would be passed to `TSSetRHSJacobianP()`
58914d0ab18SJacob Faibussowitsch 
59014d0ab18SJacob Faibussowitsch   Calling Sequence:
59114d0ab18SJacob Faibussowitsch + ts  - the `TS` context
59214d0ab18SJacob Faibussowitsch . t   - current timestep
59314d0ab18SJacob Faibussowitsch . U   - input vector (current ODE solution)
59414d0ab18SJacob Faibussowitsch . A   - output matrix
59514d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context
59614d0ab18SJacob Faibussowitsch 
59714d0ab18SJacob Faibussowitsch   Level: beginner
59814d0ab18SJacob Faibussowitsch 
599d1c5d1fcSBarry Smith   Note:
6008434afd1SBarry Smith   The deprecated `TSRHSJacobianP` still works as a replacement for `TSRHSJacobianPFn` *.
601d1c5d1fcSBarry Smith 
60214d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetRHSJacobianP()`, `TSGetRHSJacobianP()`
60314d0ab18SJacob Faibussowitsch S*/
604*2a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSRHSJacobianPFn(TS ts, PetscReal t, Vec U, Mat A, PetscCtx ctx);
60514d0ab18SJacob Faibussowitsch 
6068434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSRHSJacobianPFn *TSRHSJacobianP;
607d1c5d1fcSBarry Smith 
608*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS, Vec, TSRHSFunctionFn *, PetscCtx);
609*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS, Vec *, TSRHSFunctionFn **, PetscCtxRt);
610*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS, Mat, Mat, TSRHSJacobianFn *, PetscCtx);
611*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS, Mat *, Mat *, TSRHSJacobianFn **, PetscCtxRt);
612e1244c69SJed Brown PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS, PetscBool);
613818ad0c1SBarry Smith 
61414d0ab18SJacob Faibussowitsch /*S
6158434afd1SBarry Smith   TSSolutionFn - A prototype of a `TS` solution evaluation function that would be passed to `TSSetSolutionFunction()`
61614d0ab18SJacob Faibussowitsch 
61714d0ab18SJacob Faibussowitsch   Calling Sequence:
61814d0ab18SJacob Faibussowitsch + ts  - timestep context
61914d0ab18SJacob Faibussowitsch . t   - current time
62014d0ab18SJacob Faibussowitsch . u   - output vector
62114d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context
62214d0ab18SJacob Faibussowitsch 
62314d0ab18SJacob Faibussowitsch   Level: advanced
62414d0ab18SJacob Faibussowitsch 
625d1c5d1fcSBarry Smith   Note:
6268434afd1SBarry Smith   The deprecated `TSSolutionFunction` still works as a replacement for `TSSolutionFn` *.
627d1c5d1fcSBarry Smith 
62814d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetSolutionFunction()`, `DMTSSetSolutionFunction()`
62914d0ab18SJacob Faibussowitsch S*/
630*2a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSSolutionFn(TS ts, PetscReal t, Vec u, PetscCtx ctx);
63114d0ab18SJacob Faibussowitsch 
6328434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSSolutionFn *TSSolutionFunction;
633d1c5d1fcSBarry Smith 
634*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS, TSSolutionFn *, PetscCtx);
63514d0ab18SJacob Faibussowitsch 
63614d0ab18SJacob Faibussowitsch /*S
6378434afd1SBarry Smith   TSForcingFn - A prototype of a `TS` forcing function evaluation function that would be passed to `TSSetForcingFunction()`
63814d0ab18SJacob Faibussowitsch 
63914d0ab18SJacob Faibussowitsch   Calling Sequence:
64014d0ab18SJacob Faibussowitsch + ts  - timestep context
64114d0ab18SJacob Faibussowitsch . t   - current time
64214d0ab18SJacob Faibussowitsch . f   - output vector
64314d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context
64414d0ab18SJacob Faibussowitsch 
64514d0ab18SJacob Faibussowitsch   Level: advanced
64614d0ab18SJacob Faibussowitsch 
647d1c5d1fcSBarry Smith   Note:
6488434afd1SBarry Smith   The deprecated `TSForcingFunction` still works as a replacement for `TSForcingFn` *.
649d1c5d1fcSBarry Smith 
65014d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetForcingFunction()`, `DMTSSetForcingFunction()`
65114d0ab18SJacob Faibussowitsch S*/
652*2a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSForcingFn(TS ts, PetscReal t, Vec f, PetscCtx ctx);
65314d0ab18SJacob Faibussowitsch 
6548434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSForcingFn *TSForcingFunction;
655d1c5d1fcSBarry Smith 
656*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS, TSForcingFn *, PetscCtx);
657ef20d060SBarry Smith 
65814d0ab18SJacob Faibussowitsch /*S
6598434afd1SBarry Smith   TSIFunctionFn - A prototype of a `TS` implicit function evaluation function that would be passed to `TSSetIFunction()
66014d0ab18SJacob Faibussowitsch 
66114d0ab18SJacob Faibussowitsch   Calling Sequence:
66214d0ab18SJacob Faibussowitsch + ts  - the `TS` context obtained from `TSCreate()`
66314d0ab18SJacob Faibussowitsch . t   - time at step/stage being solved
66414d0ab18SJacob Faibussowitsch . U   - state vector
66514d0ab18SJacob Faibussowitsch . U_t - time derivative of state vector
66614d0ab18SJacob Faibussowitsch . F   - function vector
66714d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined context for function
66814d0ab18SJacob Faibussowitsch 
66914d0ab18SJacob Faibussowitsch   Level: beginner
67014d0ab18SJacob Faibussowitsch 
671d1c5d1fcSBarry Smith   Note:
6728434afd1SBarry Smith   The deprecated `TSIFunction` still works as a replacement for `TSIFunctionFn` *.
673d1c5d1fcSBarry Smith 
6748434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetIFunction()`, `DMTSSetIFunction()`, `TSIJacobianFn`, `TSRHSFunctionFn`, `TSRHSJacobianFn`
67514d0ab18SJacob Faibussowitsch S*/
676*2a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSIFunctionFn(TS ts, PetscReal t, Vec U, Vec U_t, Vec F, PetscCtx ctx);
677d1c5d1fcSBarry Smith 
6788434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSIFunctionFn *TSIFunction;
67914d0ab18SJacob Faibussowitsch 
68014d0ab18SJacob Faibussowitsch /*S
6818434afd1SBarry Smith   TSIJacobianFn - A prototype of a `TS` Jacobian evaluation function that would be passed to `TSSetIJacobian()`
68214d0ab18SJacob Faibussowitsch 
68314d0ab18SJacob Faibussowitsch   Calling Sequence:
68414d0ab18SJacob Faibussowitsch + ts   - the `TS` context obtained from `TSCreate()`
68514d0ab18SJacob Faibussowitsch . t    - time at step/stage being solved
68614d0ab18SJacob Faibussowitsch . U    - state vector
68714d0ab18SJacob Faibussowitsch . U_t  - time derivative of state vector
68814d0ab18SJacob Faibussowitsch . a    - shift
68914d0ab18SJacob Faibussowitsch . Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
69014d0ab18SJacob Faibussowitsch . Pmat - matrix used for constructing preconditioner, usually the same as `Amat`
69114d0ab18SJacob Faibussowitsch - ctx  - [optional] user-defined context for Jacobian evaluation routine
69214d0ab18SJacob Faibussowitsch 
69314d0ab18SJacob Faibussowitsch   Level: beginner
69414d0ab18SJacob Faibussowitsch 
695d1c5d1fcSBarry Smith   Note:
6968434afd1SBarry Smith   The deprecated `TSIJacobian` still works as a replacement for `TSIJacobianFn` *.
69714d0ab18SJacob Faibussowitsch 
6988434afd1SBarry Smith .seealso: [](ch_ts), `TSSetIJacobian()`, `DMTSSetIJacobian()`, `TSIFunctionFn`, `TSRHSFunctionFn`, `TSRHSJacobianFn`
699d1c5d1fcSBarry Smith S*/
700*2a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSIJacobianFn(TS ts, PetscReal t, Vec U, Vec U_t, PetscReal a, Mat Amat, Mat Pmat, PetscCtx ctx);
701d1c5d1fcSBarry Smith 
7028434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSIJacobianFn *TSIJacobian;
703d1c5d1fcSBarry Smith 
704*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSSetIFunction(TS, Vec, TSIFunctionFn *, PetscCtx);
705*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSGetIFunction(TS, Vec *, TSIFunctionFn **, PetscCtxRt);
706*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS, Mat, Mat, TSIJacobianFn *, PetscCtx);
707*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS, Mat *, Mat *, TSIJacobianFn **, PetscCtxRt);
708316643e7SJed Brown 
70914d0ab18SJacob Faibussowitsch /*S
7108434afd1SBarry Smith   TSI2FunctionFn - A prototype of a `TS` implicit function evaluation function for 2nd order systems that would be passed to `TSSetI2Function()`
71114d0ab18SJacob Faibussowitsch 
71214d0ab18SJacob Faibussowitsch   Calling Sequence:
71314d0ab18SJacob Faibussowitsch + ts   - the `TS` context obtained from `TSCreate()`
71414d0ab18SJacob Faibussowitsch . t    - time at step/stage being solved
71514d0ab18SJacob Faibussowitsch . U    - state vector
71614d0ab18SJacob Faibussowitsch . U_t  - time derivative of state vector
71714d0ab18SJacob Faibussowitsch . U_tt - second time derivative of state vector
71814d0ab18SJacob Faibussowitsch . F    - function vector
71914d0ab18SJacob Faibussowitsch - ctx  - [optional] user-defined context for matrix evaluation routine (may be `NULL`)
72014d0ab18SJacob Faibussowitsch 
72114d0ab18SJacob Faibussowitsch   Level: advanced
72214d0ab18SJacob Faibussowitsch 
723d1c5d1fcSBarry Smith   Note:
7248434afd1SBarry Smith   The deprecated `TSI2Function` still works as a replacement for `TSI2FunctionFn` *.
725d1c5d1fcSBarry Smith 
7268434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetI2Function()`, `DMTSSetI2Function()`, `TSIFunctionFn`
72714d0ab18SJacob Faibussowitsch S*/
728*2a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSI2FunctionFn(TS ts, PetscReal t, Vec U, Vec U_t, Vec U_tt, Vec F, PetscCtx ctx);
729d1c5d1fcSBarry Smith 
7308434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSI2FunctionFn *TSI2Function;
73114d0ab18SJacob Faibussowitsch 
73214d0ab18SJacob Faibussowitsch /*S
7338434afd1SBarry Smith   TSI2JacobianFn - A prototype of a `TS` implicit Jacobian evaluation function for 2nd order systems that would be passed to `TSSetI2Jacobian()`
73414d0ab18SJacob Faibussowitsch 
73514d0ab18SJacob Faibussowitsch   Calling Sequence:
73614d0ab18SJacob Faibussowitsch + ts   - the `TS` context obtained from `TSCreate()`
73714d0ab18SJacob Faibussowitsch . t    - time at step/stage being solved
73814d0ab18SJacob Faibussowitsch . U    - state vector
73914d0ab18SJacob Faibussowitsch . U_t  - time derivative of state vector
74014d0ab18SJacob Faibussowitsch . U_tt - second time derivative of state vector
74114d0ab18SJacob Faibussowitsch . v    - shift for U_t
74214d0ab18SJacob Faibussowitsch . a    - shift for U_tt
74314d0ab18SJacob Faibussowitsch . 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
74414d0ab18SJacob Faibussowitsch . jac  - matrix from which to construct the preconditioner, may be same as `J`
74514d0ab18SJacob Faibussowitsch - ctx  - [optional] user-defined context for matrix evaluation routine
74614d0ab18SJacob Faibussowitsch 
74714d0ab18SJacob Faibussowitsch   Level: advanced
74814d0ab18SJacob Faibussowitsch 
749d1c5d1fcSBarry Smith   Note:
7508434afd1SBarry Smith   The deprecated `TSI2Jacobian` still works as a replacement for `TSI2JacobianFn` *.
75114d0ab18SJacob Faibussowitsch 
7528434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetI2Jacobian()`, `DMTSSetI2Jacobian()`, `TSIFunctionFn`, `TSIJacobianFn`, `TSRHSFunctionFn`, `TSRHSJacobianFn`
753d1c5d1fcSBarry Smith S*/
754*2a8381b2SBarry Smith 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);
755d1c5d1fcSBarry Smith 
7568434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSI2JacobianFn *TSI2Jacobian;
757d1c5d1fcSBarry Smith 
758*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSSetI2Function(TS, Vec, TSI2FunctionFn *, PetscCtx);
759*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSGetI2Function(TS, Vec *, TSI2FunctionFn **, PetscCtxRt);
760*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSSetI2Jacobian(TS, Mat, Mat, TSI2JacobianFn *, PetscCtx);
761*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSGetI2Jacobian(TS, Mat *, Mat *, TSI2JacobianFn **, PetscCtxRt);
762efe9872eSLisandro Dalcin 
763545aaa6fSHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitSetIS(TS, const char[], IS);
764545aaa6fSHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetIS(TS, const char[], IS *);
765*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSRHSSplitSetRHSFunction(TS, const char[], Vec, TSRHSFunctionFn *, PetscCtx);
766*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSRHSSplitSetIFunction(TS, const char[], Vec, TSIFunctionFn *, PetscCtx);
767*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSRHSSplitSetIJacobian(TS, const char[], Mat, Mat, TSIJacobianFn *, PetscCtx);
7681d06f6b3SHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTS(TS, const char[], TS *);
7691d06f6b3SHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTSs(TS, PetscInt *, TS *[]);
7700fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSSetUseSplitRHSFunction(TS, PetscBool);
7710fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSGetUseSplitRHSFunction(TS, PetscBool *);
7724bd3aaa3SHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetSNES(TS, SNES *);
7734bd3aaa3SHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitSetSNES(TS, SNES);
774d9194312SHong Zhang 
7758434afd1SBarry Smith PETSC_EXTERN TSRHSFunctionFn TSComputeRHSFunctionLinear;
7768434afd1SBarry Smith PETSC_EXTERN TSRHSJacobianFn TSComputeRHSJacobianConstant;
777*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode  TSComputeIFunctionLinear(TS, PetscReal, Vec, Vec, Vec, PetscCtx);
778*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode  TSComputeIJacobianConstant(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscCtx);
7791c8a10a0SJed Brown PETSC_EXTERN PetscErrorCode  TSComputeSolutionFunction(TS, PetscReal, Vec);
7809b7cd975SBarry Smith PETSC_EXTERN PetscErrorCode  TSComputeForcingFunction(TS, PetscReal, Vec);
781*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode  TSComputeIJacobianDefaultColor(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscCtx);
782d7cfae9bSHong Zhang PETSC_EXTERN PetscErrorCode  TSPruneIJacobianColor(TS, Mat, Mat);
783e34be4c2SBarry Smith 
784014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
785b8123daeSJed Brown PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS, PetscReal));
7869be3e283SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS, PetscReal, PetscInt, Vec *));
787c688d042SShri Abhyankar PETSC_EXTERN PetscErrorCode TSSetPostEvaluate(TS, PetscErrorCode (*)(TS));
788014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
789*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSSetResize(TS, PetscBool, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscBool *, PetscCtx), PetscErrorCode (*)(TS, PetscInt, Vec[], Vec[], PetscCtx), PetscCtx);
790014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPreStep(TS);
791b8123daeSJed Brown PETSC_EXTERN PetscErrorCode TSPreStage(TS, PetscReal);
7920fc7ecceSBarry Smith PETSC_EXTERN PetscErrorCode TSPostStage(TS, PetscReal, PetscInt, Vec[]);
793c688d042SShri Abhyankar PETSC_EXTERN PetscErrorCode TSPostEvaluate(TS);
794014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPostStep(TS);
7956bd3a4fdSStefano Zampini PETSC_EXTERN PetscErrorCode TSResize(TS);
7966bd3a4fdSStefano Zampini PETSC_EXTERN PetscErrorCode TSResizeRetrieveVec(TS, const char *, Vec *);
7976bd3a4fdSStefano Zampini PETSC_EXTERN PetscErrorCode TSResizeRegisterVec(TS, const char *, Vec);
798c6bf8827SStefano Zampini PETSC_EXTERN PetscErrorCode TSGetStepResize(TS, PetscBool *);
7996bd3a4fdSStefano Zampini 
800014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSInterpolate(TS, PetscReal, Vec);
801014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetTolerances(TS, PetscReal, Vec, PetscReal, Vec);
802c5033834SJed Brown PETSC_EXTERN PetscErrorCode TSGetTolerances(TS, PetscReal *, Vec *, PetscReal *, Vec *);
8037453f775SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *);
8048a175baeSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm(TS, Vec, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *);
805014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS, PetscReal);
806014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS, PetscReal *);
807d183316bSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSSetFunctionDomainError(TS, PetscErrorCode (*)(TS, PetscReal, Vec, PetscBool *));
808d183316bSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSFunctionDomainError(TS, PetscReal, Vec, PetscBool *);
809000e7ae3SMatthew Knepley 
810*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS, PetscErrorCode (*)(TS, PetscReal *, PetscCtx), PetscCtx);
811*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS, PetscReal *, PetscCtx);
812014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS, PetscReal);
813*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS, PetscErrorCode (*)(TS, Vec, PetscCtx, PetscReal *, PetscBool *), PetscCtx);
814014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS, PetscReal);
815014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);
8160fe2f5c2SJames Wright PETSC_EXTERN PetscErrorCode TSPseudoComputeFunction(TS, Vec, Vec *, PetscReal *);
81721c89e3eSBarry Smith 
818014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPythonSetType(TS, const char[]);
819ebead697SStefano Zampini PETSC_EXTERN PetscErrorCode TSPythonGetType(TS, const char *[]);
8201d6018f0SLisandro Dalcin 
821014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS, PetscReal, Vec, Vec);
822d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS, PetscReal, Vec, Mat, Mat);
823014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS, PetscReal, Vec, Vec, Vec, PetscBool);
824d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscBool);
825efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSComputeI2Function(TS, PetscReal, Vec, Vec, Vec, Vec);
826efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSComputeI2Jacobian(TS, PetscReal, Vec, Vec, Vec, PetscReal, PetscReal, Mat, Mat);
827f9c1d6abSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS, PetscReal, PetscReal, PetscReal *, PetscReal *);
828818ad0c1SBarry Smith 
829014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS, Vec, Vec);
830d6ebe24aSShri Abhyankar 
831*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetBoundaryLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, PetscCtx), PetscCtx);
832*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM, TSRHSFunctionFn *, PetscCtx);
833*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM, TSRHSFunctionFn **, PetscCtxRt);
83449abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionContextDestroy(DM, PetscCtxDestroyFn *);
835*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM, TSRHSJacobianFn *, PetscCtx);
836*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM, TSRHSJacobianFn **, PetscCtxRt);
83749abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobianContextDestroy(DM, PetscCtxDestroyFn *);
838*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM, TSIFunctionFn *, PetscCtx);
839*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM, TSIFunctionFn **, PetscCtxRt);
84049abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetIFunctionContextDestroy(DM, PetscCtxDestroyFn *);
841*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM, TSIJacobianFn *, PetscCtx);
842*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM, TSIJacobianFn **, PetscCtxRt);
84349abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetIJacobianContextDestroy(DM, PetscCtxDestroyFn *);
844*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetI2Function(DM, TSI2FunctionFn *, PetscCtx);
845*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetI2Function(DM, TSI2FunctionFn **, PetscCtxRt);
84649abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetI2FunctionContextDestroy(DM, PetscCtxDestroyFn *);
847*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetI2Jacobian(DM, TSI2JacobianFn *, PetscCtx);
848*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetI2Jacobian(DM, TSI2JacobianFn **, PetscCtxRt);
84949abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetI2JacobianContextDestroy(DM, PetscCtxDestroyFn *);
850efe9872eSLisandro Dalcin 
85114d0ab18SJacob Faibussowitsch /*S
8528434afd1SBarry Smith   TSTransientVariableFn - A prototype of a function to transform from state to transient variables that would be passed to `TSSetTransientVariable()`
85314d0ab18SJacob Faibussowitsch 
85414d0ab18SJacob Faibussowitsch   Calling Sequence:
85514d0ab18SJacob Faibussowitsch + ts  - timestep context
85614d0ab18SJacob Faibussowitsch . p   - input vector (primitive form)
85714d0ab18SJacob Faibussowitsch . c   - output vector, transient variables (conservative form)
85814d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context
85914d0ab18SJacob Faibussowitsch 
86014d0ab18SJacob Faibussowitsch   Level: advanced
86114d0ab18SJacob Faibussowitsch 
862d1c5d1fcSBarry Smith   Note:
8638434afd1SBarry Smith   The deprecated `TSTransientVariable` still works as a replacement for `TSTransientVariableFn` *.
864d1c5d1fcSBarry Smith 
86514d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetTransientVariable()`, `DMTSSetTransientVariable()`
86614d0ab18SJacob Faibussowitsch S*/
867*2a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSTransientVariableFn(TS ts, Vec p, Vec c, PetscCtx ctx);
86814d0ab18SJacob Faibussowitsch 
8698434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSTransientVariableFn *TSTransientVariable;
870d1c5d1fcSBarry Smith 
871*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSSetTransientVariable(TS, TSTransientVariableFn *, PetscCtx);
872*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetTransientVariable(DM, TSTransientVariableFn *, PetscCtx);
873*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetTransientVariable(DM, TSTransientVariableFn **, PetscCtx);
874e3c11fc1SJed Brown PETSC_EXTERN PetscErrorCode TSComputeTransientVariable(TS, Vec, Vec);
875e3c11fc1SJed Brown PETSC_EXTERN PetscErrorCode TSHasTransientVariable(TS, PetscBool *);
876e3c11fc1SJed Brown 
877*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM, TSSolutionFn *, PetscCtx);
878*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM, TSSolutionFn **, PetscCtxRt);
879*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM, TSForcingFn *, PetscCtx);
880*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM, TSForcingFn **, PetscCtxRt);
881e17bd7bbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCheckResidual(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscReal *);
882e17bd7bbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCheckJacobian(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscBool *, PetscReal *);
8837f96f943SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec);
8849b7cd975SBarry Smith 
885*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetIFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, Vec, PetscCtx), PetscCtxRt);
886*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, Vec, PetscCtx), PetscCtx);
887*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetIJacobianLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscCtx), PetscCtxRt);
888*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscCtx), PetscCtx);
889*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetRHSFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, PetscCtx), PetscCtxRt);
890*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, PetscCtx), PetscCtx);
891b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrix(DM);
892b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrixLumped(DM);
893b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSDestroyRHSMassMatrix(DM);
8946c6b9e74SPeter Brune 
895*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM, PetscErrorCode (*)(PetscCtx, PetscViewer), PetscErrorCode (*)(PetscCtx *, PetscViewer));
896*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM, PetscErrorCode (*)(PetscCtx, PetscViewer), PetscErrorCode (*)(PetscCtx *, PetscViewer));
89724989b8cSPeter Brune 
89814d0ab18SJacob Faibussowitsch /*S
899dd8e379bSPierre Jolivet   DMDATSRHSFunctionLocalFn - A prototype of a local `TS` right-hand side residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetRHSFunctionLocal()`
9006c6b9e74SPeter Brune 
90114d0ab18SJacob Faibussowitsch   Calling Sequence:
90214d0ab18SJacob Faibussowitsch + info - defines the subdomain to evaluate the residual on
90314d0ab18SJacob Faibussowitsch . t    - time at which to evaluate residual
90414d0ab18SJacob Faibussowitsch . x    - array of local state information
90514d0ab18SJacob Faibussowitsch . f    - output array of local residual information
90614d0ab18SJacob Faibussowitsch - ctx  - optional user context
90714d0ab18SJacob Faibussowitsch 
90814d0ab18SJacob Faibussowitsch   Level: beginner
90914d0ab18SJacob Faibussowitsch 
910d1c5d1fcSBarry Smith   Note:
9118434afd1SBarry Smith   The deprecated `DMDATSRHSFunctionLocal` still works as a replacement for `DMDATSRHSFunctionLocalFn` *.
912d1c5d1fcSBarry Smith 
9138434afd1SBarry Smith .seealso: `DMDA`, `DMDATSSetRHSFunctionLocal()`, `TSRHSFunctionFn`, `DMDATSRHSJacobianLocalFn`, `DMDATSIJacobianLocalFn`, `DMDATSIFunctionLocalFn`
91414d0ab18SJacob Faibussowitsch S*/
915*2a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDATSRHSFunctionLocalFn(DMDALocalInfo *info, PetscReal t, void *x, void *f, PetscCtx ctx);
916d1c5d1fcSBarry Smith 
9178434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef DMDATSRHSFunctionLocalFn *DMDATSRHSFunctionLocal;
91814d0ab18SJacob Faibussowitsch 
91914d0ab18SJacob Faibussowitsch /*S
9208434afd1SBarry Smith   DMDATSRHSJacobianLocalFn - A prototype of a local residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetRHSJacobianLocal()`
92114d0ab18SJacob Faibussowitsch 
92214d0ab18SJacob Faibussowitsch   Calling Sequence:
92314d0ab18SJacob Faibussowitsch + info - defines the subdomain to evaluate the residual on
92414d0ab18SJacob Faibussowitsch . t    - time at which to evaluate residual
92514d0ab18SJacob Faibussowitsch . x    - array of local state information
92614d0ab18SJacob Faibussowitsch . J    - Jacobian matrix
92714d0ab18SJacob Faibussowitsch . B    - matrix from which to construct the preconditioner; often same as `J`
92814d0ab18SJacob Faibussowitsch - ctx  - optional context
92914d0ab18SJacob Faibussowitsch 
93014d0ab18SJacob Faibussowitsch   Level: beginner
93114d0ab18SJacob Faibussowitsch 
932d1c5d1fcSBarry Smith   Note:
9338434afd1SBarry Smith   The deprecated `DMDATSRHSJacobianLocal` still works as a replacement for `DMDATSRHSJacobianLocalFn` *.
934d1c5d1fcSBarry Smith 
9358434afd1SBarry Smith .seealso: `DMDA`, `DMDATSSetRHSJacobianLocal()`, `TSRHSJacobianFn`, `DMDATSRHSFunctionLocalFn`, `DMDATSIJacobianLocalFn`, `DMDATSIFunctionLocalFn`
93614d0ab18SJacob Faibussowitsch S*/
937*2a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDATSRHSJacobianLocalFn(DMDALocalInfo *info, PetscReal t, void *x, Mat J, Mat B, PetscCtx ctx);
938d1c5d1fcSBarry Smith 
9398434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef DMDATSRHSJacobianLocalFn *DMDATSRHSJacobianLocal;
94014d0ab18SJacob Faibussowitsch 
94114d0ab18SJacob Faibussowitsch /*S
9428434afd1SBarry Smith   DMDATSIFunctionLocalFn - A prototype of a local residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetIFunctionLocal()`
94314d0ab18SJacob Faibussowitsch 
94414d0ab18SJacob Faibussowitsch   Calling Sequence:
94514d0ab18SJacob Faibussowitsch + info  - defines the subdomain to evaluate the residual on
94614d0ab18SJacob Faibussowitsch . t     - time at which to evaluate residual
94714d0ab18SJacob Faibussowitsch . x     - array of local state information
94814d0ab18SJacob Faibussowitsch . xdot  - array of local time derivative information
94914d0ab18SJacob Faibussowitsch . imode - output array of local function evaluation information
95014d0ab18SJacob Faibussowitsch - ctx   - optional context
95114d0ab18SJacob Faibussowitsch 
95214d0ab18SJacob Faibussowitsch   Level: beginner
95314d0ab18SJacob Faibussowitsch 
954d1c5d1fcSBarry Smith   Note:
9558434afd1SBarry Smith   The deprecated `DMDATSIFunctionLocal` still works as a replacement for `DMDATSIFunctionLocalFn` *.
956d1c5d1fcSBarry Smith 
9578434afd1SBarry Smith .seealso: `DMDA`, `DMDATSSetIFunctionLocal()`, `DMDATSIJacobianLocalFn`, `TSIFunctionFn`
95814d0ab18SJacob Faibussowitsch S*/
959*2a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDATSIFunctionLocalFn(DMDALocalInfo *info, PetscReal t, void *x, void *xdot, void *imode, PetscCtx ctx);
960d1c5d1fcSBarry Smith 
9618434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef DMDATSIFunctionLocalFn *DMDATSIFunctionLocal;
96214d0ab18SJacob Faibussowitsch 
96314d0ab18SJacob Faibussowitsch /*S
9648434afd1SBarry Smith   DMDATSIJacobianLocalFn - A prototype of a local residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetIJacobianLocal()`
96514d0ab18SJacob Faibussowitsch 
96614d0ab18SJacob Faibussowitsch   Calling Sequence:
96714d0ab18SJacob Faibussowitsch + info  - defines the subdomain to evaluate the residual on
96814d0ab18SJacob Faibussowitsch . t     - time at which to evaluate the jacobian
96914d0ab18SJacob Faibussowitsch . x     - array of local state information
97014d0ab18SJacob Faibussowitsch . xdot  - time derivative at this state
97114d0ab18SJacob Faibussowitsch . shift - see `TSSetIJacobian()` for the meaning of this parameter
97214d0ab18SJacob Faibussowitsch . J     - Jacobian matrix
97314d0ab18SJacob Faibussowitsch . B     - matrix from which to construct the preconditioner; often same as `J`
97414d0ab18SJacob Faibussowitsch - ctx   - optional context
97514d0ab18SJacob Faibussowitsch 
97614d0ab18SJacob Faibussowitsch   Level: beginner
97714d0ab18SJacob Faibussowitsch 
978d1c5d1fcSBarry Smith   Note:
9798434afd1SBarry Smith   The deprecated `DMDATSIJacobianLocal` still works as a replacement for `DMDATSIJacobianLocalFn` *.
98014d0ab18SJacob Faibussowitsch 
9818434afd1SBarry Smith .seealso: `DMDA` `DMDATSSetIJacobianLocal()`, `TSIJacobianFn`, `DMDATSIFunctionLocalFn`, `DMDATSRHSFunctionLocalFn`,  `DMDATSRHSJacobianlocal()`
982d1c5d1fcSBarry Smith S*/
983*2a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDATSIJacobianLocalFn(DMDALocalInfo *info, PetscReal t, void *x, void *xdot, PetscReal shift, Mat J, Mat B, PetscCtx ctx);
984d1c5d1fcSBarry Smith 
9858434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef DMDATSIJacobianLocalFn *DMDATSIJacobianLocal;
986d1c5d1fcSBarry Smith 
9878434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM, InsertMode, DMDATSRHSFunctionLocalFn *, void *);
9888434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM, DMDATSRHSJacobianLocalFn *, void *);
9898434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM, InsertMode, DMDATSIFunctionLocalFn *, void *);
9908434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM, DMDATSIJacobianLocalFn *, void *);
9916c6b9e74SPeter Brune 
992be1b0d75SMatthew G. Knepley typedef struct _n_TSMonitorLGCtx *TSMonitorLGCtx;
993d1212d36SBarry Smith typedef struct {
994d1212d36SBarry Smith   Vec            ray;
995d1212d36SBarry Smith   VecScatter     scatter;
996d1212d36SBarry Smith   PetscViewer    viewer;
99751b4a12fSMatthew G. Knepley   TSMonitorLGCtx lgctx;
998d1212d36SBarry Smith } TSMonitorDMDARayCtx;
999*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(PetscCtxRt);
1000d1212d36SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS, PetscInt, PetscReal, Vec, void *);
100151b4a12fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS, PetscInt, PetscReal, Vec, void *);
1002d1212d36SBarry Smith 
1003bdad233fSMatthew Knepley /* Dynamic creation and loading functions */
1004140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList TSList;
100519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode    TSGetType(TS, TSType *);
100619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode    TSSetType(TS, TSType);
1007bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode    TSRegister(const char[], PetscErrorCode (*)(TS));
100830de9b25SBarry Smith 
1009014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetSNES(TS, SNES *);
1010deb2cd25SJed Brown PETSC_EXTERN PetscErrorCode TSSetSNES(TS, SNES);
1011014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetKSP(TS, KSP *);
1012818ad0c1SBarry Smith 
1013014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSView(TS, PetscViewer);
101455849f57SBarry Smith PETSC_EXTERN PetscErrorCode TSLoad(TS, PetscViewer);
1015fe2efc57SMark PETSC_EXTERN PetscErrorCode TSViewFromOptions(TS, PetscObject, const char[]);
1016fe2efc57SMark PETSC_EXTERN PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory, PetscObject, const char[]);
101755849f57SBarry Smith 
101855849f57SBarry Smith #define TS_FILE_CLASSID 1211225
101921c89e3eSBarry Smith 
1020*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS, PetscCtx);
1021*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS, PetscCtxRt);
102221c89e3eSBarry Smith 
1023a9f9c1f6SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtx *);
1024a9f9c1f6SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *);
10254f09c107SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS, PetscInt, PetscReal, Vec, void *);
10264f09c107SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS, PetscInt, PetscReal, Vec, void *);
102731152f8aSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS, const char *const *);
1028b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS, const char *const **);
1029e673d494SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx, const char *const *);
1030a66092f1SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS, const char *const *);
1031a66092f1SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx, const char *const *);
103249abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS, PetscErrorCode (*)(void *, Vec, Vec *), PetscCtxDestroyFn *, void *);
103349abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx, PetscErrorCode (*)(void *, Vec, Vec *), PetscCtxDestroyFn *, void *);
10344f09c107SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS, PetscInt, PetscReal, Vec, void *);
1035201da799SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS, PetscInt, PetscReal, Vec, void *);
1036201da799SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS, PetscInt, PetscReal, Vec, void *);
1037edbaebb3SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorError(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
1038d0c080abSJoseph Pusztay PETSC_EXTERN PetscErrorCode TSDMSwarmMonitorMoments(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
1039ef20d060SBarry Smith 
1040e669de00SBarry Smith struct _n_TSMonitorLGCtxNetwork {
1041e669de00SBarry Smith   PetscInt     nlg;
1042e669de00SBarry Smith   PetscDrawLG *lg;
1043e669de00SBarry Smith   PetscBool    semilogy;
1044e669de00SBarry Smith   PetscInt     howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
1045e669de00SBarry Smith };
1046e669de00SBarry Smith typedef struct _n_TSMonitorLGCtxNetwork *TSMonitorLGCtxNetwork;
1047e669de00SBarry Smith PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkDestroy(TSMonitorLGCtxNetwork *);
1048e669de00SBarry Smith PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkCreate(TS, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtxNetwork *);
1049e669de00SBarry Smith PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkSolution(TS, PetscInt, PetscReal, Vec, void *);
1050e669de00SBarry Smith 
1051b3d3934dSBarry Smith typedef struct _n_TSMonitorEnvelopeCtx *TSMonitorEnvelopeCtx;
1052b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeCtxCreate(TS, TSMonitorEnvelopeCtx *);
1053b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode             TSMonitorEnvelope(TS, PetscInt, PetscReal, Vec, void *);
1054b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeGetBounds(TS, Vec *, Vec *);
1055b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *);
1056b3d3934dSBarry Smith 
10578189c53fSBarry Smith typedef struct _n_TSMonitorSPEigCtx *TSMonitorSPEigCtx;
10588189c53fSBarry Smith PETSC_EXTERN PetscErrorCode          TSMonitorSPEigCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorSPEigCtx *);
10598189c53fSBarry Smith PETSC_EXTERN PetscErrorCode          TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *);
10608189c53fSBarry Smith PETSC_EXTERN PetscErrorCode          TSMonitorSPEig(TS, PetscInt, PetscReal, Vec, void *);
10613914022bSBarry Smith 
10621b575b74SJoseph Pusztay typedef struct _n_TSMonitorSPCtx *TSMonitorSPCtx;
106360e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorSPCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscBool, PetscBool, TSMonitorSPCtx *);
10641b575b74SJoseph Pusztay PETSC_EXTERN PetscErrorCode       TSMonitorSPCtxDestroy(TSMonitorSPCtx *);
10650ec8ee2bSJoseph Pusztay PETSC_EXTERN PetscErrorCode       TSMonitorSPSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
10661b575b74SJoseph Pusztay 
106760e16b1bSMatthew G. Knepley typedef struct _n_TSMonitorHGCtx *TSMonitorHGCtx;
106860e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorHGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscInt, PetscBool, TSMonitorHGCtx *);
106960e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
107060e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorHGCtxDestroy(TSMonitorHGCtx *);
107160e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
107260e16b1bSMatthew G. Knepley 
1073ca4445c7SIlya Fursov PETSC_EXTERN PetscErrorCode TSSetEventHandler(TS, PetscInt, PetscInt[], PetscBool[], PetscErrorCode (*)(TS, PetscReal, Vec, PetscReal[], void *), PetscErrorCode (*)(TS, PetscInt, PetscInt[], PetscReal, Vec, PetscBool, void *), void *);
1074ca4445c7SIlya Fursov PETSC_EXTERN PetscErrorCode TSSetPostEventStep(TS, PetscReal);
1075fe4ad979SIlya Fursov PETSC_EXTERN PetscErrorCode TSSetPostEventSecondStep(TS, PetscReal);
TSSetPostEventIntervalStep(TS ts,PetscReal dt)1076fe4ad979SIlya Fursov PETSC_DEPRECATED_FUNCTION(3, 21, 0, "TSSetPostEventSecondStep()", ) static inline PetscErrorCode TSSetPostEventIntervalStep(TS ts, PetscReal dt)
1077fe4ad979SIlya Fursov {
1078fe4ad979SIlya Fursov   return TSSetPostEventSecondStep(ts, dt);
1079fe4ad979SIlya Fursov }
10806427ac75SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS, PetscReal, PetscReal[]);
10811ea83e56SMiguel PETSC_EXTERN PetscErrorCode TSGetNumEvents(TS, PetscInt *);
10821ea83e56SMiguel 
108376bdecfbSBarry Smith /*J
1084195e9b02SBarry Smith    TSSSPType - string with the name of a `TSSSP` scheme.
1085815f1ad5SJed Brown 
1086815f1ad5SJed Brown    Level: beginner
1087815f1ad5SJed Brown 
10881cc06b55SBarry Smith .seealso: [](ch_ts), `TSSSPSetType()`, `TS`, `TSSSP`
108976bdecfbSBarry Smith J*/
109019fd82e9SBarry Smith typedef const char *TSSSPType;
1091815f1ad5SJed Brown #define TSSSPRKS2  "rks2"
1092815f1ad5SJed Brown #define TSSSPRKS3  "rks3"
1093815f1ad5SJed Brown #define TSSSPRK104 "rk104"
1094815f1ad5SJed Brown 
109519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode    TSSSPSetType(TS, TSSSPType);
109619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode    TSSSPGetType(TS, TSSSPType *);
1097014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSSPSetNumStages(TS, PetscInt);
1098014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSSPGetNumStages(TS, PetscInt *);
1099787849ffSJed Brown PETSC_EXTERN PetscErrorCode    TSSSPInitializePackage(void);
1100560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode    TSSSPFinalizePackage(void);
1101013797aaSBarry Smith PETSC_EXTERN PetscFunctionList TSSSPList;
1102815f1ad5SJed Brown 
1103ed657a08SJed Brown /*S
110484df9cb4SJed Brown    TSAdapt - Abstract object that manages time-step adaptivity
110584df9cb4SJed Brown 
110684df9cb4SJed Brown    Level: beginner
110784df9cb4SJed Brown 
1108efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TS`, `TSGetAdapt()`, `TSAdaptCreate()`, `TSAdaptType`
110984df9cb4SJed Brown S*/
111084df9cb4SJed Brown typedef struct _p_TSAdapt *TSAdapt;
111184df9cb4SJed Brown 
1112f0fc11ceSJed Brown /*J
111387497f52SBarry Smith    TSAdaptType - String with the name of `TSAdapt` scheme.
111484df9cb4SJed Brown 
111584df9cb4SJed Brown    Level: beginner
111684df9cb4SJed Brown 
1117efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSGetAdapt()`, `TSAdaptSetType()`, `TS`, `TSAdapt`
1118f0fc11ceSJed Brown J*/
1119f8963224SJed Brown typedef const char *TSAdaptType;
1120b0f836d7SJed Brown #define TSADAPTNONE    "none"
11211917a363SLisandro Dalcin #define TSADAPTBASIC   "basic"
1122cb7ab186SLisandro Dalcin #define TSADAPTDSP     "dsp"
11238d59e960SJed Brown #define TSADAPTCFL     "cfl"
11241917a363SLisandro Dalcin #define TSADAPTGLEE    "glee"
1125d949e4a4SStefano Zampini #define TSADAPTHISTORY "history"
112684df9cb4SJed Brown 
1127552698daSJed Brown PETSC_EXTERN PetscErrorCode TSGetAdapt(TS, TSAdapt *);
1128bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[], PetscErrorCode (*)(TSAdapt));
1129607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void);
1130014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
1131014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm, TSAdapt *);
113219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt, TSAdaptType);
1133d0288e4fSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetType(TSAdapt, TSAdaptType *);
1134014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt, const char[]);
1135014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
1136014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt, const char[], PetscInt, PetscInt, PetscReal, PetscReal, PetscBool);
1137014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt, PetscInt *, const PetscInt **, const PetscInt **, const PetscReal **, const PetscReal **);
1138014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt, TS, PetscReal, PetscInt *, PetscReal *, PetscBool *);
1139b295832fSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt, TS, PetscReal, Vec, PetscBool *);
1140014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt, PetscViewer);
1141ad6bc421SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt, PetscViewer);
1142ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(TSAdapt, PetscOptionItems);
1143099af0a3SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt);
1144014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt *);
1145014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt, PetscBool);
1146bf997491SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt, PetscBool);
11471917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptSetSafety(TSAdapt, PetscReal, PetscReal);
11481917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetSafety(TSAdapt, PetscReal *, PetscReal *);
114976cddca1SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt, PetscReal);
115076cddca1SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt, PetscReal *);
11511917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptSetClip(TSAdapt, PetscReal, PetscReal);
11521917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetClip(TSAdapt, PetscReal *, PetscReal *);
115362c23b28SMark PETSC_EXTERN PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt, PetscReal);
115462c23b28SMark PETSC_EXTERN PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt, PetscReal *);
1155014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt, PetscReal, PetscReal);
11561917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetStepLimits(TSAdapt, PetscReal *, PetscReal *);
1157b295832fSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt, PetscErrorCode (*)(TSAdapt, TS, PetscReal, Vec, PetscBool *));
1158ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptHistorySetHistory(TSAdapt, PetscInt, PetscReal[], PetscBool);
115975017583SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptHistorySetTrajectory(TSAdapt, TSTrajectory, PetscBool);
116075017583SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptHistoryGetStep(TSAdapt, PetscInt, PetscReal *, PetscReal *);
1161de50f1caSBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt, PetscInt);
1162cb7ab186SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptDSPSetFilter(TSAdapt, const char *);
1163cb7ab186SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptDSPSetPID(TSAdapt, PetscReal, PetscReal, PetscReal);
116484df9cb4SJed Brown 
116584df9cb4SJed Brown /*S
116687497f52SBarry Smith    TSGLLEAdapt - Abstract object that manages time-step adaptivity for `TSGLLE`
1167ed657a08SJed Brown 
1168ed657a08SJed Brown    Level: beginner
1169ed657a08SJed Brown 
1170195e9b02SBarry Smith    Developer Note:
117187497f52SBarry Smith    This functionality should be replaced by the `TSAdapt`.
117284df9cb4SJed Brown 
11731cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLEAdaptCreate()`, `TSGLLEAdaptType`
1174ed657a08SJed Brown S*/
117526d28e4eSEmil Constantinescu typedef struct _p_TSGLLEAdapt *TSGLLEAdapt;
1176ed657a08SJed Brown 
117776bdecfbSBarry Smith /*J
117887497f52SBarry Smith    TSGLLEAdaptType - String with the name of `TSGLLEAdapt` scheme
1179ed657a08SJed Brown 
1180ed657a08SJed Brown    Level: beginner
1181ed657a08SJed Brown 
1182195e9b02SBarry Smith    Developer Note:
1183195e9b02SBarry Smith    This functionality should be replaced by the `TSAdaptType`.
1184195e9b02SBarry Smith 
11851cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLLEAdaptSetType()`, `TS`
118676bdecfbSBarry Smith J*/
118726d28e4eSEmil Constantinescu typedef const char *TSGLLEAdaptType;
118826d28e4eSEmil Constantinescu #define TSGLLEADAPT_NONE "none"
118926d28e4eSEmil Constantinescu #define TSGLLEADAPT_SIZE "size"
119026d28e4eSEmil Constantinescu #define TSGLLEADAPT_BOTH "both"
1191ed657a08SJed Brown 
119226d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptRegister(const char[], PetscErrorCode (*)(TSGLLEAdapt));
119326d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptInitializePackage(void);
119426d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptFinalizePackage(void);
119526d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate(MPI_Comm, TSGLLEAdapt *);
119626d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt, TSGLLEAdaptType);
119726d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt, const char[]);
119826d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptChoose(TSGLLEAdapt, PetscInt, const PetscInt[], const PetscReal[], const PetscReal[], PetscInt, PetscReal, PetscReal, PetscInt *, PetscReal *, PetscBool *);
119926d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt, PetscViewer);
1200ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetFromOptions(TSGLLEAdapt, PetscOptionItems);
120126d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt *);
1202ed657a08SJed Brown 
120376bdecfbSBarry Smith /*J
120487497f52SBarry Smith    TSGLLEAcceptType - String with the name of `TSGLLEAccept` scheme
1205ed657a08SJed Brown 
1206ed657a08SJed Brown    Level: beginner
1207ed657a08SJed Brown 
12081cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLLESetAcceptType()`, `TS`, `TSGLLEAccept`
120976bdecfbSBarry Smith J*/
121026d28e4eSEmil Constantinescu typedef const char *TSGLLEAcceptType;
121126d28e4eSEmil Constantinescu #define TSGLLEACCEPT_ALWAYS "always"
1212ed657a08SJed Brown 
1213d1c5d1fcSBarry Smith /*S
12148434afd1SBarry Smith   TSGLLEAcceptFn - A prototype of a `TS` accept function that would be passed to `TSGLLEAcceptRegister()`
1215d1c5d1fcSBarry Smith 
1216d1c5d1fcSBarry Smith   Calling Sequence:
1217d1c5d1fcSBarry Smith + ts  - timestep context
1218d1c5d1fcSBarry Smith . nt - time to end of solution time
1219d1c5d1fcSBarry Smith . h - the proposed step-size
1220d1c5d1fcSBarry Smith . enorm - unknown
1221d1c5d1fcSBarry Smith - accept - output, if the proposal is accepted
1222d1c5d1fcSBarry Smith 
1223d1c5d1fcSBarry Smith   Level: beginner
1224d1c5d1fcSBarry Smith 
1225d1c5d1fcSBarry Smith   Note:
12268434afd1SBarry Smith   The deprecated `TSGLLEAcceptFunction` still works as a replacement for `TSGLLEAcceptFn` *
1227d1c5d1fcSBarry Smith 
12288434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetRHSFunction()`, `DMTSSetRHSFunction()`, `TSIFunctionFn`,
12298434afd1SBarry Smith `TSIJacobianFn`, `TSRHSJacobianFn`, `TSGLLEAcceptRegister()`
1230d1c5d1fcSBarry Smith S*/
12318434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSGLLEAcceptFn)(TS ts, PetscReal nt, PetscReal h, const PetscReal enorm[], PetscBool *accept);
1232d1c5d1fcSBarry Smith 
12338434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSGLLEAcceptFn *TSGLLEAcceptFunction;
1234d1c5d1fcSBarry Smith 
12358434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGLLEAcceptRegister(const char[], TSGLLEAcceptFn *);
1236ed657a08SJed Brown 
123776bdecfbSBarry Smith /*J
1238195e9b02SBarry Smith   TSGLLEType - string with the name of a General Linear `TSGLLE` type
123918b56cb9SJed Brown 
124018b56cb9SJed Brown   Level: beginner
124118b56cb9SJed Brown 
12421cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLESetType()`, `TSGLLERegister()`, `TSGLLEAccept`
124376bdecfbSBarry Smith J*/
124426d28e4eSEmil Constantinescu typedef const char *TSGLLEType;
124526d28e4eSEmil Constantinescu #define TSGLLE_IRKS "irks"
124618b56cb9SJed Brown 
124726d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLERegister(const char[], PetscErrorCode (*)(TS));
124826d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEInitializePackage(void);
124926d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEFinalizePackage(void);
125026d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLESetType(TS, TSGLLEType);
125126d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEGetAdapt(TS, TSGLLEAdapt *);
125226d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLESetAcceptType(TS, TSGLLEAcceptType);
125318b56cb9SJed Brown 
1254020d8f30SJed Brown /*J
1255195e9b02SBarry Smith    TSEIMEXType - String with the name of an Extrapolated IMEX `TSEIMEX` type
12567d5697caSHong Zhang 
12577d5697caSHong Zhang    Level: beginner
12587d5697caSHong Zhang 
12591cc06b55SBarry Smith .seealso: [](ch_ts), `TSEIMEXSetType()`, `TS`, `TSEIMEX`, `TSEIMEXRegister()`
12607d5697caSHong Zhang J*/
12617d5697caSHong Zhang #define TSEIMEXType char *
12627d5697caSHong Zhang 
1263ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS, PetscInt);
1264ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS, PetscInt, PetscInt);
12657d5697caSHong Zhang PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS, PetscBool);
12667d5697caSHong Zhang 
12677d5697caSHong Zhang /*J
1268195e9b02SBarry Smith    TSRKType - String with the name of a Runge-Kutta `TSRK` type
1269f68a32c8SEmil Constantinescu 
1270f68a32c8SEmil Constantinescu    Level: beginner
1271f68a32c8SEmil Constantinescu 
12721cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSRKSetType()`, `TS`, `TSRK`, `TSRKRegister()`
1273f68a32c8SEmil Constantinescu J*/
1274f68a32c8SEmil Constantinescu typedef const char *TSRKType;
1275f68a32c8SEmil Constantinescu #define TSRK1FE "1fe"
1276fdefd5e5SDebojyoti Ghosh #define TSRK2A  "2a"
1277e7685601SHong Zhang #define TSRK2B  "2b"
1278f68a32c8SEmil Constantinescu #define TSRK3   "3"
1279fdefd5e5SDebojyoti Ghosh #define TSRK3BS "3bs"
1280f68a32c8SEmil Constantinescu #define TSRK4   "4"
1281fdefd5e5SDebojyoti Ghosh #define TSRK5F  "5f"
1282fdefd5e5SDebojyoti Ghosh #define TSRK5DP "5dp"
128305e23783SLisandro Dalcin #define TSRK5BS "5bs"
128437acaa02SHendrik Ranocha #define TSRK6VR "6vr"
128537acaa02SHendrik Ranocha #define TSRK7VR "7vr"
128637acaa02SHendrik Ranocha #define TSRK8VR "8vr"
128705e23783SLisandro Dalcin 
12882ea87ba9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRKGetOrder(TS, PetscInt *);
12890f7a28cdSStefano Zampini PETSC_EXTERN PetscErrorCode TSRKGetType(TS, TSRKType *);
12900f7a28cdSStefano Zampini PETSC_EXTERN PetscErrorCode TSRKSetType(TS, TSRKType);
12910f7a28cdSStefano Zampini PETSC_EXTERN PetscErrorCode TSRKGetTableau(TS, PetscInt *, const PetscReal **, const PetscReal **, const PetscReal **, const PetscReal **, PetscInt *, const PetscReal **, PetscBool *);
12920fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSRKSetMultirate(TS, PetscBool);
12930fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSRKGetMultirate(TS, PetscBool *);
1294f68a32c8SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
1295f68a32c8SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void);
1296560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void);
1297f68a32c8SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void);
1298f68a32c8SEmil Constantinescu 
1299f68a32c8SEmil Constantinescu /*J
1300af27ebaaSBarry Smith    TSMPRKType - String with the name of a partitioned Runge-Kutta `TSMPRK` type
13019f537349Sluzhanghpp 
13029f537349Sluzhanghpp    Level: beginner
13039f537349Sluzhanghpp 
13041cc06b55SBarry Smith .seealso: [](ch_ts), `TSMPRKSetType()`, `TS`, `TSMPRK`, `TSMPRKRegister()`
13059f537349Sluzhanghpp J*/
13064b84eec9SHong Zhang typedef const char *TSMPRKType;
130719c2959aSHong Zhang #define TSMPRK2A22 "2a22"
130819c2959aSHong Zhang #define TSMPRK2A23 "2a23"
130919c2959aSHong Zhang #define TSMPRK2A32 "2a32"
131019c2959aSHong Zhang #define TSMPRK2A33 "2a33"
13114b84eec9SHong Zhang #define TSMPRKP2   "p2"
13124b84eec9SHong Zhang #define TSMPRKP3   "p3"
13139f537349Sluzhanghpp 
1314ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode TSMPRKGetType(TS, TSMPRKType *);
1315ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode TSMPRKSetType(TS, TSMPRKType);
131679bc8a77SHong Zhang 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[]);
13174b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKInitializePackage(void);
13184b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKFinalizePackage(void);
13194b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKRegisterDestroy(void);
13209f537349Sluzhanghpp 
13219f537349Sluzhanghpp /*J
1322195e9b02SBarry Smith    TSIRKType - String with the name of an implicit Runge-Kutta `TSIRK` type
1323d2567f34SHong Zhang 
1324d2567f34SHong Zhang    Level: beginner
1325d2567f34SHong Zhang 
13261cc06b55SBarry Smith .seealso: [](ch_ts), `TSIRKSetType()`, `TS`, `TSIRK`, `TSIRKRegister()`
1327d2567f34SHong Zhang J*/
1328d2567f34SHong Zhang typedef const char *TSIRKType;
1329d2567f34SHong Zhang #define TSIRKGAUSS "gauss"
1330d2567f34SHong Zhang 
1331d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKGetType(TS, TSIRKType *);
1332d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKSetType(TS, TSIRKType);
1333d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKGetNumStages(TS, PetscInt *);
1334d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKSetNumStages(TS, PetscInt);
1335d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKRegister(const char[], PetscErrorCode (*function)(TS));
1336d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKTableauCreate(TS, PetscInt, const PetscReal *, const PetscReal *, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, const PetscScalar *);
1337d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKInitializePackage(void);
1338d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKFinalizePackage(void);
1339d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKRegisterDestroy(void);
1340d2567f34SHong Zhang 
1341d2567f34SHong Zhang /*J
1342195e9b02SBarry Smith    TSGLEEType - String with the name of a General Linear with Error Estimation `TSGLEE` type
1343b6a60446SDebojyoti Ghosh 
1344b6a60446SDebojyoti Ghosh    Level: beginner
1345b6a60446SDebojyoti Ghosh 
13461cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEESetType()`, `TS`, `TSGLEE`, `TSGLEERegister()`
1347b6a60446SDebojyoti Ghosh J*/
1348b6a60446SDebojyoti Ghosh typedef const char *TSGLEEType;
1349ab8c5912SEmil Constantinescu #define TSGLEEi1      "BE1"
1350b6a60446SDebojyoti Ghosh #define TSGLEE23      "23"
1351b6a60446SDebojyoti Ghosh #define TSGLEE24      "24"
1352b6a60446SDebojyoti Ghosh #define TSGLEE25I     "25i"
1353b6a60446SDebojyoti Ghosh #define TSGLEE35      "35"
1354b6a60446SDebojyoti Ghosh #define TSGLEEEXRK2A  "exrk2a"
1355b6a60446SDebojyoti Ghosh #define TSGLEERK32G1  "rk32g1"
1356b6a60446SDebojyoti Ghosh #define TSGLEERK285EX "rk285ex"
135716a05f60SBarry Smith 
1358b6a60446SDebojyoti Ghosh /*J
1359195e9b02SBarry Smith    TSGLEEMode - String with the mode of error estimation for a General Linear with Error Estimation `TSGLEE` type
1360b6a60446SDebojyoti Ghosh 
1361b6a60446SDebojyoti Ghosh    Level: beginner
1362b6a60446SDebojyoti Ghosh 
13631cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEESetMode()`, `TS`, `TSGLEE`, `TSGLEERegister()`
1364b6a60446SDebojyoti Ghosh J*/
1365ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode TSGLEEGetType(TS, TSGLEEType *);
1366ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode TSGLEESetType(TS, TSGLEEType);
136757df6a1bSDebojyoti Ghosh 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[]);
13684bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSGLEERegisterAll(void);
1369b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEEFinalizePackage(void);
1370b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEEInitializePackage(void);
1371b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEERegisterDestroy(void);
1372b6a60446SDebojyoti Ghosh 
1373b6a60446SDebojyoti Ghosh /*J
1374195e9b02SBarry Smith    TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX `TSARKIMEX` type
1375020d8f30SJed Brown 
1376020d8f30SJed Brown    Level: beginner
1377020d8f30SJed Brown 
13781cc06b55SBarry Smith .seealso: [](ch_ts), `TSARKIMEXSetType()`, `TS`, `TSARKIMEX`, `TSARKIMEXRegister()`
1379020d8f30SJed Brown J*/
138019fd82e9SBarry Smith typedef const char *TSARKIMEXType;
1381e817cc15SEmil Constantinescu #define TSARKIMEX1BEE   "1bee"
13821f80e275SEmil Constantinescu #define TSARKIMEXA2     "a2"
13831f80e275SEmil Constantinescu #define TSARKIMEXL2     "l2"
13841f80e275SEmil Constantinescu #define TSARKIMEXARS122 "ars122"
13851f80e275SEmil Constantinescu #define TSARKIMEX2C     "2c"
13868a381b04SJed Brown #define TSARKIMEX2D     "2d"
1387a3a57f36SJed Brown #define TSARKIMEX2E     "2e"
13886cf0794eSJed Brown #define TSARKIMEXPRSSP2 "prssp2"
13898a381b04SJed Brown #define TSARKIMEX3      "3"
13906cf0794eSJed Brown #define TSARKIMEXBPR3   "bpr3"
13916cf0794eSJed Brown #define TSARKIMEXARS443 "ars443"
13928a381b04SJed Brown #define TSARKIMEX4      "4"
13938a381b04SJed Brown #define TSARKIMEX5      "5"
1394ce78bad3SBarry Smith 
1395ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS, TSARKIMEXType *);
1396ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS, TSARKIMEXType);
1397014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS, PetscBool);
13983a28c0e5SStefano Zampini PETSC_EXTERN PetscErrorCode TSARKIMEXGetFullyImplicit(TS, PetscBool *);
13993a2a065fSHong Zhang PETSC_EXTERN PetscErrorCode TSARKIMEXSetFastSlowSplit(TS, PetscBool);
14003a2a065fSHong Zhang PETSC_EXTERN PetscErrorCode TSARKIMEXGetFastSlowSplit(TS, PetscBool *);
140119fd82e9SBarry Smith 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[]);
1402607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void);
1403560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
1404014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);
14058a381b04SJed Brown 
1406020d8f30SJed Brown /*J
1407d27334e2SStefano Zampini    TSDIRKType - String with the name of a Diagonally Implicit Runge-Kutta `TSDIRK` type
1408d27334e2SStefano Zampini 
1409d27334e2SStefano Zampini    Level: beginner
1410d27334e2SStefano Zampini 
1411d27334e2SStefano Zampini .seealso: [](ch_ts), `TSDIRKSetType()`, `TS`, `TSDIRK`, `TSDIRKRegister()`
1412d27334e2SStefano Zampini J*/
1413d27334e2SStefano Zampini typedef const char *TSDIRKType;
1414d27334e2SStefano Zampini #define TSDIRKS212      "s212"
14153405e92cSStefano Zampini #define TSDIRKES122SAL  "es122sal"
1416d27334e2SStefano Zampini #define TSDIRKES213SAL  "es213sal"
1417d27334e2SStefano Zampini #define TSDIRKES324SAL  "es324sal"
1418d27334e2SStefano Zampini #define TSDIRKES325SAL  "es325sal"
1419d27334e2SStefano Zampini #define TSDIRK657A      "657a"
1420d27334e2SStefano Zampini #define TSDIRKES648SA   "es648sa"
1421d27334e2SStefano Zampini #define TSDIRK658A      "658a"
1422d27334e2SStefano Zampini #define TSDIRKS659A     "s659a"
1423d27334e2SStefano Zampini #define TSDIRK7510SAL   "7510sal"
1424d27334e2SStefano Zampini #define TSDIRKES7510SA  "es7510sa"
1425d27334e2SStefano Zampini #define TSDIRK759A      "759a"
1426d27334e2SStefano Zampini #define TSDIRKS7511SAL  "s7511sal"
1427d27334e2SStefano Zampini #define TSDIRK8614A     "8614a"
1428d27334e2SStefano Zampini #define TSDIRK8616SAL   "8616sal"
1429d27334e2SStefano Zampini #define TSDIRKES8516SAL "es8516sal"
1430ce78bad3SBarry Smith 
1431ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode TSDIRKGetType(TS, TSDIRKType *);
1432ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode TSDIRKSetType(TS, TSDIRKType);
1433d27334e2SStefano Zampini PETSC_EXTERN PetscErrorCode TSDIRKRegister(TSDIRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
1434d27334e2SStefano Zampini 
1435d27334e2SStefano Zampini /*J
1436195e9b02SBarry Smith    TSRosWType - String with the name of a Rosenbrock-W `TSROSW` type
1437020d8f30SJed Brown 
1438020d8f30SJed Brown    Level: beginner
1439020d8f30SJed Brown 
14401cc06b55SBarry Smith .seealso: [](ch_ts), `TSRosWSetType()`, `TS`, `TSROSW`, `TSRosWRegister()`
1441020d8f30SJed Brown J*/
144219fd82e9SBarry Smith typedef const char *TSRosWType;
144361692a83SJed Brown #define TSROSW2M          "2m"
144461692a83SJed Brown #define TSROSW2P          "2p"
1445fe7e6d57SJed Brown #define TSROSWRA3PW       "ra3pw"
1446fe7e6d57SJed Brown #define TSROSWRA34PW2     "ra34pw2"
1447337ef93bSJed Brown #define TSROSWR34PRW      "r34prw"
1448337ef93bSJed Brown #define TSROSWR3PRL2      "r3prl2"
1449ef3c5b88SJed Brown #define TSROSWRODAS3      "rodas3"
145048665b53SJed Brown #define TSROSWRODASPR     "rodaspr"
1451337ef93bSJed Brown #define TSROSWRODASPR2    "rodaspr2"
1452ef3c5b88SJed Brown #define TSROSWSANDU3      "sandu3"
1453961f28d0SJed Brown #define TSROSWASSP3P3S1C  "assp3p3s1c"
1454961f28d0SJed Brown #define TSROSWLASSP3P4S2C "lassp3p4s2c"
145543b21953SEmil Constantinescu #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
1456753f8adbSEmil Constantinescu #define TSROSWARK3        "ark3"
14573606a31eSEmil Constantinescu #define TSROSWTHETA1      "theta1"
14583606a31eSEmil Constantinescu #define TSROSWTHETA2      "theta2"
145942faf41dSJed Brown #define TSROSWGRK4T       "grk4t"
146042faf41dSJed Brown #define TSROSWSHAMP4      "shamp4"
146142faf41dSJed Brown #define TSROSWVELDD4      "veldd4"
146242faf41dSJed Brown #define TSROSW4L          "4l"
1463b1c69cc3SEmil Constantinescu 
14646bd7aeb5SHong Zhang PETSC_EXTERN PetscErrorCode TSRosWGetType(TS, TSRosWType *);
14656bd7aeb5SHong Zhang PETSC_EXTERN PetscErrorCode TSRosWSetType(TS, TSRosWType);
1466014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS, PetscBool);
146719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
146819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal);
1469607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void);
1470560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
1471014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);
1472e27a552bSJed Brown 
1473211a84d6SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSBDFSetOrder(TS, PetscInt);
1474211a84d6SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSBDFGetOrder(TS, PetscInt *);
1475211a84d6SLisandro Dalcin 
14766bd7aeb5SHong Zhang /*J
1477195e9b02SBarry Smith   TSBasicSymplecticType - String with the name of a basic symplectic integration `TSBASICSYMPLECTIC` type
14786bd7aeb5SHong Zhang 
14796bd7aeb5SHong Zhang   Level: beginner
14806bd7aeb5SHong Zhang 
14811cc06b55SBarry Smith .seealso: [](ch_ts), `TSBasicSymplecticSetType()`, `TS`, `TSBASICSYMPLECTIC`, `TSBasicSymplecticRegister()`
14826bd7aeb5SHong Zhang J*/
1483e49d4f37SHong Zhang typedef const char *TSBasicSymplecticType;
1484e49d4f37SHong Zhang #define TSBASICSYMPLECTICSIEULER   "1"
1485e49d4f37SHong Zhang #define TSBASICSYMPLECTICVELVERLET "2"
1486e49d4f37SHong Zhang #define TSBASICSYMPLECTIC3         "3"
1487e49d4f37SHong Zhang #define TSBASICSYMPLECTIC4         "4"
1488ce78bad3SBarry Smith 
1489e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticSetType(TS, TSBasicSymplecticType);
1490e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticGetType(TS, TSBasicSymplecticType *);
1491e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegister(TSBasicSymplecticType, PetscInt, PetscInt, PetscReal[], PetscReal[]);
14924bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterAll(void);
1493e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticInitializePackage(void);
1494e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticFinalizePackage(void);
1495e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterDestroy(void);
14966bd7aeb5SHong Zhang 
1497db4653c2SDaniel Finn /*J
1498195e9b02SBarry Smith   TSDISCGRAD - The Discrete Gradient integrator is a timestepper for Hamiltonian systems designed to conserve the first integral (energy),
149987497f52SBarry Smith   but also has the property for some systems of monotonicity in a functional.
1500db4653c2SDaniel Finn 
1501db4653c2SDaniel Finn   Level: beginner
1502db4653c2SDaniel Finn 
1503f940b0e3Sdanofinn .seealso: [](ch_ts), `TS`, TSDiscGradSetFormulation()`, `TSDiscGradGetFormulation()`, `TSDiscGradSetType()`, `TSDiscGradGetType()`
1504db4653c2SDaniel Finn J*/
1505f940b0e3Sdanofinn typedef enum {
1506f940b0e3Sdanofinn   TS_DG_GONZALEZ,
1507f940b0e3Sdanofinn   TS_DG_AVERAGE,
1508f940b0e3Sdanofinn   TS_DG_NONE
1509f940b0e3Sdanofinn } TSDGType;
151040be0ff1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSDiscGradSetFormulation(TS, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (*)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (*)(TS, PetscReal, Vec, Vec, void *), void *);
15114bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSDiscGradGetFormulation(TS, PetscErrorCode (**)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (**)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (**)(TS, PetscReal, Vec, Vec, void *), void *);
1512f940b0e3Sdanofinn PETSC_EXTERN PetscErrorCode TSDiscGradSetType(TS, TSDGType);
1513f940b0e3Sdanofinn PETSC_EXTERN PetscErrorCode TSDiscGradGetType(TS, TSDGType *);
151440be0ff1SMatthew G. Knepley 
151583e2fdc7SBarry Smith /*
15160f3b3ca1SHong Zhang        PETSc interface to Sundials
151783e2fdc7SBarry Smith */
1518e808b789SPatrick Sanan #ifdef PETSC_HAVE_SUNDIALS2
15199371c9d4SSatish Balay typedef enum {
15209371c9d4SSatish Balay   SUNDIALS_ADAMS = 1,
15219371c9d4SSatish Balay   SUNDIALS_BDF   = 2
15229371c9d4SSatish Balay } TSSundialsLmmType;
15236a6fc655SJed Brown PETSC_EXTERN const char *const TSSundialsLmmTypes[];
15249371c9d4SSatish Balay typedef enum {
15259371c9d4SSatish Balay   SUNDIALS_MODIFIED_GS  = 1,
15269371c9d4SSatish Balay   SUNDIALS_CLASSICAL_GS = 2
15279371c9d4SSatish Balay } TSSundialsGramSchmidtType;
15286a6fc655SJed Brown PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
15294bf303faSJacob Faibussowitsch 
1530014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS, TSSundialsLmmType);
1531014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS, PC *);
1532014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS, PetscReal, PetscReal);
1533014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS, PetscReal);
1534014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS, PetscReal);
1535014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS, PetscInt *, PetscInt *);
1536014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS, TSSundialsGramSchmidtType);
1537014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS, PetscInt);
1538014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS, PetscReal);
1539014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS, PetscBool);
1540014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS, PetscInt);
1541c4e80e11SFlorian PETSC_EXTERN PetscErrorCode TSSundialsSetMaxord(TS, PetscInt);
1542918687b8SPatrick Sanan PETSC_EXTERN PetscErrorCode TSSundialsSetUseDense(TS, PetscBool);
15436dc17ff5SMatthew Knepley #endif
154483e2fdc7SBarry Smith 
1545014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS, PetscReal);
1546014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS, PetscReal *);
1547014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS, PetscBool *);
1548014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS, PetscBool);
15490de4c49aSJed Brown 
1550014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS, PetscReal);
1551014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS, PetscReal, PetscReal, PetscReal);
1552014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS, PetscReal *, PetscReal *, PetscReal *);
155388df8a41SLisandro Dalcin 
1554818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAlpha2SetRadius(TS, PetscReal);
1555818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAlpha2SetParams(TS, PetscReal, PetscReal, PetscReal, PetscReal);
1556818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAlpha2GetParams(TS, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
1557818efac9SLisandro Dalcin 
1558220f924aSDavid Kamensky /*S
15598434afd1SBarry Smith   TSAlpha2PredictorFn - A callback to set the predictor (i.e., the initial guess for the nonlinear solver) in
1560220f924aSDavid Kamensky   a second-order generalized-alpha time integrator.
1561220f924aSDavid Kamensky 
1562220f924aSDavid Kamensky   Calling Sequence:
1563220f924aSDavid Kamensky + ts   - the `TS` context obtained from `TSCreate()`
1564220f924aSDavid Kamensky . X0   - the previous time step's state vector
1565220f924aSDavid Kamensky . V0   - the previous time step's first derivative of the state vector
1566220f924aSDavid Kamensky . A0   - the previous time step's second derivative of the state vector
1567220f924aSDavid Kamensky . X1   - the vector into which the initial guess for the current time step will be written
1568220f924aSDavid Kamensky - ctx  - [optional] user-defined context for the predictor evaluation routine (may be `NULL`)
1569220f924aSDavid Kamensky 
1570220f924aSDavid Kamensky   Level: intermediate
1571220f924aSDavid Kamensky 
1572d1c5d1fcSBarry Smith   Note:
15738434afd1SBarry Smith   The deprecated `TSAlpha2Predictor` still works as a replacement for `TSAlpha2PredictorFn` *.
1574d1c5d1fcSBarry Smith 
1575220f924aSDavid Kamensky .seealso: [](ch_ts), `TS`, `TSAlpha2SetPredictor()`
1576220f924aSDavid Kamensky S*/
1577*2a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode TSAlpha2PredictorFn(TS ts, Vec X0, Vec V0, Vec A0, Vec X1, PetscCtx ctx);
1578d1c5d1fcSBarry Smith 
15798434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSAlpha2PredictorFn *TSAlpha2Predictor;
1580d1c5d1fcSBarry Smith 
1581ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode TSAlpha2SetPredictor(TS, TSAlpha2PredictorFn *, void *);
1582220f924aSDavid Kamensky 
1583014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetDM(TS, DM);
1584014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetDM(TS, DM *);
15856c699258SBarry Smith 
1586014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES, Vec, Vec, void *);
1587d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES, Vec, Mat, Mat, void *);
15880f5c6efeSJed Brown 
1589f3b1f45cSBarry Smith PETSC_EXTERN PetscErrorCode TSRHSJacobianTest(TS, PetscBool *);
1590f3b1f45cSBarry Smith PETSC_EXTERN PetscErrorCode TSRHSJacobianTestTranspose(TS, PetscBool *);
1591aad739acSMatthew G. Knepley 
15922e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSGetComputeInitialCondition(TS, PetscErrorCode (**)(TS, Vec));
15932e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSSetComputeInitialCondition(TS, PetscErrorCode (*)(TS, Vec));
15942e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSComputeInitialCondition(TS, Vec);
15952e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSGetComputeExactError(TS, PetscErrorCode (**)(TS, Vec, Vec));
1596aad739acSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSSetComputeExactError(TS, PetscErrorCode (*)(TS, Vec, Vec));
1597aad739acSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSComputeExactError(TS, Vec, Vec);
1598f2ed2dc7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscConvEstUseTS(PetscConvEst, PetscBool);
1599d60b7d5cSBarry Smith 
1600d60b7d5cSBarry Smith PETSC_EXTERN PetscErrorCode TSSetMatStructure(TS, MatStructure);
1601