xref: /petsc/include/petscts.h (revision 49abdd8a111d9c2ef7fc48ade253ef64e07f9b37)
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);
255*49abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorSet(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, void *), void *, PetscCtxDestroyFn *);
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 *);
378*49abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointMonitorSet(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *), void *, PetscCtxDestroyFn *);
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 *);
4387f27e910SStefano Zampini 
4397f27e910SStefano Zampini typedef struct _n_TSMonitorVTKCtx *TSMonitorVTKCtx;
4407f27e910SStefano Zampini PETSC_EXTERN PetscErrorCode        TSMonitorSolutionVTK(TS, PetscInt, PetscReal, Vec, TSMonitorVTKCtx);
4417f27e910SStefano Zampini PETSC_EXTERN PetscErrorCode        TSMonitorSolutionVTKDestroy(TSMonitorVTKCtx *);
4427f27e910SStefano Zampini PETSC_EXTERN PetscErrorCode        TSMonitorSolutionVTKCtxCreate(const char *, TSMonitorVTKCtx *);
443fb1732b5SBarry Smith 
444014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSStep(TS);
4457cbde773SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSEvaluateWLTE(TS, NormType, PetscInt *, PetscReal *);
446014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS, PetscInt, Vec, PetscBool *);
447cc708dedSBarry Smith PETSC_EXTERN PetscErrorCode TSSolve(TS, Vec);
448e817cc15SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetEquationType(TS, TSEquationType *);
449e817cc15SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSSetEquationType(TS, TSEquationType);
450014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS, TSConvergedReason *);
451d6ad946cSShri Abhyankar PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS, TSConvergedReason);
452487e0bb9SJed Brown PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS, PetscReal *);
4535ef26d82SJed Brown PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS, PetscInt *);
4545ef26d82SJed Brown PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS, PetscInt *);
455cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS, PetscInt *);
456cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS, PetscInt);
457cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS, PetscInt *);
458cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS, PetscInt);
459cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS, PetscBool);
460dcb233daSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRestartStep(TS);
46124655328SShri PETSC_EXTERN PetscErrorCode TSRollBack(TS);
4629dc50cb5SStefano Zampini PETSC_EXTERN PetscErrorCode TSGetStepRollBack(TS, PetscBool *);
463d2daff3dSHong Zhang 
464cc4f23bcSHong Zhang PETSC_EXTERN PetscErrorCode TSGetStages(TS, PetscInt *, Vec *[]);
465818ad0c1SBarry Smith 
466014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetTime(TS, PetscReal *);
467014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetTime(TS, PetscReal);
468e5e524a1SHong Zhang PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS, PetscReal *);
46980275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS, PetscReal *);
47080275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS, PetscReal);
47180275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetStepNumber(TS, PetscInt *);
47280275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetStepNumber(TS, PetscInt);
473818ad0c1SBarry Smith 
47414d0ab18SJacob Faibussowitsch /*S
4758434afd1SBarry Smith   TSRHSFunctionFn - A prototype of a `TS` right-hand-side evaluation function that would be passed to `TSSetRHSFunction()`
47614d0ab18SJacob Faibussowitsch 
47714d0ab18SJacob Faibussowitsch   Calling Sequence:
47814d0ab18SJacob Faibussowitsch + ts  - timestep context
47914d0ab18SJacob Faibussowitsch . t   - current time
48014d0ab18SJacob Faibussowitsch . u   - input vector
48114d0ab18SJacob Faibussowitsch . F   - function vector
48214d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context
48314d0ab18SJacob Faibussowitsch 
48414d0ab18SJacob Faibussowitsch   Level: beginner
48514d0ab18SJacob Faibussowitsch 
486d1c5d1fcSBarry Smith   Note:
4878434afd1SBarry Smith   The deprecated `TSRHSFunction` still works as a replacement for `TSRHSFunctionFn` *.
488d1c5d1fcSBarry Smith 
4898434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetRHSFunction()`, `DMTSSetRHSFunction()`, `TSIFunctionFn`,
4908434afd1SBarry Smith `TSIJacobianFn`, `TSRHSJacobianFn`
49114d0ab18SJacob Faibussowitsch S*/
4928434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSRHSFunctionFn)(TS ts, PetscReal t, Vec u, Vec F, void *ctx);
493d1c5d1fcSBarry Smith 
4948434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSRHSFunctionFn *TSRHSFunction;
49514d0ab18SJacob Faibussowitsch 
49614d0ab18SJacob Faibussowitsch /*S
4978434afd1SBarry Smith   TSRHSJacobianFn - A prototype of a `TS` right-hand-side Jacobian evaluation function that would be passed to `TSSetRHSJacobian()`
49814d0ab18SJacob Faibussowitsch 
49914d0ab18SJacob Faibussowitsch   Calling Sequence:
50014d0ab18SJacob Faibussowitsch + ts   - the `TS` context obtained from `TSCreate()`
50114d0ab18SJacob Faibussowitsch . t    - current time
50214d0ab18SJacob Faibussowitsch . u    - input vector
50314d0ab18SJacob Faibussowitsch . Amat - (approximate) Jacobian matrix
504af27ebaaSBarry Smith . Pmat - matrix from which preconditioner is to be constructed (usually the same as `Amat`)
50514d0ab18SJacob Faibussowitsch - ctx  - [optional] user-defined context for matrix evaluation routine
50614d0ab18SJacob Faibussowitsch 
50714d0ab18SJacob Faibussowitsch   Level: beginner
50814d0ab18SJacob Faibussowitsch 
509d1c5d1fcSBarry Smith   Note:
5108434afd1SBarry Smith   The deprecated `TSRHSJacobian` still works as a replacement for `TSRHSJacobianFn` *.
511d1c5d1fcSBarry Smith 
5128434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetRHSJacobian()`, `DMTSSetRHSJacobian()`, `TSRHSFunctionFn`,
5138434afd1SBarry Smith `TSIFunctionFn`, `TSIJacobianFn`
51414d0ab18SJacob Faibussowitsch S*/
5158434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSRHSJacobianFn)(TS ts, PetscReal t, Vec u, Mat Amat, Mat Pmat, void *ctx);
516d1c5d1fcSBarry Smith 
5178434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSRHSJacobianFn *TSRHSJacobian;
51814d0ab18SJacob Faibussowitsch 
51914d0ab18SJacob Faibussowitsch /*S
5208434afd1SBarry Smith   TSRHSJacobianPFn - A prototype of a function that computes the Jacobian of G w.r.t. the parameters P where
521af27ebaaSBarry Smith   U_t = G(U,P,t), as well as the location to store the matrix that would be passed to `TSSetRHSJacobianP()`
52214d0ab18SJacob Faibussowitsch 
52314d0ab18SJacob Faibussowitsch   Calling Sequence:
52414d0ab18SJacob Faibussowitsch + ts  - the `TS` context
52514d0ab18SJacob Faibussowitsch . t   - current timestep
52614d0ab18SJacob Faibussowitsch . U   - input vector (current ODE solution)
52714d0ab18SJacob Faibussowitsch . A   - output matrix
52814d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context
52914d0ab18SJacob Faibussowitsch 
53014d0ab18SJacob Faibussowitsch   Level: beginner
53114d0ab18SJacob Faibussowitsch 
532d1c5d1fcSBarry Smith   Note:
5338434afd1SBarry Smith   The deprecated `TSRHSJacobianP` still works as a replacement for `TSRHSJacobianPFn` *.
534d1c5d1fcSBarry Smith 
53514d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetRHSJacobianP()`, `TSGetRHSJacobianP()`
53614d0ab18SJacob Faibussowitsch S*/
5378434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSRHSJacobianPFn)(TS ts, PetscReal t, Vec U, Mat A, void *ctx);
53814d0ab18SJacob Faibussowitsch 
5398434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSRHSJacobianPFn *TSRHSJacobianP;
540d1c5d1fcSBarry Smith 
5418434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS, Vec, TSRHSFunctionFn *, void *);
5428434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS, Vec *, TSRHSFunctionFn **, void **);
5438434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS, Mat, Mat, TSRHSJacobianFn *, void *);
5448434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS, Mat *, Mat *, TSRHSJacobianFn **, void **);
545e1244c69SJed Brown PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS, PetscBool);
546818ad0c1SBarry Smith 
54714d0ab18SJacob Faibussowitsch /*S
5488434afd1SBarry Smith   TSSolutionFn - A prototype of a `TS` solution evaluation function that would be passed to `TSSetSolutionFunction()`
54914d0ab18SJacob Faibussowitsch 
55014d0ab18SJacob Faibussowitsch   Calling Sequence:
55114d0ab18SJacob Faibussowitsch + ts  - timestep context
55214d0ab18SJacob Faibussowitsch . t   - current time
55314d0ab18SJacob Faibussowitsch . u   - output vector
55414d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context
55514d0ab18SJacob Faibussowitsch 
55614d0ab18SJacob Faibussowitsch   Level: advanced
55714d0ab18SJacob Faibussowitsch 
558d1c5d1fcSBarry Smith   Note:
5598434afd1SBarry Smith   The deprecated `TSSolutionFunction` still works as a replacement for `TSSolutionFn` *.
560d1c5d1fcSBarry Smith 
56114d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetSolutionFunction()`, `DMTSSetSolutionFunction()`
56214d0ab18SJacob Faibussowitsch S*/
5638434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSSolutionFn)(TS ts, PetscReal t, Vec u, void *ctx);
56414d0ab18SJacob Faibussowitsch 
5658434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSSolutionFn *TSSolutionFunction;
566d1c5d1fcSBarry Smith 
5678434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS, TSSolutionFn *, void *);
56814d0ab18SJacob Faibussowitsch 
56914d0ab18SJacob Faibussowitsch /*S
5708434afd1SBarry Smith   TSForcingFn - A prototype of a `TS` forcing function evaluation function that would be passed to `TSSetForcingFunction()`
57114d0ab18SJacob Faibussowitsch 
57214d0ab18SJacob Faibussowitsch   Calling Sequence:
57314d0ab18SJacob Faibussowitsch + ts  - timestep context
57414d0ab18SJacob Faibussowitsch . t   - current time
57514d0ab18SJacob Faibussowitsch . f   - output vector
57614d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context
57714d0ab18SJacob Faibussowitsch 
57814d0ab18SJacob Faibussowitsch   Level: advanced
57914d0ab18SJacob Faibussowitsch 
580d1c5d1fcSBarry Smith   Note:
5818434afd1SBarry Smith   The deprecated `TSForcingFunction` still works as a replacement for `TSForcingFn` *.
582d1c5d1fcSBarry Smith 
58314d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetForcingFunction()`, `DMTSSetForcingFunction()`
58414d0ab18SJacob Faibussowitsch S*/
5858434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSForcingFn)(TS ts, PetscReal t, Vec f, void *ctx);
58614d0ab18SJacob Faibussowitsch 
5878434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSForcingFn *TSForcingFunction;
588d1c5d1fcSBarry Smith 
5898434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS, TSForcingFn *, void *);
590ef20d060SBarry Smith 
59114d0ab18SJacob Faibussowitsch /*S
5928434afd1SBarry Smith   TSIFunctionFn - A prototype of a `TS` implicit function evaluation function that would be passed to `TSSetIFunction()
59314d0ab18SJacob Faibussowitsch 
59414d0ab18SJacob Faibussowitsch   Calling Sequence:
59514d0ab18SJacob Faibussowitsch + ts  - the `TS` context obtained from `TSCreate()`
59614d0ab18SJacob Faibussowitsch . t   - time at step/stage being solved
59714d0ab18SJacob Faibussowitsch . U   - state vector
59814d0ab18SJacob Faibussowitsch . U_t - time derivative of state vector
59914d0ab18SJacob Faibussowitsch . F   - function vector
60014d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined context for function
60114d0ab18SJacob Faibussowitsch 
60214d0ab18SJacob Faibussowitsch   Level: beginner
60314d0ab18SJacob Faibussowitsch 
604d1c5d1fcSBarry Smith   Note:
6058434afd1SBarry Smith   The deprecated `TSIFunction` still works as a replacement for `TSIFunctionFn` *.
606d1c5d1fcSBarry Smith 
6078434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetIFunction()`, `DMTSSetIFunction()`, `TSIJacobianFn`, `TSRHSFunctionFn`, `TSRHSJacobianFn`
60814d0ab18SJacob Faibussowitsch S*/
6098434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSIFunctionFn)(TS ts, PetscReal t, Vec U, Vec U_t, Vec F, void *ctx);
610d1c5d1fcSBarry Smith 
6118434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSIFunctionFn *TSIFunction;
61214d0ab18SJacob Faibussowitsch 
61314d0ab18SJacob Faibussowitsch /*S
6148434afd1SBarry Smith   TSIJacobianFn - A prototype of a `TS` Jacobian evaluation function that would be passed to `TSSetIJacobian()`
61514d0ab18SJacob Faibussowitsch 
61614d0ab18SJacob Faibussowitsch   Calling Sequence:
61714d0ab18SJacob Faibussowitsch + ts   - the `TS` context obtained from `TSCreate()`
61814d0ab18SJacob Faibussowitsch . t    - time at step/stage being solved
61914d0ab18SJacob Faibussowitsch . U    - state vector
62014d0ab18SJacob Faibussowitsch . U_t  - time derivative of state vector
62114d0ab18SJacob Faibussowitsch . a    - shift
62214d0ab18SJacob Faibussowitsch . Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
62314d0ab18SJacob Faibussowitsch . Pmat - matrix used for constructing preconditioner, usually the same as `Amat`
62414d0ab18SJacob Faibussowitsch - ctx  - [optional] user-defined context for Jacobian evaluation routine
62514d0ab18SJacob Faibussowitsch 
62614d0ab18SJacob Faibussowitsch   Level: beginner
62714d0ab18SJacob Faibussowitsch 
628d1c5d1fcSBarry Smith   Note:
6298434afd1SBarry Smith   The deprecated `TSIJacobian` still works as a replacement for `TSIJacobianFn` *.
63014d0ab18SJacob Faibussowitsch 
6318434afd1SBarry Smith .seealso: [](ch_ts), `TSSetIJacobian()`, `DMTSSetIJacobian()`, `TSIFunctionFn`, `TSRHSFunctionFn`, `TSRHSJacobianFn`
632d1c5d1fcSBarry Smith S*/
6338434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSIJacobianFn)(TS ts, PetscReal t, Vec U, Vec U_t, PetscReal a, Mat Amat, Mat Pmat, void *ctx);
634d1c5d1fcSBarry Smith 
6358434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSIJacobianFn *TSIJacobian;
636d1c5d1fcSBarry Smith 
6378434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetIFunction(TS, Vec, TSIFunctionFn *, void *);
6388434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGetIFunction(TS, Vec *, TSIFunctionFn **, void **);
6398434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS, Mat, Mat, TSIJacobianFn *, void *);
6408434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS, Mat *, Mat *, TSIJacobianFn **, void **);
641316643e7SJed Brown 
64214d0ab18SJacob Faibussowitsch /*S
6438434afd1SBarry Smith   TSI2FunctionFn - A prototype of a `TS` implicit function evaluation function for 2nd order systems that would be passed to `TSSetI2Function()`
64414d0ab18SJacob Faibussowitsch 
64514d0ab18SJacob Faibussowitsch   Calling Sequence:
64614d0ab18SJacob Faibussowitsch + ts   - the `TS` context obtained from `TSCreate()`
64714d0ab18SJacob Faibussowitsch . t    - time at step/stage being solved
64814d0ab18SJacob Faibussowitsch . U    - state vector
64914d0ab18SJacob Faibussowitsch . U_t  - time derivative of state vector
65014d0ab18SJacob Faibussowitsch . U_tt - second time derivative of state vector
65114d0ab18SJacob Faibussowitsch . F    - function vector
65214d0ab18SJacob Faibussowitsch - ctx  - [optional] user-defined context for matrix evaluation routine (may be `NULL`)
65314d0ab18SJacob Faibussowitsch 
65414d0ab18SJacob Faibussowitsch   Level: advanced
65514d0ab18SJacob Faibussowitsch 
656d1c5d1fcSBarry Smith   Note:
6578434afd1SBarry Smith   The deprecated `TSI2Function` still works as a replacement for `TSI2FunctionFn` *.
658d1c5d1fcSBarry Smith 
6598434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetI2Function()`, `DMTSSetI2Function()`, `TSIFunctionFn`
66014d0ab18SJacob Faibussowitsch S*/
6618434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSI2FunctionFn)(TS ts, PetscReal t, Vec U, Vec U_t, Vec U_tt, Vec F, void *ctx);
662d1c5d1fcSBarry Smith 
6638434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSI2FunctionFn *TSI2Function;
66414d0ab18SJacob Faibussowitsch 
66514d0ab18SJacob Faibussowitsch /*S
6668434afd1SBarry Smith   TSI2JacobianFn - A prototype of a `TS` implicit Jacobian evaluation function for 2nd order systems that would be passed to `TSSetI2Jacobian()`
66714d0ab18SJacob Faibussowitsch 
66814d0ab18SJacob Faibussowitsch   Calling Sequence:
66914d0ab18SJacob Faibussowitsch + ts   - the `TS` context obtained from `TSCreate()`
67014d0ab18SJacob Faibussowitsch . t    - time at step/stage being solved
67114d0ab18SJacob Faibussowitsch . U    - state vector
67214d0ab18SJacob Faibussowitsch . U_t  - time derivative of state vector
67314d0ab18SJacob Faibussowitsch . U_tt - second time derivative of state vector
67414d0ab18SJacob Faibussowitsch . v    - shift for U_t
67514d0ab18SJacob Faibussowitsch . a    - shift for U_tt
67614d0ab18SJacob 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
67714d0ab18SJacob Faibussowitsch . jac  - matrix from which to construct the preconditioner, may be same as `J`
67814d0ab18SJacob Faibussowitsch - ctx  - [optional] user-defined context for matrix evaluation routine
67914d0ab18SJacob Faibussowitsch 
68014d0ab18SJacob Faibussowitsch   Level: advanced
68114d0ab18SJacob Faibussowitsch 
682d1c5d1fcSBarry Smith   Note:
6838434afd1SBarry Smith   The deprecated `TSI2Jacobian` still works as a replacement for `TSI2JacobianFn` *.
68414d0ab18SJacob Faibussowitsch 
6858434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetI2Jacobian()`, `DMTSSetI2Jacobian()`, `TSIFunctionFn`, `TSIJacobianFn`, `TSRHSFunctionFn`, `TSRHSJacobianFn`
686d1c5d1fcSBarry Smith S*/
6878434afd1SBarry 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);
688d1c5d1fcSBarry Smith 
6898434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSI2JacobianFn *TSI2Jacobian;
690d1c5d1fcSBarry Smith 
6918434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetI2Function(TS, Vec, TSI2FunctionFn *, void *);
6928434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGetI2Function(TS, Vec *, TSI2FunctionFn **, void **);
6938434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetI2Jacobian(TS, Mat, Mat, TSI2JacobianFn *, void *);
6948434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGetI2Jacobian(TS, Mat *, Mat *, TSI2JacobianFn **, void **);
695efe9872eSLisandro Dalcin 
696545aaa6fSHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitSetIS(TS, const char[], IS);
697545aaa6fSHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetIS(TS, const char[], IS *);
6988434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSRHSSplitSetRHSFunction(TS, const char[], Vec, TSRHSFunctionFn *, void *);
6993a2a065fSHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitSetIFunction(TS, const char[], Vec, TSIFunctionFn *, void *);
7003a2a065fSHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitSetIJacobian(TS, const char[], Mat, Mat, TSIJacobianFn *, void *);
7011d06f6b3SHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTS(TS, const char[], TS *);
7021d06f6b3SHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTSs(TS, PetscInt *, TS *[]);
7030fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSSetUseSplitRHSFunction(TS, PetscBool);
7040fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSGetUseSplitRHSFunction(TS, PetscBool *);
7054bd3aaa3SHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetSNES(TS, SNES *);
7064bd3aaa3SHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitSetSNES(TS, SNES);
707d9194312SHong Zhang 
7088434afd1SBarry Smith PETSC_EXTERN TSRHSFunctionFn TSComputeRHSFunctionLinear;
7098434afd1SBarry Smith PETSC_EXTERN TSRHSJacobianFn TSComputeRHSJacobianConstant;
710014dd563SJed Brown PETSC_EXTERN PetscErrorCode  TSComputeIFunctionLinear(TS, PetscReal, Vec, Vec, Vec, void *);
711d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode  TSComputeIJacobianConstant(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
7121c8a10a0SJed Brown PETSC_EXTERN PetscErrorCode  TSComputeSolutionFunction(TS, PetscReal, Vec);
7139b7cd975SBarry Smith PETSC_EXTERN PetscErrorCode  TSComputeForcingFunction(TS, PetscReal, Vec);
714847ff0e1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode  TSComputeIJacobianDefaultColor(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
715d7cfae9bSHong Zhang PETSC_EXTERN PetscErrorCode  TSPruneIJacobianColor(TS, Mat, Mat);
716e34be4c2SBarry Smith 
717014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
718b8123daeSJed Brown PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS, PetscReal));
7199be3e283SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS, PetscReal, PetscInt, Vec *));
720c688d042SShri Abhyankar PETSC_EXTERN PetscErrorCode TSSetPostEvaluate(TS, PetscErrorCode (*)(TS));
721014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
722ecc87898SStefano Zampini PETSC_EXTERN PetscErrorCode TSSetResize(TS, PetscBool, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscBool *, void *), PetscErrorCode (*)(TS, PetscInt, Vec[], Vec[], void *), void *);
723014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPreStep(TS);
724b8123daeSJed Brown PETSC_EXTERN PetscErrorCode TSPreStage(TS, PetscReal);
7259be3e283SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSPostStage(TS, PetscReal, PetscInt, Vec *);
726c688d042SShri Abhyankar PETSC_EXTERN PetscErrorCode TSPostEvaluate(TS);
727014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPostStep(TS);
7286bd3a4fdSStefano Zampini PETSC_EXTERN PetscErrorCode TSResize(TS);
7296bd3a4fdSStefano Zampini PETSC_EXTERN PetscErrorCode TSResizeRetrieveVec(TS, const char *, Vec *);
7306bd3a4fdSStefano Zampini PETSC_EXTERN PetscErrorCode TSResizeRegisterVec(TS, const char *, Vec);
731c6bf8827SStefano Zampini PETSC_EXTERN PetscErrorCode TSGetStepResize(TS, PetscBool *);
7326bd3a4fdSStefano Zampini 
733014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSInterpolate(TS, PetscReal, Vec);
734014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetTolerances(TS, PetscReal, Vec, PetscReal, Vec);
735c5033834SJed Brown PETSC_EXTERN PetscErrorCode TSGetTolerances(TS, PetscReal *, Vec *, PetscReal *, Vec *);
7367453f775SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *);
7378a175baeSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm(TS, Vec, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *);
738014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS, PetscReal);
739014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS, PetscReal *);
740d183316bSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSSetFunctionDomainError(TS, PetscErrorCode (*)(TS, PetscReal, Vec, PetscBool *));
741d183316bSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSFunctionDomainError(TS, PetscReal, Vec, PetscBool *);
742000e7ae3SMatthew Knepley 
743014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS, PetscErrorCode (*)(TS, PetscReal *, void *), void *);
7448d359177SBarry Smith PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS, PetscReal *, void *);
745014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS, PetscReal *);
746014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS, PetscReal);
747014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS, PetscErrorCode (*)(TS, Vec, void *, PetscReal *, PetscBool *), void *);
7488d359177SBarry Smith PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS, Vec, void *, PetscReal *, PetscBool *);
749014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS, Vec, PetscReal *, PetscBool *);
750014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS, PetscReal);
751014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);
75221c89e3eSBarry Smith 
753014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPythonSetType(TS, const char[]);
754ebead697SStefano Zampini PETSC_EXTERN PetscErrorCode TSPythonGetType(TS, const char *[]);
7551d6018f0SLisandro Dalcin 
756014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS, PetscReal, Vec, Vec);
757d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS, PetscReal, Vec, Mat, Mat);
758014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS, PetscReal, Vec, Vec, Vec, PetscBool);
759d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscBool);
760efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSComputeI2Function(TS, PetscReal, Vec, Vec, Vec, Vec);
761efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSComputeI2Jacobian(TS, PetscReal, Vec, Vec, Vec, PetscReal, PetscReal, Mat, Mat);
762f9c1d6abSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS, PetscReal, PetscReal, PetscReal *, PetscReal *);
763818ad0c1SBarry Smith 
764014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS, Vec, Vec);
765d6ebe24aSShri Abhyankar 
766967bb25aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetBoundaryLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *);
7678434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM, TSRHSFunctionFn *, void *);
7688434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM, TSRHSFunctionFn **, void **);
769*49abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionContextDestroy(DM, PetscCtxDestroyFn *);
7708434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM, TSRHSJacobianFn *, void *);
7718434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM, TSRHSJacobianFn **, void **);
772*49abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobianContextDestroy(DM, PetscCtxDestroyFn *);
7738434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM, TSIFunctionFn *, void *);
7748434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM, TSIFunctionFn **, void **);
775*49abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetIFunctionContextDestroy(DM, PetscCtxDestroyFn *);
7768434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM, TSIJacobianFn *, void *);
7778434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM, TSIJacobianFn **, void **);
778*49abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetIJacobianContextDestroy(DM, PetscCtxDestroyFn *);
7798434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetI2Function(DM, TSI2FunctionFn *, void *);
7808434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetI2Function(DM, TSI2FunctionFn **, void **);
781*49abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetI2FunctionContextDestroy(DM, PetscCtxDestroyFn *);
7828434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetI2Jacobian(DM, TSI2JacobianFn *, void *);
7838434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetI2Jacobian(DM, TSI2JacobianFn **, void **);
784*49abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetI2JacobianContextDestroy(DM, PetscCtxDestroyFn *);
785efe9872eSLisandro Dalcin 
78614d0ab18SJacob Faibussowitsch /*S
7878434afd1SBarry Smith   TSTransientVariableFn - A prototype of a function to transform from state to transient variables that would be passed to `TSSetTransientVariable()`
78814d0ab18SJacob Faibussowitsch 
78914d0ab18SJacob Faibussowitsch   Calling Sequence:
79014d0ab18SJacob Faibussowitsch + ts  - timestep context
79114d0ab18SJacob Faibussowitsch . p   - input vector (primitive form)
79214d0ab18SJacob Faibussowitsch . c   - output vector, transient variables (conservative form)
79314d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context
79414d0ab18SJacob Faibussowitsch 
79514d0ab18SJacob Faibussowitsch   Level: advanced
79614d0ab18SJacob Faibussowitsch 
797d1c5d1fcSBarry Smith   Note:
7988434afd1SBarry Smith   The deprecated `TSTransientVariable` still works as a replacement for `TSTransientVariableFn` *.
799d1c5d1fcSBarry Smith 
80014d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetTransientVariable()`, `DMTSSetTransientVariable()`
80114d0ab18SJacob Faibussowitsch S*/
8028434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSTransientVariableFn)(TS ts, Vec p, Vec c, void *ctx);
80314d0ab18SJacob Faibussowitsch 
8048434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSTransientVariableFn *TSTransientVariable;
805d1c5d1fcSBarry Smith 
8068434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetTransientVariable(TS, TSTransientVariableFn *, void *);
8078434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetTransientVariable(DM, TSTransientVariableFn *, void *);
8088434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetTransientVariable(DM, TSTransientVariableFn **, void *);
809e3c11fc1SJed Brown PETSC_EXTERN PetscErrorCode TSComputeTransientVariable(TS, Vec, Vec);
810e3c11fc1SJed Brown PETSC_EXTERN PetscErrorCode TSHasTransientVariable(TS, PetscBool *);
811e3c11fc1SJed Brown 
8128434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM, TSSolutionFn *, void *);
8138434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM, TSSolutionFn **, void **);
8148434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM, TSForcingFn *, void *);
8158434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM, TSForcingFn **, void **);
816e17bd7bbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCheckResidual(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscReal *);
817e17bd7bbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCheckJacobian(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscBool *, PetscReal *);
8187f96f943SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec);
8199b7cd975SBarry Smith 
8200c9409e4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSGetIFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, Vec, void *), void **);
821d433e6cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, Vec, void *), void *);
8220c9409e4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSGetIJacobianLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void **);
823d433e6cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void *);
8240c9409e4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSGetRHSFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, void *), void **);
825d433e6cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *);
826b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrix(DM);
827b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrixLumped(DM);
828b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSDestroyRHSMassMatrix(DM);
8296c6b9e74SPeter Brune 
8306c6b9e74SPeter Brune PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM, PetscErrorCode (*)(void *, PetscViewer), PetscErrorCode (*)(void **, PetscViewer));
8316c6b9e74SPeter Brune PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM, PetscErrorCode (*)(void *, PetscViewer), PetscErrorCode (*)(void **, PetscViewer));
8326c6b9e74SPeter Brune 
83314d0ab18SJacob Faibussowitsch /*S
834dd8e379bSPierre Jolivet   DMDATSRHSFunctionLocalFn - A prototype of a local `TS` right-hand side residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetRHSFunctionLocal()`
8356c6b9e74SPeter Brune 
83614d0ab18SJacob Faibussowitsch   Calling Sequence:
83714d0ab18SJacob Faibussowitsch + info - defines the subdomain to evaluate the residual on
83814d0ab18SJacob Faibussowitsch . t    - time at which to evaluate residual
83914d0ab18SJacob Faibussowitsch . x    - array of local state information
84014d0ab18SJacob Faibussowitsch . f    - output array of local residual information
84114d0ab18SJacob Faibussowitsch - ctx  - optional user context
84214d0ab18SJacob Faibussowitsch 
84314d0ab18SJacob Faibussowitsch   Level: beginner
84414d0ab18SJacob Faibussowitsch 
845d1c5d1fcSBarry Smith   Note:
8468434afd1SBarry Smith   The deprecated `DMDATSRHSFunctionLocal` still works as a replacement for `DMDATSRHSFunctionLocalFn` *.
847d1c5d1fcSBarry Smith 
8488434afd1SBarry Smith .seealso: `DMDA`, `DMDATSSetRHSFunctionLocal()`, `TSRHSFunctionFn`, `DMDATSRHSJacobianLocalFn`, `DMDATSIJacobianLocalFn`, `DMDATSIFunctionLocalFn`
84914d0ab18SJacob Faibussowitsch S*/
8508434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDATSRHSFunctionLocalFn)(DMDALocalInfo *info, PetscReal t, void *x, void *f, void *ctx);
851d1c5d1fcSBarry Smith 
8528434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef DMDATSRHSFunctionLocalFn *DMDATSRHSFunctionLocal;
85314d0ab18SJacob Faibussowitsch 
85414d0ab18SJacob Faibussowitsch /*S
8558434afd1SBarry Smith   DMDATSRHSJacobianLocalFn - A prototype of a local residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetRHSJacobianLocal()`
85614d0ab18SJacob Faibussowitsch 
85714d0ab18SJacob Faibussowitsch   Calling Sequence:
85814d0ab18SJacob Faibussowitsch + info - defines the subdomain to evaluate the residual on
85914d0ab18SJacob Faibussowitsch . t    - time at which to evaluate residual
86014d0ab18SJacob Faibussowitsch . x    - array of local state information
86114d0ab18SJacob Faibussowitsch . J    - Jacobian matrix
86214d0ab18SJacob Faibussowitsch . B    - matrix from which to construct the preconditioner; often same as `J`
86314d0ab18SJacob Faibussowitsch - ctx  - optional context
86414d0ab18SJacob Faibussowitsch 
86514d0ab18SJacob Faibussowitsch   Level: beginner
86614d0ab18SJacob Faibussowitsch 
867d1c5d1fcSBarry Smith   Note:
8688434afd1SBarry Smith   The deprecated `DMDATSRHSJacobianLocal` still works as a replacement for `DMDATSRHSJacobianLocalFn` *.
869d1c5d1fcSBarry Smith 
8708434afd1SBarry Smith .seealso: `DMDA`, `DMDATSSetRHSJacobianLocal()`, `TSRHSJacobianFn`, `DMDATSRHSFunctionLocalFn`, `DMDATSIJacobianLocalFn`, `DMDATSIFunctionLocalFn`
87114d0ab18SJacob Faibussowitsch S*/
8728434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDATSRHSJacobianLocalFn)(DMDALocalInfo *info, PetscReal t, void *x, Mat J, Mat B, void *ctx);
873d1c5d1fcSBarry Smith 
8748434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef DMDATSRHSJacobianLocalFn *DMDATSRHSJacobianLocal;
87514d0ab18SJacob Faibussowitsch 
87614d0ab18SJacob Faibussowitsch /*S
8778434afd1SBarry Smith   DMDATSIFunctionLocalFn - A prototype of a local residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetIFunctionLocal()`
87814d0ab18SJacob Faibussowitsch 
87914d0ab18SJacob Faibussowitsch   Calling Sequence:
88014d0ab18SJacob Faibussowitsch + info  - defines the subdomain to evaluate the residual on
88114d0ab18SJacob Faibussowitsch . t     - time at which to evaluate residual
88214d0ab18SJacob Faibussowitsch . x     - array of local state information
88314d0ab18SJacob Faibussowitsch . xdot  - array of local time derivative information
88414d0ab18SJacob Faibussowitsch . imode - output array of local function evaluation information
88514d0ab18SJacob Faibussowitsch - ctx   - optional context
88614d0ab18SJacob Faibussowitsch 
88714d0ab18SJacob Faibussowitsch   Level: beginner
88814d0ab18SJacob Faibussowitsch 
889d1c5d1fcSBarry Smith   Note:
8908434afd1SBarry Smith   The deprecated `DMDATSIFunctionLocal` still works as a replacement for `DMDATSIFunctionLocalFn` *.
891d1c5d1fcSBarry Smith 
8928434afd1SBarry Smith .seealso: `DMDA`, `DMDATSSetIFunctionLocal()`, `DMDATSIJacobianLocalFn`, `TSIFunctionFn`
89314d0ab18SJacob Faibussowitsch S*/
8948434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDATSIFunctionLocalFn)(DMDALocalInfo *info, PetscReal t, void *x, void *xdot, void *imode, void *ctx);
895d1c5d1fcSBarry Smith 
8968434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef DMDATSIFunctionLocalFn *DMDATSIFunctionLocal;
89714d0ab18SJacob Faibussowitsch 
89814d0ab18SJacob Faibussowitsch /*S
8998434afd1SBarry Smith   DMDATSIJacobianLocalFn - A prototype of a local residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetIJacobianLocal()`
90014d0ab18SJacob Faibussowitsch 
90114d0ab18SJacob Faibussowitsch   Calling Sequence:
90214d0ab18SJacob Faibussowitsch + info  - defines the subdomain to evaluate the residual on
90314d0ab18SJacob Faibussowitsch . t     - time at which to evaluate the jacobian
90414d0ab18SJacob Faibussowitsch . x     - array of local state information
90514d0ab18SJacob Faibussowitsch . xdot  - time derivative at this state
90614d0ab18SJacob Faibussowitsch . shift - see `TSSetIJacobian()` for the meaning of this parameter
90714d0ab18SJacob Faibussowitsch . J     - Jacobian matrix
90814d0ab18SJacob Faibussowitsch . B     - matrix from which to construct the preconditioner; often same as `J`
90914d0ab18SJacob Faibussowitsch - ctx   - optional context
91014d0ab18SJacob Faibussowitsch 
91114d0ab18SJacob Faibussowitsch   Level: beginner
91214d0ab18SJacob Faibussowitsch 
913d1c5d1fcSBarry Smith   Note:
9148434afd1SBarry Smith   The deprecated `DMDATSIJacobianLocal` still works as a replacement for `DMDATSIJacobianLocalFn` *.
91514d0ab18SJacob Faibussowitsch 
9168434afd1SBarry Smith .seealso: `DMDA` `DMDATSSetIJacobianLocal()`, `TSIJacobianFn`, `DMDATSIFunctionLocalFn`, `DMDATSRHSFunctionLocalFn`,  `DMDATSRHSJacobianlocal()`
917d1c5d1fcSBarry Smith S*/
9188434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDATSIJacobianLocalFn)(DMDALocalInfo *info, PetscReal t, void *x, void *xdot, PetscReal shift, Mat J, Mat B, void *ctx);
919d1c5d1fcSBarry Smith 
9208434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef DMDATSIJacobianLocalFn *DMDATSIJacobianLocal;
921d1c5d1fcSBarry Smith 
9228434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM, InsertMode, DMDATSRHSFunctionLocalFn *, void *);
9238434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM, DMDATSRHSJacobianLocalFn *, void *);
9248434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM, InsertMode, DMDATSIFunctionLocalFn *, void *);
9258434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM, DMDATSIJacobianLocalFn *, void *);
9266c6b9e74SPeter Brune 
927be1b0d75SMatthew G. Knepley typedef struct _n_TSMonitorLGCtx *TSMonitorLGCtx;
928d1212d36SBarry Smith typedef struct {
929d1212d36SBarry Smith   Vec            ray;
930d1212d36SBarry Smith   VecScatter     scatter;
931d1212d36SBarry Smith   PetscViewer    viewer;
93251b4a12fSMatthew G. Knepley   TSMonitorLGCtx lgctx;
933d1212d36SBarry Smith } TSMonitorDMDARayCtx;
934d1212d36SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void **);
935d1212d36SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS, PetscInt, PetscReal, Vec, void *);
93651b4a12fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS, PetscInt, PetscReal, Vec, void *);
937d1212d36SBarry Smith 
938bdad233fSMatthew Knepley /* Dynamic creation and loading functions */
939140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList TSList;
94019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode    TSGetType(TS, TSType *);
94119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode    TSSetType(TS, TSType);
942bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode    TSRegister(const char[], PetscErrorCode (*)(TS));
94330de9b25SBarry Smith 
944014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetSNES(TS, SNES *);
945deb2cd25SJed Brown PETSC_EXTERN PetscErrorCode TSSetSNES(TS, SNES);
946014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetKSP(TS, KSP *);
947818ad0c1SBarry Smith 
948014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSView(TS, PetscViewer);
94955849f57SBarry Smith PETSC_EXTERN PetscErrorCode TSLoad(TS, PetscViewer);
950fe2efc57SMark PETSC_EXTERN PetscErrorCode TSViewFromOptions(TS, PetscObject, const char[]);
951fe2efc57SMark PETSC_EXTERN PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory, PetscObject, const char[]);
95255849f57SBarry Smith 
95355849f57SBarry Smith #define TS_FILE_CLASSID 1211225
95421c89e3eSBarry Smith 
955014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS, void *);
956014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS, void *);
95721c89e3eSBarry Smith 
958a9f9c1f6SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtx *);
959a9f9c1f6SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *);
9604f09c107SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS, PetscInt, PetscReal, Vec, void *);
9614f09c107SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS, PetscInt, PetscReal, Vec, void *);
96231152f8aSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS, const char *const *);
963b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS, const char *const **);
964e673d494SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx, const char *const *);
965a66092f1SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS, const char *const *);
966a66092f1SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx, const char *const *);
967*49abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS, PetscErrorCode (*)(void *, Vec, Vec *), PetscCtxDestroyFn *, void *);
968*49abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx, PetscErrorCode (*)(void *, Vec, Vec *), PetscCtxDestroyFn *, void *);
9694f09c107SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS, PetscInt, PetscReal, Vec, void *);
970201da799SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS, PetscInt, PetscReal, Vec, void *);
971201da799SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS, PetscInt, PetscReal, Vec, void *);
972edbaebb3SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorError(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
973d0c080abSJoseph Pusztay PETSC_EXTERN PetscErrorCode TSDMSwarmMonitorMoments(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
974ef20d060SBarry Smith 
975e669de00SBarry Smith struct _n_TSMonitorLGCtxNetwork {
976e669de00SBarry Smith   PetscInt     nlg;
977e669de00SBarry Smith   PetscDrawLG *lg;
978e669de00SBarry Smith   PetscBool    semilogy;
979e669de00SBarry Smith   PetscInt     howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
980e669de00SBarry Smith };
981e669de00SBarry Smith typedef struct _n_TSMonitorLGCtxNetwork *TSMonitorLGCtxNetwork;
982e669de00SBarry Smith PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkDestroy(TSMonitorLGCtxNetwork *);
983e669de00SBarry Smith PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkCreate(TS, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtxNetwork *);
984e669de00SBarry Smith PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkSolution(TS, PetscInt, PetscReal, Vec, void *);
985e669de00SBarry Smith 
986b3d3934dSBarry Smith typedef struct _n_TSMonitorEnvelopeCtx *TSMonitorEnvelopeCtx;
987b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeCtxCreate(TS, TSMonitorEnvelopeCtx *);
988b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode             TSMonitorEnvelope(TS, PetscInt, PetscReal, Vec, void *);
989b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeGetBounds(TS, Vec *, Vec *);
990b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *);
991b3d3934dSBarry Smith 
9928189c53fSBarry Smith typedef struct _n_TSMonitorSPEigCtx *TSMonitorSPEigCtx;
9938189c53fSBarry Smith PETSC_EXTERN PetscErrorCode          TSMonitorSPEigCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorSPEigCtx *);
9948189c53fSBarry Smith PETSC_EXTERN PetscErrorCode          TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *);
9958189c53fSBarry Smith PETSC_EXTERN PetscErrorCode          TSMonitorSPEig(TS, PetscInt, PetscReal, Vec, void *);
9963914022bSBarry Smith 
9971b575b74SJoseph Pusztay typedef struct _n_TSMonitorSPCtx *TSMonitorSPCtx;
99860e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorSPCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscBool, PetscBool, TSMonitorSPCtx *);
9991b575b74SJoseph Pusztay PETSC_EXTERN PetscErrorCode       TSMonitorSPCtxDestroy(TSMonitorSPCtx *);
10000ec8ee2bSJoseph Pusztay PETSC_EXTERN PetscErrorCode       TSMonitorSPSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
10011b575b74SJoseph Pusztay 
100260e16b1bSMatthew G. Knepley typedef struct _n_TSMonitorHGCtx *TSMonitorHGCtx;
100360e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorHGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscInt, PetscBool, TSMonitorHGCtx *);
100460e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
100560e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorHGCtxDestroy(TSMonitorHGCtx *);
100660e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
100760e16b1bSMatthew G. Knepley 
1008ca4445c7SIlya Fursov PETSC_EXTERN PetscErrorCode TSSetEventHandler(TS, PetscInt, PetscInt[], PetscBool[], PetscErrorCode (*)(TS, PetscReal, Vec, PetscReal[], void *), PetscErrorCode (*)(TS, PetscInt, PetscInt[], PetscReal, Vec, PetscBool, void *), void *);
1009ca4445c7SIlya Fursov PETSC_EXTERN PetscErrorCode TSSetPostEventStep(TS, PetscReal);
1010fe4ad979SIlya Fursov PETSC_EXTERN PetscErrorCode TSSetPostEventSecondStep(TS, PetscReal);
1011fe4ad979SIlya Fursov PETSC_DEPRECATED_FUNCTION(3, 21, 0, "TSSetPostEventSecondStep()", ) static inline PetscErrorCode TSSetPostEventIntervalStep(TS ts, PetscReal dt)
1012fe4ad979SIlya Fursov {
1013fe4ad979SIlya Fursov   return TSSetPostEventSecondStep(ts, dt);
1014fe4ad979SIlya Fursov }
10156427ac75SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS, PetscReal, PetscReal[]);
10161ea83e56SMiguel PETSC_EXTERN PetscErrorCode TSGetNumEvents(TS, PetscInt *);
10171ea83e56SMiguel 
101876bdecfbSBarry Smith /*J
1019195e9b02SBarry Smith    TSSSPType - string with the name of a `TSSSP` scheme.
1020815f1ad5SJed Brown 
1021815f1ad5SJed Brown    Level: beginner
1022815f1ad5SJed Brown 
10231cc06b55SBarry Smith .seealso: [](ch_ts), `TSSSPSetType()`, `TS`, `TSSSP`
102476bdecfbSBarry Smith J*/
102519fd82e9SBarry Smith typedef const char *TSSSPType;
1026815f1ad5SJed Brown #define TSSSPRKS2  "rks2"
1027815f1ad5SJed Brown #define TSSSPRKS3  "rks3"
1028815f1ad5SJed Brown #define TSSSPRK104 "rk104"
1029815f1ad5SJed Brown 
103019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode    TSSSPSetType(TS, TSSSPType);
103119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode    TSSSPGetType(TS, TSSSPType *);
1032014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSSPSetNumStages(TS, PetscInt);
1033014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSSPGetNumStages(TS, PetscInt *);
1034787849ffSJed Brown PETSC_EXTERN PetscErrorCode    TSSSPInitializePackage(void);
1035560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode    TSSSPFinalizePackage(void);
1036013797aaSBarry Smith PETSC_EXTERN PetscFunctionList TSSSPList;
1037815f1ad5SJed Brown 
1038ed657a08SJed Brown /*S
103984df9cb4SJed Brown    TSAdapt - Abstract object that manages time-step adaptivity
104084df9cb4SJed Brown 
104184df9cb4SJed Brown    Level: beginner
104284df9cb4SJed Brown 
1043af27ebaaSBarry Smith .seealso: [](ch_ts), `TS`, `TSGetAdapt()`, `TSAdaptCreate()`, `TSAdaptType`
104484df9cb4SJed Brown S*/
104584df9cb4SJed Brown typedef struct _p_TSAdapt *TSAdapt;
104684df9cb4SJed Brown 
1047f0fc11ceSJed Brown /*J
104887497f52SBarry Smith    TSAdaptType - String with the name of `TSAdapt` scheme.
104984df9cb4SJed Brown 
105084df9cb4SJed Brown    Level: beginner
105184df9cb4SJed Brown 
1052af27ebaaSBarry Smith .seealso: [](ch_ts), `TSGetAdapt()`, `TSAdaptSetType()`, `TS`, `TSAdapt`
1053f0fc11ceSJed Brown J*/
1054f8963224SJed Brown typedef const char *TSAdaptType;
1055b0f836d7SJed Brown #define TSADAPTNONE    "none"
10561917a363SLisandro Dalcin #define TSADAPTBASIC   "basic"
1057cb7ab186SLisandro Dalcin #define TSADAPTDSP     "dsp"
10588d59e960SJed Brown #define TSADAPTCFL     "cfl"
10591917a363SLisandro Dalcin #define TSADAPTGLEE    "glee"
1060d949e4a4SStefano Zampini #define TSADAPTHISTORY "history"
106184df9cb4SJed Brown 
1062552698daSJed Brown PETSC_EXTERN PetscErrorCode TSGetAdapt(TS, TSAdapt *);
1063bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[], PetscErrorCode (*)(TSAdapt));
1064607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void);
1065014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
1066014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm, TSAdapt *);
106719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt, TSAdaptType);
1068d0288e4fSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetType(TSAdapt, TSAdaptType *);
1069014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt, const char[]);
1070014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
1071014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt, const char[], PetscInt, PetscInt, PetscReal, PetscReal, PetscBool);
1072014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt, PetscInt *, const PetscInt **, const PetscInt **, const PetscReal **, const PetscReal **);
1073014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt, TS, PetscReal, PetscInt *, PetscReal *, PetscBool *);
1074b295832fSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt, TS, PetscReal, Vec, PetscBool *);
1075014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt, PetscViewer);
1076ad6bc421SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt, PetscViewer);
1077dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(TSAdapt, PetscOptionItems *);
1078099af0a3SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt);
1079014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt *);
1080014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt, PetscBool);
1081bf997491SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt, PetscBool);
10821917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptSetSafety(TSAdapt, PetscReal, PetscReal);
10831917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetSafety(TSAdapt, PetscReal *, PetscReal *);
108476cddca1SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt, PetscReal);
108576cddca1SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt, PetscReal *);
10861917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptSetClip(TSAdapt, PetscReal, PetscReal);
10871917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetClip(TSAdapt, PetscReal *, PetscReal *);
108862c23b28SMark PETSC_EXTERN PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt, PetscReal);
108962c23b28SMark PETSC_EXTERN PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt, PetscReal *);
1090014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt, PetscReal, PetscReal);
10911917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetStepLimits(TSAdapt, PetscReal *, PetscReal *);
1092b295832fSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt, PetscErrorCode (*)(TSAdapt, TS, PetscReal, Vec, PetscBool *));
1093d949e4a4SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptHistorySetHistory(TSAdapt, PetscInt n, PetscReal hist[], PetscBool);
109475017583SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptHistorySetTrajectory(TSAdapt, TSTrajectory, PetscBool);
109575017583SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptHistoryGetStep(TSAdapt, PetscInt, PetscReal *, PetscReal *);
1096de50f1caSBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt, PetscInt);
1097cb7ab186SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptDSPSetFilter(TSAdapt, const char *);
1098cb7ab186SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptDSPSetPID(TSAdapt, PetscReal, PetscReal, PetscReal);
109984df9cb4SJed Brown 
110084df9cb4SJed Brown /*S
110187497f52SBarry Smith    TSGLLEAdapt - Abstract object that manages time-step adaptivity for `TSGLLE`
1102ed657a08SJed Brown 
1103ed657a08SJed Brown    Level: beginner
1104ed657a08SJed Brown 
1105195e9b02SBarry Smith    Developer Note:
110687497f52SBarry Smith    This functionality should be replaced by the `TSAdapt`.
110784df9cb4SJed Brown 
11081cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLEAdaptCreate()`, `TSGLLEAdaptType`
1109ed657a08SJed Brown S*/
111026d28e4eSEmil Constantinescu typedef struct _p_TSGLLEAdapt *TSGLLEAdapt;
1111ed657a08SJed Brown 
111276bdecfbSBarry Smith /*J
111387497f52SBarry Smith    TSGLLEAdaptType - String with the name of `TSGLLEAdapt` scheme
1114ed657a08SJed Brown 
1115ed657a08SJed Brown    Level: beginner
1116ed657a08SJed Brown 
1117195e9b02SBarry Smith    Developer Note:
1118195e9b02SBarry Smith    This functionality should be replaced by the `TSAdaptType`.
1119195e9b02SBarry Smith 
11201cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLLEAdaptSetType()`, `TS`
112176bdecfbSBarry Smith J*/
112226d28e4eSEmil Constantinescu typedef const char *TSGLLEAdaptType;
112326d28e4eSEmil Constantinescu #define TSGLLEADAPT_NONE "none"
112426d28e4eSEmil Constantinescu #define TSGLLEADAPT_SIZE "size"
112526d28e4eSEmil Constantinescu #define TSGLLEADAPT_BOTH "both"
1126ed657a08SJed Brown 
112726d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptRegister(const char[], PetscErrorCode (*)(TSGLLEAdapt));
112826d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptInitializePackage(void);
112926d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptFinalizePackage(void);
113026d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate(MPI_Comm, TSGLLEAdapt *);
113126d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt, TSGLLEAdaptType);
113226d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt, const char[]);
113326d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptChoose(TSGLLEAdapt, PetscInt, const PetscInt[], const PetscReal[], const PetscReal[], PetscInt, PetscReal, PetscReal, PetscInt *, PetscReal *, PetscBool *);
113426d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt, PetscViewer);
1135dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetFromOptions(TSGLLEAdapt, PetscOptionItems *);
113626d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt *);
1137ed657a08SJed Brown 
113876bdecfbSBarry Smith /*J
113987497f52SBarry Smith    TSGLLEAcceptType - String with the name of `TSGLLEAccept` scheme
1140ed657a08SJed Brown 
1141ed657a08SJed Brown    Level: beginner
1142ed657a08SJed Brown 
11431cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLLESetAcceptType()`, `TS`, `TSGLLEAccept`
114476bdecfbSBarry Smith J*/
114526d28e4eSEmil Constantinescu typedef const char *TSGLLEAcceptType;
114626d28e4eSEmil Constantinescu #define TSGLLEACCEPT_ALWAYS "always"
1147ed657a08SJed Brown 
1148d1c5d1fcSBarry Smith /*S
11498434afd1SBarry Smith   TSGLLEAcceptFn - A prototype of a `TS` accept function that would be passed to `TSGLLEAcceptRegister()`
1150d1c5d1fcSBarry Smith 
1151d1c5d1fcSBarry Smith   Calling Sequence:
1152d1c5d1fcSBarry Smith + ts  - timestep context
1153d1c5d1fcSBarry Smith . nt - time to end of solution time
1154d1c5d1fcSBarry Smith . h - the proposed step-size
1155d1c5d1fcSBarry Smith . enorm - unknown
1156d1c5d1fcSBarry Smith - accept - output, if the proposal is accepted
1157d1c5d1fcSBarry Smith 
1158d1c5d1fcSBarry Smith   Level: beginner
1159d1c5d1fcSBarry Smith 
1160d1c5d1fcSBarry Smith   Note:
11618434afd1SBarry Smith   The deprecated `TSGLLEAcceptFunction` still works as a replacement for `TSGLLEAcceptFn` *
1162d1c5d1fcSBarry Smith 
11638434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetRHSFunction()`, `DMTSSetRHSFunction()`, `TSIFunctionFn`,
11648434afd1SBarry Smith `TSIJacobianFn`, `TSRHSJacobianFn`, `TSGLLEAcceptRegister()`
1165d1c5d1fcSBarry Smith S*/
11668434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSGLLEAcceptFn)(TS ts, PetscReal nt, PetscReal h, const PetscReal enorm[], PetscBool *accept);
1167d1c5d1fcSBarry Smith 
11688434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSGLLEAcceptFn *TSGLLEAcceptFunction;
1169d1c5d1fcSBarry Smith 
11708434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGLLEAcceptRegister(const char[], TSGLLEAcceptFn *);
1171ed657a08SJed Brown 
117276bdecfbSBarry Smith /*J
1173195e9b02SBarry Smith   TSGLLEType - string with the name of a General Linear `TSGLLE` type
117418b56cb9SJed Brown 
117518b56cb9SJed Brown   Level: beginner
117618b56cb9SJed Brown 
11771cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLESetType()`, `TSGLLERegister()`, `TSGLLEAccept`
117876bdecfbSBarry Smith J*/
117926d28e4eSEmil Constantinescu typedef const char *TSGLLEType;
118026d28e4eSEmil Constantinescu #define TSGLLE_IRKS "irks"
118118b56cb9SJed Brown 
118226d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLERegister(const char[], PetscErrorCode (*)(TS));
118326d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEInitializePackage(void);
118426d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEFinalizePackage(void);
118526d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLESetType(TS, TSGLLEType);
118626d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEGetAdapt(TS, TSGLLEAdapt *);
118726d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLESetAcceptType(TS, TSGLLEAcceptType);
118818b56cb9SJed Brown 
1189020d8f30SJed Brown /*J
1190195e9b02SBarry Smith    TSEIMEXType - String with the name of an Extrapolated IMEX `TSEIMEX` type
11917d5697caSHong Zhang 
11927d5697caSHong Zhang    Level: beginner
11937d5697caSHong Zhang 
11941cc06b55SBarry Smith .seealso: [](ch_ts), `TSEIMEXSetType()`, `TS`, `TSEIMEX`, `TSEIMEXRegister()`
11957d5697caSHong Zhang J*/
11967d5697caSHong Zhang #define TSEIMEXType char *
11977d5697caSHong Zhang 
11987d5697caSHong Zhang PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts, PetscInt);
11997d5697caSHong Zhang PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts, PetscInt, PetscInt);
12007d5697caSHong Zhang PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS, PetscBool);
12017d5697caSHong Zhang 
12027d5697caSHong Zhang /*J
1203195e9b02SBarry Smith    TSRKType - String with the name of a Runge-Kutta `TSRK` type
1204f68a32c8SEmil Constantinescu 
1205f68a32c8SEmil Constantinescu    Level: beginner
1206f68a32c8SEmil Constantinescu 
12071cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSRKSetType()`, `TS`, `TSRK`, `TSRKRegister()`
1208f68a32c8SEmil Constantinescu J*/
1209f68a32c8SEmil Constantinescu typedef const char *TSRKType;
1210f68a32c8SEmil Constantinescu #define TSRK1FE "1fe"
1211fdefd5e5SDebojyoti Ghosh #define TSRK2A  "2a"
1212e7685601SHong Zhang #define TSRK2B  "2b"
1213f68a32c8SEmil Constantinescu #define TSRK3   "3"
1214fdefd5e5SDebojyoti Ghosh #define TSRK3BS "3bs"
1215f68a32c8SEmil Constantinescu #define TSRK4   "4"
1216fdefd5e5SDebojyoti Ghosh #define TSRK5F  "5f"
1217fdefd5e5SDebojyoti Ghosh #define TSRK5DP "5dp"
121805e23783SLisandro Dalcin #define TSRK5BS "5bs"
121937acaa02SHendrik Ranocha #define TSRK6VR "6vr"
122037acaa02SHendrik Ranocha #define TSRK7VR "7vr"
122137acaa02SHendrik Ranocha #define TSRK8VR "8vr"
122205e23783SLisandro Dalcin 
12232ea87ba9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRKGetOrder(TS, PetscInt *);
12240f7a28cdSStefano Zampini PETSC_EXTERN PetscErrorCode TSRKGetType(TS, TSRKType *);
12250f7a28cdSStefano Zampini PETSC_EXTERN PetscErrorCode TSRKSetType(TS, TSRKType);
12260f7a28cdSStefano Zampini PETSC_EXTERN PetscErrorCode TSRKGetTableau(TS, PetscInt *, const PetscReal **, const PetscReal **, const PetscReal **, const PetscReal **, PetscInt *, const PetscReal **, PetscBool *);
12270fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSRKSetMultirate(TS, PetscBool);
12280fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSRKGetMultirate(TS, PetscBool *);
1229f68a32c8SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
1230f68a32c8SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void);
1231560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void);
1232f68a32c8SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void);
1233f68a32c8SEmil Constantinescu 
1234f68a32c8SEmil Constantinescu /*J
1235af27ebaaSBarry Smith    TSMPRKType - String with the name of a partitioned Runge-Kutta `TSMPRK` type
12369f537349Sluzhanghpp 
12379f537349Sluzhanghpp    Level: beginner
12389f537349Sluzhanghpp 
12391cc06b55SBarry Smith .seealso: [](ch_ts), `TSMPRKSetType()`, `TS`, `TSMPRK`, `TSMPRKRegister()`
12409f537349Sluzhanghpp J*/
12414b84eec9SHong Zhang typedef const char *TSMPRKType;
124219c2959aSHong Zhang #define TSMPRK2A22 "2a22"
124319c2959aSHong Zhang #define TSMPRK2A23 "2a23"
124419c2959aSHong Zhang #define TSMPRK2A32 "2a32"
124519c2959aSHong Zhang #define TSMPRK2A33 "2a33"
12464b84eec9SHong Zhang #define TSMPRKP2   "p2"
12474b84eec9SHong Zhang #define TSMPRKP3   "p3"
12489f537349Sluzhanghpp 
12494b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKGetType(TS ts, TSMPRKType *);
12504b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKSetType(TS ts, TSMPRKType);
125179bc8a77SHong 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[]);
12524b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKInitializePackage(void);
12534b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKFinalizePackage(void);
12544b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKRegisterDestroy(void);
12559f537349Sluzhanghpp 
12569f537349Sluzhanghpp /*J
1257195e9b02SBarry Smith    TSIRKType - String with the name of an implicit Runge-Kutta `TSIRK` type
1258d2567f34SHong Zhang 
1259d2567f34SHong Zhang    Level: beginner
1260d2567f34SHong Zhang 
12611cc06b55SBarry Smith .seealso: [](ch_ts), `TSIRKSetType()`, `TS`, `TSIRK`, `TSIRKRegister()`
1262d2567f34SHong Zhang J*/
1263d2567f34SHong Zhang typedef const char *TSIRKType;
1264d2567f34SHong Zhang #define TSIRKGAUSS "gauss"
1265d2567f34SHong Zhang 
1266d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKGetType(TS, TSIRKType *);
1267d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKSetType(TS, TSIRKType);
1268d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKGetNumStages(TS, PetscInt *);
1269d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKSetNumStages(TS, PetscInt);
1270d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKRegister(const char[], PetscErrorCode (*function)(TS));
1271d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKTableauCreate(TS, PetscInt, const PetscReal *, const PetscReal *, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, const PetscScalar *);
1272d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKInitializePackage(void);
1273d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKFinalizePackage(void);
1274d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKRegisterDestroy(void);
1275d2567f34SHong Zhang 
1276d2567f34SHong Zhang /*J
1277195e9b02SBarry Smith    TSGLEEType - String with the name of a General Linear with Error Estimation `TSGLEE` type
1278b6a60446SDebojyoti Ghosh 
1279b6a60446SDebojyoti Ghosh    Level: beginner
1280b6a60446SDebojyoti Ghosh 
12811cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEESetType()`, `TS`, `TSGLEE`, `TSGLEERegister()`
1282b6a60446SDebojyoti Ghosh J*/
1283b6a60446SDebojyoti Ghosh typedef const char *TSGLEEType;
1284ab8c5912SEmil Constantinescu #define TSGLEEi1      "BE1"
1285b6a60446SDebojyoti Ghosh #define TSGLEE23      "23"
1286b6a60446SDebojyoti Ghosh #define TSGLEE24      "24"
1287b6a60446SDebojyoti Ghosh #define TSGLEE25I     "25i"
1288b6a60446SDebojyoti Ghosh #define TSGLEE35      "35"
1289b6a60446SDebojyoti Ghosh #define TSGLEEEXRK2A  "exrk2a"
1290b6a60446SDebojyoti Ghosh #define TSGLEERK32G1  "rk32g1"
1291b6a60446SDebojyoti Ghosh #define TSGLEERK285EX "rk285ex"
129216a05f60SBarry Smith 
1293b6a60446SDebojyoti Ghosh /*J
1294195e9b02SBarry Smith    TSGLEEMode - String with the mode of error estimation for a General Linear with Error Estimation `TSGLEE` type
1295b6a60446SDebojyoti Ghosh 
1296b6a60446SDebojyoti Ghosh    Level: beginner
1297b6a60446SDebojyoti Ghosh 
12981cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEESetMode()`, `TS`, `TSGLEE`, `TSGLEERegister()`
1299b6a60446SDebojyoti Ghosh J*/
1300b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEEGetType(TS ts, TSGLEEType *);
1301b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEESetType(TS ts, TSGLEEType);
130257df6a1bSDebojyoti 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[]);
13034bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSGLEERegisterAll(void);
1304b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEEFinalizePackage(void);
1305b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEEInitializePackage(void);
1306b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEERegisterDestroy(void);
1307b6a60446SDebojyoti Ghosh 
1308b6a60446SDebojyoti Ghosh /*J
1309195e9b02SBarry Smith    TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX `TSARKIMEX` type
1310020d8f30SJed Brown 
1311020d8f30SJed Brown    Level: beginner
1312020d8f30SJed Brown 
13131cc06b55SBarry Smith .seealso: [](ch_ts), `TSARKIMEXSetType()`, `TS`, `TSARKIMEX`, `TSARKIMEXRegister()`
1314020d8f30SJed Brown J*/
131519fd82e9SBarry Smith typedef const char *TSARKIMEXType;
1316e817cc15SEmil Constantinescu #define TSARKIMEX1BEE   "1bee"
13171f80e275SEmil Constantinescu #define TSARKIMEXA2     "a2"
13181f80e275SEmil Constantinescu #define TSARKIMEXL2     "l2"
13191f80e275SEmil Constantinescu #define TSARKIMEXARS122 "ars122"
13201f80e275SEmil Constantinescu #define TSARKIMEX2C     "2c"
13218a381b04SJed Brown #define TSARKIMEX2D     "2d"
1322a3a57f36SJed Brown #define TSARKIMEX2E     "2e"
13236cf0794eSJed Brown #define TSARKIMEXPRSSP2 "prssp2"
13248a381b04SJed Brown #define TSARKIMEX3      "3"
13256cf0794eSJed Brown #define TSARKIMEXBPR3   "bpr3"
13266cf0794eSJed Brown #define TSARKIMEXARS443 "ars443"
13278a381b04SJed Brown #define TSARKIMEX4      "4"
13288a381b04SJed Brown #define TSARKIMEX5      "5"
132919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts, TSARKIMEXType *);
133019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts, TSARKIMEXType);
1331014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS, PetscBool);
13323a28c0e5SStefano Zampini PETSC_EXTERN PetscErrorCode TSARKIMEXGetFullyImplicit(TS, PetscBool *);
13333a2a065fSHong Zhang PETSC_EXTERN PetscErrorCode TSARKIMEXSetFastSlowSplit(TS, PetscBool);
13343a2a065fSHong Zhang PETSC_EXTERN PetscErrorCode TSARKIMEXGetFastSlowSplit(TS, PetscBool *);
133519fd82e9SBarry 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[]);
1336607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void);
1337560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
1338014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);
13398a381b04SJed Brown 
1340020d8f30SJed Brown /*J
1341d27334e2SStefano Zampini    TSDIRKType - String with the name of a Diagonally Implicit Runge-Kutta `TSDIRK` type
1342d27334e2SStefano Zampini 
1343d27334e2SStefano Zampini    Level: beginner
1344d27334e2SStefano Zampini 
1345d27334e2SStefano Zampini .seealso: [](ch_ts), `TSDIRKSetType()`, `TS`, `TSDIRK`, `TSDIRKRegister()`
1346d27334e2SStefano Zampini J*/
1347d27334e2SStefano Zampini typedef const char *TSDIRKType;
1348d27334e2SStefano Zampini #define TSDIRKS212      "s212"
13493405e92cSStefano Zampini #define TSDIRKES122SAL  "es122sal"
1350d27334e2SStefano Zampini #define TSDIRKES213SAL  "es213sal"
1351d27334e2SStefano Zampini #define TSDIRKES324SAL  "es324sal"
1352d27334e2SStefano Zampini #define TSDIRKES325SAL  "es325sal"
1353d27334e2SStefano Zampini #define TSDIRK657A      "657a"
1354d27334e2SStefano Zampini #define TSDIRKES648SA   "es648sa"
1355d27334e2SStefano Zampini #define TSDIRK658A      "658a"
1356d27334e2SStefano Zampini #define TSDIRKS659A     "s659a"
1357d27334e2SStefano Zampini #define TSDIRK7510SAL   "7510sal"
1358d27334e2SStefano Zampini #define TSDIRKES7510SA  "es7510sa"
1359d27334e2SStefano Zampini #define TSDIRK759A      "759a"
1360d27334e2SStefano Zampini #define TSDIRKS7511SAL  "s7511sal"
1361d27334e2SStefano Zampini #define TSDIRK8614A     "8614a"
1362d27334e2SStefano Zampini #define TSDIRK8616SAL   "8616sal"
1363d27334e2SStefano Zampini #define TSDIRKES8516SAL "es8516sal"
1364d27334e2SStefano Zampini PETSC_EXTERN PetscErrorCode TSDIRKGetType(TS ts, TSDIRKType *);
1365d27334e2SStefano Zampini PETSC_EXTERN PetscErrorCode TSDIRKSetType(TS ts, TSDIRKType);
1366d27334e2SStefano Zampini PETSC_EXTERN PetscErrorCode TSDIRKRegister(TSDIRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
1367d27334e2SStefano Zampini 
1368d27334e2SStefano Zampini /*J
1369195e9b02SBarry Smith    TSRosWType - String with the name of a Rosenbrock-W `TSROSW` type
1370020d8f30SJed Brown 
1371020d8f30SJed Brown    Level: beginner
1372020d8f30SJed Brown 
13731cc06b55SBarry Smith .seealso: [](ch_ts), `TSRosWSetType()`, `TS`, `TSROSW`, `TSRosWRegister()`
1374020d8f30SJed Brown J*/
137519fd82e9SBarry Smith typedef const char *TSRosWType;
137661692a83SJed Brown #define TSROSW2M          "2m"
137761692a83SJed Brown #define TSROSW2P          "2p"
1378fe7e6d57SJed Brown #define TSROSWRA3PW       "ra3pw"
1379fe7e6d57SJed Brown #define TSROSWRA34PW2     "ra34pw2"
1380337ef93bSJed Brown #define TSROSWR34PRW      "r34prw"
1381337ef93bSJed Brown #define TSROSWR3PRL2      "r3prl2"
1382ef3c5b88SJed Brown #define TSROSWRODAS3      "rodas3"
138348665b53SJed Brown #define TSROSWRODASPR     "rodaspr"
1384337ef93bSJed Brown #define TSROSWRODASPR2    "rodaspr2"
1385ef3c5b88SJed Brown #define TSROSWSANDU3      "sandu3"
1386961f28d0SJed Brown #define TSROSWASSP3P3S1C  "assp3p3s1c"
1387961f28d0SJed Brown #define TSROSWLASSP3P4S2C "lassp3p4s2c"
138843b21953SEmil Constantinescu #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
1389753f8adbSEmil Constantinescu #define TSROSWARK3        "ark3"
13903606a31eSEmil Constantinescu #define TSROSWTHETA1      "theta1"
13913606a31eSEmil Constantinescu #define TSROSWTHETA2      "theta2"
139242faf41dSJed Brown #define TSROSWGRK4T       "grk4t"
139342faf41dSJed Brown #define TSROSWSHAMP4      "shamp4"
139442faf41dSJed Brown #define TSROSWVELDD4      "veldd4"
139542faf41dSJed Brown #define TSROSW4L          "4l"
1396b1c69cc3SEmil Constantinescu 
13976bd7aeb5SHong Zhang PETSC_EXTERN PetscErrorCode TSRosWGetType(TS, TSRosWType *);
13986bd7aeb5SHong Zhang PETSC_EXTERN PetscErrorCode TSRosWSetType(TS, TSRosWType);
1399014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS, PetscBool);
140019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
140119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal);
1402607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void);
1403560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
1404014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);
1405e27a552bSJed Brown 
1406211a84d6SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSBDFSetOrder(TS, PetscInt);
1407211a84d6SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSBDFGetOrder(TS, PetscInt *);
1408211a84d6SLisandro Dalcin 
14096bd7aeb5SHong Zhang /*J
1410195e9b02SBarry Smith   TSBasicSymplecticType - String with the name of a basic symplectic integration `TSBASICSYMPLECTIC` type
14116bd7aeb5SHong Zhang 
14126bd7aeb5SHong Zhang   Level: beginner
14136bd7aeb5SHong Zhang 
14141cc06b55SBarry Smith .seealso: [](ch_ts), `TSBasicSymplecticSetType()`, `TS`, `TSBASICSYMPLECTIC`, `TSBasicSymplecticRegister()`
14156bd7aeb5SHong Zhang J*/
1416e49d4f37SHong Zhang typedef const char *TSBasicSymplecticType;
1417e49d4f37SHong Zhang #define TSBASICSYMPLECTICSIEULER   "1"
1418e49d4f37SHong Zhang #define TSBASICSYMPLECTICVELVERLET "2"
1419e49d4f37SHong Zhang #define TSBASICSYMPLECTIC3         "3"
1420e49d4f37SHong Zhang #define TSBASICSYMPLECTIC4         "4"
1421e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticSetType(TS, TSBasicSymplecticType);
1422e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticGetType(TS, TSBasicSymplecticType *);
1423e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegister(TSBasicSymplecticType, PetscInt, PetscInt, PetscReal[], PetscReal[]);
14244bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterAll(void);
1425e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticInitializePackage(void);
1426e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticFinalizePackage(void);
1427e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterDestroy(void);
14286bd7aeb5SHong Zhang 
1429db4653c2SDaniel Finn /*J
1430195e9b02SBarry Smith   TSDISCGRAD - The Discrete Gradient integrator is a timestepper for Hamiltonian systems designed to conserve the first integral (energy),
143187497f52SBarry Smith   but also has the property for some systems of monotonicity in a functional.
1432db4653c2SDaniel Finn 
1433db4653c2SDaniel Finn   Level: beginner
1434db4653c2SDaniel Finn 
14351cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, TSDiscGradSetFormulation()`, `TSDiscGradGetFormulation()`
1436db4653c2SDaniel Finn J*/
143740be0ff1SMatthew 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 *);
14384bf303faSJacob 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 *);
1439db4653c2SDaniel Finn PETSC_EXTERN PetscErrorCode TSDiscGradIsGonzalez(TS, PetscBool *);
1440db4653c2SDaniel Finn PETSC_EXTERN PetscErrorCode TSDiscGradUseGonzalez(TS, PetscBool);
144140be0ff1SMatthew G. Knepley 
144283e2fdc7SBarry Smith /*
14430f3b3ca1SHong Zhang        PETSc interface to Sundials
144483e2fdc7SBarry Smith */
1445e808b789SPatrick Sanan #ifdef PETSC_HAVE_SUNDIALS2
14469371c9d4SSatish Balay typedef enum {
14479371c9d4SSatish Balay   SUNDIALS_ADAMS = 1,
14489371c9d4SSatish Balay   SUNDIALS_BDF   = 2
14499371c9d4SSatish Balay } TSSundialsLmmType;
14506a6fc655SJed Brown PETSC_EXTERN const char *const TSSundialsLmmTypes[];
14519371c9d4SSatish Balay typedef enum {
14529371c9d4SSatish Balay   SUNDIALS_MODIFIED_GS  = 1,
14539371c9d4SSatish Balay   SUNDIALS_CLASSICAL_GS = 2
14549371c9d4SSatish Balay } TSSundialsGramSchmidtType;
14556a6fc655SJed Brown PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
14564bf303faSJacob Faibussowitsch 
1457014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS, TSSundialsLmmType);
1458014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS, PC *);
1459014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS, PetscReal, PetscReal);
1460014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS, PetscReal);
1461014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS, PetscReal);
1462014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS, PetscInt *, PetscInt *);
1463014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS, TSSundialsGramSchmidtType);
1464014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS, PetscInt);
1465014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS, PetscReal);
1466014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS, PetscBool);
1467014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS, PetscInt);
1468c4e80e11SFlorian PETSC_EXTERN PetscErrorCode TSSundialsSetMaxord(TS, PetscInt);
1469918687b8SPatrick Sanan PETSC_EXTERN PetscErrorCode TSSundialsSetUseDense(TS, PetscBool);
14706dc17ff5SMatthew Knepley #endif
147183e2fdc7SBarry Smith 
1472014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS, PetscReal);
1473014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS, PetscReal *);
1474014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS, PetscBool *);
1475014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS, PetscBool);
14760de4c49aSJed Brown 
1477014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS, PetscReal);
1478014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS, PetscReal, PetscReal, PetscReal);
1479014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS, PetscReal *, PetscReal *, PetscReal *);
148088df8a41SLisandro Dalcin 
1481818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAlpha2SetRadius(TS, PetscReal);
1482818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAlpha2SetParams(TS, PetscReal, PetscReal, PetscReal, PetscReal);
1483818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAlpha2GetParams(TS, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
1484818efac9SLisandro Dalcin 
1485220f924aSDavid Kamensky /*S
14868434afd1SBarry Smith   TSAlpha2PredictorFn - A callback to set the predictor (i.e., the initial guess for the nonlinear solver) in
1487220f924aSDavid Kamensky   a second-order generalized-alpha time integrator.
1488220f924aSDavid Kamensky 
1489220f924aSDavid Kamensky   Calling Sequence:
1490220f924aSDavid Kamensky + ts   - the `TS` context obtained from `TSCreate()`
1491220f924aSDavid Kamensky . X0   - the previous time step's state vector
1492220f924aSDavid Kamensky . V0   - the previous time step's first derivative of the state vector
1493220f924aSDavid Kamensky . A0   - the previous time step's second derivative of the state vector
1494220f924aSDavid Kamensky . X1   - the vector into which the initial guess for the current time step will be written
1495220f924aSDavid Kamensky - ctx  - [optional] user-defined context for the predictor evaluation routine (may be `NULL`)
1496220f924aSDavid Kamensky 
1497220f924aSDavid Kamensky   Level: intermediate
1498220f924aSDavid Kamensky 
1499d1c5d1fcSBarry Smith   Note:
15008434afd1SBarry Smith   The deprecated `TSAlpha2Predictor` still works as a replacement for `TSAlpha2PredictorFn` *.
1501d1c5d1fcSBarry Smith 
1502220f924aSDavid Kamensky .seealso: [](ch_ts), `TS`, `TSAlpha2SetPredictor()`
1503220f924aSDavid Kamensky S*/
15048434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSAlpha2PredictorFn)(TS ts, Vec X0, Vec V0, Vec A0, Vec X1, void *ctx);
1505d1c5d1fcSBarry Smith 
15068434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSAlpha2PredictorFn *TSAlpha2Predictor;
1507d1c5d1fcSBarry Smith 
15088434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSAlpha2SetPredictor(TS, TSAlpha2PredictorFn *, void *ctx);
1509220f924aSDavid Kamensky 
1510014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetDM(TS, DM);
1511014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetDM(TS, DM *);
15126c699258SBarry Smith 
1513014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES, Vec, Vec, void *);
1514d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES, Vec, Mat, Mat, void *);
15150f5c6efeSJed Brown 
1516f3b1f45cSBarry Smith PETSC_EXTERN PetscErrorCode TSRHSJacobianTest(TS, PetscBool *);
1517f3b1f45cSBarry Smith PETSC_EXTERN PetscErrorCode TSRHSJacobianTestTranspose(TS, PetscBool *);
1518aad739acSMatthew G. Knepley 
15192e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSGetComputeInitialCondition(TS, PetscErrorCode (**)(TS, Vec));
15202e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSSetComputeInitialCondition(TS, PetscErrorCode (*)(TS, Vec));
15212e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSComputeInitialCondition(TS, Vec);
15222e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSGetComputeExactError(TS, PetscErrorCode (**)(TS, Vec, Vec));
1523aad739acSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSSetComputeExactError(TS, PetscErrorCode (*)(TS, Vec, Vec));
1524aad739acSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSComputeExactError(TS, Vec, Vec);
1525f2ed2dc7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscConvEstUseTS(PetscConvEst, PetscBool);
1526d60b7d5cSBarry Smith 
1527d60b7d5cSBarry Smith PETSC_EXTERN PetscErrorCode TSSetMatStructure(TS, MatStructure);
1528