xref: /petsc/include/petscts.h (revision 9dc50cb58a67b42abbba6aecdb4cbe9393c4ef0e)
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
15435da068SBarry Smith    TS - Abstract PETSc object that manages all time-steppers (ODE integrators)
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
2487497f52SBarry Smith    TSType - String with the name of a PETSc `TS` method.
25435da068SBarry Smith 
26435da068SBarry Smith    Level: beginner
27435da068SBarry Smith 
281cc06b55SBarry Smith .seealso: [](integrator_table), [](ch_ts), `TSSetType()`, `TS`, `TSRegister()`
2976bdecfbSBarry Smith J*/
3019fd82e9SBarry Smith typedef const char *TSType;
319596e0b4SJed Brown #define TSEULER           "euler"
329596e0b4SJed Brown #define TSBEULER          "beuler"
33e49d4f37SHong Zhang #define TSBASICSYMPLECTIC "basicsymplectic"
349596e0b4SJed Brown #define TSPSEUDO          "pseudo"
354d91e141SJed Brown #define TSCN              "cn"
369596e0b4SJed Brown #define TSSUNDIALS        "sundials"
374d91e141SJed Brown #define TSRK              "rk"
389596e0b4SJed Brown #define TSPYTHON          "python"
399596e0b4SJed Brown #define TSTHETA           "theta"
4088df8a41SLisandro Dalcin #define TSALPHA           "alpha"
41818efac9SLisandro Dalcin #define TSALPHA2          "alpha2"
4226d28e4eSEmil Constantinescu #define TSGLLE            "glle"
43b6a60446SDebojyoti Ghosh #define TSGLEE            "glee"
44ef7bb5aaSJed Brown #define TSSSP             "ssp"
458a381b04SJed Brown #define TSARKIMEX         "arkimex"
46e27a552bSJed Brown #define TSROSW            "rosw"
477d5697caSHong Zhang #define TSEIMEX           "eimex"
48abadfa0fSMatthew G. Knepley #define TSMIMEX           "mimex"
49211a84d6SLisandro Dalcin #define TSBDF             "bdf"
50d249a78fSBarry Smith #define TSRADAU5          "radau5"
514b84eec9SHong Zhang #define TSMPRK            "mprk"
5240be0ff1SMatthew G. Knepley #define TSDISCGRAD        "discgrad"
53d2567f34SHong Zhang #define TSIRK             "irk"
54d27334e2SStefano Zampini #define TSDIRK            "dirk"
5540be0ff1SMatthew G. Knepley 
56435da068SBarry Smith /*E
5787497f52SBarry Smith    TSProblemType - Determines the type of problem this `TS` object is to be used to solve
58435da068SBarry Smith 
59195e9b02SBarry Smith    Values:
60195e9b02SBarry Smith  + `TS_LINEAR`    - a linear ODE or DAE
61195e9b02SBarry Smith  - `TS_NONLINEAR` - a nonlinear ODE or DAE
62195e9b02SBarry Smith 
63435da068SBarry Smith    Level: beginner
64435da068SBarry Smith 
651cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSCreate()`
66435da068SBarry Smith E*/
679371c9d4SSatish Balay typedef enum {
689371c9d4SSatish Balay   TS_LINEAR,
699371c9d4SSatish Balay   TS_NONLINEAR
709371c9d4SSatish Balay } TSProblemType;
71818ad0c1SBarry Smith 
72b93fea0eSJed Brown /*E
7387497f52SBarry Smith    TSEquationType - type of `TS` problem that is solved
74e817cc15SEmil Constantinescu 
75195e9b02SBarry Smith    Values:
76195e9b02SBarry Smith +  `TS_EQ_UNSPECIFIED` - (default)
7716a05f60SBarry 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
7816a05f60SBarry Smith -  `TS_EQ_IMPLICIT`    - {ODE and DAE index 1, 2, 3, HI} F(t,U,U_t) = 0
79e817cc15SEmil Constantinescu 
8095bd0b28SBarry Smith    Level: beginner
8195bd0b28SBarry Smith 
821cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGetEquationType()`, `TSSetEquationType()`
83e817cc15SEmil Constantinescu E*/
84e817cc15SEmil Constantinescu typedef enum {
85e817cc15SEmil Constantinescu   TS_EQ_UNSPECIFIED               = -1,
86e817cc15SEmil Constantinescu   TS_EQ_EXPLICIT                  = 0,
87e817cc15SEmil Constantinescu   TS_EQ_ODE_EXPLICIT              = 1,
88e817cc15SEmil Constantinescu   TS_EQ_DAE_SEMI_EXPLICIT_INDEX1  = 100,
89e817cc15SEmil Constantinescu   TS_EQ_DAE_SEMI_EXPLICIT_INDEX2  = 200,
90e817cc15SEmil Constantinescu   TS_EQ_DAE_SEMI_EXPLICIT_INDEX3  = 300,
91e817cc15SEmil Constantinescu   TS_EQ_DAE_SEMI_EXPLICIT_INDEXHI = 500,
92e817cc15SEmil Constantinescu   TS_EQ_IMPLICIT                  = 1000,
93e817cc15SEmil Constantinescu   TS_EQ_ODE_IMPLICIT              = 1001,
94e817cc15SEmil Constantinescu   TS_EQ_DAE_IMPLICIT_INDEX1       = 1100,
95e817cc15SEmil Constantinescu   TS_EQ_DAE_IMPLICIT_INDEX2       = 1200,
96e817cc15SEmil Constantinescu   TS_EQ_DAE_IMPLICIT_INDEX3       = 1300,
9719436ca2SJed Brown   TS_EQ_DAE_IMPLICIT_INDEXHI      = 1500
98e817cc15SEmil Constantinescu } TSEquationType;
99e817cc15SEmil Constantinescu PETSC_EXTERN const char *const *TSEquationTypes;
100e817cc15SEmil Constantinescu 
101e817cc15SEmil Constantinescu /*E
10216a05f60SBarry Smith    TSConvergedReason - reason a `TS` method has converged (integrated to the requested time) or not
103b93fea0eSJed Brown 
104195e9b02SBarry Smith    Values:
105195e9b02SBarry Smith +  `TS_CONVERGED_ITERATING`          - this only occurs if `TSGetConvergedReason()` is called during the `TSSolve()`
106195e9b02SBarry Smith .  `TS_CONVERGED_TIME`               - the final time was reached
107195e9b02SBarry Smith .  `TS_CONVERGED_ITS`                - the maximum number of iterations (time-steps) was reached prior to the final time
108195e9b02SBarry Smith .  `TS_CONVERGED_USER`               - user requested termination
109195e9b02SBarry Smith .  `TS_CONVERGED_EVENT`              - user requested termination on event detection
110195e9b02SBarry Smith .  `TS_CONVERGED_PSEUDO_FATOL`       - stops when function norm decreased by a set amount, used only for `TSPSEUDO`
111195e9b02SBarry Smith .  `TS_CONVERGED_PSEUDO_FRTOL`       - stops when function norm decreases below a set amount, used only for `TSPSEUDO`
112195e9b02SBarry Smith .  `TS_DIVERGED_NONLINEAR_SOLVE`     - too many nonlinear solve failures have occurred
113195e9b02SBarry Smith .  `TS_DIVERGED_STEP_REJECTED`       - too many steps were rejected
114195e9b02SBarry Smith .  `TSFORWARD_DIVERGED_LINEAR_SOLVE` - tangent linear solve failed
115195e9b02SBarry Smith -  `TSADJOINT_DIVERGED_LINEAR_SOLVE` - transposed linear solve failed
116195e9b02SBarry Smith 
117b93fea0eSJed Brown    Level: beginner
118b93fea0eSJed Brown 
1191cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGetConvergedReason()`
120b93fea0eSJed Brown E*/
121193ac0bcSJed Brown typedef enum {
122193ac0bcSJed Brown   TS_CONVERGED_ITERATING          = 0,
123193ac0bcSJed Brown   TS_CONVERGED_TIME               = 1,
124193ac0bcSJed Brown   TS_CONVERGED_ITS                = 2,
125d6ad946cSShri Abhyankar   TS_CONVERGED_USER               = 3,
1262b7db910SShri Abhyankar   TS_CONVERGED_EVENT              = 4,
1273118ae5eSBarry Smith   TS_CONVERGED_PSEUDO_FATOL       = 5,
1283118ae5eSBarry Smith   TS_CONVERGED_PSEUDO_FRTOL       = 6,
129193ac0bcSJed Brown   TS_DIVERGED_NONLINEAR_SOLVE     = -1,
130b5b4017aSHong Zhang   TS_DIVERGED_STEP_REJECTED       = -2,
13105c9950eSHong Zhang   TSFORWARD_DIVERGED_LINEAR_SOLVE = -3,
13205c9950eSHong Zhang   TSADJOINT_DIVERGED_LINEAR_SOLVE = -4
133193ac0bcSJed Brown } TSConvergedReason;
134014dd563SJed Brown PETSC_EXTERN const char *const *TSConvergedReasons;
1351957e957SBarry Smith 
136b93fea0eSJed Brown /*MC
13787497f52SBarry Smith    TS_CONVERGED_ITERATING - this only occurs if `TSGetConvergedReason()` is called during the `TSSolve()`
138b93fea0eSJed Brown 
139b93fea0eSJed Brown    Level: beginner
140b93fea0eSJed Brown 
1411cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`
142b93fea0eSJed Brown M*/
143b93fea0eSJed Brown 
144b93fea0eSJed Brown /*MC
145b93fea0eSJed Brown    TS_CONVERGED_TIME - the final time was reached
146b93fea0eSJed Brown 
147b93fea0eSJed Brown    Level: beginner
148b93fea0eSJed Brown 
1491cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxTime()`, `TSGetMaxTime()`, `TSGetSolveTime()`
150b93fea0eSJed Brown M*/
151b93fea0eSJed Brown 
152b93fea0eSJed Brown /*MC
1531957e957SBarry Smith    TS_CONVERGED_ITS - the maximum number of iterations (time-steps) was reached prior to the final time
154b93fea0eSJed Brown 
155b93fea0eSJed Brown    Level: beginner
156b93fea0eSJed Brown 
1571cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxSteps()`, `TSGetMaxSteps()`
158b93fea0eSJed Brown M*/
1591957e957SBarry Smith 
160d6ad946cSShri Abhyankar /*MC
161d6ad946cSShri Abhyankar    TS_CONVERGED_USER - user requested termination
162d6ad946cSShri Abhyankar 
163d6ad946cSShri Abhyankar    Level: beginner
164d6ad946cSShri Abhyankar 
1651cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`
166b93fea0eSJed Brown M*/
167b93fea0eSJed Brown 
168b93fea0eSJed Brown /*MC
169d76fc68cSShri Abhyankar    TS_CONVERGED_EVENT - user requested termination on event detection
170d76fc68cSShri Abhyankar 
171d76fc68cSShri Abhyankar    Level: beginner
172d76fc68cSShri Abhyankar 
1731cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`
174d76fc68cSShri Abhyankar M*/
175d76fc68cSShri Abhyankar 
176d76fc68cSShri Abhyankar /*MC
17787497f52SBarry Smith    TS_CONVERGED_PSEUDO_FRTOL - stops when function norm decreased by a set amount, used only for `TSPSEUDO`
1783118ae5eSBarry Smith 
1793c7db156SBarry Smith    Options Database Key:
180d5b43468SJose E. Roman .   -ts_pseudo_frtol <rtol> - use specified rtol
1813118ae5eSBarry Smith 
182195e9b02SBarry Smith    Level: beginner
183195e9b02SBarry Smith 
1841cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`, `TS_CONVERGED_PSEUDO_FATOL`
1853118ae5eSBarry Smith M*/
1863118ae5eSBarry Smith 
1873118ae5eSBarry Smith /*MC
18887497f52SBarry Smith    TS_CONVERGED_PSEUDO_FATOL - stops when function norm decreases below a set amount, used only for `TSPSEUDO`
1893118ae5eSBarry Smith 
1903c7db156SBarry Smith    Options Database Key:
19167b8a455SSatish Balay .   -ts_pseudo_fatol <atol> - use specified atol
1923118ae5eSBarry Smith 
193195e9b02SBarry Smith    Level: beginner
194195e9b02SBarry Smith 
1951cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`, `TS_CONVERGED_PSEUDO_FRTOL`
1963118ae5eSBarry Smith M*/
1973118ae5eSBarry Smith 
1983118ae5eSBarry Smith /*MC
199b93fea0eSJed Brown    TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed
200b93fea0eSJed Brown 
201b93fea0eSJed Brown    Level: beginner
202b93fea0eSJed Brown 
203195e9b02SBarry Smith    Note:
204195e9b02SBarry Smith    See `TSSetMaxSNESFailures()` for how to allow more nonlinear solver failures.
2051957e957SBarry Smith 
2061cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSGetSNES()`, `SNESGetConvergedReason()`, `TSSetMaxSNESFailures()`
207b93fea0eSJed Brown M*/
208b93fea0eSJed Brown 
209b93fea0eSJed Brown /*MC
210b93fea0eSJed Brown    TS_DIVERGED_STEP_REJECTED - too many steps were rejected
211b93fea0eSJed Brown 
212b93fea0eSJed Brown    Level: beginner
213b93fea0eSJed Brown 
21495452b02SPatrick Sanan    Notes:
215195e9b02SBarry Smith    See `TSSetMaxStepRejections()` for how to allow more step rejections.
2161957e957SBarry Smith 
2171cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxStepRejections()`
218b93fea0eSJed Brown M*/
219b93fea0eSJed Brown 
220d6ad946cSShri Abhyankar /*E
221d6ad946cSShri Abhyankar    TSExactFinalTimeOption - option for handling of final time step
222d6ad946cSShri Abhyankar 
223195e9b02SBarry Smith    Values:
22416a05f60SBarry Smith +  `TS_EXACTFINALTIME_STEPOVER`    - Don't do anything if requested final time is exceeded
225195e9b02SBarry Smith .  `TS_EXACTFINALTIME_INTERPOLATE` - Interpolate back to final time
22616a05f60SBarry Smith -  `TS_EXACTFINALTIME_MATCHSTEP`   - Adapt final time step to match the final time requested
227195e9b02SBarry Smith 
228d6ad946cSShri Abhyankar    Level: beginner
229d6ad946cSShri Abhyankar 
2301cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGetConvergedReason()`, `TSSetExactFinalTime()`, `TSGetExactFinalTime()`
231d6ad946cSShri Abhyankar E*/
2329371c9d4SSatish Balay typedef enum {
2339371c9d4SSatish Balay   TS_EXACTFINALTIME_UNSPECIFIED = 0,
2349371c9d4SSatish Balay   TS_EXACTFINALTIME_STEPOVER    = 1,
2359371c9d4SSatish Balay   TS_EXACTFINALTIME_INTERPOLATE = 2,
2369371c9d4SSatish Balay   TS_EXACTFINALTIME_MATCHSTEP   = 3
2379371c9d4SSatish Balay } TSExactFinalTimeOption;
238d6ad946cSShri Abhyankar PETSC_EXTERN const char *const TSExactFinalTimeOptions[];
239d6ad946cSShri Abhyankar 
240000e7ae3SMatthew Knepley /* Logging support */
241014dd563SJed Brown PETSC_EXTERN PetscClassId TS_CLASSID;
242d74926cbSBarry Smith PETSC_EXTERN PetscClassId DMTS_CLASSID;
2431b9b13dfSLisandro Dalcin PETSC_EXTERN PetscClassId TSADAPT_CLASSID;
244000e7ae3SMatthew Knepley 
245607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSInitializePackage(void);
246560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSFinalizePackage(void);
2478ba1e511SMatthew Knepley 
248014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm, TS *);
249baa10174SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSClone(TS, TS *);
250014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSDestroy(TS *);
251818ad0c1SBarry Smith 
252014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetProblemType(TS, TSProblemType);
253014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetProblemType(TS, TSProblemType *);
254014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSMonitor(TS, PetscInt, PetscReal, Vec);
255014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSMonitorSet(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, void *), void *, PetscErrorCode (*)(void **));
256721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorSetFromOptions(TS, const char[], const char[], const char[], PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *), PetscErrorCode (*)(TS, PetscViewerAndFormat *));
257014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS);
258818ad0c1SBarry Smith 
259014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS, const char[]);
260014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS, const char[]);
261014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS, const char *[]);
262014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS);
263014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetUp(TS);
264014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSReset(TS);
265818ad0c1SBarry Smith 
266014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetSolution(TS, Vec);
267014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetSolution(TS, Vec *);
268818ad0c1SBarry Smith 
269efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TS2SetSolution(TS, Vec, Vec);
270efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TS2GetSolution(TS, Vec *, Vec *);
27157df6a1bSDebojyoti Ghosh 
272642f3702SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetSolutionComponents(TS, PetscInt *, Vec *);
273642f3702SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetAuxSolution(TS, Vec *);
2740a01e1b2SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetTimeError(TS, PetscInt, Vec *);
27557df6a1bSDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSSetTimeError(TS, Vec);
2766e2e232bSDebojyoti Ghosh 
277a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSSetRHSJacobianP(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), void *);
278cd4cee2dSHong Zhang PETSC_EXTERN PetscErrorCode TSGetRHSJacobianP(TS, Mat *, PetscErrorCode (**)(TS, PetscReal, Vec, Mat, void *), void **);
279a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianP(TS, PetscReal, Vec, Mat);
28033f72c5dSHong Zhang PETSC_EXTERN PetscErrorCode TSSetIJacobianP(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Vec, PetscReal, Mat, void *), void *);
28133f72c5dSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIJacobianP(TS, PetscReal, Vec, Vec, PetscReal, Mat, PetscBool);
282edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS() then TSComputeRHSJacobianP()", ) PetscErrorCode TSComputeDRDPFunction(TS, PetscReal, Vec, Vec *);
283edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS() then TSComputeRHSJacobian()", ) PetscErrorCode TSComputeDRDUFunction(TS, PetscReal, Vec, Vec *);
2849371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TSSetIHessianProduct(TS, Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), void *);
285b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUU(TS, PetscReal, Vec, Vec *, Vec, Vec *);
286b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUP(TS, PetscReal, Vec, Vec *, Vec, Vec *);
287b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPU(TS, PetscReal, Vec, Vec *, Vec, Vec *);
288b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPP(TS, PetscReal, Vec, Vec *, Vec, Vec *);
2899371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TSSetRHSHessianProduct(TS, Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), void *);
290b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS, PetscReal, Vec, Vec *, Vec, Vec *);
291b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS, PetscReal, Vec, Vec *, Vec, Vec *);
292b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS, PetscReal, Vec, Vec *, Vec, Vec *);
293b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS, PetscReal, Vec, Vec *, Vec, Vec *);
294b5b4017aSHong Zhang PETSC_EXTERN PetscErrorCode TSSetCostHessianProducts(TS, PetscInt, Vec *, Vec *, Vec);
295b5b4017aSHong Zhang PETSC_EXTERN PetscErrorCode TSGetCostHessianProducts(TS, PetscInt *, Vec **, Vec **, Vec *);
296ba0675f6SHong Zhang PETSC_EXTERN PetscErrorCode TSComputeSNESJacobian(TS, Vec, Mat, Mat);
297a05bf03eSHong Zhang 
298bc952696SBarry Smith /*S
299e91eccc2SStefano Zampini    TSTrajectory - Abstract PETSc object that stores the trajectory (solution of ODE/DAE at each time step)
300bc952696SBarry Smith 
301bc952696SBarry Smith    Level: advanced
302bc952696SBarry Smith 
3031cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectorySetType()`, `TSTrajectoryDestroy()`, `TSTrajectoryReset()`
304bc952696SBarry Smith S*/
305bc952696SBarry Smith typedef struct _p_TSTrajectory *TSTrajectory;
306bc952696SBarry Smith 
307bc952696SBarry Smith /*J
308195e9b02SBarry Smith    TSTrajectoryType - String with the name of a PETSc `TS` trajectory storage method
309bc952696SBarry Smith 
310bc952696SBarry Smith    Level: intermediate
311bc952696SBarry Smith 
3121cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectoryDestroy()`
313bc952696SBarry Smith J*/
314bc952696SBarry Smith typedef const char *TSTrajectoryType;
315bc952696SBarry Smith #define TSTRAJECTORYBASIC         "basic"
3161c8c567eSBarry Smith #define TSTRAJECTORYSINGLEFILE    "singlefile"
3179a53571cSHong Zhang #define TSTRAJECTORYMEMORY        "memory"
3182b043167SHong Zhang #define TSTRAJECTORYVISUALIZATION "visualization"
319bc952696SBarry Smith 
320bc952696SBarry Smith PETSC_EXTERN PetscFunctionList TSTrajectoryList;
321bc952696SBarry Smith PETSC_EXTERN PetscClassId      TSTRAJECTORY_CLASSID;
322bc952696SBarry Smith PETSC_EXTERN PetscBool         TSTrajectoryRegisterAllCalled;
323bc952696SBarry Smith 
324bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS);
3252d29f1f2SStefano Zampini PETSC_EXTERN PetscErrorCode TSResetTrajectory(TS);
32667a3cfb0SHong Zhang PETSC_EXTERN PetscErrorCode TSRemoveTrajectory(TS);
327bc952696SBarry Smith 
328bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectoryCreate(MPI_Comm, TSTrajectory *);
3299a992471SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryReset(TSTrajectory);
330bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy(TSTrajectory *);
331560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSTrajectoryView(TSTrajectory, PetscViewer);
332fd9d3c67SJed Brown PETSC_EXTERN PetscErrorCode TSTrajectorySetType(TSTrajectory, TS, TSTrajectoryType);
333881c1a9bSHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryGetType(TSTrajectory, TS, TSTrajectoryType *);
334bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySet(TSTrajectory, TS, PetscInt, PetscReal, Vec);
335c679fc15SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryGet(TSTrajectory, TS, PetscInt, PetscReal *);
336fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetVecs(TSTrajectory, TS, PetscInt, PetscReal *, Vec, Vec);
337fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetUpdatedHistoryVecs(TSTrajectory, TS, PetscReal, Vec *, Vec *);
338fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetNumSteps(TSTrajectory, PetscInt *);
339fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryRestoreUpdatedHistoryVecs(TSTrajectory, Vec *, Vec *);
340972caf09SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory, TS);
341560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSTrajectoryRegister(const char[], PetscErrorCode (*)(TSTrajectory, TS));
342bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectoryRegisterAll(void);
34368bece0bSHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetUp(TSTrajectory, TS);
3449ffb3502SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetUseHistory(TSTrajectory, PetscBool);
3452bee684fSHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetMonitor(TSTrajectory, PetscBool);
34678fbdcc8SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySetVariableNames(TSTrajectory, const char *const *);
34708347785SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySetTransform(TSTrajectory, PetscErrorCode (*)(void *, Vec, Vec *), PetscErrorCode (*)(void *), void *);
348fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectorySetSolutionOnly(TSTrajectory, PetscBool);
349fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetSolutionOnly(TSTrajectory, PetscBool *);
35064fc91eeSBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySetKeepFiles(TSTrajectory, PetscBool);
35164e38db7SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetDirname(TSTrajectory, const char[]);
35264e38db7SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetFiletemplate(TSTrajectory, const char[]);
35378fbdcc8SBarry Smith PETSC_EXTERN PetscErrorCode TSGetTrajectory(TS, TSTrajectory *);
354bc952696SBarry Smith 
3554bf303faSJacob Faibussowitsch typedef enum {
3564bf303faSJacob Faibussowitsch   TJ_REVOLVE,
3574bf303faSJacob Faibussowitsch   TJ_CAMS,
3584bf303faSJacob Faibussowitsch   TJ_PETSC
3594bf303faSJacob Faibussowitsch } TSTrajectoryMemoryType;
3604bf303faSJacob Faibussowitsch PETSC_EXTERN const char *const TSTrajectoryMemoryTypes[];
3614bf303faSJacob Faibussowitsch 
3624bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectoryMemorySetType(TSTrajectory, TSTrajectoryMemoryType);
3634bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxCpsRAM(TSTrajectory, PetscInt);
3644bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxCpsDisk(TSTrajectory, PetscInt);
3654bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxUnitsRAM(TSTrajectory, PetscInt);
3664bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxUnitsDisk(TSTrajectory, PetscInt);
3674bf303faSJacob Faibussowitsch 
368dfb21088SHong Zhang PETSC_EXTERN PetscErrorCode TSSetCostGradients(TS, PetscInt, Vec *, Vec *);
369dfb21088SHong Zhang PETSC_EXTERN PetscErrorCode TSGetCostGradients(TS, PetscInt *, Vec **, Vec **);
370edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 12, 0, "TSCreateQuadratureTS() and TSForwardSetSensitivities()", ) PetscErrorCode TSSetCostIntegrand(TS, PetscInt, Vec, PetscErrorCode (*)(TS, PetscReal, Vec, Vec, void *), PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, void *), PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, void *), PetscBool, void *);
371dfb21088SHong Zhang PETSC_EXTERN PetscErrorCode TSGetCostIntegral(TS, Vec *);
372022c081aSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeCostIntegrand(TS, PetscReal, Vec, Vec);
373cd4cee2dSHong Zhang PETSC_EXTERN PetscErrorCode TSCreateQuadratureTS(TS, PetscBool, TS *);
374cd4cee2dSHong Zhang PETSC_EXTERN PetscErrorCode TSGetQuadratureTS(TS, PetscBool *, TS *);
375d4aa0a58SBarry Smith 
376dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSetFromOptions(TS, PetscOptionItems *);
377a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitor(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *);
378a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitorSet(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *), void *, PetscErrorCode (*)(void **));
379a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitorCancel(TS);
380a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitorSetFromOptions(TS, const char[], const char[], const char[], PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, PetscViewerAndFormat *), PetscErrorCode (*)(TS, PetscViewerAndFormat *));
381a05bf03eSHong Zhang 
382edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSSetRHSJacobianP()", ) PetscErrorCode TSAdjointSetRHSJacobian(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), void *);
383edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSComputeRHSJacobianP()", ) PetscErrorCode TSAdjointComputeRHSJacobian(TS, PetscReal, Vec, Mat);
384edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS()", ) PetscErrorCode TSAdjointComputeDRDPFunction(TS, PetscReal, Vec, Vec *);
385edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS()", ) PetscErrorCode TSAdjointComputeDRDYFunction(TS, PetscReal, Vec, Vec *);
3862c39e106SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS);
38773b18844SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSetSteps(TS, PetscInt);
388d4aa0a58SBarry Smith 
38973b18844SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointStep(TS);
390f6a906c0SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS);
391ece46509SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointReset(TS);
392999a3721SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointCostIntegral(TS);
393ecf68647SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointSetForward(TS, Mat);
394ecf68647SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointResetForward(TS);
395715f1b00SHong Zhang 
39613b1b0ffSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardSetSensitivities(TS, PetscInt, Mat);
39713b1b0ffSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardGetSensitivities(TS, PetscInt *, Mat *);
398edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 12, 0, "TSCreateQuadratureTS()", ) PetscErrorCode TSForwardSetIntegralGradients(TS, PetscInt, Vec *);
399edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 12, 0, "TSForwardGetSensitivities()", ) PetscErrorCode TSForwardGetIntegralGradients(TS, PetscInt *, Vec **);
400715f1b00SHong Zhang PETSC_EXTERN PetscErrorCode TSForwardSetUp(TS);
4017adebddeSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardReset(TS);
402999a3721SHong Zhang PETSC_EXTERN PetscErrorCode TSForwardCostIntegral(TS);
403715f1b00SHong Zhang PETSC_EXTERN PetscErrorCode TSForwardStep(TS);
404b5b4017aSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardSetInitialSensitivities(TS, Mat);
405cc4f23bcSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardGetStages(TS, PetscInt *, Mat *[]);
406c235aa8dSHong Zhang 
407618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetMaxSteps(TS, PetscInt);
408618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetMaxSteps(TS, PetscInt *);
409618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetMaxTime(TS, PetscReal);
410618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetMaxTime(TS, PetscReal *);
41149354f04SShri Abhyankar PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS, TSExactFinalTimeOption);
412f6953c82SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetExactFinalTime(TS, TSExactFinalTimeOption *);
4134a658b32SHong Zhang PETSC_EXTERN PetscErrorCode TSSetTimeSpan(TS, PetscInt, PetscReal *);
4144a658b32SHong Zhang PETSC_EXTERN PetscErrorCode TSGetTimeSpan(TS, PetscInt *, const PetscReal **);
4154a658b32SHong Zhang PETSC_EXTERN PetscErrorCode TSGetTimeSpanSolutions(TS, PetscInt *, Vec **);
416818ad0c1SBarry Smith 
417edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSSetTime()", ) PetscErrorCode TSSetInitialTimeStep(TS, PetscReal, PetscReal);
418edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSSetMax()", ) PetscErrorCode TSSetDuration(TS, PetscInt, PetscReal);
419edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetMax()", ) PetscErrorCode TSGetDuration(TS, PetscInt *, PetscReal *);
420edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetStepNumber()", ) PetscErrorCode TSGetTimeStepNumber(TS, PetscInt *);
421edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetStepNumber()", ) PetscErrorCode TSGetTotalSteps(TS, PetscInt *);
42219eac22cSLisandro Dalcin 
423721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
424cc9c3a59SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorExtreme(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
42583a4ac43SBarry Smith 
42683a4ac43SBarry Smith typedef struct _n_TSMonitorDrawCtx *TSMonitorDrawCtx;
42783a4ac43SBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorDrawCtx *);
42883a4ac43SBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *);
42983a4ac43SBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawSolution(TS, PetscInt, PetscReal, Vec, void *);
4302d5ee99bSBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawSolutionPhase(TS, PetscInt, PetscReal, Vec, void *);
43183a4ac43SBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawError(TS, PetscInt, PetscReal, Vec, void *);
4320ed3bfb6SBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawSolutionFunction(TS, PetscInt, PetscReal, Vec, void *);
43383a4ac43SBarry Smith 
434721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointMonitorDefault(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, PetscViewerAndFormat *);
4359110b2e7SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitorDrawSensi(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *);
4369110b2e7SHong Zhang 
437721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorSolution(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
438014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS, PetscInt, PetscReal, Vec, void *);
439014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void *);
440fb1732b5SBarry Smith 
441014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSStep(TS);
4427cbde773SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSEvaluateWLTE(TS, NormType, PetscInt *, PetscReal *);
443014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS, PetscInt, Vec, PetscBool *);
444cc708dedSBarry Smith PETSC_EXTERN PetscErrorCode TSSolve(TS, Vec);
445e817cc15SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetEquationType(TS, TSEquationType *);
446e817cc15SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSSetEquationType(TS, TSEquationType);
447014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS, TSConvergedReason *);
448d6ad946cSShri Abhyankar PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS, TSConvergedReason);
449487e0bb9SJed Brown PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS, PetscReal *);
4505ef26d82SJed Brown PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS, PetscInt *);
4515ef26d82SJed Brown PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS, PetscInt *);
452cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS, PetscInt *);
453cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS, PetscInt);
454cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS, PetscInt *);
455cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS, PetscInt);
456cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS, PetscBool);
457dcb233daSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRestartStep(TS);
45824655328SShri PETSC_EXTERN PetscErrorCode TSRollBack(TS);
459*9dc50cb5SStefano Zampini PETSC_EXTERN PetscErrorCode TSGetStepRollBack(TS, PetscBool *);
460d2daff3dSHong Zhang 
461cc4f23bcSHong Zhang PETSC_EXTERN PetscErrorCode TSGetStages(TS, PetscInt *, Vec *[]);
462818ad0c1SBarry Smith 
463014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetTime(TS, PetscReal *);
464014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetTime(TS, PetscReal);
465e5e524a1SHong Zhang PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS, PetscReal *);
46680275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS, PetscReal *);
46780275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS, PetscReal);
46880275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetStepNumber(TS, PetscInt *);
46980275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetStepNumber(TS, PetscInt);
470818ad0c1SBarry Smith 
47114d0ab18SJacob Faibussowitsch /*S
4728434afd1SBarry Smith   TSRHSFunctionFn - A prototype of a `TS` right-hand-side evaluation function that would be passed to `TSSetRHSFunction()`
47314d0ab18SJacob Faibussowitsch 
47414d0ab18SJacob Faibussowitsch   Calling Sequence:
47514d0ab18SJacob Faibussowitsch + ts  - timestep context
47614d0ab18SJacob Faibussowitsch . t   - current time
47714d0ab18SJacob Faibussowitsch . u   - input vector
47814d0ab18SJacob Faibussowitsch . F   - function vector
47914d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context
48014d0ab18SJacob Faibussowitsch 
48114d0ab18SJacob Faibussowitsch   Level: beginner
48214d0ab18SJacob Faibussowitsch 
483d1c5d1fcSBarry Smith   Note:
4848434afd1SBarry Smith   The deprecated `TSRHSFunction` still works as a replacement for `TSRHSFunctionFn` *.
485d1c5d1fcSBarry Smith 
4868434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetRHSFunction()`, `DMTSSetRHSFunction()`, `TSIFunctionFn`,
4878434afd1SBarry Smith `TSIJacobianFn`, `TSRHSJacobianFn`
48814d0ab18SJacob Faibussowitsch S*/
4898434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSRHSFunctionFn)(TS ts, PetscReal t, Vec u, Vec F, void *ctx);
490d1c5d1fcSBarry Smith 
4918434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSRHSFunctionFn *TSRHSFunction;
49214d0ab18SJacob Faibussowitsch 
49314d0ab18SJacob Faibussowitsch /*S
4948434afd1SBarry Smith   TSRHSJacobianFn - A prototype of a `TS` right-hand-side Jacobian evaluation function that would be passed to `TSSetRHSJacobian()`
49514d0ab18SJacob Faibussowitsch 
49614d0ab18SJacob Faibussowitsch   Calling Sequence:
49714d0ab18SJacob Faibussowitsch + ts   - the `TS` context obtained from `TSCreate()`
49814d0ab18SJacob Faibussowitsch . t    - current time
49914d0ab18SJacob Faibussowitsch . u    - input vector
50014d0ab18SJacob Faibussowitsch . Amat - (approximate) Jacobian matrix
501af27ebaaSBarry Smith . Pmat - matrix from which preconditioner is to be constructed (usually the same as `Amat`)
50214d0ab18SJacob Faibussowitsch - ctx  - [optional] user-defined context for matrix evaluation routine
50314d0ab18SJacob Faibussowitsch 
50414d0ab18SJacob Faibussowitsch   Level: beginner
50514d0ab18SJacob Faibussowitsch 
506d1c5d1fcSBarry Smith   Note:
5078434afd1SBarry Smith   The deprecated `TSRHSJacobian` still works as a replacement for `TSRHSJacobianFn` *.
508d1c5d1fcSBarry Smith 
5098434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetRHSJacobian()`, `DMTSSetRHSJacobian()`, `TSRHSFunctionFn`,
5108434afd1SBarry Smith `TSIFunctionFn`, `TSIJacobianFn`
51114d0ab18SJacob Faibussowitsch S*/
5128434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSRHSJacobianFn)(TS ts, PetscReal t, Vec u, Mat Amat, Mat Pmat, void *ctx);
513d1c5d1fcSBarry Smith 
5148434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSRHSJacobianFn *TSRHSJacobian;
51514d0ab18SJacob Faibussowitsch 
51614d0ab18SJacob Faibussowitsch /*S
5178434afd1SBarry Smith   TSRHSJacobianPFn - A prototype of a function that computes the Jacobian of G w.r.t. the parameters P where
518af27ebaaSBarry Smith   U_t = G(U,P,t), as well as the location to store the matrix that would be passed to `TSSetRHSJacobianP()`
51914d0ab18SJacob Faibussowitsch 
52014d0ab18SJacob Faibussowitsch   Calling Sequence:
52114d0ab18SJacob Faibussowitsch + ts  - the `TS` context
52214d0ab18SJacob Faibussowitsch . t   - current timestep
52314d0ab18SJacob Faibussowitsch . U   - input vector (current ODE solution)
52414d0ab18SJacob Faibussowitsch . A   - output matrix
52514d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context
52614d0ab18SJacob Faibussowitsch 
52714d0ab18SJacob Faibussowitsch   Level: beginner
52814d0ab18SJacob Faibussowitsch 
529d1c5d1fcSBarry Smith   Note:
5308434afd1SBarry Smith   The deprecated `TSRHSJacobianP` still works as a replacement for `TSRHSJacobianPFn` *.
531d1c5d1fcSBarry Smith 
53214d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetRHSJacobianP()`, `TSGetRHSJacobianP()`
53314d0ab18SJacob Faibussowitsch S*/
5348434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSRHSJacobianPFn)(TS ts, PetscReal t, Vec U, Mat A, void *ctx);
53514d0ab18SJacob Faibussowitsch 
5368434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSRHSJacobianPFn *TSRHSJacobianP;
537d1c5d1fcSBarry Smith 
5388434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS, Vec, TSRHSFunctionFn *, void *);
5398434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS, Vec *, TSRHSFunctionFn **, void **);
5408434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS, Mat, Mat, TSRHSJacobianFn *, void *);
5418434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS, Mat *, Mat *, TSRHSJacobianFn **, void **);
542e1244c69SJed Brown PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS, PetscBool);
543818ad0c1SBarry Smith 
54414d0ab18SJacob Faibussowitsch /*S
5458434afd1SBarry Smith   TSSolutionFn - A prototype of a `TS` solution evaluation function that would be passed to `TSSetSolutionFunction()`
54614d0ab18SJacob Faibussowitsch 
54714d0ab18SJacob Faibussowitsch   Calling Sequence:
54814d0ab18SJacob Faibussowitsch + ts  - timestep context
54914d0ab18SJacob Faibussowitsch . t   - current time
55014d0ab18SJacob Faibussowitsch . u   - output vector
55114d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context
55214d0ab18SJacob Faibussowitsch 
55314d0ab18SJacob Faibussowitsch   Level: advanced
55414d0ab18SJacob Faibussowitsch 
555d1c5d1fcSBarry Smith   Note:
5568434afd1SBarry Smith   The deprecated `TSSolutionFunction` still works as a replacement for `TSSolutionFn` *.
557d1c5d1fcSBarry Smith 
55814d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetSolutionFunction()`, `DMTSSetSolutionFunction()`
55914d0ab18SJacob Faibussowitsch S*/
5608434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSSolutionFn)(TS ts, PetscReal t, Vec u, void *ctx);
56114d0ab18SJacob Faibussowitsch 
5628434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSSolutionFn *TSSolutionFunction;
563d1c5d1fcSBarry Smith 
5648434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS, TSSolutionFn *, void *);
56514d0ab18SJacob Faibussowitsch 
56614d0ab18SJacob Faibussowitsch /*S
5678434afd1SBarry Smith   TSForcingFn - A prototype of a `TS` forcing function evaluation function that would be passed to `TSSetForcingFunction()`
56814d0ab18SJacob Faibussowitsch 
56914d0ab18SJacob Faibussowitsch   Calling Sequence:
57014d0ab18SJacob Faibussowitsch + ts  - timestep context
57114d0ab18SJacob Faibussowitsch . t   - current time
57214d0ab18SJacob Faibussowitsch . f   - output vector
57314d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context
57414d0ab18SJacob Faibussowitsch 
57514d0ab18SJacob Faibussowitsch   Level: advanced
57614d0ab18SJacob Faibussowitsch 
577d1c5d1fcSBarry Smith   Note:
5788434afd1SBarry Smith   The deprecated `TSForcingFunction` still works as a replacement for `TSForcingFn` *.
579d1c5d1fcSBarry Smith 
58014d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetForcingFunction()`, `DMTSSetForcingFunction()`
58114d0ab18SJacob Faibussowitsch S*/
5828434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSForcingFn)(TS ts, PetscReal t, Vec f, void *ctx);
58314d0ab18SJacob Faibussowitsch 
5848434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSForcingFn *TSForcingFunction;
585d1c5d1fcSBarry Smith 
5868434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS, TSForcingFn *, void *);
587ef20d060SBarry Smith 
58814d0ab18SJacob Faibussowitsch /*S
5898434afd1SBarry Smith   TSIFunctionFn - A prototype of a `TS` implicit function evaluation function that would be passed to `TSSetIFunction()
59014d0ab18SJacob Faibussowitsch 
59114d0ab18SJacob Faibussowitsch   Calling Sequence:
59214d0ab18SJacob Faibussowitsch + ts  - the `TS` context obtained from `TSCreate()`
59314d0ab18SJacob Faibussowitsch . t   - time at step/stage being solved
59414d0ab18SJacob Faibussowitsch . U   - state vector
59514d0ab18SJacob Faibussowitsch . U_t - time derivative of state vector
59614d0ab18SJacob Faibussowitsch . F   - function vector
59714d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined context for function
59814d0ab18SJacob Faibussowitsch 
59914d0ab18SJacob Faibussowitsch   Level: beginner
60014d0ab18SJacob Faibussowitsch 
601d1c5d1fcSBarry Smith   Note:
6028434afd1SBarry Smith   The deprecated `TSIFunction` still works as a replacement for `TSIFunctionFn` *.
603d1c5d1fcSBarry Smith 
6048434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetIFunction()`, `DMTSSetIFunction()`, `TSIJacobianFn`, `TSRHSFunctionFn`, `TSRHSJacobianFn`
60514d0ab18SJacob Faibussowitsch S*/
6068434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSIFunctionFn)(TS ts, PetscReal t, Vec U, Vec U_t, Vec F, void *ctx);
607d1c5d1fcSBarry Smith 
6088434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSIFunctionFn *TSIFunction;
60914d0ab18SJacob Faibussowitsch 
61014d0ab18SJacob Faibussowitsch /*S
6118434afd1SBarry Smith   TSIJacobianFn - A prototype of a `TS` Jacobian evaluation function that would be passed to `TSSetIJacobian()`
61214d0ab18SJacob Faibussowitsch 
61314d0ab18SJacob Faibussowitsch   Calling Sequence:
61414d0ab18SJacob Faibussowitsch + ts   - the `TS` context obtained from `TSCreate()`
61514d0ab18SJacob Faibussowitsch . t    - time at step/stage being solved
61614d0ab18SJacob Faibussowitsch . U    - state vector
61714d0ab18SJacob Faibussowitsch . U_t  - time derivative of state vector
61814d0ab18SJacob Faibussowitsch . a    - shift
61914d0ab18SJacob Faibussowitsch . Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
62014d0ab18SJacob Faibussowitsch . Pmat - matrix used for constructing preconditioner, usually the same as `Amat`
62114d0ab18SJacob Faibussowitsch - ctx  - [optional] user-defined context for Jacobian evaluation routine
62214d0ab18SJacob Faibussowitsch 
62314d0ab18SJacob Faibussowitsch   Level: beginner
62414d0ab18SJacob Faibussowitsch 
625d1c5d1fcSBarry Smith   Note:
6268434afd1SBarry Smith   The deprecated `TSIJacobian` still works as a replacement for `TSIJacobianFn` *.
62714d0ab18SJacob Faibussowitsch 
6288434afd1SBarry Smith .seealso: [](ch_ts), `TSSetIJacobian()`, `DMTSSetIJacobian()`, `TSIFunctionFn`, `TSRHSFunctionFn`, `TSRHSJacobianFn`
629d1c5d1fcSBarry Smith S*/
6308434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSIJacobianFn)(TS ts, PetscReal t, Vec U, Vec U_t, PetscReal a, Mat Amat, Mat Pmat, void *ctx);
631d1c5d1fcSBarry Smith 
6328434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSIJacobianFn *TSIJacobian;
633d1c5d1fcSBarry Smith 
6348434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetIFunction(TS, Vec, TSIFunctionFn *, void *);
6358434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGetIFunction(TS, Vec *, TSIFunctionFn **, void **);
6368434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS, Mat, Mat, TSIJacobianFn *, void *);
6378434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS, Mat *, Mat *, TSIJacobianFn **, void **);
638316643e7SJed Brown 
63914d0ab18SJacob Faibussowitsch /*S
6408434afd1SBarry Smith   TSI2FunctionFn - A prototype of a `TS` implicit function evaluation function for 2nd order systems that would be passed to `TSSetI2Function()`
64114d0ab18SJacob Faibussowitsch 
64214d0ab18SJacob Faibussowitsch   Calling Sequence:
64314d0ab18SJacob Faibussowitsch + ts   - the `TS` context obtained from `TSCreate()`
64414d0ab18SJacob Faibussowitsch . t    - time at step/stage being solved
64514d0ab18SJacob Faibussowitsch . U    - state vector
64614d0ab18SJacob Faibussowitsch . U_t  - time derivative of state vector
64714d0ab18SJacob Faibussowitsch . U_tt - second time derivative of state vector
64814d0ab18SJacob Faibussowitsch . F    - function vector
64914d0ab18SJacob Faibussowitsch - ctx  - [optional] user-defined context for matrix evaluation routine (may be `NULL`)
65014d0ab18SJacob Faibussowitsch 
65114d0ab18SJacob Faibussowitsch   Level: advanced
65214d0ab18SJacob Faibussowitsch 
653d1c5d1fcSBarry Smith   Note:
6548434afd1SBarry Smith   The deprecated `TSI2Function` still works as a replacement for `TSI2FunctionFn` *.
655d1c5d1fcSBarry Smith 
6568434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetI2Function()`, `DMTSSetI2Function()`, `TSIFunctionFn`
65714d0ab18SJacob Faibussowitsch S*/
6588434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSI2FunctionFn)(TS ts, PetscReal t, Vec U, Vec U_t, Vec U_tt, Vec F, void *ctx);
659d1c5d1fcSBarry Smith 
6608434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSI2FunctionFn *TSI2Function;
66114d0ab18SJacob Faibussowitsch 
66214d0ab18SJacob Faibussowitsch /*S
6638434afd1SBarry Smith   TSI2JacobianFn - A prototype of a `TS` implicit Jacobian evaluation function for 2nd order systems that would be passed to `TSSetI2Jacobian()`
66414d0ab18SJacob Faibussowitsch 
66514d0ab18SJacob Faibussowitsch   Calling Sequence:
66614d0ab18SJacob Faibussowitsch + ts   - the `TS` context obtained from `TSCreate()`
66714d0ab18SJacob Faibussowitsch . t    - time at step/stage being solved
66814d0ab18SJacob Faibussowitsch . U    - state vector
66914d0ab18SJacob Faibussowitsch . U_t  - time derivative of state vector
67014d0ab18SJacob Faibussowitsch . U_tt - second time derivative of state vector
67114d0ab18SJacob Faibussowitsch . v    - shift for U_t
67214d0ab18SJacob Faibussowitsch . a    - shift for U_tt
67314d0ab18SJacob 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
67414d0ab18SJacob Faibussowitsch . jac  - matrix from which to construct the preconditioner, may be same as `J`
67514d0ab18SJacob Faibussowitsch - ctx  - [optional] user-defined context for matrix evaluation routine
67614d0ab18SJacob Faibussowitsch 
67714d0ab18SJacob Faibussowitsch   Level: advanced
67814d0ab18SJacob Faibussowitsch 
679d1c5d1fcSBarry Smith   Note:
6808434afd1SBarry Smith   The deprecated `TSI2Jacobian` still works as a replacement for `TSI2JacobianFn` *.
68114d0ab18SJacob Faibussowitsch 
6828434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetI2Jacobian()`, `DMTSSetI2Jacobian()`, `TSIFunctionFn`, `TSIJacobianFn`, `TSRHSFunctionFn`, `TSRHSJacobianFn`
683d1c5d1fcSBarry Smith S*/
6848434afd1SBarry 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, void *ctx);
685d1c5d1fcSBarry Smith 
6868434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSI2JacobianFn *TSI2Jacobian;
687d1c5d1fcSBarry Smith 
6888434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetI2Function(TS, Vec, TSI2FunctionFn *, void *);
6898434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGetI2Function(TS, Vec *, TSI2FunctionFn **, void **);
6908434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetI2Jacobian(TS, Mat, Mat, TSI2JacobianFn *, void *);
6918434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGetI2Jacobian(TS, Mat *, Mat *, TSI2JacobianFn **, void **);
692efe9872eSLisandro Dalcin 
693545aaa6fSHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitSetIS(TS, const char[], IS);
694545aaa6fSHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetIS(TS, const char[], IS *);
6958434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSRHSSplitSetRHSFunction(TS, const char[], Vec, TSRHSFunctionFn *, void *);
6961d06f6b3SHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTS(TS, const char[], TS *);
6971d06f6b3SHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTSs(TS, PetscInt *, TS *[]);
6980fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSSetUseSplitRHSFunction(TS, PetscBool);
6990fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSGetUseSplitRHSFunction(TS, PetscBool *);
700d9194312SHong Zhang 
7018434afd1SBarry Smith PETSC_EXTERN TSRHSFunctionFn TSComputeRHSFunctionLinear;
7028434afd1SBarry Smith PETSC_EXTERN TSRHSJacobianFn TSComputeRHSJacobianConstant;
703014dd563SJed Brown PETSC_EXTERN PetscErrorCode  TSComputeIFunctionLinear(TS, PetscReal, Vec, Vec, Vec, void *);
704d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode  TSComputeIJacobianConstant(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
7051c8a10a0SJed Brown PETSC_EXTERN PetscErrorCode  TSComputeSolutionFunction(TS, PetscReal, Vec);
7069b7cd975SBarry Smith PETSC_EXTERN PetscErrorCode  TSComputeForcingFunction(TS, PetscReal, Vec);
707847ff0e1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode  TSComputeIJacobianDefaultColor(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
708d7cfae9bSHong Zhang PETSC_EXTERN PetscErrorCode  TSPruneIJacobianColor(TS, Mat, Mat);
709e34be4c2SBarry Smith 
710014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
711b8123daeSJed Brown PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS, PetscReal));
7129be3e283SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS, PetscReal, PetscInt, Vec *));
713c688d042SShri Abhyankar PETSC_EXTERN PetscErrorCode TSSetPostEvaluate(TS, PetscErrorCode (*)(TS));
714014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
7156bd3a4fdSStefano Zampini PETSC_EXTERN PetscErrorCode TSSetResize(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscBool *, void *), PetscErrorCode (*)(TS, PetscInt, Vec[], Vec[], void *), void *);
716014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPreStep(TS);
717b8123daeSJed Brown PETSC_EXTERN PetscErrorCode TSPreStage(TS, PetscReal);
7189be3e283SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSPostStage(TS, PetscReal, PetscInt, Vec *);
719c688d042SShri Abhyankar PETSC_EXTERN PetscErrorCode TSPostEvaluate(TS);
720014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPostStep(TS);
7216bd3a4fdSStefano Zampini PETSC_EXTERN PetscErrorCode TSResize(TS);
7226bd3a4fdSStefano Zampini PETSC_EXTERN PetscErrorCode TSResizeRetrieveVec(TS, const char *, Vec *);
7236bd3a4fdSStefano Zampini PETSC_EXTERN PetscErrorCode TSResizeRegisterVec(TS, const char *, Vec);
7246bd3a4fdSStefano Zampini 
725014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSInterpolate(TS, PetscReal, Vec);
726014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetTolerances(TS, PetscReal, Vec, PetscReal, Vec);
727c5033834SJed Brown PETSC_EXTERN PetscErrorCode TSGetTolerances(TS, PetscReal *, Vec *, PetscReal *, Vec *);
7287453f775SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *);
7298a175baeSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm(TS, Vec, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *);
730014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS, PetscReal);
731014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS, PetscReal *);
732d183316bSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSSetFunctionDomainError(TS, PetscErrorCode (*)(TS, PetscReal, Vec, PetscBool *));
733d183316bSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSFunctionDomainError(TS, PetscReal, Vec, PetscBool *);
734000e7ae3SMatthew Knepley 
735014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS, PetscErrorCode (*)(TS, PetscReal *, void *), void *);
7368d359177SBarry Smith PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS, PetscReal *, void *);
737014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS, PetscReal *);
738014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS, PetscReal);
739014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS, PetscErrorCode (*)(TS, Vec, void *, PetscReal *, PetscBool *), void *);
7408d359177SBarry Smith PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS, Vec, void *, PetscReal *, PetscBool *);
741014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS, Vec, PetscReal *, PetscBool *);
742014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS, PetscReal);
743014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);
74421c89e3eSBarry Smith 
745014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPythonSetType(TS, const char[]);
746ebead697SStefano Zampini PETSC_EXTERN PetscErrorCode TSPythonGetType(TS, const char *[]);
7471d6018f0SLisandro Dalcin 
748014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS, PetscReal, Vec, Vec);
749d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS, PetscReal, Vec, Mat, Mat);
750014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS, PetscReal, Vec, Vec, Vec, PetscBool);
751d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscBool);
752efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSComputeI2Function(TS, PetscReal, Vec, Vec, Vec, Vec);
753efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSComputeI2Jacobian(TS, PetscReal, Vec, Vec, Vec, PetscReal, PetscReal, Mat, Mat);
754f9c1d6abSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS, PetscReal, PetscReal, PetscReal *, PetscReal *);
755818ad0c1SBarry Smith 
756014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS, Vec, Vec);
757d6ebe24aSShri Abhyankar 
758967bb25aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetBoundaryLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *);
7598434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM, TSRHSFunctionFn *, void *);
7608434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM, TSRHSFunctionFn **, void **);
761800f99ffSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionContextDestroy(DM, PetscErrorCode (*)(void *));
7628434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM, TSRHSJacobianFn *, void *);
7638434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM, TSRHSJacobianFn **, void **);
764800f99ffSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobianContextDestroy(DM, PetscErrorCode (*)(void *));
7658434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM, TSIFunctionFn *, void *);
7668434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM, TSIFunctionFn **, void **);
767800f99ffSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMTSSetIFunctionContextDestroy(DM, PetscErrorCode (*)(void *));
7688434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM, TSIJacobianFn *, void *);
7698434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM, TSIJacobianFn **, void **);
770800f99ffSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMTSSetIJacobianContextDestroy(DM, PetscErrorCode (*)(void *));
7718434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetI2Function(DM, TSI2FunctionFn *, void *);
7728434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetI2Function(DM, TSI2FunctionFn **, void **);
773800f99ffSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMTSSetI2FunctionContextDestroy(DM, PetscErrorCode (*)(void *));
7748434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetI2Jacobian(DM, TSI2JacobianFn *, void *);
7758434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetI2Jacobian(DM, TSI2JacobianFn **, void **);
776800f99ffSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMTSSetI2JacobianContextDestroy(DM, PetscErrorCode (*)(void *));
777efe9872eSLisandro Dalcin 
77814d0ab18SJacob Faibussowitsch /*S
7798434afd1SBarry Smith   TSTransientVariableFn - A prototype of a function to transform from state to transient variables that would be passed to `TSSetTransientVariable()`
78014d0ab18SJacob Faibussowitsch 
78114d0ab18SJacob Faibussowitsch   Calling Sequence:
78214d0ab18SJacob Faibussowitsch + ts  - timestep context
78314d0ab18SJacob Faibussowitsch . p   - input vector (primitive form)
78414d0ab18SJacob Faibussowitsch . c   - output vector, transient variables (conservative form)
78514d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context
78614d0ab18SJacob Faibussowitsch 
78714d0ab18SJacob Faibussowitsch   Level: advanced
78814d0ab18SJacob Faibussowitsch 
789d1c5d1fcSBarry Smith   Note:
7908434afd1SBarry Smith   The deprecated `TSTransientVariable` still works as a replacement for `TSTransientVariableFn` *.
791d1c5d1fcSBarry Smith 
79214d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetTransientVariable()`, `DMTSSetTransientVariable()`
79314d0ab18SJacob Faibussowitsch S*/
7948434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSTransientVariableFn)(TS ts, Vec p, Vec c, void *ctx);
79514d0ab18SJacob Faibussowitsch 
7968434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSTransientVariableFn *TSTransientVariable;
797d1c5d1fcSBarry Smith 
7988434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetTransientVariable(TS, TSTransientVariableFn *, void *);
7998434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetTransientVariable(DM, TSTransientVariableFn *, void *);
8008434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetTransientVariable(DM, TSTransientVariableFn **, void *);
801e3c11fc1SJed Brown PETSC_EXTERN PetscErrorCode TSComputeTransientVariable(TS, Vec, Vec);
802e3c11fc1SJed Brown PETSC_EXTERN PetscErrorCode TSHasTransientVariable(TS, PetscBool *);
803e3c11fc1SJed Brown 
8048434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM, TSSolutionFn *, void *);
8058434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM, TSSolutionFn **, void **);
8068434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM, TSForcingFn *, void *);
8078434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM, TSForcingFn **, void **);
808e17bd7bbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCheckResidual(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscReal *);
809e17bd7bbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCheckJacobian(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscBool *, PetscReal *);
8107f96f943SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec);
8119b7cd975SBarry Smith 
8120c9409e4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSGetIFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, Vec, void *), void **);
813d433e6cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, Vec, void *), void *);
8140c9409e4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSGetIJacobianLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void **);
815d433e6cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void *);
8160c9409e4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSGetRHSFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, void *), void **);
817d433e6cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *);
818b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrix(DM);
819b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrixLumped(DM);
820b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSDestroyRHSMassMatrix(DM);
8216c6b9e74SPeter Brune 
8226c6b9e74SPeter Brune PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM, PetscErrorCode (*)(void *, PetscViewer), PetscErrorCode (*)(void **, PetscViewer));
8236c6b9e74SPeter Brune PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM, PetscErrorCode (*)(void *, PetscViewer), PetscErrorCode (*)(void **, PetscViewer));
8246c6b9e74SPeter Brune 
82514d0ab18SJacob Faibussowitsch /*S
8268434afd1SBarry Smith   DMDATSRHSFunctionLocalFn - A prototype of a local `TS` right hand side residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetRHSFunctionLocal()`
8276c6b9e74SPeter Brune 
82814d0ab18SJacob Faibussowitsch   Calling Sequence:
82914d0ab18SJacob Faibussowitsch + info - defines the subdomain to evaluate the residual on
83014d0ab18SJacob Faibussowitsch . t    - time at which to evaluate residual
83114d0ab18SJacob Faibussowitsch . x    - array of local state information
83214d0ab18SJacob Faibussowitsch . f    - output array of local residual information
83314d0ab18SJacob Faibussowitsch - ctx  - optional user context
83414d0ab18SJacob Faibussowitsch 
83514d0ab18SJacob Faibussowitsch   Level: beginner
83614d0ab18SJacob Faibussowitsch 
837d1c5d1fcSBarry Smith   Note:
8388434afd1SBarry Smith   The deprecated `DMDATSRHSFunctionLocal` still works as a replacement for `DMDATSRHSFunctionLocalFn` *.
839d1c5d1fcSBarry Smith 
8408434afd1SBarry Smith .seealso: `DMDA`, `DMDATSSetRHSFunctionLocal()`, `TSRHSFunctionFn`, `DMDATSRHSJacobianLocalFn`, `DMDATSIJacobianLocalFn`, `DMDATSIFunctionLocalFn`
84114d0ab18SJacob Faibussowitsch S*/
8428434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDATSRHSFunctionLocalFn)(DMDALocalInfo *info, PetscReal t, void *x, void *f, void *ctx);
843d1c5d1fcSBarry Smith 
8448434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef DMDATSRHSFunctionLocalFn *DMDATSRHSFunctionLocal;
84514d0ab18SJacob Faibussowitsch 
84614d0ab18SJacob Faibussowitsch /*S
8478434afd1SBarry Smith   DMDATSRHSJacobianLocalFn - A prototype of a local residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetRHSJacobianLocal()`
84814d0ab18SJacob Faibussowitsch 
84914d0ab18SJacob Faibussowitsch   Calling Sequence:
85014d0ab18SJacob Faibussowitsch + info - defines the subdomain to evaluate the residual on
85114d0ab18SJacob Faibussowitsch . t    - time at which to evaluate residual
85214d0ab18SJacob Faibussowitsch . x    - array of local state information
85314d0ab18SJacob Faibussowitsch . J    - Jacobian matrix
85414d0ab18SJacob Faibussowitsch . B    - matrix from which to construct the preconditioner; often same as `J`
85514d0ab18SJacob Faibussowitsch - ctx  - optional context
85614d0ab18SJacob Faibussowitsch 
85714d0ab18SJacob Faibussowitsch   Level: beginner
85814d0ab18SJacob Faibussowitsch 
859d1c5d1fcSBarry Smith   Note:
8608434afd1SBarry Smith   The deprecated `DMDATSRHSJacobianLocal` still works as a replacement for `DMDATSRHSJacobianLocalFn` *.
861d1c5d1fcSBarry Smith 
8628434afd1SBarry Smith .seealso: `DMDA`, `DMDATSSetRHSJacobianLocal()`, `TSRHSJacobianFn`, `DMDATSRHSFunctionLocalFn`, `DMDATSIJacobianLocalFn`, `DMDATSIFunctionLocalFn`
86314d0ab18SJacob Faibussowitsch S*/
8648434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDATSRHSJacobianLocalFn)(DMDALocalInfo *info, PetscReal t, void *x, Mat J, Mat B, void *ctx);
865d1c5d1fcSBarry Smith 
8668434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef DMDATSRHSJacobianLocalFn *DMDATSRHSJacobianLocal;
86714d0ab18SJacob Faibussowitsch 
86814d0ab18SJacob Faibussowitsch /*S
8698434afd1SBarry Smith   DMDATSIFunctionLocalFn - A prototype of a local residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetIFunctionLocal()`
87014d0ab18SJacob Faibussowitsch 
87114d0ab18SJacob Faibussowitsch   Calling Sequence:
87214d0ab18SJacob Faibussowitsch + info  - defines the subdomain to evaluate the residual on
87314d0ab18SJacob Faibussowitsch . t     - time at which to evaluate residual
87414d0ab18SJacob Faibussowitsch . x     - array of local state information
87514d0ab18SJacob Faibussowitsch . xdot  - array of local time derivative information
87614d0ab18SJacob Faibussowitsch . imode - output array of local function evaluation information
87714d0ab18SJacob Faibussowitsch - ctx   - optional context
87814d0ab18SJacob Faibussowitsch 
87914d0ab18SJacob Faibussowitsch   Level: beginner
88014d0ab18SJacob Faibussowitsch 
881d1c5d1fcSBarry Smith   Note:
8828434afd1SBarry Smith   The deprecated `DMDATSIFunctionLocal` still works as a replacement for `DMDATSIFunctionLocalFn` *.
883d1c5d1fcSBarry Smith 
8848434afd1SBarry Smith .seealso: `DMDA`, `DMDATSSetIFunctionLocal()`, `DMDATSIJacobianLocalFn`, `TSIFunctionFn`
88514d0ab18SJacob Faibussowitsch S*/
8868434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDATSIFunctionLocalFn)(DMDALocalInfo *info, PetscReal t, void *x, void *xdot, void *imode, void *ctx);
887d1c5d1fcSBarry Smith 
8888434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef DMDATSIFunctionLocalFn *DMDATSIFunctionLocal;
88914d0ab18SJacob Faibussowitsch 
89014d0ab18SJacob Faibussowitsch /*S
8918434afd1SBarry Smith   DMDATSIJacobianLocalFn - A prototype of a local residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetIJacobianLocal()`
89214d0ab18SJacob Faibussowitsch 
89314d0ab18SJacob Faibussowitsch   Calling Sequence:
89414d0ab18SJacob Faibussowitsch + info  - defines the subdomain to evaluate the residual on
89514d0ab18SJacob Faibussowitsch . t     - time at which to evaluate the jacobian
89614d0ab18SJacob Faibussowitsch . x     - array of local state information
89714d0ab18SJacob Faibussowitsch . xdot  - time derivative at this state
89814d0ab18SJacob Faibussowitsch . shift - see `TSSetIJacobian()` for the meaning of this parameter
89914d0ab18SJacob Faibussowitsch . J     - Jacobian matrix
90014d0ab18SJacob Faibussowitsch . B     - matrix from which to construct the preconditioner; often same as `J`
90114d0ab18SJacob Faibussowitsch - ctx   - optional context
90214d0ab18SJacob Faibussowitsch 
90314d0ab18SJacob Faibussowitsch   Level: beginner
90414d0ab18SJacob Faibussowitsch 
905d1c5d1fcSBarry Smith   Note:
9068434afd1SBarry Smith   The deprecated `DMDATSIJacobianLocal` still works as a replacement for `DMDATSIJacobianLocalFn` *.
90714d0ab18SJacob Faibussowitsch 
9088434afd1SBarry Smith .seealso: `DMDA` `DMDATSSetIJacobianLocal()`, `TSIJacobianFn`, `DMDATSIFunctionLocalFn`, `DMDATSRHSFunctionLocalFn`,  `DMDATSRHSJacobianlocal()`
909d1c5d1fcSBarry Smith S*/
9108434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDATSIJacobianLocalFn)(DMDALocalInfo *info, PetscReal t, void *x, void *xdot, PetscReal shift, Mat J, Mat B, void *ctx);
911d1c5d1fcSBarry Smith 
9128434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef DMDATSIJacobianLocalFn *DMDATSIJacobianLocal;
913d1c5d1fcSBarry Smith 
9148434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM, InsertMode, DMDATSRHSFunctionLocalFn *, void *);
9158434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM, DMDATSRHSJacobianLocalFn *, void *);
9168434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM, InsertMode, DMDATSIFunctionLocalFn *, void *);
9178434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM, DMDATSIJacobianLocalFn *, void *);
9186c6b9e74SPeter Brune 
919be1b0d75SMatthew G. Knepley typedef struct _n_TSMonitorLGCtx *TSMonitorLGCtx;
920d1212d36SBarry Smith typedef struct {
921d1212d36SBarry Smith   Vec            ray;
922d1212d36SBarry Smith   VecScatter     scatter;
923d1212d36SBarry Smith   PetscViewer    viewer;
92451b4a12fSMatthew G. Knepley   TSMonitorLGCtx lgctx;
925d1212d36SBarry Smith } TSMonitorDMDARayCtx;
926d1212d36SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void **);
927d1212d36SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS, PetscInt, PetscReal, Vec, void *);
92851b4a12fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS, PetscInt, PetscReal, Vec, void *);
929d1212d36SBarry Smith 
930bdad233fSMatthew Knepley /* Dynamic creation and loading functions */
931140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList TSList;
93219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode    TSGetType(TS, TSType *);
93319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode    TSSetType(TS, TSType);
934bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode    TSRegister(const char[], PetscErrorCode (*)(TS));
93530de9b25SBarry Smith 
936014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetSNES(TS, SNES *);
937deb2cd25SJed Brown PETSC_EXTERN PetscErrorCode TSSetSNES(TS, SNES);
938014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetKSP(TS, KSP *);
939818ad0c1SBarry Smith 
940014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSView(TS, PetscViewer);
94155849f57SBarry Smith PETSC_EXTERN PetscErrorCode TSLoad(TS, PetscViewer);
942fe2efc57SMark PETSC_EXTERN PetscErrorCode TSViewFromOptions(TS, PetscObject, const char[]);
943fe2efc57SMark PETSC_EXTERN PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory, PetscObject, const char[]);
94455849f57SBarry Smith 
94555849f57SBarry Smith #define TS_FILE_CLASSID 1211225
94621c89e3eSBarry Smith 
947014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS, void *);
948014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS, void *);
94921c89e3eSBarry Smith 
950a9f9c1f6SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtx *);
951a9f9c1f6SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *);
9524f09c107SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS, PetscInt, PetscReal, Vec, void *);
9534f09c107SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS, PetscInt, PetscReal, Vec, void *);
95431152f8aSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS, const char *const *);
955b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS, const char *const **);
956e673d494SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx, const char *const *);
957a66092f1SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS, const char *const *);
958a66092f1SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx, const char *const *);
9597684fa3eSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS, PetscErrorCode (*)(void *, Vec, Vec *), PetscErrorCode (*)(void *), void *);
9607684fa3eSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx, PetscErrorCode (*)(void *, Vec, Vec *), PetscErrorCode (*)(void *), void *);
9614f09c107SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS, PetscInt, PetscReal, Vec, void *);
962201da799SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS, PetscInt, PetscReal, Vec, void *);
963201da799SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS, PetscInt, PetscReal, Vec, void *);
964edbaebb3SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorError(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
965d0c080abSJoseph Pusztay PETSC_EXTERN PetscErrorCode TSDMSwarmMonitorMoments(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
966ef20d060SBarry Smith 
967e669de00SBarry Smith struct _n_TSMonitorLGCtxNetwork {
968e669de00SBarry Smith   PetscInt     nlg;
969e669de00SBarry Smith   PetscDrawLG *lg;
970e669de00SBarry Smith   PetscBool    semilogy;
971e669de00SBarry Smith   PetscInt     howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
972e669de00SBarry Smith };
973e669de00SBarry Smith typedef struct _n_TSMonitorLGCtxNetwork *TSMonitorLGCtxNetwork;
974e669de00SBarry Smith PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkDestroy(TSMonitorLGCtxNetwork *);
975e669de00SBarry Smith PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkCreate(TS, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtxNetwork *);
976e669de00SBarry Smith PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkSolution(TS, PetscInt, PetscReal, Vec, void *);
977e669de00SBarry Smith 
978b3d3934dSBarry Smith typedef struct _n_TSMonitorEnvelopeCtx *TSMonitorEnvelopeCtx;
979b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeCtxCreate(TS, TSMonitorEnvelopeCtx *);
980b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode             TSMonitorEnvelope(TS, PetscInt, PetscReal, Vec, void *);
981b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeGetBounds(TS, Vec *, Vec *);
982b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *);
983b3d3934dSBarry Smith 
9848189c53fSBarry Smith typedef struct _n_TSMonitorSPEigCtx *TSMonitorSPEigCtx;
9858189c53fSBarry Smith PETSC_EXTERN PetscErrorCode          TSMonitorSPEigCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorSPEigCtx *);
9868189c53fSBarry Smith PETSC_EXTERN PetscErrorCode          TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *);
9878189c53fSBarry Smith PETSC_EXTERN PetscErrorCode          TSMonitorSPEig(TS, PetscInt, PetscReal, Vec, void *);
9883914022bSBarry Smith 
9891b575b74SJoseph Pusztay typedef struct _n_TSMonitorSPCtx *TSMonitorSPCtx;
99060e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorSPCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscBool, PetscBool, TSMonitorSPCtx *);
9911b575b74SJoseph Pusztay PETSC_EXTERN PetscErrorCode       TSMonitorSPCtxDestroy(TSMonitorSPCtx *);
9920ec8ee2bSJoseph Pusztay PETSC_EXTERN PetscErrorCode       TSMonitorSPSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
9931b575b74SJoseph Pusztay 
99460e16b1bSMatthew G. Knepley typedef struct _n_TSMonitorHGCtx *TSMonitorHGCtx;
99560e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorHGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscInt, PetscBool, TSMonitorHGCtx *);
99660e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
99760e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorHGCtxDestroy(TSMonitorHGCtx *);
99860e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
99960e16b1bSMatthew G. Knepley 
1000ca4445c7SIlya Fursov PETSC_EXTERN PetscErrorCode TSSetEventHandler(TS, PetscInt, PetscInt[], PetscBool[], PetscErrorCode (*)(TS, PetscReal, Vec, PetscReal[], void *), PetscErrorCode (*)(TS, PetscInt, PetscInt[], PetscReal, Vec, PetscBool, void *), void *);
1001ca4445c7SIlya Fursov PETSC_EXTERN PetscErrorCode TSSetPostEventStep(TS, PetscReal);
1002fe4ad979SIlya Fursov PETSC_EXTERN PetscErrorCode TSSetPostEventSecondStep(TS, PetscReal);
1003fe4ad979SIlya Fursov PETSC_DEPRECATED_FUNCTION(3, 21, 0, "TSSetPostEventSecondStep()", ) static inline PetscErrorCode TSSetPostEventIntervalStep(TS ts, PetscReal dt)
1004fe4ad979SIlya Fursov {
1005fe4ad979SIlya Fursov   return TSSetPostEventSecondStep(ts, dt);
1006fe4ad979SIlya Fursov }
10076427ac75SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS, PetscReal, PetscReal[]);
10081ea83e56SMiguel PETSC_EXTERN PetscErrorCode TSGetNumEvents(TS, PetscInt *);
10091ea83e56SMiguel 
101076bdecfbSBarry Smith /*J
1011195e9b02SBarry Smith    TSSSPType - string with the name of a `TSSSP` scheme.
1012815f1ad5SJed Brown 
1013815f1ad5SJed Brown    Level: beginner
1014815f1ad5SJed Brown 
10151cc06b55SBarry Smith .seealso: [](ch_ts), `TSSSPSetType()`, `TS`, `TSSSP`
101676bdecfbSBarry Smith J*/
101719fd82e9SBarry Smith typedef const char *TSSSPType;
1018815f1ad5SJed Brown #define TSSSPRKS2  "rks2"
1019815f1ad5SJed Brown #define TSSSPRKS3  "rks3"
1020815f1ad5SJed Brown #define TSSSPRK104 "rk104"
1021815f1ad5SJed Brown 
102219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode    TSSSPSetType(TS, TSSSPType);
102319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode    TSSSPGetType(TS, TSSSPType *);
1024014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSSPSetNumStages(TS, PetscInt);
1025014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSSPGetNumStages(TS, PetscInt *);
1026787849ffSJed Brown PETSC_EXTERN PetscErrorCode    TSSSPInitializePackage(void);
1027560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode    TSSSPFinalizePackage(void);
1028013797aaSBarry Smith PETSC_EXTERN PetscFunctionList TSSSPList;
1029815f1ad5SJed Brown 
1030ed657a08SJed Brown /*S
103184df9cb4SJed Brown    TSAdapt - Abstract object that manages time-step adaptivity
103284df9cb4SJed Brown 
103384df9cb4SJed Brown    Level: beginner
103484df9cb4SJed Brown 
1035af27ebaaSBarry Smith .seealso: [](ch_ts), `TS`, `TSGetAdapt()`, `TSAdaptCreate()`, `TSAdaptType`
103684df9cb4SJed Brown S*/
103784df9cb4SJed Brown typedef struct _p_TSAdapt *TSAdapt;
103884df9cb4SJed Brown 
1039f0fc11ceSJed Brown /*J
104087497f52SBarry Smith    TSAdaptType - String with the name of `TSAdapt` scheme.
104184df9cb4SJed Brown 
104284df9cb4SJed Brown    Level: beginner
104384df9cb4SJed Brown 
1044af27ebaaSBarry Smith .seealso: [](ch_ts), `TSGetAdapt()`, `TSAdaptSetType()`, `TS`, `TSAdapt`
1045f0fc11ceSJed Brown J*/
1046f8963224SJed Brown typedef const char *TSAdaptType;
1047b0f836d7SJed Brown #define TSADAPTNONE    "none"
10481917a363SLisandro Dalcin #define TSADAPTBASIC   "basic"
1049cb7ab186SLisandro Dalcin #define TSADAPTDSP     "dsp"
10508d59e960SJed Brown #define TSADAPTCFL     "cfl"
10511917a363SLisandro Dalcin #define TSADAPTGLEE    "glee"
1052d949e4a4SStefano Zampini #define TSADAPTHISTORY "history"
105384df9cb4SJed Brown 
1054552698daSJed Brown PETSC_EXTERN PetscErrorCode TSGetAdapt(TS, TSAdapt *);
1055bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[], PetscErrorCode (*)(TSAdapt));
1056607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void);
1057014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
1058014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm, TSAdapt *);
105919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt, TSAdaptType);
1060d0288e4fSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetType(TSAdapt, TSAdaptType *);
1061014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt, const char[]);
1062014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
1063014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt, const char[], PetscInt, PetscInt, PetscReal, PetscReal, PetscBool);
1064014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt, PetscInt *, const PetscInt **, const PetscInt **, const PetscReal **, const PetscReal **);
1065014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt, TS, PetscReal, PetscInt *, PetscReal *, PetscBool *);
1066b295832fSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt, TS, PetscReal, Vec, PetscBool *);
1067014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt, PetscViewer);
1068ad6bc421SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt, PetscViewer);
1069dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(TSAdapt, PetscOptionItems *);
1070099af0a3SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt);
1071014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt *);
1072014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt, PetscBool);
1073bf997491SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt, PetscBool);
10741917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptSetSafety(TSAdapt, PetscReal, PetscReal);
10751917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetSafety(TSAdapt, PetscReal *, PetscReal *);
107676cddca1SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt, PetscReal);
107776cddca1SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt, PetscReal *);
10781917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptSetClip(TSAdapt, PetscReal, PetscReal);
10791917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetClip(TSAdapt, PetscReal *, PetscReal *);
108062c23b28SMark PETSC_EXTERN PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt, PetscReal);
108162c23b28SMark PETSC_EXTERN PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt, PetscReal *);
1082014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt, PetscReal, PetscReal);
10831917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetStepLimits(TSAdapt, PetscReal *, PetscReal *);
1084b295832fSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt, PetscErrorCode (*)(TSAdapt, TS, PetscReal, Vec, PetscBool *));
1085d949e4a4SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptHistorySetHistory(TSAdapt, PetscInt n, PetscReal hist[], PetscBool);
108675017583SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptHistorySetTrajectory(TSAdapt, TSTrajectory, PetscBool);
108775017583SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptHistoryGetStep(TSAdapt, PetscInt, PetscReal *, PetscReal *);
1088de50f1caSBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt, PetscInt);
1089cb7ab186SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptDSPSetFilter(TSAdapt, const char *);
1090cb7ab186SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptDSPSetPID(TSAdapt, PetscReal, PetscReal, PetscReal);
109184df9cb4SJed Brown 
109284df9cb4SJed Brown /*S
109387497f52SBarry Smith    TSGLLEAdapt - Abstract object that manages time-step adaptivity for `TSGLLE`
1094ed657a08SJed Brown 
1095ed657a08SJed Brown    Level: beginner
1096ed657a08SJed Brown 
1097195e9b02SBarry Smith    Developer Note:
109887497f52SBarry Smith    This functionality should be replaced by the `TSAdapt`.
109984df9cb4SJed Brown 
11001cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLEAdaptCreate()`, `TSGLLEAdaptType`
1101ed657a08SJed Brown S*/
110226d28e4eSEmil Constantinescu typedef struct _p_TSGLLEAdapt *TSGLLEAdapt;
1103ed657a08SJed Brown 
110476bdecfbSBarry Smith /*J
110587497f52SBarry Smith    TSGLLEAdaptType - String with the name of `TSGLLEAdapt` scheme
1106ed657a08SJed Brown 
1107ed657a08SJed Brown    Level: beginner
1108ed657a08SJed Brown 
1109195e9b02SBarry Smith    Developer Note:
1110195e9b02SBarry Smith    This functionality should be replaced by the `TSAdaptType`.
1111195e9b02SBarry Smith 
11121cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLLEAdaptSetType()`, `TS`
111376bdecfbSBarry Smith J*/
111426d28e4eSEmil Constantinescu typedef const char *TSGLLEAdaptType;
111526d28e4eSEmil Constantinescu #define TSGLLEADAPT_NONE "none"
111626d28e4eSEmil Constantinescu #define TSGLLEADAPT_SIZE "size"
111726d28e4eSEmil Constantinescu #define TSGLLEADAPT_BOTH "both"
1118ed657a08SJed Brown 
111926d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptRegister(const char[], PetscErrorCode (*)(TSGLLEAdapt));
112026d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptInitializePackage(void);
112126d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptFinalizePackage(void);
112226d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate(MPI_Comm, TSGLLEAdapt *);
112326d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt, TSGLLEAdaptType);
112426d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt, const char[]);
112526d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptChoose(TSGLLEAdapt, PetscInt, const PetscInt[], const PetscReal[], const PetscReal[], PetscInt, PetscReal, PetscReal, PetscInt *, PetscReal *, PetscBool *);
112626d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt, PetscViewer);
1127dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetFromOptions(TSGLLEAdapt, PetscOptionItems *);
112826d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt *);
1129ed657a08SJed Brown 
113076bdecfbSBarry Smith /*J
113187497f52SBarry Smith    TSGLLEAcceptType - String with the name of `TSGLLEAccept` scheme
1132ed657a08SJed Brown 
1133ed657a08SJed Brown    Level: beginner
1134ed657a08SJed Brown 
11351cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLLESetAcceptType()`, `TS`, `TSGLLEAccept`
113676bdecfbSBarry Smith J*/
113726d28e4eSEmil Constantinescu typedef const char *TSGLLEAcceptType;
113826d28e4eSEmil Constantinescu #define TSGLLEACCEPT_ALWAYS "always"
1139ed657a08SJed Brown 
1140d1c5d1fcSBarry Smith /*S
11418434afd1SBarry Smith   TSGLLEAcceptFn - A prototype of a `TS` accept function that would be passed to `TSGLLEAcceptRegister()`
1142d1c5d1fcSBarry Smith 
1143d1c5d1fcSBarry Smith   Calling Sequence:
1144d1c5d1fcSBarry Smith + ts  - timestep context
1145d1c5d1fcSBarry Smith . nt - time to end of solution time
1146d1c5d1fcSBarry Smith . h - the proposed step-size
1147d1c5d1fcSBarry Smith . enorm - unknown
1148d1c5d1fcSBarry Smith - accept - output, if the proposal is accepted
1149d1c5d1fcSBarry Smith 
1150d1c5d1fcSBarry Smith   Level: beginner
1151d1c5d1fcSBarry Smith 
1152d1c5d1fcSBarry Smith   Note:
11538434afd1SBarry Smith   The deprecated `TSGLLEAcceptFunction` still works as a replacement for `TSGLLEAcceptFn` *
1154d1c5d1fcSBarry Smith 
11558434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetRHSFunction()`, `DMTSSetRHSFunction()`, `TSIFunctionFn`,
11568434afd1SBarry Smith `TSIJacobianFn`, `TSRHSJacobianFn`, `TSGLLEAcceptRegister()`
1157d1c5d1fcSBarry Smith S*/
11588434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSGLLEAcceptFn)(TS ts, PetscReal nt, PetscReal h, const PetscReal enorm[], PetscBool *accept);
1159d1c5d1fcSBarry Smith 
11608434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSGLLEAcceptFn *TSGLLEAcceptFunction;
1161d1c5d1fcSBarry Smith 
11628434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGLLEAcceptRegister(const char[], TSGLLEAcceptFn *);
1163ed657a08SJed Brown 
116476bdecfbSBarry Smith /*J
1165195e9b02SBarry Smith   TSGLLEType - string with the name of a General Linear `TSGLLE` type
116618b56cb9SJed Brown 
116718b56cb9SJed Brown   Level: beginner
116818b56cb9SJed Brown 
11691cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLESetType()`, `TSGLLERegister()`, `TSGLLEAccept`
117076bdecfbSBarry Smith J*/
117126d28e4eSEmil Constantinescu typedef const char *TSGLLEType;
117226d28e4eSEmil Constantinescu #define TSGLLE_IRKS "irks"
117318b56cb9SJed Brown 
117426d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLERegister(const char[], PetscErrorCode (*)(TS));
117526d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEInitializePackage(void);
117626d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEFinalizePackage(void);
117726d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLESetType(TS, TSGLLEType);
117826d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEGetAdapt(TS, TSGLLEAdapt *);
117926d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLESetAcceptType(TS, TSGLLEAcceptType);
118018b56cb9SJed Brown 
1181020d8f30SJed Brown /*J
1182195e9b02SBarry Smith    TSEIMEXType - String with the name of an Extrapolated IMEX `TSEIMEX` type
11837d5697caSHong Zhang 
11847d5697caSHong Zhang    Level: beginner
11857d5697caSHong Zhang 
11861cc06b55SBarry Smith .seealso: [](ch_ts), `TSEIMEXSetType()`, `TS`, `TSEIMEX`, `TSEIMEXRegister()`
11877d5697caSHong Zhang J*/
11887d5697caSHong Zhang #define TSEIMEXType char *
11897d5697caSHong Zhang 
11907d5697caSHong Zhang PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts, PetscInt);
11917d5697caSHong Zhang PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts, PetscInt, PetscInt);
11927d5697caSHong Zhang PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS, PetscBool);
11937d5697caSHong Zhang 
11947d5697caSHong Zhang /*J
1195195e9b02SBarry Smith    TSRKType - String with the name of a Runge-Kutta `TSRK` type
1196f68a32c8SEmil Constantinescu 
1197f68a32c8SEmil Constantinescu    Level: beginner
1198f68a32c8SEmil Constantinescu 
11991cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSRKSetType()`, `TS`, `TSRK`, `TSRKRegister()`
1200f68a32c8SEmil Constantinescu J*/
1201f68a32c8SEmil Constantinescu typedef const char *TSRKType;
1202f68a32c8SEmil Constantinescu #define TSRK1FE "1fe"
1203fdefd5e5SDebojyoti Ghosh #define TSRK2A  "2a"
1204e7685601SHong Zhang #define TSRK2B  "2b"
1205f68a32c8SEmil Constantinescu #define TSRK3   "3"
1206fdefd5e5SDebojyoti Ghosh #define TSRK3BS "3bs"
1207f68a32c8SEmil Constantinescu #define TSRK4   "4"
1208fdefd5e5SDebojyoti Ghosh #define TSRK5F  "5f"
1209fdefd5e5SDebojyoti Ghosh #define TSRK5DP "5dp"
121005e23783SLisandro Dalcin #define TSRK5BS "5bs"
121137acaa02SHendrik Ranocha #define TSRK6VR "6vr"
121237acaa02SHendrik Ranocha #define TSRK7VR "7vr"
121337acaa02SHendrik Ranocha #define TSRK8VR "8vr"
121405e23783SLisandro Dalcin 
12152ea87ba9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRKGetOrder(TS, PetscInt *);
12160f7a28cdSStefano Zampini PETSC_EXTERN PetscErrorCode TSRKGetType(TS, TSRKType *);
12170f7a28cdSStefano Zampini PETSC_EXTERN PetscErrorCode TSRKSetType(TS, TSRKType);
12180f7a28cdSStefano Zampini PETSC_EXTERN PetscErrorCode TSRKGetTableau(TS, PetscInt *, const PetscReal **, const PetscReal **, const PetscReal **, const PetscReal **, PetscInt *, const PetscReal **, PetscBool *);
12190fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSRKSetMultirate(TS, PetscBool);
12200fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSRKGetMultirate(TS, PetscBool *);
1221f68a32c8SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
1222f68a32c8SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void);
1223560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void);
1224f68a32c8SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void);
1225f68a32c8SEmil Constantinescu 
1226f68a32c8SEmil Constantinescu /*J
1227af27ebaaSBarry Smith    TSMPRKType - String with the name of a partitioned Runge-Kutta `TSMPRK` type
12289f537349Sluzhanghpp 
12299f537349Sluzhanghpp    Level: beginner
12309f537349Sluzhanghpp 
12311cc06b55SBarry Smith .seealso: [](ch_ts), `TSMPRKSetType()`, `TS`, `TSMPRK`, `TSMPRKRegister()`
12329f537349Sluzhanghpp J*/
12334b84eec9SHong Zhang typedef const char *TSMPRKType;
123419c2959aSHong Zhang #define TSMPRK2A22 "2a22"
123519c2959aSHong Zhang #define TSMPRK2A23 "2a23"
123619c2959aSHong Zhang #define TSMPRK2A32 "2a32"
123719c2959aSHong Zhang #define TSMPRK2A33 "2a33"
12384b84eec9SHong Zhang #define TSMPRKP2   "p2"
12394b84eec9SHong Zhang #define TSMPRKP3   "p3"
12409f537349Sluzhanghpp 
12414b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKGetType(TS ts, TSMPRKType *);
12424b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKSetType(TS ts, TSMPRKType);
124379bc8a77SHong 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[]);
12444b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKInitializePackage(void);
12454b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKFinalizePackage(void);
12464b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKRegisterDestroy(void);
12479f537349Sluzhanghpp 
12489f537349Sluzhanghpp /*J
1249195e9b02SBarry Smith    TSIRKType - String with the name of an implicit Runge-Kutta `TSIRK` type
1250d2567f34SHong Zhang 
1251d2567f34SHong Zhang    Level: beginner
1252d2567f34SHong Zhang 
12531cc06b55SBarry Smith .seealso: [](ch_ts), `TSIRKSetType()`, `TS`, `TSIRK`, `TSIRKRegister()`
1254d2567f34SHong Zhang J*/
1255d2567f34SHong Zhang typedef const char *TSIRKType;
1256d2567f34SHong Zhang #define TSIRKGAUSS "gauss"
1257d2567f34SHong Zhang 
1258d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKGetType(TS, TSIRKType *);
1259d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKSetType(TS, TSIRKType);
1260d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKGetNumStages(TS, PetscInt *);
1261d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKSetNumStages(TS, PetscInt);
1262d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKRegister(const char[], PetscErrorCode (*function)(TS));
1263d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKTableauCreate(TS, PetscInt, const PetscReal *, const PetscReal *, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, const PetscScalar *);
1264d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKInitializePackage(void);
1265d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKFinalizePackage(void);
1266d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKRegisterDestroy(void);
1267d2567f34SHong Zhang 
1268d2567f34SHong Zhang /*J
1269195e9b02SBarry Smith    TSGLEEType - String with the name of a General Linear with Error Estimation `TSGLEE` type
1270b6a60446SDebojyoti Ghosh 
1271b6a60446SDebojyoti Ghosh    Level: beginner
1272b6a60446SDebojyoti Ghosh 
12731cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEESetType()`, `TS`, `TSGLEE`, `TSGLEERegister()`
1274b6a60446SDebojyoti Ghosh J*/
1275b6a60446SDebojyoti Ghosh typedef const char *TSGLEEType;
1276ab8c5912SEmil Constantinescu #define TSGLEEi1      "BE1"
1277b6a60446SDebojyoti Ghosh #define TSGLEE23      "23"
1278b6a60446SDebojyoti Ghosh #define TSGLEE24      "24"
1279b6a60446SDebojyoti Ghosh #define TSGLEE25I     "25i"
1280b6a60446SDebojyoti Ghosh #define TSGLEE35      "35"
1281b6a60446SDebojyoti Ghosh #define TSGLEEEXRK2A  "exrk2a"
1282b6a60446SDebojyoti Ghosh #define TSGLEERK32G1  "rk32g1"
1283b6a60446SDebojyoti Ghosh #define TSGLEERK285EX "rk285ex"
128416a05f60SBarry Smith 
1285b6a60446SDebojyoti Ghosh /*J
1286195e9b02SBarry Smith    TSGLEEMode - String with the mode of error estimation for a General Linear with Error Estimation `TSGLEE` type
1287b6a60446SDebojyoti Ghosh 
1288b6a60446SDebojyoti Ghosh    Level: beginner
1289b6a60446SDebojyoti Ghosh 
12901cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEESetMode()`, `TS`, `TSGLEE`, `TSGLEERegister()`
1291b6a60446SDebojyoti Ghosh J*/
1292b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEEGetType(TS ts, TSGLEEType *);
1293b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEESetType(TS ts, TSGLEEType);
129457df6a1bSDebojyoti 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[]);
12954bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSGLEERegisterAll(void);
1296b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEEFinalizePackage(void);
1297b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEEInitializePackage(void);
1298b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEERegisterDestroy(void);
1299b6a60446SDebojyoti Ghosh 
1300b6a60446SDebojyoti Ghosh /*J
1301195e9b02SBarry Smith    TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX `TSARKIMEX` type
1302020d8f30SJed Brown 
1303020d8f30SJed Brown    Level: beginner
1304020d8f30SJed Brown 
13051cc06b55SBarry Smith .seealso: [](ch_ts), `TSARKIMEXSetType()`, `TS`, `TSARKIMEX`, `TSARKIMEXRegister()`
1306020d8f30SJed Brown J*/
130719fd82e9SBarry Smith typedef const char *TSARKIMEXType;
1308e817cc15SEmil Constantinescu #define TSARKIMEX1BEE   "1bee"
13091f80e275SEmil Constantinescu #define TSARKIMEXA2     "a2"
13101f80e275SEmil Constantinescu #define TSARKIMEXL2     "l2"
13111f80e275SEmil Constantinescu #define TSARKIMEXARS122 "ars122"
13121f80e275SEmil Constantinescu #define TSARKIMEX2C     "2c"
13138a381b04SJed Brown #define TSARKIMEX2D     "2d"
1314a3a57f36SJed Brown #define TSARKIMEX2E     "2e"
13156cf0794eSJed Brown #define TSARKIMEXPRSSP2 "prssp2"
13168a381b04SJed Brown #define TSARKIMEX3      "3"
13176cf0794eSJed Brown #define TSARKIMEXBPR3   "bpr3"
13186cf0794eSJed Brown #define TSARKIMEXARS443 "ars443"
13198a381b04SJed Brown #define TSARKIMEX4      "4"
13208a381b04SJed Brown #define TSARKIMEX5      "5"
132119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts, TSARKIMEXType *);
132219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts, TSARKIMEXType);
1323014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS, PetscBool);
13243a28c0e5SStefano Zampini PETSC_EXTERN PetscErrorCode TSARKIMEXGetFullyImplicit(TS, PetscBool *);
132519fd82e9SBarry 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[]);
1326607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void);
1327560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
1328014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);
13298a381b04SJed Brown 
1330020d8f30SJed Brown /*J
1331d27334e2SStefano Zampini    TSDIRKType - String with the name of a Diagonally Implicit Runge-Kutta `TSDIRK` type
1332d27334e2SStefano Zampini 
1333d27334e2SStefano Zampini    Level: beginner
1334d27334e2SStefano Zampini 
1335d27334e2SStefano Zampini .seealso: [](ch_ts), `TSDIRKSetType()`, `TS`, `TSDIRK`, `TSDIRKRegister()`
1336d27334e2SStefano Zampini J*/
1337d27334e2SStefano Zampini typedef const char *TSDIRKType;
1338d27334e2SStefano Zampini #define TSDIRKS212      "s212"
13393405e92cSStefano Zampini #define TSDIRKES122SAL  "es122sal"
1340d27334e2SStefano Zampini #define TSDIRKES213SAL  "es213sal"
1341d27334e2SStefano Zampini #define TSDIRKES324SAL  "es324sal"
1342d27334e2SStefano Zampini #define TSDIRKES325SAL  "es325sal"
1343d27334e2SStefano Zampini #define TSDIRK657A      "657a"
1344d27334e2SStefano Zampini #define TSDIRKES648SA   "es648sa"
1345d27334e2SStefano Zampini #define TSDIRK658A      "658a"
1346d27334e2SStefano Zampini #define TSDIRKS659A     "s659a"
1347d27334e2SStefano Zampini #define TSDIRK7510SAL   "7510sal"
1348d27334e2SStefano Zampini #define TSDIRKES7510SA  "es7510sa"
1349d27334e2SStefano Zampini #define TSDIRK759A      "759a"
1350d27334e2SStefano Zampini #define TSDIRKS7511SAL  "s7511sal"
1351d27334e2SStefano Zampini #define TSDIRK8614A     "8614a"
1352d27334e2SStefano Zampini #define TSDIRK8616SAL   "8616sal"
1353d27334e2SStefano Zampini #define TSDIRKES8516SAL "es8516sal"
1354d27334e2SStefano Zampini PETSC_EXTERN PetscErrorCode TSDIRKGetType(TS ts, TSDIRKType *);
1355d27334e2SStefano Zampini PETSC_EXTERN PetscErrorCode TSDIRKSetType(TS ts, TSDIRKType);
1356d27334e2SStefano Zampini PETSC_EXTERN PetscErrorCode TSDIRKRegister(TSDIRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
1357d27334e2SStefano Zampini 
1358d27334e2SStefano Zampini /*J
1359195e9b02SBarry Smith    TSRosWType - String with the name of a Rosenbrock-W `TSROSW` type
1360020d8f30SJed Brown 
1361020d8f30SJed Brown    Level: beginner
1362020d8f30SJed Brown 
13631cc06b55SBarry Smith .seealso: [](ch_ts), `TSRosWSetType()`, `TS`, `TSROSW`, `TSRosWRegister()`
1364020d8f30SJed Brown J*/
136519fd82e9SBarry Smith typedef const char *TSRosWType;
136661692a83SJed Brown #define TSROSW2M          "2m"
136761692a83SJed Brown #define TSROSW2P          "2p"
1368fe7e6d57SJed Brown #define TSROSWRA3PW       "ra3pw"
1369fe7e6d57SJed Brown #define TSROSWRA34PW2     "ra34pw2"
1370ef3c5b88SJed Brown #define TSROSWRODAS3      "rodas3"
1371ef3c5b88SJed Brown #define TSROSWSANDU3      "sandu3"
1372961f28d0SJed Brown #define TSROSWASSP3P3S1C  "assp3p3s1c"
1373961f28d0SJed Brown #define TSROSWLASSP3P4S2C "lassp3p4s2c"
137443b21953SEmil Constantinescu #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
1375753f8adbSEmil Constantinescu #define TSROSWARK3        "ark3"
13763606a31eSEmil Constantinescu #define TSROSWTHETA1      "theta1"
13773606a31eSEmil Constantinescu #define TSROSWTHETA2      "theta2"
137842faf41dSJed Brown #define TSROSWGRK4T       "grk4t"
137942faf41dSJed Brown #define TSROSWSHAMP4      "shamp4"
138042faf41dSJed Brown #define TSROSWVELDD4      "veldd4"
138142faf41dSJed Brown #define TSROSW4L          "4l"
1382b1c69cc3SEmil Constantinescu 
13836bd7aeb5SHong Zhang PETSC_EXTERN PetscErrorCode TSRosWGetType(TS, TSRosWType *);
13846bd7aeb5SHong Zhang PETSC_EXTERN PetscErrorCode TSRosWSetType(TS, TSRosWType);
1385014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS, PetscBool);
138619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
138719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal);
1388607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void);
1389560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
1390014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);
1391e27a552bSJed Brown 
1392211a84d6SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSBDFSetOrder(TS, PetscInt);
1393211a84d6SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSBDFGetOrder(TS, PetscInt *);
1394211a84d6SLisandro Dalcin 
13956bd7aeb5SHong Zhang /*J
1396195e9b02SBarry Smith   TSBasicSymplecticType - String with the name of a basic symplectic integration `TSBASICSYMPLECTIC` type
13976bd7aeb5SHong Zhang 
13986bd7aeb5SHong Zhang   Level: beginner
13996bd7aeb5SHong Zhang 
14001cc06b55SBarry Smith .seealso: [](ch_ts), `TSBasicSymplecticSetType()`, `TS`, `TSBASICSYMPLECTIC`, `TSBasicSymplecticRegister()`
14016bd7aeb5SHong Zhang J*/
1402e49d4f37SHong Zhang typedef const char *TSBasicSymplecticType;
1403e49d4f37SHong Zhang #define TSBASICSYMPLECTICSIEULER   "1"
1404e49d4f37SHong Zhang #define TSBASICSYMPLECTICVELVERLET "2"
1405e49d4f37SHong Zhang #define TSBASICSYMPLECTIC3         "3"
1406e49d4f37SHong Zhang #define TSBASICSYMPLECTIC4         "4"
1407e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticSetType(TS, TSBasicSymplecticType);
1408e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticGetType(TS, TSBasicSymplecticType *);
1409e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegister(TSBasicSymplecticType, PetscInt, PetscInt, PetscReal[], PetscReal[]);
14104bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterAll(void);
1411e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticInitializePackage(void);
1412e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticFinalizePackage(void);
1413e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterDestroy(void);
14146bd7aeb5SHong Zhang 
1415db4653c2SDaniel Finn /*J
1416195e9b02SBarry Smith   TSDISCGRAD - The Discrete Gradient integrator is a timestepper for Hamiltonian systems designed to conserve the first integral (energy),
141787497f52SBarry Smith   but also has the property for some systems of monotonicity in a functional.
1418db4653c2SDaniel Finn 
1419db4653c2SDaniel Finn   Level: beginner
1420db4653c2SDaniel Finn 
14211cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, TSDiscGradSetFormulation()`, `TSDiscGradGetFormulation()`
1422db4653c2SDaniel Finn J*/
142340be0ff1SMatthew 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 *);
14244bf303faSJacob 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 *);
1425db4653c2SDaniel Finn PETSC_EXTERN PetscErrorCode TSDiscGradIsGonzalez(TS, PetscBool *);
1426db4653c2SDaniel Finn PETSC_EXTERN PetscErrorCode TSDiscGradUseGonzalez(TS, PetscBool);
142740be0ff1SMatthew G. Knepley 
142883e2fdc7SBarry Smith /*
14290f3b3ca1SHong Zhang        PETSc interface to Sundials
143083e2fdc7SBarry Smith */
1431e808b789SPatrick Sanan #ifdef PETSC_HAVE_SUNDIALS2
14329371c9d4SSatish Balay typedef enum {
14339371c9d4SSatish Balay   SUNDIALS_ADAMS = 1,
14349371c9d4SSatish Balay   SUNDIALS_BDF   = 2
14359371c9d4SSatish Balay } TSSundialsLmmType;
14366a6fc655SJed Brown PETSC_EXTERN const char *const TSSundialsLmmTypes[];
14379371c9d4SSatish Balay typedef enum {
14389371c9d4SSatish Balay   SUNDIALS_MODIFIED_GS  = 1,
14399371c9d4SSatish Balay   SUNDIALS_CLASSICAL_GS = 2
14409371c9d4SSatish Balay } TSSundialsGramSchmidtType;
14416a6fc655SJed Brown PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
14424bf303faSJacob Faibussowitsch 
1443014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS, TSSundialsLmmType);
1444014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS, PC *);
1445014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS, PetscReal, PetscReal);
1446014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS, PetscReal);
1447014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS, PetscReal);
1448014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS, PetscInt *, PetscInt *);
1449014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS, TSSundialsGramSchmidtType);
1450014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS, PetscInt);
1451014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS, PetscReal);
1452014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS, PetscBool);
1453014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS, PetscInt);
1454c4e80e11SFlorian PETSC_EXTERN PetscErrorCode TSSundialsSetMaxord(TS, PetscInt);
1455918687b8SPatrick Sanan PETSC_EXTERN PetscErrorCode TSSundialsSetUseDense(TS, PetscBool);
14566dc17ff5SMatthew Knepley #endif
145783e2fdc7SBarry Smith 
1458014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS, PetscReal);
1459014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS, PetscReal *);
1460014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS, PetscBool *);
1461014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS, PetscBool);
14620de4c49aSJed Brown 
1463014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS, PetscReal);
1464014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS, PetscReal, PetscReal, PetscReal);
1465014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS, PetscReal *, PetscReal *, PetscReal *);
146688df8a41SLisandro Dalcin 
1467818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAlpha2SetRadius(TS, PetscReal);
1468818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAlpha2SetParams(TS, PetscReal, PetscReal, PetscReal, PetscReal);
1469818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAlpha2GetParams(TS, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
1470818efac9SLisandro Dalcin 
1471220f924aSDavid Kamensky /*S
14728434afd1SBarry Smith   TSAlpha2PredictorFn - A callback to set the predictor (i.e., the initial guess for the nonlinear solver) in
1473220f924aSDavid Kamensky   a second-order generalized-alpha time integrator.
1474220f924aSDavid Kamensky 
1475220f924aSDavid Kamensky   Calling Sequence:
1476220f924aSDavid Kamensky + ts   - the `TS` context obtained from `TSCreate()`
1477220f924aSDavid Kamensky . X0   - the previous time step's state vector
1478220f924aSDavid Kamensky . V0   - the previous time step's first derivative of the state vector
1479220f924aSDavid Kamensky . A0   - the previous time step's second derivative of the state vector
1480220f924aSDavid Kamensky . X1   - the vector into which the initial guess for the current time step will be written
1481220f924aSDavid Kamensky - ctx  - [optional] user-defined context for the predictor evaluation routine (may be `NULL`)
1482220f924aSDavid Kamensky 
1483220f924aSDavid Kamensky   Level: intermediate
1484220f924aSDavid Kamensky 
1485d1c5d1fcSBarry Smith   Note:
14868434afd1SBarry Smith   The deprecated `TSAlpha2Predictor` still works as a replacement for `TSAlpha2PredictorFn` *.
1487d1c5d1fcSBarry Smith 
1488220f924aSDavid Kamensky .seealso: [](ch_ts), `TS`, `TSAlpha2SetPredictor()`
1489220f924aSDavid Kamensky S*/
14908434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSAlpha2PredictorFn)(TS ts, Vec X0, Vec V0, Vec A0, Vec X1, void *ctx);
1491d1c5d1fcSBarry Smith 
14928434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSAlpha2PredictorFn *TSAlpha2Predictor;
1493d1c5d1fcSBarry Smith 
14948434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSAlpha2SetPredictor(TS, TSAlpha2PredictorFn *, void *ctx);
1495220f924aSDavid Kamensky 
1496014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetDM(TS, DM);
1497014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetDM(TS, DM *);
14986c699258SBarry Smith 
1499014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES, Vec, Vec, void *);
1500d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES, Vec, Mat, Mat, void *);
15010f5c6efeSJed Brown 
1502f3b1f45cSBarry Smith PETSC_EXTERN PetscErrorCode TSRHSJacobianTest(TS, PetscBool *);
1503f3b1f45cSBarry Smith PETSC_EXTERN PetscErrorCode TSRHSJacobianTestTranspose(TS, PetscBool *);
1504aad739acSMatthew G. Knepley 
15052e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSGetComputeInitialCondition(TS, PetscErrorCode (**)(TS, Vec));
15062e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSSetComputeInitialCondition(TS, PetscErrorCode (*)(TS, Vec));
15072e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSComputeInitialCondition(TS, Vec);
15082e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSGetComputeExactError(TS, PetscErrorCode (**)(TS, Vec, Vec));
1509aad739acSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSSetComputeExactError(TS, PetscErrorCode (*)(TS, Vec, Vec));
1510aad739acSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSComputeExactError(TS, Vec, Vec);
1511f2ed2dc7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscConvEstUseTS(PetscConvEst, PetscBool);
1512d60b7d5cSBarry Smith 
1513d60b7d5cSBarry Smith PETSC_EXTERN PetscErrorCode TSSetMatStructure(TS, MatStructure);
1514