xref: /petsc/include/petscts.h (revision c80d56d95e45277402b921e5bfcd988e56b8d46e)
1818ad0c1SBarry Smith /*
2316643e7SJed Brown    User interface for the timestepping package. This package
3f64a0f93SLois Curfman McInnes    is for use in solving time-dependent PDEs.
4818ad0c1SBarry Smith */
5a4963045SJacob Faibussowitsch #pragma once
6ac09b921SBarry Smith 
72c8e378dSBarry Smith #include <petscsnes.h>
8900f6b5bSMatthew G. Knepley #include <petscconvest.h>
9818ad0c1SBarry Smith 
1014d0ab18SJacob Faibussowitsch /*I <petscts.h> I*/
1114d0ab18SJacob Faibussowitsch 
12ac09b921SBarry Smith /* SUBMANSEC = TS */
13ac09b921SBarry Smith 
14435da068SBarry Smith /*S
150b4b7b1cSBarry Smith    TS - Abstract PETSc object that manages integrating an ODE.
16435da068SBarry Smith 
17435da068SBarry Smith    Level: beginner
18435da068SBarry Smith 
191cc06b55SBarry Smith .seealso: [](integrator_table), [](ch_ts), `TSCreate()`, `TSSetType()`, `TSType`, `SNES`, `KSP`, `PC`, `TSDestroy()`
20435da068SBarry Smith S*/
21f09e8eb9SSatish Balay typedef struct _p_TS *TS;
22435da068SBarry Smith 
2376bdecfbSBarry Smith /*J
240b4b7b1cSBarry Smith    TSType - String with the name of a PETSc `TS` method. These are all the time/ODE integrators that PETSc provides.
25435da068SBarry Smith 
26435da068SBarry Smith    Level: beginner
27435da068SBarry Smith 
280b4b7b1cSBarry Smith    Note:
290b4b7b1cSBarry Smith    Use `TSSetType()` or the options database key `-ts_type` to set the ODE integrator method to use with a given `TS` object
300b4b7b1cSBarry Smith 
311cc06b55SBarry Smith .seealso: [](integrator_table), [](ch_ts), `TSSetType()`, `TS`, `TSRegister()`
3276bdecfbSBarry Smith J*/
3319fd82e9SBarry Smith typedef const char *TSType;
349596e0b4SJed Brown #define TSEULER           "euler"
359596e0b4SJed Brown #define TSBEULER          "beuler"
36e49d4f37SHong Zhang #define TSBASICSYMPLECTIC "basicsymplectic"
379596e0b4SJed Brown #define TSPSEUDO          "pseudo"
384d91e141SJed Brown #define TSCN              "cn"
399596e0b4SJed Brown #define TSSUNDIALS        "sundials"
404d91e141SJed Brown #define TSRK              "rk"
419596e0b4SJed Brown #define TSPYTHON          "python"
429596e0b4SJed Brown #define TSTHETA           "theta"
4388df8a41SLisandro Dalcin #define TSALPHA           "alpha"
44818efac9SLisandro Dalcin #define TSALPHA2          "alpha2"
4526d28e4eSEmil Constantinescu #define TSGLLE            "glle"
46b6a60446SDebojyoti Ghosh #define TSGLEE            "glee"
47ef7bb5aaSJed Brown #define TSSSP             "ssp"
488a381b04SJed Brown #define TSARKIMEX         "arkimex"
49e27a552bSJed Brown #define TSROSW            "rosw"
507d5697caSHong Zhang #define TSEIMEX           "eimex"
51abadfa0fSMatthew G. Knepley #define TSMIMEX           "mimex"
52211a84d6SLisandro Dalcin #define TSBDF             "bdf"
53d249a78fSBarry Smith #define TSRADAU5          "radau5"
544b84eec9SHong Zhang #define TSMPRK            "mprk"
5540be0ff1SMatthew G. Knepley #define TSDISCGRAD        "discgrad"
56d2567f34SHong Zhang #define TSIRK             "irk"
57d27334e2SStefano Zampini #define TSDIRK            "dirk"
5840be0ff1SMatthew G. Knepley 
59435da068SBarry Smith /*E
6087497f52SBarry Smith    TSProblemType - Determines the type of problem this `TS` object is to be used to solve
61435da068SBarry Smith 
62195e9b02SBarry Smith    Values:
63195e9b02SBarry Smith  + `TS_LINEAR`    - a linear ODE or DAE
64195e9b02SBarry Smith  - `TS_NONLINEAR` - a nonlinear ODE or DAE
65195e9b02SBarry Smith 
66435da068SBarry Smith    Level: beginner
67435da068SBarry Smith 
681cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSCreate()`
69435da068SBarry Smith E*/
709371c9d4SSatish Balay typedef enum {
719371c9d4SSatish Balay   TS_LINEAR,
729371c9d4SSatish Balay   TS_NONLINEAR
739371c9d4SSatish Balay } TSProblemType;
74818ad0c1SBarry Smith 
75b93fea0eSJed Brown /*E
7687497f52SBarry Smith    TSEquationType - type of `TS` problem that is solved
77e817cc15SEmil Constantinescu 
78195e9b02SBarry Smith    Values:
79195e9b02SBarry Smith +  `TS_EQ_UNSPECIFIED` - (default)
8016a05f60SBarry Smith .  `TS_EQ_EXPLICIT`    - {ODE and DAE index 1, 2, 3, HI} F(t,U,U_t) := M(t) U_t - G(U,t) = 0
8116a05f60SBarry Smith -  `TS_EQ_IMPLICIT`    - {ODE and DAE index 1, 2, 3, HI} F(t,U,U_t) = 0
82e817cc15SEmil Constantinescu 
8395bd0b28SBarry Smith    Level: beginner
8495bd0b28SBarry Smith 
851cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGetEquationType()`, `TSSetEquationType()`
86e817cc15SEmil Constantinescu E*/
87e817cc15SEmil Constantinescu typedef enum {
88e817cc15SEmil Constantinescu   TS_EQ_UNSPECIFIED               = -1,
89e817cc15SEmil Constantinescu   TS_EQ_EXPLICIT                  = 0,
90e817cc15SEmil Constantinescu   TS_EQ_ODE_EXPLICIT              = 1,
91e817cc15SEmil Constantinescu   TS_EQ_DAE_SEMI_EXPLICIT_INDEX1  = 100,
92e817cc15SEmil Constantinescu   TS_EQ_DAE_SEMI_EXPLICIT_INDEX2  = 200,
93e817cc15SEmil Constantinescu   TS_EQ_DAE_SEMI_EXPLICIT_INDEX3  = 300,
94e817cc15SEmil Constantinescu   TS_EQ_DAE_SEMI_EXPLICIT_INDEXHI = 500,
95e817cc15SEmil Constantinescu   TS_EQ_IMPLICIT                  = 1000,
96e817cc15SEmil Constantinescu   TS_EQ_ODE_IMPLICIT              = 1001,
97e817cc15SEmil Constantinescu   TS_EQ_DAE_IMPLICIT_INDEX1       = 1100,
98e817cc15SEmil Constantinescu   TS_EQ_DAE_IMPLICIT_INDEX2       = 1200,
99e817cc15SEmil Constantinescu   TS_EQ_DAE_IMPLICIT_INDEX3       = 1300,
10019436ca2SJed Brown   TS_EQ_DAE_IMPLICIT_INDEXHI      = 1500
101e817cc15SEmil Constantinescu } TSEquationType;
102e817cc15SEmil Constantinescu PETSC_EXTERN const char *const *TSEquationTypes;
103e817cc15SEmil Constantinescu 
104e817cc15SEmil Constantinescu /*E
10516a05f60SBarry Smith    TSConvergedReason - reason a `TS` method has converged (integrated to the requested time) or not
106b93fea0eSJed Brown 
107195e9b02SBarry Smith    Values:
108195e9b02SBarry Smith +  `TS_CONVERGED_ITERATING`          - this only occurs if `TSGetConvergedReason()` is called during the `TSSolve()`
109195e9b02SBarry Smith .  `TS_CONVERGED_TIME`               - the final time was reached
110195e9b02SBarry Smith .  `TS_CONVERGED_ITS`                - the maximum number of iterations (time-steps) was reached prior to the final time
111195e9b02SBarry Smith .  `TS_CONVERGED_USER`               - user requested termination
112195e9b02SBarry Smith .  `TS_CONVERGED_EVENT`              - user requested termination on event detection
113195e9b02SBarry Smith .  `TS_CONVERGED_PSEUDO_FATOL`       - stops when function norm decreased by a set amount, used only for `TSPSEUDO`
114195e9b02SBarry Smith .  `TS_CONVERGED_PSEUDO_FRTOL`       - stops when function norm decreases below a set amount, used only for `TSPSEUDO`
115195e9b02SBarry Smith .  `TS_DIVERGED_NONLINEAR_SOLVE`     - too many nonlinear solve failures have occurred
116195e9b02SBarry Smith .  `TS_DIVERGED_STEP_REJECTED`       - too many steps were rejected
117195e9b02SBarry Smith .  `TSFORWARD_DIVERGED_LINEAR_SOLVE` - tangent linear solve failed
118195e9b02SBarry Smith -  `TSADJOINT_DIVERGED_LINEAR_SOLVE` - transposed linear solve failed
119195e9b02SBarry Smith 
120b93fea0eSJed Brown    Level: beginner
121b93fea0eSJed Brown 
1221cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGetConvergedReason()`
123b93fea0eSJed Brown E*/
124193ac0bcSJed Brown typedef enum {
125193ac0bcSJed Brown   TS_CONVERGED_ITERATING          = 0,
126193ac0bcSJed Brown   TS_CONVERGED_TIME               = 1,
127193ac0bcSJed Brown   TS_CONVERGED_ITS                = 2,
128d6ad946cSShri Abhyankar   TS_CONVERGED_USER               = 3,
1292b7db910SShri Abhyankar   TS_CONVERGED_EVENT              = 4,
1303118ae5eSBarry Smith   TS_CONVERGED_PSEUDO_FATOL       = 5,
1313118ae5eSBarry Smith   TS_CONVERGED_PSEUDO_FRTOL       = 6,
132193ac0bcSJed Brown   TS_DIVERGED_NONLINEAR_SOLVE     = -1,
133b5b4017aSHong Zhang   TS_DIVERGED_STEP_REJECTED       = -2,
13405c9950eSHong Zhang   TSFORWARD_DIVERGED_LINEAR_SOLVE = -3,
13505c9950eSHong Zhang   TSADJOINT_DIVERGED_LINEAR_SOLVE = -4
136193ac0bcSJed Brown } TSConvergedReason;
137014dd563SJed Brown PETSC_EXTERN const char *const *TSConvergedReasons;
1381957e957SBarry Smith 
139b93fea0eSJed Brown /*MC
14087497f52SBarry Smith    TS_CONVERGED_ITERATING - this only occurs if `TSGetConvergedReason()` is called during the `TSSolve()`
141b93fea0eSJed Brown 
142b93fea0eSJed Brown    Level: beginner
143b93fea0eSJed Brown 
1441cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`
145b93fea0eSJed Brown M*/
146b93fea0eSJed Brown 
147b93fea0eSJed Brown /*MC
148b93fea0eSJed Brown    TS_CONVERGED_TIME - the final time was reached
149b93fea0eSJed Brown 
150b93fea0eSJed Brown    Level: beginner
151b93fea0eSJed Brown 
1521cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxTime()`, `TSGetMaxTime()`, `TSGetSolveTime()`
153b93fea0eSJed Brown M*/
154b93fea0eSJed Brown 
155b93fea0eSJed Brown /*MC
1561957e957SBarry Smith    TS_CONVERGED_ITS - the maximum number of iterations (time-steps) was reached prior to the final time
157b93fea0eSJed Brown 
158b93fea0eSJed Brown    Level: beginner
159b93fea0eSJed Brown 
1601cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxSteps()`, `TSGetMaxSteps()`
161b93fea0eSJed Brown M*/
1621957e957SBarry Smith 
163d6ad946cSShri Abhyankar /*MC
164d6ad946cSShri Abhyankar    TS_CONVERGED_USER - user requested termination
165d6ad946cSShri Abhyankar 
166d6ad946cSShri Abhyankar    Level: beginner
167d6ad946cSShri Abhyankar 
1681cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`
169b93fea0eSJed Brown M*/
170b93fea0eSJed Brown 
171b93fea0eSJed Brown /*MC
172d76fc68cSShri Abhyankar    TS_CONVERGED_EVENT - user requested termination on event detection
173d76fc68cSShri Abhyankar 
174d76fc68cSShri Abhyankar    Level: beginner
175d76fc68cSShri Abhyankar 
1761cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`
177d76fc68cSShri Abhyankar M*/
178d76fc68cSShri Abhyankar 
179d76fc68cSShri Abhyankar /*MC
18087497f52SBarry Smith    TS_CONVERGED_PSEUDO_FRTOL - stops when function norm decreased by a set amount, used only for `TSPSEUDO`
1813118ae5eSBarry Smith 
1823c7db156SBarry Smith    Options Database Key:
183d5b43468SJose E. Roman .   -ts_pseudo_frtol <rtol> - use specified rtol
1843118ae5eSBarry Smith 
185195e9b02SBarry Smith    Level: beginner
186195e9b02SBarry Smith 
1871cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`, `TS_CONVERGED_PSEUDO_FATOL`
1883118ae5eSBarry Smith M*/
1893118ae5eSBarry Smith 
1903118ae5eSBarry Smith /*MC
19187497f52SBarry Smith    TS_CONVERGED_PSEUDO_FATOL - stops when function norm decreases below a set amount, used only for `TSPSEUDO`
1923118ae5eSBarry Smith 
1933c7db156SBarry Smith    Options Database Key:
19467b8a455SSatish Balay .   -ts_pseudo_fatol <atol> - use specified atol
1953118ae5eSBarry Smith 
196195e9b02SBarry Smith    Level: beginner
197195e9b02SBarry Smith 
1981cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`, `TS_CONVERGED_PSEUDO_FRTOL`
1993118ae5eSBarry Smith M*/
2003118ae5eSBarry Smith 
2013118ae5eSBarry Smith /*MC
202b93fea0eSJed Brown    TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed
203b93fea0eSJed Brown 
204b93fea0eSJed Brown    Level: beginner
205b93fea0eSJed Brown 
206195e9b02SBarry Smith    Note:
207195e9b02SBarry Smith    See `TSSetMaxSNESFailures()` for how to allow more nonlinear solver failures.
2081957e957SBarry Smith 
2091cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSGetSNES()`, `SNESGetConvergedReason()`, `TSSetMaxSNESFailures()`
210b93fea0eSJed Brown M*/
211b93fea0eSJed Brown 
212b93fea0eSJed Brown /*MC
213b93fea0eSJed Brown    TS_DIVERGED_STEP_REJECTED - too many steps were rejected
214b93fea0eSJed Brown 
215b93fea0eSJed Brown    Level: beginner
216b93fea0eSJed Brown 
21795452b02SPatrick Sanan    Notes:
218195e9b02SBarry Smith    See `TSSetMaxStepRejections()` for how to allow more step rejections.
2191957e957SBarry Smith 
2201cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxStepRejections()`
221b93fea0eSJed Brown M*/
222b93fea0eSJed Brown 
223d6ad946cSShri Abhyankar /*E
224d6ad946cSShri Abhyankar    TSExactFinalTimeOption - option for handling of final time step
225d6ad946cSShri Abhyankar 
226195e9b02SBarry Smith    Values:
22716a05f60SBarry Smith +  `TS_EXACTFINALTIME_STEPOVER`    - Don't do anything if requested final time is exceeded
228195e9b02SBarry Smith .  `TS_EXACTFINALTIME_INTERPOLATE` - Interpolate back to final time
22916a05f60SBarry Smith -  `TS_EXACTFINALTIME_MATCHSTEP`   - Adapt final time step to match the final time requested
230195e9b02SBarry Smith 
231d6ad946cSShri Abhyankar    Level: beginner
232d6ad946cSShri Abhyankar 
2331cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGetConvergedReason()`, `TSSetExactFinalTime()`, `TSGetExactFinalTime()`
234d6ad946cSShri Abhyankar E*/
2359371c9d4SSatish Balay typedef enum {
2369371c9d4SSatish Balay   TS_EXACTFINALTIME_UNSPECIFIED = 0,
2379371c9d4SSatish Balay   TS_EXACTFINALTIME_STEPOVER    = 1,
2389371c9d4SSatish Balay   TS_EXACTFINALTIME_INTERPOLATE = 2,
2399371c9d4SSatish Balay   TS_EXACTFINALTIME_MATCHSTEP   = 3
2409371c9d4SSatish Balay } TSExactFinalTimeOption;
241d6ad946cSShri Abhyankar PETSC_EXTERN const char *const TSExactFinalTimeOptions[];
242d6ad946cSShri Abhyankar 
243000e7ae3SMatthew Knepley /* Logging support */
244014dd563SJed Brown PETSC_EXTERN PetscClassId TS_CLASSID;
245d74926cbSBarry Smith PETSC_EXTERN PetscClassId DMTS_CLASSID;
2461b9b13dfSLisandro Dalcin PETSC_EXTERN PetscClassId TSADAPT_CLASSID;
247000e7ae3SMatthew Knepley 
248607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSInitializePackage(void);
249560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSFinalizePackage(void);
2508ba1e511SMatthew Knepley 
251014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm, TS *);
252baa10174SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSClone(TS, TS *);
253014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSDestroy(TS *);
254818ad0c1SBarry Smith 
255014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetProblemType(TS, TSProblemType);
256014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetProblemType(TS, TSProblemType *);
257014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSMonitor(TS, PetscInt, PetscReal, Vec);
25849abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorSet(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, void *), void *, PetscCtxDestroyFn *);
259721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorSetFromOptions(TS, const char[], const char[], const char[], PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *), PetscErrorCode (*)(TS, PetscViewerAndFormat *));
260014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS);
261818ad0c1SBarry Smith 
262014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS, const char[]);
263014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS, const char[]);
264014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS, const char *[]);
265014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS);
266014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetUp(TS);
267014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSReset(TS);
268818ad0c1SBarry Smith 
269014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetSolution(TS, Vec);
270014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetSolution(TS, Vec *);
271818ad0c1SBarry Smith 
272efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TS2SetSolution(TS, Vec, Vec);
273efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TS2GetSolution(TS, Vec *, Vec *);
27457df6a1bSDebojyoti Ghosh 
275642f3702SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetSolutionComponents(TS, PetscInt *, Vec *);
276642f3702SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetAuxSolution(TS, Vec *);
2770a01e1b2SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetTimeError(TS, PetscInt, Vec *);
27857df6a1bSDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSSetTimeError(TS, Vec);
2796e2e232bSDebojyoti Ghosh 
280a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSSetRHSJacobianP(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), void *);
281cd4cee2dSHong Zhang PETSC_EXTERN PetscErrorCode TSGetRHSJacobianP(TS, Mat *, PetscErrorCode (**)(TS, PetscReal, Vec, Mat, void *), void **);
282a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianP(TS, PetscReal, Vec, Mat);
28333f72c5dSHong Zhang PETSC_EXTERN PetscErrorCode TSSetIJacobianP(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Vec, PetscReal, Mat, void *), void *);
28433f72c5dSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIJacobianP(TS, PetscReal, Vec, Vec, PetscReal, Mat, PetscBool);
285edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS() then TSComputeRHSJacobianP()", ) PetscErrorCode TSComputeDRDPFunction(TS, PetscReal, Vec, Vec *);
286edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS() then TSComputeRHSJacobian()", ) PetscErrorCode TSComputeDRDUFunction(TS, PetscReal, Vec, Vec *);
2879371c9d4SSatish 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 *);
288b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUU(TS, PetscReal, Vec, Vec *, Vec, Vec *);
289b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUP(TS, PetscReal, Vec, Vec *, Vec, Vec *);
290b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPU(TS, PetscReal, Vec, Vec *, Vec, Vec *);
291b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPP(TS, PetscReal, Vec, Vec *, Vec, Vec *);
2929371c9d4SSatish 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 *);
293b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS, PetscReal, Vec, Vec *, Vec, Vec *);
294b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS, PetscReal, Vec, Vec *, Vec, Vec *);
295b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS, PetscReal, Vec, Vec *, Vec, Vec *);
296b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS, PetscReal, Vec, Vec *, Vec, Vec *);
297b5b4017aSHong Zhang PETSC_EXTERN PetscErrorCode TSSetCostHessianProducts(TS, PetscInt, Vec *, Vec *, Vec);
298b5b4017aSHong Zhang PETSC_EXTERN PetscErrorCode TSGetCostHessianProducts(TS, PetscInt *, Vec **, Vec **, Vec *);
299ba0675f6SHong Zhang PETSC_EXTERN PetscErrorCode TSComputeSNESJacobian(TS, Vec, Mat, Mat);
300a05bf03eSHong Zhang 
301bc952696SBarry Smith /*S
302e91eccc2SStefano Zampini    TSTrajectory - Abstract PETSc object that stores the trajectory (solution of ODE/DAE at each time step)
303bc952696SBarry Smith 
304bc952696SBarry Smith    Level: advanced
305bc952696SBarry Smith 
3061cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectorySetType()`, `TSTrajectoryDestroy()`, `TSTrajectoryReset()`
307bc952696SBarry Smith S*/
308bc952696SBarry Smith typedef struct _p_TSTrajectory *TSTrajectory;
309bc952696SBarry Smith 
310bc952696SBarry Smith /*J
311195e9b02SBarry Smith    TSTrajectoryType - String with the name of a PETSc `TS` trajectory storage method
312bc952696SBarry Smith 
313bc952696SBarry Smith    Level: intermediate
314bc952696SBarry Smith 
3151cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectoryDestroy()`
316bc952696SBarry Smith J*/
317bc952696SBarry Smith typedef const char *TSTrajectoryType;
318bc952696SBarry Smith #define TSTRAJECTORYBASIC         "basic"
3191c8c567eSBarry Smith #define TSTRAJECTORYSINGLEFILE    "singlefile"
3209a53571cSHong Zhang #define TSTRAJECTORYMEMORY        "memory"
3212b043167SHong Zhang #define TSTRAJECTORYVISUALIZATION "visualization"
322bc952696SBarry Smith 
323bc952696SBarry Smith PETSC_EXTERN PetscFunctionList TSTrajectoryList;
324bc952696SBarry Smith PETSC_EXTERN PetscClassId      TSTRAJECTORY_CLASSID;
325bc952696SBarry Smith PETSC_EXTERN PetscBool         TSTrajectoryRegisterAllCalled;
326bc952696SBarry Smith 
327bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS);
3282d29f1f2SStefano Zampini PETSC_EXTERN PetscErrorCode TSResetTrajectory(TS);
32967a3cfb0SHong Zhang PETSC_EXTERN PetscErrorCode TSRemoveTrajectory(TS);
330bc952696SBarry Smith 
331bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectoryCreate(MPI_Comm, TSTrajectory *);
3329a992471SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryReset(TSTrajectory);
333bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy(TSTrajectory *);
334560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSTrajectoryView(TSTrajectory, PetscViewer);
335fd9d3c67SJed Brown PETSC_EXTERN PetscErrorCode TSTrajectorySetType(TSTrajectory, TS, TSTrajectoryType);
336881c1a9bSHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryGetType(TSTrajectory, TS, TSTrajectoryType *);
337bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySet(TSTrajectory, TS, PetscInt, PetscReal, Vec);
338c679fc15SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryGet(TSTrajectory, TS, PetscInt, PetscReal *);
339fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetVecs(TSTrajectory, TS, PetscInt, PetscReal *, Vec, Vec);
340fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetUpdatedHistoryVecs(TSTrajectory, TS, PetscReal, Vec *, Vec *);
341fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetNumSteps(TSTrajectory, PetscInt *);
342fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryRestoreUpdatedHistoryVecs(TSTrajectory, Vec *, Vec *);
343972caf09SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory, TS);
344560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSTrajectoryRegister(const char[], PetscErrorCode (*)(TSTrajectory, TS));
345bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectoryRegisterAll(void);
34668bece0bSHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetUp(TSTrajectory, TS);
3479ffb3502SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetUseHistory(TSTrajectory, PetscBool);
3482bee684fSHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetMonitor(TSTrajectory, PetscBool);
34978fbdcc8SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySetVariableNames(TSTrajectory, const char *const *);
35008347785SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySetTransform(TSTrajectory, PetscErrorCode (*)(void *, Vec, Vec *), PetscErrorCode (*)(void *), void *);
351fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectorySetSolutionOnly(TSTrajectory, PetscBool);
352fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetSolutionOnly(TSTrajectory, PetscBool *);
35364fc91eeSBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySetKeepFiles(TSTrajectory, PetscBool);
35464e38db7SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetDirname(TSTrajectory, const char[]);
35564e38db7SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetFiletemplate(TSTrajectory, const char[]);
35678fbdcc8SBarry Smith PETSC_EXTERN PetscErrorCode TSGetTrajectory(TS, TSTrajectory *);
357bc952696SBarry Smith 
3584bf303faSJacob Faibussowitsch typedef enum {
3594bf303faSJacob Faibussowitsch   TJ_REVOLVE,
3604bf303faSJacob Faibussowitsch   TJ_CAMS,
3614bf303faSJacob Faibussowitsch   TJ_PETSC
3624bf303faSJacob Faibussowitsch } TSTrajectoryMemoryType;
3634bf303faSJacob Faibussowitsch PETSC_EXTERN const char *const TSTrajectoryMemoryTypes[];
3644bf303faSJacob Faibussowitsch 
3654bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectoryMemorySetType(TSTrajectory, TSTrajectoryMemoryType);
3664bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxCpsRAM(TSTrajectory, PetscInt);
3674bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxCpsDisk(TSTrajectory, PetscInt);
3684bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxUnitsRAM(TSTrajectory, PetscInt);
3694bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxUnitsDisk(TSTrajectory, PetscInt);
3704bf303faSJacob Faibussowitsch 
371dfb21088SHong Zhang PETSC_EXTERN PetscErrorCode TSSetCostGradients(TS, PetscInt, Vec *, Vec *);
372dfb21088SHong Zhang PETSC_EXTERN PetscErrorCode TSGetCostGradients(TS, PetscInt *, Vec **, Vec **);
373edd03b47SJacob 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 *);
374dfb21088SHong Zhang PETSC_EXTERN PetscErrorCode TSGetCostIntegral(TS, Vec *);
375022c081aSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeCostIntegrand(TS, PetscReal, Vec, Vec);
376cd4cee2dSHong Zhang PETSC_EXTERN PetscErrorCode TSCreateQuadratureTS(TS, PetscBool, TS *);
377cd4cee2dSHong Zhang PETSC_EXTERN PetscErrorCode TSGetQuadratureTS(TS, PetscBool *, TS *);
378d4aa0a58SBarry Smith 
379dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSetFromOptions(TS, PetscOptionItems *);
380a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitor(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *);
38149abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointMonitorSet(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *), void *, PetscCtxDestroyFn *);
382a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitorCancel(TS);
383a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitorSetFromOptions(TS, const char[], const char[], const char[], PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, PetscViewerAndFormat *), PetscErrorCode (*)(TS, PetscViewerAndFormat *));
384a05bf03eSHong Zhang 
385edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSSetRHSJacobianP()", ) PetscErrorCode TSAdjointSetRHSJacobian(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), void *);
386edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSComputeRHSJacobianP()", ) PetscErrorCode TSAdjointComputeRHSJacobian(TS, PetscReal, Vec, Mat);
387edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS()", ) PetscErrorCode TSAdjointComputeDRDPFunction(TS, PetscReal, Vec, Vec *);
388edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS()", ) PetscErrorCode TSAdjointComputeDRDYFunction(TS, PetscReal, Vec, Vec *);
3892c39e106SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS);
39073b18844SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSetSteps(TS, PetscInt);
391d4aa0a58SBarry Smith 
39273b18844SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointStep(TS);
393f6a906c0SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS);
394ece46509SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointReset(TS);
395999a3721SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointCostIntegral(TS);
396ecf68647SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointSetForward(TS, Mat);
397ecf68647SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointResetForward(TS);
398715f1b00SHong Zhang 
39913b1b0ffSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardSetSensitivities(TS, PetscInt, Mat);
40013b1b0ffSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardGetSensitivities(TS, PetscInt *, Mat *);
401edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 12, 0, "TSCreateQuadratureTS()", ) PetscErrorCode TSForwardSetIntegralGradients(TS, PetscInt, Vec *);
402edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 12, 0, "TSForwardGetSensitivities()", ) PetscErrorCode TSForwardGetIntegralGradients(TS, PetscInt *, Vec **);
403715f1b00SHong Zhang PETSC_EXTERN PetscErrorCode TSForwardSetUp(TS);
4047adebddeSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardReset(TS);
405999a3721SHong Zhang PETSC_EXTERN PetscErrorCode TSForwardCostIntegral(TS);
406715f1b00SHong Zhang PETSC_EXTERN PetscErrorCode TSForwardStep(TS);
407b5b4017aSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardSetInitialSensitivities(TS, Mat);
408cc4f23bcSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardGetStages(TS, PetscInt *, Mat *[]);
409c235aa8dSHong Zhang 
410618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetMaxSteps(TS, PetscInt);
411618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetMaxSteps(TS, PetscInt *);
412618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetMaxTime(TS, PetscReal);
413618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetMaxTime(TS, PetscReal *);
41449354f04SShri Abhyankar PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS, TSExactFinalTimeOption);
415f6953c82SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetExactFinalTime(TS, TSExactFinalTimeOption *);
4168343f784SJames Wright PETSC_EXTERN PetscErrorCode TSSetEvaluationTimes(TS, PetscInt, PetscReal *);
4178343f784SJames Wright PETSC_EXTERN PetscErrorCode TSGetEvaluationTimes(TS, PetscInt *, const PetscReal **);
418*c80d56d9SJames Wright PETSC_EXTERN PetscErrorCode TSGetEvaluationSolutions(TS, PetscInt *, const PetscReal **, Vec **);
419*c80d56d9SJames Wright PETSC_EXTERN PetscErrorCode TSSetTimeSpan(TS, PetscInt, PetscReal *);
420*c80d56d9SJames Wright 
421*c80d56d9SJames Wright /*@C
422*c80d56d9SJames Wright   TSGetTimeSpan - gets the time span set with `TSSetTimeSpan()`
423*c80d56d9SJames Wright 
424*c80d56d9SJames Wright   Not Collective
425*c80d56d9SJames Wright 
426*c80d56d9SJames Wright   Input Parameter:
427*c80d56d9SJames Wright . ts - the time-stepper
428*c80d56d9SJames Wright 
429*c80d56d9SJames Wright   Output Parameters:
430*c80d56d9SJames Wright + n          - number of the time points (>=2)
431*c80d56d9SJames Wright - span_times - array of the time points. The first element and the last element are the initial time and the final time respectively.
432*c80d56d9SJames Wright 
433*c80d56d9SJames Wright   Level: deprecated
434*c80d56d9SJames Wright 
435*c80d56d9SJames Wright   Note:
436*c80d56d9SJames Wright   Deprecated, use `TSGetEvaluationTimes()`.
437*c80d56d9SJames Wright 
438*c80d56d9SJames Wright   The values obtained are valid until the `TS` object is destroyed.
439*c80d56d9SJames Wright 
440*c80d56d9SJames Wright   Both `n` and `span_times` can be `NULL`.
441*c80d56d9SJames Wright 
442*c80d56d9SJames Wright .seealso: [](ch_ts), `TS`, `TSGetEvaluationTimes()`, `TSSetTimeSpan()`, `TSSetEvaluationTimes()`, `TSGetEvaluationSolutions()`
443*c80d56d9SJames Wright  @*/
444*c80d56d9SJames Wright PETSC_DEPRECATED_FUNCTION(3, 23, 0, "TSGetEvaluationTimes()", ) static inline PetscErrorCode TSGetTimeSpan(TS ts, PetscInt *n, const PetscReal *span_times[])
445*c80d56d9SJames Wright {
446*c80d56d9SJames Wright   return TSGetEvaluationTimes(ts, n, span_times);
447*c80d56d9SJames Wright }
448*c80d56d9SJames Wright 
449*c80d56d9SJames Wright /*@C
450*c80d56d9SJames Wright   TSGetTimeSpanSolutions - Get the number of solutions and the solutions at the time points specified by the time span.
451*c80d56d9SJames Wright 
452*c80d56d9SJames Wright   Input Parameter:
453*c80d56d9SJames Wright . ts - the `TS` context obtained from `TSCreate()`
454*c80d56d9SJames Wright 
455*c80d56d9SJames Wright   Output Parameters:
456*c80d56d9SJames Wright + nsol - the number of solutions
457*c80d56d9SJames Wright - Sols - the solution vectors
458*c80d56d9SJames Wright 
459*c80d56d9SJames Wright   Level: deprecated
460*c80d56d9SJames Wright 
461*c80d56d9SJames Wright   Notes:
462*c80d56d9SJames Wright   Deprecated, use `TSGetEvaluationSolutions()`.
463*c80d56d9SJames Wright 
464*c80d56d9SJames Wright   Both `nsol` and `Sols` can be `NULL`.
465*c80d56d9SJames Wright 
466*c80d56d9SJames Wright   Some time points in the time span may be skipped by `TS` so that `nsol` is less than the number of points specified by `TSSetTimeSpan()`.
467*c80d56d9SJames Wright   For example, manipulating the step size, especially with a reduced precision, may cause `TS` to step over certain points in the span.
468*c80d56d9SJames Wright   This issue is alleviated in `TSGetEvaluationSolutions()` by returning the solution times that `Sols` were recorded at.
469*c80d56d9SJames Wright 
470*c80d56d9SJames Wright .seealso: [](ch_ts), `TS`, `TSGetEvaluationSolutions()`, `TSSetTimeSpan()`, `TSGetEvaluationTimes()`, `TSSetEvaluationTimes()`
471*c80d56d9SJames Wright  @*/
472*c80d56d9SJames Wright PETSC_DEPRECATED_FUNCTION(3, 23, 0, "TSGetEvaluationSolutions()", ) static inline PetscErrorCode TSGetTimeSpanSolutions(TS ts, PetscInt *nsol, Vec **Sols)
473*c80d56d9SJames Wright {
474*c80d56d9SJames Wright   return TSGetEvaluationSolutions(ts, nsol, NULL, Sols);
475*c80d56d9SJames Wright }
476818ad0c1SBarry Smith 
477edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSSetTime()", ) PetscErrorCode TSSetInitialTimeStep(TS, PetscReal, PetscReal);
478edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSSetMax()", ) PetscErrorCode TSSetDuration(TS, PetscInt, PetscReal);
479edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetMax()", ) PetscErrorCode TSGetDuration(TS, PetscInt *, PetscReal *);
480edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetStepNumber()", ) PetscErrorCode TSGetTimeStepNumber(TS, PetscInt *);
481edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetStepNumber()", ) PetscErrorCode TSGetTotalSteps(TS, PetscInt *);
48219eac22cSLisandro Dalcin 
483721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
484cc9c3a59SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorExtreme(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
48583a4ac43SBarry Smith 
48683a4ac43SBarry Smith typedef struct _n_TSMonitorDrawCtx *TSMonitorDrawCtx;
48783a4ac43SBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorDrawCtx *);
48883a4ac43SBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *);
48983a4ac43SBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawSolution(TS, PetscInt, PetscReal, Vec, void *);
4902d5ee99bSBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawSolutionPhase(TS, PetscInt, PetscReal, Vec, void *);
49183a4ac43SBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawError(TS, PetscInt, PetscReal, Vec, void *);
4920ed3bfb6SBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawSolutionFunction(TS, PetscInt, PetscReal, Vec, void *);
49383a4ac43SBarry Smith 
494721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointMonitorDefault(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, PetscViewerAndFormat *);
4959110b2e7SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitorDrawSensi(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *);
4969110b2e7SHong Zhang 
497721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorSolution(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
4987f27e910SStefano Zampini 
4997f27e910SStefano Zampini typedef struct _n_TSMonitorVTKCtx *TSMonitorVTKCtx;
5007f27e910SStefano Zampini PETSC_EXTERN PetscErrorCode        TSMonitorSolutionVTK(TS, PetscInt, PetscReal, Vec, TSMonitorVTKCtx);
5017f27e910SStefano Zampini PETSC_EXTERN PetscErrorCode        TSMonitorSolutionVTKDestroy(TSMonitorVTKCtx *);
5027f27e910SStefano Zampini PETSC_EXTERN PetscErrorCode        TSMonitorSolutionVTKCtxCreate(const char *, TSMonitorVTKCtx *);
503fb1732b5SBarry Smith 
504014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSStep(TS);
5057cbde773SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSEvaluateWLTE(TS, NormType, PetscInt *, PetscReal *);
506014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS, PetscInt, Vec, PetscBool *);
507cc708dedSBarry Smith PETSC_EXTERN PetscErrorCode TSSolve(TS, Vec);
508e817cc15SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetEquationType(TS, TSEquationType *);
509e817cc15SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSSetEquationType(TS, TSEquationType);
510014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS, TSConvergedReason *);
511d6ad946cSShri Abhyankar PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS, TSConvergedReason);
512487e0bb9SJed Brown PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS, PetscReal *);
5135ef26d82SJed Brown PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS, PetscInt *);
5145ef26d82SJed Brown PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS, PetscInt *);
515cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS, PetscInt *);
516cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS, PetscInt);
517cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS, PetscInt *);
518cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS, PetscInt);
519cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS, PetscBool);
520dcb233daSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRestartStep(TS);
52124655328SShri PETSC_EXTERN PetscErrorCode TSRollBack(TS);
5229dc50cb5SStefano Zampini PETSC_EXTERN PetscErrorCode TSGetStepRollBack(TS, PetscBool *);
523d2daff3dSHong Zhang 
524cc4f23bcSHong Zhang PETSC_EXTERN PetscErrorCode TSGetStages(TS, PetscInt *, Vec *[]);
525818ad0c1SBarry Smith 
526014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetTime(TS, PetscReal *);
527014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetTime(TS, PetscReal);
528e5e524a1SHong Zhang PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS, PetscReal *);
52980275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS, PetscReal *);
53080275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS, PetscReal);
53180275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetStepNumber(TS, PetscInt *);
53280275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetStepNumber(TS, PetscInt);
533818ad0c1SBarry Smith 
53414d0ab18SJacob Faibussowitsch /*S
5358434afd1SBarry Smith   TSRHSFunctionFn - A prototype of a `TS` right-hand-side evaluation function that would be passed to `TSSetRHSFunction()`
53614d0ab18SJacob Faibussowitsch 
53714d0ab18SJacob Faibussowitsch   Calling Sequence:
53814d0ab18SJacob Faibussowitsch + ts  - timestep context
53914d0ab18SJacob Faibussowitsch . t   - current time
54014d0ab18SJacob Faibussowitsch . u   - input vector
54114d0ab18SJacob Faibussowitsch . F   - function vector
54214d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context
54314d0ab18SJacob Faibussowitsch 
54414d0ab18SJacob Faibussowitsch   Level: beginner
54514d0ab18SJacob Faibussowitsch 
546d1c5d1fcSBarry Smith   Note:
5478434afd1SBarry Smith   The deprecated `TSRHSFunction` still works as a replacement for `TSRHSFunctionFn` *.
548d1c5d1fcSBarry Smith 
5498434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetRHSFunction()`, `DMTSSetRHSFunction()`, `TSIFunctionFn`,
5508434afd1SBarry Smith `TSIJacobianFn`, `TSRHSJacobianFn`
55114d0ab18SJacob Faibussowitsch S*/
5528434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSRHSFunctionFn)(TS ts, PetscReal t, Vec u, Vec F, void *ctx);
553d1c5d1fcSBarry Smith 
5548434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSRHSFunctionFn *TSRHSFunction;
55514d0ab18SJacob Faibussowitsch 
55614d0ab18SJacob Faibussowitsch /*S
5578434afd1SBarry Smith   TSRHSJacobianFn - A prototype of a `TS` right-hand-side Jacobian evaluation function that would be passed to `TSSetRHSJacobian()`
55814d0ab18SJacob Faibussowitsch 
55914d0ab18SJacob Faibussowitsch   Calling Sequence:
56014d0ab18SJacob Faibussowitsch + ts   - the `TS` context obtained from `TSCreate()`
56114d0ab18SJacob Faibussowitsch . t    - current time
56214d0ab18SJacob Faibussowitsch . u    - input vector
56314d0ab18SJacob Faibussowitsch . Amat - (approximate) Jacobian matrix
564af27ebaaSBarry Smith . Pmat - matrix from which preconditioner is to be constructed (usually the same as `Amat`)
56514d0ab18SJacob Faibussowitsch - ctx  - [optional] user-defined context for matrix evaluation routine
56614d0ab18SJacob Faibussowitsch 
56714d0ab18SJacob Faibussowitsch   Level: beginner
56814d0ab18SJacob Faibussowitsch 
569d1c5d1fcSBarry Smith   Note:
5708434afd1SBarry Smith   The deprecated `TSRHSJacobian` still works as a replacement for `TSRHSJacobianFn` *.
571d1c5d1fcSBarry Smith 
5728434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetRHSJacobian()`, `DMTSSetRHSJacobian()`, `TSRHSFunctionFn`,
5738434afd1SBarry Smith `TSIFunctionFn`, `TSIJacobianFn`
57414d0ab18SJacob Faibussowitsch S*/
5758434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSRHSJacobianFn)(TS ts, PetscReal t, Vec u, Mat Amat, Mat Pmat, void *ctx);
576d1c5d1fcSBarry Smith 
5778434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSRHSJacobianFn *TSRHSJacobian;
57814d0ab18SJacob Faibussowitsch 
57914d0ab18SJacob Faibussowitsch /*S
5808434afd1SBarry Smith   TSRHSJacobianPFn - A prototype of a function that computes the Jacobian of G w.r.t. the parameters P where
581af27ebaaSBarry Smith   U_t = G(U,P,t), as well as the location to store the matrix that would be passed to `TSSetRHSJacobianP()`
58214d0ab18SJacob Faibussowitsch 
58314d0ab18SJacob Faibussowitsch   Calling Sequence:
58414d0ab18SJacob Faibussowitsch + ts  - the `TS` context
58514d0ab18SJacob Faibussowitsch . t   - current timestep
58614d0ab18SJacob Faibussowitsch . U   - input vector (current ODE solution)
58714d0ab18SJacob Faibussowitsch . A   - output matrix
58814d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context
58914d0ab18SJacob Faibussowitsch 
59014d0ab18SJacob Faibussowitsch   Level: beginner
59114d0ab18SJacob Faibussowitsch 
592d1c5d1fcSBarry Smith   Note:
5938434afd1SBarry Smith   The deprecated `TSRHSJacobianP` still works as a replacement for `TSRHSJacobianPFn` *.
594d1c5d1fcSBarry Smith 
59514d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetRHSJacobianP()`, `TSGetRHSJacobianP()`
59614d0ab18SJacob Faibussowitsch S*/
5978434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSRHSJacobianPFn)(TS ts, PetscReal t, Vec U, Mat A, void *ctx);
59814d0ab18SJacob Faibussowitsch 
5998434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSRHSJacobianPFn *TSRHSJacobianP;
600d1c5d1fcSBarry Smith 
6018434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS, Vec, TSRHSFunctionFn *, void *);
6028434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS, Vec *, TSRHSFunctionFn **, void **);
6038434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS, Mat, Mat, TSRHSJacobianFn *, void *);
6048434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS, Mat *, Mat *, TSRHSJacobianFn **, void **);
605e1244c69SJed Brown PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS, PetscBool);
606818ad0c1SBarry Smith 
60714d0ab18SJacob Faibussowitsch /*S
6088434afd1SBarry Smith   TSSolutionFn - A prototype of a `TS` solution evaluation function that would be passed to `TSSetSolutionFunction()`
60914d0ab18SJacob Faibussowitsch 
61014d0ab18SJacob Faibussowitsch   Calling Sequence:
61114d0ab18SJacob Faibussowitsch + ts  - timestep context
61214d0ab18SJacob Faibussowitsch . t   - current time
61314d0ab18SJacob Faibussowitsch . u   - output vector
61414d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context
61514d0ab18SJacob Faibussowitsch 
61614d0ab18SJacob Faibussowitsch   Level: advanced
61714d0ab18SJacob Faibussowitsch 
618d1c5d1fcSBarry Smith   Note:
6198434afd1SBarry Smith   The deprecated `TSSolutionFunction` still works as a replacement for `TSSolutionFn` *.
620d1c5d1fcSBarry Smith 
62114d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetSolutionFunction()`, `DMTSSetSolutionFunction()`
62214d0ab18SJacob Faibussowitsch S*/
6238434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSSolutionFn)(TS ts, PetscReal t, Vec u, void *ctx);
62414d0ab18SJacob Faibussowitsch 
6258434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSSolutionFn *TSSolutionFunction;
626d1c5d1fcSBarry Smith 
6278434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS, TSSolutionFn *, void *);
62814d0ab18SJacob Faibussowitsch 
62914d0ab18SJacob Faibussowitsch /*S
6308434afd1SBarry Smith   TSForcingFn - A prototype of a `TS` forcing function evaluation function that would be passed to `TSSetForcingFunction()`
63114d0ab18SJacob Faibussowitsch 
63214d0ab18SJacob Faibussowitsch   Calling Sequence:
63314d0ab18SJacob Faibussowitsch + ts  - timestep context
63414d0ab18SJacob Faibussowitsch . t   - current time
63514d0ab18SJacob Faibussowitsch . f   - output vector
63614d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context
63714d0ab18SJacob Faibussowitsch 
63814d0ab18SJacob Faibussowitsch   Level: advanced
63914d0ab18SJacob Faibussowitsch 
640d1c5d1fcSBarry Smith   Note:
6418434afd1SBarry Smith   The deprecated `TSForcingFunction` still works as a replacement for `TSForcingFn` *.
642d1c5d1fcSBarry Smith 
64314d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetForcingFunction()`, `DMTSSetForcingFunction()`
64414d0ab18SJacob Faibussowitsch S*/
6458434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSForcingFn)(TS ts, PetscReal t, Vec f, void *ctx);
64614d0ab18SJacob Faibussowitsch 
6478434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSForcingFn *TSForcingFunction;
648d1c5d1fcSBarry Smith 
6498434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS, TSForcingFn *, void *);
650ef20d060SBarry Smith 
65114d0ab18SJacob Faibussowitsch /*S
6528434afd1SBarry Smith   TSIFunctionFn - A prototype of a `TS` implicit function evaluation function that would be passed to `TSSetIFunction()
65314d0ab18SJacob Faibussowitsch 
65414d0ab18SJacob Faibussowitsch   Calling Sequence:
65514d0ab18SJacob Faibussowitsch + ts  - the `TS` context obtained from `TSCreate()`
65614d0ab18SJacob Faibussowitsch . t   - time at step/stage being solved
65714d0ab18SJacob Faibussowitsch . U   - state vector
65814d0ab18SJacob Faibussowitsch . U_t - time derivative of state vector
65914d0ab18SJacob Faibussowitsch . F   - function vector
66014d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined context for function
66114d0ab18SJacob Faibussowitsch 
66214d0ab18SJacob Faibussowitsch   Level: beginner
66314d0ab18SJacob Faibussowitsch 
664d1c5d1fcSBarry Smith   Note:
6658434afd1SBarry Smith   The deprecated `TSIFunction` still works as a replacement for `TSIFunctionFn` *.
666d1c5d1fcSBarry Smith 
6678434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetIFunction()`, `DMTSSetIFunction()`, `TSIJacobianFn`, `TSRHSFunctionFn`, `TSRHSJacobianFn`
66814d0ab18SJacob Faibussowitsch S*/
6698434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSIFunctionFn)(TS ts, PetscReal t, Vec U, Vec U_t, Vec F, void *ctx);
670d1c5d1fcSBarry Smith 
6718434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSIFunctionFn *TSIFunction;
67214d0ab18SJacob Faibussowitsch 
67314d0ab18SJacob Faibussowitsch /*S
6748434afd1SBarry Smith   TSIJacobianFn - A prototype of a `TS` Jacobian evaluation function that would be passed to `TSSetIJacobian()`
67514d0ab18SJacob Faibussowitsch 
67614d0ab18SJacob Faibussowitsch   Calling Sequence:
67714d0ab18SJacob Faibussowitsch + ts   - the `TS` context obtained from `TSCreate()`
67814d0ab18SJacob Faibussowitsch . t    - time at step/stage being solved
67914d0ab18SJacob Faibussowitsch . U    - state vector
68014d0ab18SJacob Faibussowitsch . U_t  - time derivative of state vector
68114d0ab18SJacob Faibussowitsch . a    - shift
68214d0ab18SJacob Faibussowitsch . Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
68314d0ab18SJacob Faibussowitsch . Pmat - matrix used for constructing preconditioner, usually the same as `Amat`
68414d0ab18SJacob Faibussowitsch - ctx  - [optional] user-defined context for Jacobian evaluation routine
68514d0ab18SJacob Faibussowitsch 
68614d0ab18SJacob Faibussowitsch   Level: beginner
68714d0ab18SJacob Faibussowitsch 
688d1c5d1fcSBarry Smith   Note:
6898434afd1SBarry Smith   The deprecated `TSIJacobian` still works as a replacement for `TSIJacobianFn` *.
69014d0ab18SJacob Faibussowitsch 
6918434afd1SBarry Smith .seealso: [](ch_ts), `TSSetIJacobian()`, `DMTSSetIJacobian()`, `TSIFunctionFn`, `TSRHSFunctionFn`, `TSRHSJacobianFn`
692d1c5d1fcSBarry Smith S*/
6938434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSIJacobianFn)(TS ts, PetscReal t, Vec U, Vec U_t, PetscReal a, Mat Amat, Mat Pmat, void *ctx);
694d1c5d1fcSBarry Smith 
6958434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSIJacobianFn *TSIJacobian;
696d1c5d1fcSBarry Smith 
6978434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetIFunction(TS, Vec, TSIFunctionFn *, void *);
6988434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGetIFunction(TS, Vec *, TSIFunctionFn **, void **);
6998434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS, Mat, Mat, TSIJacobianFn *, void *);
7008434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS, Mat *, Mat *, TSIJacobianFn **, void **);
701316643e7SJed Brown 
70214d0ab18SJacob Faibussowitsch /*S
7038434afd1SBarry Smith   TSI2FunctionFn - A prototype of a `TS` implicit function evaluation function for 2nd order systems that would be passed to `TSSetI2Function()`
70414d0ab18SJacob Faibussowitsch 
70514d0ab18SJacob Faibussowitsch   Calling Sequence:
70614d0ab18SJacob Faibussowitsch + ts   - the `TS` context obtained from `TSCreate()`
70714d0ab18SJacob Faibussowitsch . t    - time at step/stage being solved
70814d0ab18SJacob Faibussowitsch . U    - state vector
70914d0ab18SJacob Faibussowitsch . U_t  - time derivative of state vector
71014d0ab18SJacob Faibussowitsch . U_tt - second time derivative of state vector
71114d0ab18SJacob Faibussowitsch . F    - function vector
71214d0ab18SJacob Faibussowitsch - ctx  - [optional] user-defined context for matrix evaluation routine (may be `NULL`)
71314d0ab18SJacob Faibussowitsch 
71414d0ab18SJacob Faibussowitsch   Level: advanced
71514d0ab18SJacob Faibussowitsch 
716d1c5d1fcSBarry Smith   Note:
7178434afd1SBarry Smith   The deprecated `TSI2Function` still works as a replacement for `TSI2FunctionFn` *.
718d1c5d1fcSBarry Smith 
7198434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetI2Function()`, `DMTSSetI2Function()`, `TSIFunctionFn`
72014d0ab18SJacob Faibussowitsch S*/
7218434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSI2FunctionFn)(TS ts, PetscReal t, Vec U, Vec U_t, Vec U_tt, Vec F, void *ctx);
722d1c5d1fcSBarry Smith 
7238434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSI2FunctionFn *TSI2Function;
72414d0ab18SJacob Faibussowitsch 
72514d0ab18SJacob Faibussowitsch /*S
7268434afd1SBarry Smith   TSI2JacobianFn - A prototype of a `TS` implicit Jacobian evaluation function for 2nd order systems that would be passed to `TSSetI2Jacobian()`
72714d0ab18SJacob Faibussowitsch 
72814d0ab18SJacob Faibussowitsch   Calling Sequence:
72914d0ab18SJacob Faibussowitsch + ts   - the `TS` context obtained from `TSCreate()`
73014d0ab18SJacob Faibussowitsch . t    - time at step/stage being solved
73114d0ab18SJacob Faibussowitsch . U    - state vector
73214d0ab18SJacob Faibussowitsch . U_t  - time derivative of state vector
73314d0ab18SJacob Faibussowitsch . U_tt - second time derivative of state vector
73414d0ab18SJacob Faibussowitsch . v    - shift for U_t
73514d0ab18SJacob Faibussowitsch . a    - shift for U_tt
73614d0ab18SJacob 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
73714d0ab18SJacob Faibussowitsch . jac  - matrix from which to construct the preconditioner, may be same as `J`
73814d0ab18SJacob Faibussowitsch - ctx  - [optional] user-defined context for matrix evaluation routine
73914d0ab18SJacob Faibussowitsch 
74014d0ab18SJacob Faibussowitsch   Level: advanced
74114d0ab18SJacob Faibussowitsch 
742d1c5d1fcSBarry Smith   Note:
7438434afd1SBarry Smith   The deprecated `TSI2Jacobian` still works as a replacement for `TSI2JacobianFn` *.
74414d0ab18SJacob Faibussowitsch 
7458434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetI2Jacobian()`, `DMTSSetI2Jacobian()`, `TSIFunctionFn`, `TSIJacobianFn`, `TSRHSFunctionFn`, `TSRHSJacobianFn`
746d1c5d1fcSBarry Smith S*/
7478434afd1SBarry 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);
748d1c5d1fcSBarry Smith 
7498434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSI2JacobianFn *TSI2Jacobian;
750d1c5d1fcSBarry Smith 
7518434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetI2Function(TS, Vec, TSI2FunctionFn *, void *);
7528434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGetI2Function(TS, Vec *, TSI2FunctionFn **, void **);
7538434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetI2Jacobian(TS, Mat, Mat, TSI2JacobianFn *, void *);
7548434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGetI2Jacobian(TS, Mat *, Mat *, TSI2JacobianFn **, void **);
755efe9872eSLisandro Dalcin 
756545aaa6fSHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitSetIS(TS, const char[], IS);
757545aaa6fSHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetIS(TS, const char[], IS *);
7588434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSRHSSplitSetRHSFunction(TS, const char[], Vec, TSRHSFunctionFn *, void *);
7593a2a065fSHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitSetIFunction(TS, const char[], Vec, TSIFunctionFn *, void *);
7603a2a065fSHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitSetIJacobian(TS, const char[], Mat, Mat, TSIJacobianFn *, void *);
7611d06f6b3SHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTS(TS, const char[], TS *);
7621d06f6b3SHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTSs(TS, PetscInt *, TS *[]);
7630fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSSetUseSplitRHSFunction(TS, PetscBool);
7640fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSGetUseSplitRHSFunction(TS, PetscBool *);
7654bd3aaa3SHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetSNES(TS, SNES *);
7664bd3aaa3SHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitSetSNES(TS, SNES);
767d9194312SHong Zhang 
7688434afd1SBarry Smith PETSC_EXTERN TSRHSFunctionFn TSComputeRHSFunctionLinear;
7698434afd1SBarry Smith PETSC_EXTERN TSRHSJacobianFn TSComputeRHSJacobianConstant;
770014dd563SJed Brown PETSC_EXTERN PetscErrorCode  TSComputeIFunctionLinear(TS, PetscReal, Vec, Vec, Vec, void *);
771d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode  TSComputeIJacobianConstant(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
7721c8a10a0SJed Brown PETSC_EXTERN PetscErrorCode  TSComputeSolutionFunction(TS, PetscReal, Vec);
7739b7cd975SBarry Smith PETSC_EXTERN PetscErrorCode  TSComputeForcingFunction(TS, PetscReal, Vec);
774847ff0e1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode  TSComputeIJacobianDefaultColor(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
775d7cfae9bSHong Zhang PETSC_EXTERN PetscErrorCode  TSPruneIJacobianColor(TS, Mat, Mat);
776e34be4c2SBarry Smith 
777014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
778b8123daeSJed Brown PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS, PetscReal));
7799be3e283SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS, PetscReal, PetscInt, Vec *));
780c688d042SShri Abhyankar PETSC_EXTERN PetscErrorCode TSSetPostEvaluate(TS, PetscErrorCode (*)(TS));
781014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
782ecc87898SStefano Zampini PETSC_EXTERN PetscErrorCode TSSetResize(TS, PetscBool, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscBool *, void *), PetscErrorCode (*)(TS, PetscInt, Vec[], Vec[], void *), void *);
783014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPreStep(TS);
784b8123daeSJed Brown PETSC_EXTERN PetscErrorCode TSPreStage(TS, PetscReal);
7859be3e283SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSPostStage(TS, PetscReal, PetscInt, Vec *);
786c688d042SShri Abhyankar PETSC_EXTERN PetscErrorCode TSPostEvaluate(TS);
787014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPostStep(TS);
7886bd3a4fdSStefano Zampini PETSC_EXTERN PetscErrorCode TSResize(TS);
7896bd3a4fdSStefano Zampini PETSC_EXTERN PetscErrorCode TSResizeRetrieveVec(TS, const char *, Vec *);
7906bd3a4fdSStefano Zampini PETSC_EXTERN PetscErrorCode TSResizeRegisterVec(TS, const char *, Vec);
791c6bf8827SStefano Zampini PETSC_EXTERN PetscErrorCode TSGetStepResize(TS, PetscBool *);
7926bd3a4fdSStefano Zampini 
793014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSInterpolate(TS, PetscReal, Vec);
794014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetTolerances(TS, PetscReal, Vec, PetscReal, Vec);
795c5033834SJed Brown PETSC_EXTERN PetscErrorCode TSGetTolerances(TS, PetscReal *, Vec *, PetscReal *, Vec *);
7967453f775SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *);
7978a175baeSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm(TS, Vec, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *);
798014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS, PetscReal);
799014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS, PetscReal *);
800d183316bSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSSetFunctionDomainError(TS, PetscErrorCode (*)(TS, PetscReal, Vec, PetscBool *));
801d183316bSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSFunctionDomainError(TS, PetscReal, Vec, PetscBool *);
802000e7ae3SMatthew Knepley 
803014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS, PetscErrorCode (*)(TS, PetscReal *, void *), void *);
8048d359177SBarry Smith PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS, PetscReal *, void *);
805014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS, PetscReal *);
806014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS, PetscReal);
807014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS, PetscErrorCode (*)(TS, Vec, void *, PetscReal *, PetscBool *), void *);
8088d359177SBarry Smith PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS, Vec, void *, PetscReal *, PetscBool *);
809014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS, Vec, PetscReal *, PetscBool *);
810014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS, PetscReal);
811014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);
81221c89e3eSBarry Smith 
813014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPythonSetType(TS, const char[]);
814ebead697SStefano Zampini PETSC_EXTERN PetscErrorCode TSPythonGetType(TS, const char *[]);
8151d6018f0SLisandro Dalcin 
816014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS, PetscReal, Vec, Vec);
817d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS, PetscReal, Vec, Mat, Mat);
818014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS, PetscReal, Vec, Vec, Vec, PetscBool);
819d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscBool);
820efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSComputeI2Function(TS, PetscReal, Vec, Vec, Vec, Vec);
821efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSComputeI2Jacobian(TS, PetscReal, Vec, Vec, Vec, PetscReal, PetscReal, Mat, Mat);
822f9c1d6abSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS, PetscReal, PetscReal, PetscReal *, PetscReal *);
823818ad0c1SBarry Smith 
824014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS, Vec, Vec);
825d6ebe24aSShri Abhyankar 
826967bb25aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetBoundaryLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *);
8278434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM, TSRHSFunctionFn *, void *);
8288434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM, TSRHSFunctionFn **, void **);
82949abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionContextDestroy(DM, PetscCtxDestroyFn *);
8308434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM, TSRHSJacobianFn *, void *);
8318434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM, TSRHSJacobianFn **, void **);
83249abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobianContextDestroy(DM, PetscCtxDestroyFn *);
8338434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM, TSIFunctionFn *, void *);
8348434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM, TSIFunctionFn **, void **);
83549abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetIFunctionContextDestroy(DM, PetscCtxDestroyFn *);
8368434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM, TSIJacobianFn *, void *);
8378434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM, TSIJacobianFn **, void **);
83849abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetIJacobianContextDestroy(DM, PetscCtxDestroyFn *);
8398434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetI2Function(DM, TSI2FunctionFn *, void *);
8408434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetI2Function(DM, TSI2FunctionFn **, void **);
84149abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetI2FunctionContextDestroy(DM, PetscCtxDestroyFn *);
8428434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetI2Jacobian(DM, TSI2JacobianFn *, void *);
8438434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetI2Jacobian(DM, TSI2JacobianFn **, void **);
84449abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetI2JacobianContextDestroy(DM, PetscCtxDestroyFn *);
845efe9872eSLisandro Dalcin 
84614d0ab18SJacob Faibussowitsch /*S
8478434afd1SBarry Smith   TSTransientVariableFn - A prototype of a function to transform from state to transient variables that would be passed to `TSSetTransientVariable()`
84814d0ab18SJacob Faibussowitsch 
84914d0ab18SJacob Faibussowitsch   Calling Sequence:
85014d0ab18SJacob Faibussowitsch + ts  - timestep context
85114d0ab18SJacob Faibussowitsch . p   - input vector (primitive form)
85214d0ab18SJacob Faibussowitsch . c   - output vector, transient variables (conservative form)
85314d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context
85414d0ab18SJacob Faibussowitsch 
85514d0ab18SJacob Faibussowitsch   Level: advanced
85614d0ab18SJacob Faibussowitsch 
857d1c5d1fcSBarry Smith   Note:
8588434afd1SBarry Smith   The deprecated `TSTransientVariable` still works as a replacement for `TSTransientVariableFn` *.
859d1c5d1fcSBarry Smith 
86014d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetTransientVariable()`, `DMTSSetTransientVariable()`
86114d0ab18SJacob Faibussowitsch S*/
8628434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSTransientVariableFn)(TS ts, Vec p, Vec c, void *ctx);
86314d0ab18SJacob Faibussowitsch 
8648434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSTransientVariableFn *TSTransientVariable;
865d1c5d1fcSBarry Smith 
8668434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetTransientVariable(TS, TSTransientVariableFn *, void *);
8678434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetTransientVariable(DM, TSTransientVariableFn *, void *);
8688434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetTransientVariable(DM, TSTransientVariableFn **, void *);
869e3c11fc1SJed Brown PETSC_EXTERN PetscErrorCode TSComputeTransientVariable(TS, Vec, Vec);
870e3c11fc1SJed Brown PETSC_EXTERN PetscErrorCode TSHasTransientVariable(TS, PetscBool *);
871e3c11fc1SJed Brown 
8728434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM, TSSolutionFn *, void *);
8738434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM, TSSolutionFn **, void **);
8748434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM, TSForcingFn *, void *);
8758434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM, TSForcingFn **, void **);
876e17bd7bbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCheckResidual(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscReal *);
877e17bd7bbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCheckJacobian(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscBool *, PetscReal *);
8787f96f943SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec);
8799b7cd975SBarry Smith 
8800c9409e4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSGetIFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, Vec, void *), void **);
881d433e6cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, Vec, void *), void *);
8820c9409e4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSGetIJacobianLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void **);
883d433e6cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void *);
8840c9409e4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSGetRHSFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, void *), void **);
885d433e6cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *);
886b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrix(DM);
887b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrixLumped(DM);
888b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSDestroyRHSMassMatrix(DM);
8896c6b9e74SPeter Brune 
8906c6b9e74SPeter Brune PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM, PetscErrorCode (*)(void *, PetscViewer), PetscErrorCode (*)(void **, PetscViewer));
8916c6b9e74SPeter Brune PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM, PetscErrorCode (*)(void *, PetscViewer), PetscErrorCode (*)(void **, PetscViewer));
8926c6b9e74SPeter Brune 
89314d0ab18SJacob Faibussowitsch /*S
894dd8e379bSPierre Jolivet   DMDATSRHSFunctionLocalFn - A prototype of a local `TS` right-hand side residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetRHSFunctionLocal()`
8956c6b9e74SPeter Brune 
89614d0ab18SJacob Faibussowitsch   Calling Sequence:
89714d0ab18SJacob Faibussowitsch + info - defines the subdomain to evaluate the residual on
89814d0ab18SJacob Faibussowitsch . t    - time at which to evaluate residual
89914d0ab18SJacob Faibussowitsch . x    - array of local state information
90014d0ab18SJacob Faibussowitsch . f    - output array of local residual information
90114d0ab18SJacob Faibussowitsch - ctx  - optional user context
90214d0ab18SJacob Faibussowitsch 
90314d0ab18SJacob Faibussowitsch   Level: beginner
90414d0ab18SJacob Faibussowitsch 
905d1c5d1fcSBarry Smith   Note:
9068434afd1SBarry Smith   The deprecated `DMDATSRHSFunctionLocal` still works as a replacement for `DMDATSRHSFunctionLocalFn` *.
907d1c5d1fcSBarry Smith 
9088434afd1SBarry Smith .seealso: `DMDA`, `DMDATSSetRHSFunctionLocal()`, `TSRHSFunctionFn`, `DMDATSRHSJacobianLocalFn`, `DMDATSIJacobianLocalFn`, `DMDATSIFunctionLocalFn`
90914d0ab18SJacob Faibussowitsch S*/
9108434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDATSRHSFunctionLocalFn)(DMDALocalInfo *info, PetscReal t, void *x, void *f, void *ctx);
911d1c5d1fcSBarry Smith 
9128434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef DMDATSRHSFunctionLocalFn *DMDATSRHSFunctionLocal;
91314d0ab18SJacob Faibussowitsch 
91414d0ab18SJacob Faibussowitsch /*S
9158434afd1SBarry Smith   DMDATSRHSJacobianLocalFn - A prototype of a local residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetRHSJacobianLocal()`
91614d0ab18SJacob Faibussowitsch 
91714d0ab18SJacob Faibussowitsch   Calling Sequence:
91814d0ab18SJacob Faibussowitsch + info - defines the subdomain to evaluate the residual on
91914d0ab18SJacob Faibussowitsch . t    - time at which to evaluate residual
92014d0ab18SJacob Faibussowitsch . x    - array of local state information
92114d0ab18SJacob Faibussowitsch . J    - Jacobian matrix
92214d0ab18SJacob Faibussowitsch . B    - matrix from which to construct the preconditioner; often same as `J`
92314d0ab18SJacob Faibussowitsch - ctx  - optional context
92414d0ab18SJacob Faibussowitsch 
92514d0ab18SJacob Faibussowitsch   Level: beginner
92614d0ab18SJacob Faibussowitsch 
927d1c5d1fcSBarry Smith   Note:
9288434afd1SBarry Smith   The deprecated `DMDATSRHSJacobianLocal` still works as a replacement for `DMDATSRHSJacobianLocalFn` *.
929d1c5d1fcSBarry Smith 
9308434afd1SBarry Smith .seealso: `DMDA`, `DMDATSSetRHSJacobianLocal()`, `TSRHSJacobianFn`, `DMDATSRHSFunctionLocalFn`, `DMDATSIJacobianLocalFn`, `DMDATSIFunctionLocalFn`
93114d0ab18SJacob Faibussowitsch S*/
9328434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDATSRHSJacobianLocalFn)(DMDALocalInfo *info, PetscReal t, void *x, Mat J, Mat B, void *ctx);
933d1c5d1fcSBarry Smith 
9348434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef DMDATSRHSJacobianLocalFn *DMDATSRHSJacobianLocal;
93514d0ab18SJacob Faibussowitsch 
93614d0ab18SJacob Faibussowitsch /*S
9378434afd1SBarry Smith   DMDATSIFunctionLocalFn - A prototype of a local residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetIFunctionLocal()`
93814d0ab18SJacob Faibussowitsch 
93914d0ab18SJacob Faibussowitsch   Calling Sequence:
94014d0ab18SJacob Faibussowitsch + info  - defines the subdomain to evaluate the residual on
94114d0ab18SJacob Faibussowitsch . t     - time at which to evaluate residual
94214d0ab18SJacob Faibussowitsch . x     - array of local state information
94314d0ab18SJacob Faibussowitsch . xdot  - array of local time derivative information
94414d0ab18SJacob Faibussowitsch . imode - output array of local function evaluation information
94514d0ab18SJacob Faibussowitsch - ctx   - optional context
94614d0ab18SJacob Faibussowitsch 
94714d0ab18SJacob Faibussowitsch   Level: beginner
94814d0ab18SJacob Faibussowitsch 
949d1c5d1fcSBarry Smith   Note:
9508434afd1SBarry Smith   The deprecated `DMDATSIFunctionLocal` still works as a replacement for `DMDATSIFunctionLocalFn` *.
951d1c5d1fcSBarry Smith 
9528434afd1SBarry Smith .seealso: `DMDA`, `DMDATSSetIFunctionLocal()`, `DMDATSIJacobianLocalFn`, `TSIFunctionFn`
95314d0ab18SJacob Faibussowitsch S*/
9548434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDATSIFunctionLocalFn)(DMDALocalInfo *info, PetscReal t, void *x, void *xdot, void *imode, void *ctx);
955d1c5d1fcSBarry Smith 
9568434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef DMDATSIFunctionLocalFn *DMDATSIFunctionLocal;
95714d0ab18SJacob Faibussowitsch 
95814d0ab18SJacob Faibussowitsch /*S
9598434afd1SBarry Smith   DMDATSIJacobianLocalFn - A prototype of a local residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetIJacobianLocal()`
96014d0ab18SJacob Faibussowitsch 
96114d0ab18SJacob Faibussowitsch   Calling Sequence:
96214d0ab18SJacob Faibussowitsch + info  - defines the subdomain to evaluate the residual on
96314d0ab18SJacob Faibussowitsch . t     - time at which to evaluate the jacobian
96414d0ab18SJacob Faibussowitsch . x     - array of local state information
96514d0ab18SJacob Faibussowitsch . xdot  - time derivative at this state
96614d0ab18SJacob Faibussowitsch . shift - see `TSSetIJacobian()` for the meaning of this parameter
96714d0ab18SJacob Faibussowitsch . J     - Jacobian matrix
96814d0ab18SJacob Faibussowitsch . B     - matrix from which to construct the preconditioner; often same as `J`
96914d0ab18SJacob Faibussowitsch - ctx   - optional context
97014d0ab18SJacob Faibussowitsch 
97114d0ab18SJacob Faibussowitsch   Level: beginner
97214d0ab18SJacob Faibussowitsch 
973d1c5d1fcSBarry Smith   Note:
9748434afd1SBarry Smith   The deprecated `DMDATSIJacobianLocal` still works as a replacement for `DMDATSIJacobianLocalFn` *.
97514d0ab18SJacob Faibussowitsch 
9768434afd1SBarry Smith .seealso: `DMDA` `DMDATSSetIJacobianLocal()`, `TSIJacobianFn`, `DMDATSIFunctionLocalFn`, `DMDATSRHSFunctionLocalFn`,  `DMDATSRHSJacobianlocal()`
977d1c5d1fcSBarry Smith S*/
9788434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDATSIJacobianLocalFn)(DMDALocalInfo *info, PetscReal t, void *x, void *xdot, PetscReal shift, Mat J, Mat B, void *ctx);
979d1c5d1fcSBarry Smith 
9808434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef DMDATSIJacobianLocalFn *DMDATSIJacobianLocal;
981d1c5d1fcSBarry Smith 
9828434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM, InsertMode, DMDATSRHSFunctionLocalFn *, void *);
9838434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM, DMDATSRHSJacobianLocalFn *, void *);
9848434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM, InsertMode, DMDATSIFunctionLocalFn *, void *);
9858434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM, DMDATSIJacobianLocalFn *, void *);
9866c6b9e74SPeter Brune 
987be1b0d75SMatthew G. Knepley typedef struct _n_TSMonitorLGCtx *TSMonitorLGCtx;
988d1212d36SBarry Smith typedef struct {
989d1212d36SBarry Smith   Vec            ray;
990d1212d36SBarry Smith   VecScatter     scatter;
991d1212d36SBarry Smith   PetscViewer    viewer;
99251b4a12fSMatthew G. Knepley   TSMonitorLGCtx lgctx;
993d1212d36SBarry Smith } TSMonitorDMDARayCtx;
994d1212d36SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void **);
995d1212d36SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS, PetscInt, PetscReal, Vec, void *);
99651b4a12fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS, PetscInt, PetscReal, Vec, void *);
997d1212d36SBarry Smith 
998bdad233fSMatthew Knepley /* Dynamic creation and loading functions */
999140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList TSList;
100019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode    TSGetType(TS, TSType *);
100119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode    TSSetType(TS, TSType);
1002bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode    TSRegister(const char[], PetscErrorCode (*)(TS));
100330de9b25SBarry Smith 
1004014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetSNES(TS, SNES *);
1005deb2cd25SJed Brown PETSC_EXTERN PetscErrorCode TSSetSNES(TS, SNES);
1006014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetKSP(TS, KSP *);
1007818ad0c1SBarry Smith 
1008014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSView(TS, PetscViewer);
100955849f57SBarry Smith PETSC_EXTERN PetscErrorCode TSLoad(TS, PetscViewer);
1010fe2efc57SMark PETSC_EXTERN PetscErrorCode TSViewFromOptions(TS, PetscObject, const char[]);
1011fe2efc57SMark PETSC_EXTERN PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory, PetscObject, const char[]);
101255849f57SBarry Smith 
101355849f57SBarry Smith #define TS_FILE_CLASSID 1211225
101421c89e3eSBarry Smith 
1015014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS, void *);
1016014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS, void *);
101721c89e3eSBarry Smith 
1018a9f9c1f6SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtx *);
1019a9f9c1f6SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *);
10204f09c107SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS, PetscInt, PetscReal, Vec, void *);
10214f09c107SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS, PetscInt, PetscReal, Vec, void *);
102231152f8aSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS, const char *const *);
1023b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS, const char *const **);
1024e673d494SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx, const char *const *);
1025a66092f1SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS, const char *const *);
1026a66092f1SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx, const char *const *);
102749abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS, PetscErrorCode (*)(void *, Vec, Vec *), PetscCtxDestroyFn *, void *);
102849abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx, PetscErrorCode (*)(void *, Vec, Vec *), PetscCtxDestroyFn *, void *);
10294f09c107SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS, PetscInt, PetscReal, Vec, void *);
1030201da799SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS, PetscInt, PetscReal, Vec, void *);
1031201da799SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS, PetscInt, PetscReal, Vec, void *);
1032edbaebb3SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorError(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
1033d0c080abSJoseph Pusztay PETSC_EXTERN PetscErrorCode TSDMSwarmMonitorMoments(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
1034ef20d060SBarry Smith 
1035e669de00SBarry Smith struct _n_TSMonitorLGCtxNetwork {
1036e669de00SBarry Smith   PetscInt     nlg;
1037e669de00SBarry Smith   PetscDrawLG *lg;
1038e669de00SBarry Smith   PetscBool    semilogy;
1039e669de00SBarry Smith   PetscInt     howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
1040e669de00SBarry Smith };
1041e669de00SBarry Smith typedef struct _n_TSMonitorLGCtxNetwork *TSMonitorLGCtxNetwork;
1042e669de00SBarry Smith PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkDestroy(TSMonitorLGCtxNetwork *);
1043e669de00SBarry Smith PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkCreate(TS, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtxNetwork *);
1044e669de00SBarry Smith PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkSolution(TS, PetscInt, PetscReal, Vec, void *);
1045e669de00SBarry Smith 
1046b3d3934dSBarry Smith typedef struct _n_TSMonitorEnvelopeCtx *TSMonitorEnvelopeCtx;
1047b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeCtxCreate(TS, TSMonitorEnvelopeCtx *);
1048b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode             TSMonitorEnvelope(TS, PetscInt, PetscReal, Vec, void *);
1049b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeGetBounds(TS, Vec *, Vec *);
1050b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *);
1051b3d3934dSBarry Smith 
10528189c53fSBarry Smith typedef struct _n_TSMonitorSPEigCtx *TSMonitorSPEigCtx;
10538189c53fSBarry Smith PETSC_EXTERN PetscErrorCode          TSMonitorSPEigCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorSPEigCtx *);
10548189c53fSBarry Smith PETSC_EXTERN PetscErrorCode          TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *);
10558189c53fSBarry Smith PETSC_EXTERN PetscErrorCode          TSMonitorSPEig(TS, PetscInt, PetscReal, Vec, void *);
10563914022bSBarry Smith 
10571b575b74SJoseph Pusztay typedef struct _n_TSMonitorSPCtx *TSMonitorSPCtx;
105860e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorSPCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscBool, PetscBool, TSMonitorSPCtx *);
10591b575b74SJoseph Pusztay PETSC_EXTERN PetscErrorCode       TSMonitorSPCtxDestroy(TSMonitorSPCtx *);
10600ec8ee2bSJoseph Pusztay PETSC_EXTERN PetscErrorCode       TSMonitorSPSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
10611b575b74SJoseph Pusztay 
106260e16b1bSMatthew G. Knepley typedef struct _n_TSMonitorHGCtx *TSMonitorHGCtx;
106360e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorHGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscInt, PetscBool, TSMonitorHGCtx *);
106460e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
106560e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorHGCtxDestroy(TSMonitorHGCtx *);
106660e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
106760e16b1bSMatthew G. Knepley 
1068ca4445c7SIlya Fursov PETSC_EXTERN PetscErrorCode TSSetEventHandler(TS, PetscInt, PetscInt[], PetscBool[], PetscErrorCode (*)(TS, PetscReal, Vec, PetscReal[], void *), PetscErrorCode (*)(TS, PetscInt, PetscInt[], PetscReal, Vec, PetscBool, void *), void *);
1069ca4445c7SIlya Fursov PETSC_EXTERN PetscErrorCode TSSetPostEventStep(TS, PetscReal);
1070fe4ad979SIlya Fursov PETSC_EXTERN PetscErrorCode TSSetPostEventSecondStep(TS, PetscReal);
1071fe4ad979SIlya Fursov PETSC_DEPRECATED_FUNCTION(3, 21, 0, "TSSetPostEventSecondStep()", ) static inline PetscErrorCode TSSetPostEventIntervalStep(TS ts, PetscReal dt)
1072fe4ad979SIlya Fursov {
1073fe4ad979SIlya Fursov   return TSSetPostEventSecondStep(ts, dt);
1074fe4ad979SIlya Fursov }
10756427ac75SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS, PetscReal, PetscReal[]);
10761ea83e56SMiguel PETSC_EXTERN PetscErrorCode TSGetNumEvents(TS, PetscInt *);
10771ea83e56SMiguel 
107876bdecfbSBarry Smith /*J
1079195e9b02SBarry Smith    TSSSPType - string with the name of a `TSSSP` scheme.
1080815f1ad5SJed Brown 
1081815f1ad5SJed Brown    Level: beginner
1082815f1ad5SJed Brown 
10831cc06b55SBarry Smith .seealso: [](ch_ts), `TSSSPSetType()`, `TS`, `TSSSP`
108476bdecfbSBarry Smith J*/
108519fd82e9SBarry Smith typedef const char *TSSSPType;
1086815f1ad5SJed Brown #define TSSSPRKS2  "rks2"
1087815f1ad5SJed Brown #define TSSSPRKS3  "rks3"
1088815f1ad5SJed Brown #define TSSSPRK104 "rk104"
1089815f1ad5SJed Brown 
109019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode    TSSSPSetType(TS, TSSSPType);
109119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode    TSSSPGetType(TS, TSSSPType *);
1092014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSSPSetNumStages(TS, PetscInt);
1093014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSSPGetNumStages(TS, PetscInt *);
1094787849ffSJed Brown PETSC_EXTERN PetscErrorCode    TSSSPInitializePackage(void);
1095560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode    TSSSPFinalizePackage(void);
1096013797aaSBarry Smith PETSC_EXTERN PetscFunctionList TSSSPList;
1097815f1ad5SJed Brown 
1098ed657a08SJed Brown /*S
109984df9cb4SJed Brown    TSAdapt - Abstract object that manages time-step adaptivity
110084df9cb4SJed Brown 
110184df9cb4SJed Brown    Level: beginner
110284df9cb4SJed Brown 
1103af27ebaaSBarry Smith .seealso: [](ch_ts), `TS`, `TSGetAdapt()`, `TSAdaptCreate()`, `TSAdaptType`
110484df9cb4SJed Brown S*/
110584df9cb4SJed Brown typedef struct _p_TSAdapt *TSAdapt;
110684df9cb4SJed Brown 
1107f0fc11ceSJed Brown /*J
110887497f52SBarry Smith    TSAdaptType - String with the name of `TSAdapt` scheme.
110984df9cb4SJed Brown 
111084df9cb4SJed Brown    Level: beginner
111184df9cb4SJed Brown 
1112af27ebaaSBarry Smith .seealso: [](ch_ts), `TSGetAdapt()`, `TSAdaptSetType()`, `TS`, `TSAdapt`
1113f0fc11ceSJed Brown J*/
1114f8963224SJed Brown typedef const char *TSAdaptType;
1115b0f836d7SJed Brown #define TSADAPTNONE    "none"
11161917a363SLisandro Dalcin #define TSADAPTBASIC   "basic"
1117cb7ab186SLisandro Dalcin #define TSADAPTDSP     "dsp"
11188d59e960SJed Brown #define TSADAPTCFL     "cfl"
11191917a363SLisandro Dalcin #define TSADAPTGLEE    "glee"
1120d949e4a4SStefano Zampini #define TSADAPTHISTORY "history"
112184df9cb4SJed Brown 
1122552698daSJed Brown PETSC_EXTERN PetscErrorCode TSGetAdapt(TS, TSAdapt *);
1123bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[], PetscErrorCode (*)(TSAdapt));
1124607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void);
1125014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
1126014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm, TSAdapt *);
112719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt, TSAdaptType);
1128d0288e4fSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetType(TSAdapt, TSAdaptType *);
1129014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt, const char[]);
1130014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
1131014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt, const char[], PetscInt, PetscInt, PetscReal, PetscReal, PetscBool);
1132014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt, PetscInt *, const PetscInt **, const PetscInt **, const PetscReal **, const PetscReal **);
1133014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt, TS, PetscReal, PetscInt *, PetscReal *, PetscBool *);
1134b295832fSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt, TS, PetscReal, Vec, PetscBool *);
1135014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt, PetscViewer);
1136ad6bc421SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt, PetscViewer);
1137dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(TSAdapt, PetscOptionItems *);
1138099af0a3SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt);
1139014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt *);
1140014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt, PetscBool);
1141bf997491SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt, PetscBool);
11421917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptSetSafety(TSAdapt, PetscReal, PetscReal);
11431917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetSafety(TSAdapt, PetscReal *, PetscReal *);
114476cddca1SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt, PetscReal);
114576cddca1SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt, PetscReal *);
11461917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptSetClip(TSAdapt, PetscReal, PetscReal);
11471917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetClip(TSAdapt, PetscReal *, PetscReal *);
114862c23b28SMark PETSC_EXTERN PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt, PetscReal);
114962c23b28SMark PETSC_EXTERN PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt, PetscReal *);
1150014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt, PetscReal, PetscReal);
11511917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetStepLimits(TSAdapt, PetscReal *, PetscReal *);
1152b295832fSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt, PetscErrorCode (*)(TSAdapt, TS, PetscReal, Vec, PetscBool *));
1153d949e4a4SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptHistorySetHistory(TSAdapt, PetscInt n, PetscReal hist[], PetscBool);
115475017583SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptHistorySetTrajectory(TSAdapt, TSTrajectory, PetscBool);
115575017583SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptHistoryGetStep(TSAdapt, PetscInt, PetscReal *, PetscReal *);
1156de50f1caSBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt, PetscInt);
1157cb7ab186SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptDSPSetFilter(TSAdapt, const char *);
1158cb7ab186SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptDSPSetPID(TSAdapt, PetscReal, PetscReal, PetscReal);
115984df9cb4SJed Brown 
116084df9cb4SJed Brown /*S
116187497f52SBarry Smith    TSGLLEAdapt - Abstract object that manages time-step adaptivity for `TSGLLE`
1162ed657a08SJed Brown 
1163ed657a08SJed Brown    Level: beginner
1164ed657a08SJed Brown 
1165195e9b02SBarry Smith    Developer Note:
116687497f52SBarry Smith    This functionality should be replaced by the `TSAdapt`.
116784df9cb4SJed Brown 
11681cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLEAdaptCreate()`, `TSGLLEAdaptType`
1169ed657a08SJed Brown S*/
117026d28e4eSEmil Constantinescu typedef struct _p_TSGLLEAdapt *TSGLLEAdapt;
1171ed657a08SJed Brown 
117276bdecfbSBarry Smith /*J
117387497f52SBarry Smith    TSGLLEAdaptType - String with the name of `TSGLLEAdapt` scheme
1174ed657a08SJed Brown 
1175ed657a08SJed Brown    Level: beginner
1176ed657a08SJed Brown 
1177195e9b02SBarry Smith    Developer Note:
1178195e9b02SBarry Smith    This functionality should be replaced by the `TSAdaptType`.
1179195e9b02SBarry Smith 
11801cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLLEAdaptSetType()`, `TS`
118176bdecfbSBarry Smith J*/
118226d28e4eSEmil Constantinescu typedef const char *TSGLLEAdaptType;
118326d28e4eSEmil Constantinescu #define TSGLLEADAPT_NONE "none"
118426d28e4eSEmil Constantinescu #define TSGLLEADAPT_SIZE "size"
118526d28e4eSEmil Constantinescu #define TSGLLEADAPT_BOTH "both"
1186ed657a08SJed Brown 
118726d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptRegister(const char[], PetscErrorCode (*)(TSGLLEAdapt));
118826d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptInitializePackage(void);
118926d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptFinalizePackage(void);
119026d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate(MPI_Comm, TSGLLEAdapt *);
119126d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt, TSGLLEAdaptType);
119226d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt, const char[]);
119326d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptChoose(TSGLLEAdapt, PetscInt, const PetscInt[], const PetscReal[], const PetscReal[], PetscInt, PetscReal, PetscReal, PetscInt *, PetscReal *, PetscBool *);
119426d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt, PetscViewer);
1195dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetFromOptions(TSGLLEAdapt, PetscOptionItems *);
119626d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt *);
1197ed657a08SJed Brown 
119876bdecfbSBarry Smith /*J
119987497f52SBarry Smith    TSGLLEAcceptType - String with the name of `TSGLLEAccept` scheme
1200ed657a08SJed Brown 
1201ed657a08SJed Brown    Level: beginner
1202ed657a08SJed Brown 
12031cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLLESetAcceptType()`, `TS`, `TSGLLEAccept`
120476bdecfbSBarry Smith J*/
120526d28e4eSEmil Constantinescu typedef const char *TSGLLEAcceptType;
120626d28e4eSEmil Constantinescu #define TSGLLEACCEPT_ALWAYS "always"
1207ed657a08SJed Brown 
1208d1c5d1fcSBarry Smith /*S
12098434afd1SBarry Smith   TSGLLEAcceptFn - A prototype of a `TS` accept function that would be passed to `TSGLLEAcceptRegister()`
1210d1c5d1fcSBarry Smith 
1211d1c5d1fcSBarry Smith   Calling Sequence:
1212d1c5d1fcSBarry Smith + ts  - timestep context
1213d1c5d1fcSBarry Smith . nt - time to end of solution time
1214d1c5d1fcSBarry Smith . h - the proposed step-size
1215d1c5d1fcSBarry Smith . enorm - unknown
1216d1c5d1fcSBarry Smith - accept - output, if the proposal is accepted
1217d1c5d1fcSBarry Smith 
1218d1c5d1fcSBarry Smith   Level: beginner
1219d1c5d1fcSBarry Smith 
1220d1c5d1fcSBarry Smith   Note:
12218434afd1SBarry Smith   The deprecated `TSGLLEAcceptFunction` still works as a replacement for `TSGLLEAcceptFn` *
1222d1c5d1fcSBarry Smith 
12238434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetRHSFunction()`, `DMTSSetRHSFunction()`, `TSIFunctionFn`,
12248434afd1SBarry Smith `TSIJacobianFn`, `TSRHSJacobianFn`, `TSGLLEAcceptRegister()`
1225d1c5d1fcSBarry Smith S*/
12268434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSGLLEAcceptFn)(TS ts, PetscReal nt, PetscReal h, const PetscReal enorm[], PetscBool *accept);
1227d1c5d1fcSBarry Smith 
12288434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSGLLEAcceptFn *TSGLLEAcceptFunction;
1229d1c5d1fcSBarry Smith 
12308434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGLLEAcceptRegister(const char[], TSGLLEAcceptFn *);
1231ed657a08SJed Brown 
123276bdecfbSBarry Smith /*J
1233195e9b02SBarry Smith   TSGLLEType - string with the name of a General Linear `TSGLLE` type
123418b56cb9SJed Brown 
123518b56cb9SJed Brown   Level: beginner
123618b56cb9SJed Brown 
12371cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLESetType()`, `TSGLLERegister()`, `TSGLLEAccept`
123876bdecfbSBarry Smith J*/
123926d28e4eSEmil Constantinescu typedef const char *TSGLLEType;
124026d28e4eSEmil Constantinescu #define TSGLLE_IRKS "irks"
124118b56cb9SJed Brown 
124226d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLERegister(const char[], PetscErrorCode (*)(TS));
124326d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEInitializePackage(void);
124426d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEFinalizePackage(void);
124526d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLESetType(TS, TSGLLEType);
124626d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEGetAdapt(TS, TSGLLEAdapt *);
124726d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLESetAcceptType(TS, TSGLLEAcceptType);
124818b56cb9SJed Brown 
1249020d8f30SJed Brown /*J
1250195e9b02SBarry Smith    TSEIMEXType - String with the name of an Extrapolated IMEX `TSEIMEX` type
12517d5697caSHong Zhang 
12527d5697caSHong Zhang    Level: beginner
12537d5697caSHong Zhang 
12541cc06b55SBarry Smith .seealso: [](ch_ts), `TSEIMEXSetType()`, `TS`, `TSEIMEX`, `TSEIMEXRegister()`
12557d5697caSHong Zhang J*/
12567d5697caSHong Zhang #define TSEIMEXType char *
12577d5697caSHong Zhang 
12587d5697caSHong Zhang PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts, PetscInt);
12597d5697caSHong Zhang PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts, PetscInt, PetscInt);
12607d5697caSHong Zhang PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS, PetscBool);
12617d5697caSHong Zhang 
12627d5697caSHong Zhang /*J
1263195e9b02SBarry Smith    TSRKType - String with the name of a Runge-Kutta `TSRK` type
1264f68a32c8SEmil Constantinescu 
1265f68a32c8SEmil Constantinescu    Level: beginner
1266f68a32c8SEmil Constantinescu 
12671cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSRKSetType()`, `TS`, `TSRK`, `TSRKRegister()`
1268f68a32c8SEmil Constantinescu J*/
1269f68a32c8SEmil Constantinescu typedef const char *TSRKType;
1270f68a32c8SEmil Constantinescu #define TSRK1FE "1fe"
1271fdefd5e5SDebojyoti Ghosh #define TSRK2A  "2a"
1272e7685601SHong Zhang #define TSRK2B  "2b"
1273f68a32c8SEmil Constantinescu #define TSRK3   "3"
1274fdefd5e5SDebojyoti Ghosh #define TSRK3BS "3bs"
1275f68a32c8SEmil Constantinescu #define TSRK4   "4"
1276fdefd5e5SDebojyoti Ghosh #define TSRK5F  "5f"
1277fdefd5e5SDebojyoti Ghosh #define TSRK5DP "5dp"
127805e23783SLisandro Dalcin #define TSRK5BS "5bs"
127937acaa02SHendrik Ranocha #define TSRK6VR "6vr"
128037acaa02SHendrik Ranocha #define TSRK7VR "7vr"
128137acaa02SHendrik Ranocha #define TSRK8VR "8vr"
128205e23783SLisandro Dalcin 
12832ea87ba9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRKGetOrder(TS, PetscInt *);
12840f7a28cdSStefano Zampini PETSC_EXTERN PetscErrorCode TSRKGetType(TS, TSRKType *);
12850f7a28cdSStefano Zampini PETSC_EXTERN PetscErrorCode TSRKSetType(TS, TSRKType);
12860f7a28cdSStefano Zampini PETSC_EXTERN PetscErrorCode TSRKGetTableau(TS, PetscInt *, const PetscReal **, const PetscReal **, const PetscReal **, const PetscReal **, PetscInt *, const PetscReal **, PetscBool *);
12870fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSRKSetMultirate(TS, PetscBool);
12880fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSRKGetMultirate(TS, PetscBool *);
1289f68a32c8SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
1290f68a32c8SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void);
1291560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void);
1292f68a32c8SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void);
1293f68a32c8SEmil Constantinescu 
1294f68a32c8SEmil Constantinescu /*J
1295af27ebaaSBarry Smith    TSMPRKType - String with the name of a partitioned Runge-Kutta `TSMPRK` type
12969f537349Sluzhanghpp 
12979f537349Sluzhanghpp    Level: beginner
12989f537349Sluzhanghpp 
12991cc06b55SBarry Smith .seealso: [](ch_ts), `TSMPRKSetType()`, `TS`, `TSMPRK`, `TSMPRKRegister()`
13009f537349Sluzhanghpp J*/
13014b84eec9SHong Zhang typedef const char *TSMPRKType;
130219c2959aSHong Zhang #define TSMPRK2A22 "2a22"
130319c2959aSHong Zhang #define TSMPRK2A23 "2a23"
130419c2959aSHong Zhang #define TSMPRK2A32 "2a32"
130519c2959aSHong Zhang #define TSMPRK2A33 "2a33"
13064b84eec9SHong Zhang #define TSMPRKP2   "p2"
13074b84eec9SHong Zhang #define TSMPRKP3   "p3"
13089f537349Sluzhanghpp 
13094b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKGetType(TS ts, TSMPRKType *);
13104b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKSetType(TS ts, TSMPRKType);
131179bc8a77SHong 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[]);
13124b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKInitializePackage(void);
13134b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKFinalizePackage(void);
13144b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKRegisterDestroy(void);
13159f537349Sluzhanghpp 
13169f537349Sluzhanghpp /*J
1317195e9b02SBarry Smith    TSIRKType - String with the name of an implicit Runge-Kutta `TSIRK` type
1318d2567f34SHong Zhang 
1319d2567f34SHong Zhang    Level: beginner
1320d2567f34SHong Zhang 
13211cc06b55SBarry Smith .seealso: [](ch_ts), `TSIRKSetType()`, `TS`, `TSIRK`, `TSIRKRegister()`
1322d2567f34SHong Zhang J*/
1323d2567f34SHong Zhang typedef const char *TSIRKType;
1324d2567f34SHong Zhang #define TSIRKGAUSS "gauss"
1325d2567f34SHong Zhang 
1326d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKGetType(TS, TSIRKType *);
1327d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKSetType(TS, TSIRKType);
1328d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKGetNumStages(TS, PetscInt *);
1329d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKSetNumStages(TS, PetscInt);
1330d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKRegister(const char[], PetscErrorCode (*function)(TS));
1331d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKTableauCreate(TS, PetscInt, const PetscReal *, const PetscReal *, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, const PetscScalar *);
1332d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKInitializePackage(void);
1333d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKFinalizePackage(void);
1334d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKRegisterDestroy(void);
1335d2567f34SHong Zhang 
1336d2567f34SHong Zhang /*J
1337195e9b02SBarry Smith    TSGLEEType - String with the name of a General Linear with Error Estimation `TSGLEE` type
1338b6a60446SDebojyoti Ghosh 
1339b6a60446SDebojyoti Ghosh    Level: beginner
1340b6a60446SDebojyoti Ghosh 
13411cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEESetType()`, `TS`, `TSGLEE`, `TSGLEERegister()`
1342b6a60446SDebojyoti Ghosh J*/
1343b6a60446SDebojyoti Ghosh typedef const char *TSGLEEType;
1344ab8c5912SEmil Constantinescu #define TSGLEEi1      "BE1"
1345b6a60446SDebojyoti Ghosh #define TSGLEE23      "23"
1346b6a60446SDebojyoti Ghosh #define TSGLEE24      "24"
1347b6a60446SDebojyoti Ghosh #define TSGLEE25I     "25i"
1348b6a60446SDebojyoti Ghosh #define TSGLEE35      "35"
1349b6a60446SDebojyoti Ghosh #define TSGLEEEXRK2A  "exrk2a"
1350b6a60446SDebojyoti Ghosh #define TSGLEERK32G1  "rk32g1"
1351b6a60446SDebojyoti Ghosh #define TSGLEERK285EX "rk285ex"
135216a05f60SBarry Smith 
1353b6a60446SDebojyoti Ghosh /*J
1354195e9b02SBarry Smith    TSGLEEMode - String with the mode of error estimation for a General Linear with Error Estimation `TSGLEE` type
1355b6a60446SDebojyoti Ghosh 
1356b6a60446SDebojyoti Ghosh    Level: beginner
1357b6a60446SDebojyoti Ghosh 
13581cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEESetMode()`, `TS`, `TSGLEE`, `TSGLEERegister()`
1359b6a60446SDebojyoti Ghosh J*/
1360b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEEGetType(TS ts, TSGLEEType *);
1361b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEESetType(TS ts, TSGLEEType);
136257df6a1bSDebojyoti 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[]);
13634bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSGLEERegisterAll(void);
1364b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEEFinalizePackage(void);
1365b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEEInitializePackage(void);
1366b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEERegisterDestroy(void);
1367b6a60446SDebojyoti Ghosh 
1368b6a60446SDebojyoti Ghosh /*J
1369195e9b02SBarry Smith    TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX `TSARKIMEX` type
1370020d8f30SJed Brown 
1371020d8f30SJed Brown    Level: beginner
1372020d8f30SJed Brown 
13731cc06b55SBarry Smith .seealso: [](ch_ts), `TSARKIMEXSetType()`, `TS`, `TSARKIMEX`, `TSARKIMEXRegister()`
1374020d8f30SJed Brown J*/
137519fd82e9SBarry Smith typedef const char *TSARKIMEXType;
1376e817cc15SEmil Constantinescu #define TSARKIMEX1BEE   "1bee"
13771f80e275SEmil Constantinescu #define TSARKIMEXA2     "a2"
13781f80e275SEmil Constantinescu #define TSARKIMEXL2     "l2"
13791f80e275SEmil Constantinescu #define TSARKIMEXARS122 "ars122"
13801f80e275SEmil Constantinescu #define TSARKIMEX2C     "2c"
13818a381b04SJed Brown #define TSARKIMEX2D     "2d"
1382a3a57f36SJed Brown #define TSARKIMEX2E     "2e"
13836cf0794eSJed Brown #define TSARKIMEXPRSSP2 "prssp2"
13848a381b04SJed Brown #define TSARKIMEX3      "3"
13856cf0794eSJed Brown #define TSARKIMEXBPR3   "bpr3"
13866cf0794eSJed Brown #define TSARKIMEXARS443 "ars443"
13878a381b04SJed Brown #define TSARKIMEX4      "4"
13888a381b04SJed Brown #define TSARKIMEX5      "5"
138919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts, TSARKIMEXType *);
139019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts, TSARKIMEXType);
1391014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS, PetscBool);
13923a28c0e5SStefano Zampini PETSC_EXTERN PetscErrorCode TSARKIMEXGetFullyImplicit(TS, PetscBool *);
13933a2a065fSHong Zhang PETSC_EXTERN PetscErrorCode TSARKIMEXSetFastSlowSplit(TS, PetscBool);
13943a2a065fSHong Zhang PETSC_EXTERN PetscErrorCode TSARKIMEXGetFastSlowSplit(TS, PetscBool *);
139519fd82e9SBarry 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[]);
1396607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void);
1397560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
1398014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);
13998a381b04SJed Brown 
1400020d8f30SJed Brown /*J
1401d27334e2SStefano Zampini    TSDIRKType - String with the name of a Diagonally Implicit Runge-Kutta `TSDIRK` type
1402d27334e2SStefano Zampini 
1403d27334e2SStefano Zampini    Level: beginner
1404d27334e2SStefano Zampini 
1405d27334e2SStefano Zampini .seealso: [](ch_ts), `TSDIRKSetType()`, `TS`, `TSDIRK`, `TSDIRKRegister()`
1406d27334e2SStefano Zampini J*/
1407d27334e2SStefano Zampini typedef const char *TSDIRKType;
1408d27334e2SStefano Zampini #define TSDIRKS212      "s212"
14093405e92cSStefano Zampini #define TSDIRKES122SAL  "es122sal"
1410d27334e2SStefano Zampini #define TSDIRKES213SAL  "es213sal"
1411d27334e2SStefano Zampini #define TSDIRKES324SAL  "es324sal"
1412d27334e2SStefano Zampini #define TSDIRKES325SAL  "es325sal"
1413d27334e2SStefano Zampini #define TSDIRK657A      "657a"
1414d27334e2SStefano Zampini #define TSDIRKES648SA   "es648sa"
1415d27334e2SStefano Zampini #define TSDIRK658A      "658a"
1416d27334e2SStefano Zampini #define TSDIRKS659A     "s659a"
1417d27334e2SStefano Zampini #define TSDIRK7510SAL   "7510sal"
1418d27334e2SStefano Zampini #define TSDIRKES7510SA  "es7510sa"
1419d27334e2SStefano Zampini #define TSDIRK759A      "759a"
1420d27334e2SStefano Zampini #define TSDIRKS7511SAL  "s7511sal"
1421d27334e2SStefano Zampini #define TSDIRK8614A     "8614a"
1422d27334e2SStefano Zampini #define TSDIRK8616SAL   "8616sal"
1423d27334e2SStefano Zampini #define TSDIRKES8516SAL "es8516sal"
1424d27334e2SStefano Zampini PETSC_EXTERN PetscErrorCode TSDIRKGetType(TS ts, TSDIRKType *);
1425d27334e2SStefano Zampini PETSC_EXTERN PetscErrorCode TSDIRKSetType(TS ts, TSDIRKType);
1426d27334e2SStefano Zampini PETSC_EXTERN PetscErrorCode TSDIRKRegister(TSDIRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
1427d27334e2SStefano Zampini 
1428d27334e2SStefano Zampini /*J
1429195e9b02SBarry Smith    TSRosWType - String with the name of a Rosenbrock-W `TSROSW` type
1430020d8f30SJed Brown 
1431020d8f30SJed Brown    Level: beginner
1432020d8f30SJed Brown 
14331cc06b55SBarry Smith .seealso: [](ch_ts), `TSRosWSetType()`, `TS`, `TSROSW`, `TSRosWRegister()`
1434020d8f30SJed Brown J*/
143519fd82e9SBarry Smith typedef const char *TSRosWType;
143661692a83SJed Brown #define TSROSW2M          "2m"
143761692a83SJed Brown #define TSROSW2P          "2p"
1438fe7e6d57SJed Brown #define TSROSWRA3PW       "ra3pw"
1439fe7e6d57SJed Brown #define TSROSWRA34PW2     "ra34pw2"
1440337ef93bSJed Brown #define TSROSWR34PRW      "r34prw"
1441337ef93bSJed Brown #define TSROSWR3PRL2      "r3prl2"
1442ef3c5b88SJed Brown #define TSROSWRODAS3      "rodas3"
144348665b53SJed Brown #define TSROSWRODASPR     "rodaspr"
1444337ef93bSJed Brown #define TSROSWRODASPR2    "rodaspr2"
1445ef3c5b88SJed Brown #define TSROSWSANDU3      "sandu3"
1446961f28d0SJed Brown #define TSROSWASSP3P3S1C  "assp3p3s1c"
1447961f28d0SJed Brown #define TSROSWLASSP3P4S2C "lassp3p4s2c"
144843b21953SEmil Constantinescu #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
1449753f8adbSEmil Constantinescu #define TSROSWARK3        "ark3"
14503606a31eSEmil Constantinescu #define TSROSWTHETA1      "theta1"
14513606a31eSEmil Constantinescu #define TSROSWTHETA2      "theta2"
145242faf41dSJed Brown #define TSROSWGRK4T       "grk4t"
145342faf41dSJed Brown #define TSROSWSHAMP4      "shamp4"
145442faf41dSJed Brown #define TSROSWVELDD4      "veldd4"
145542faf41dSJed Brown #define TSROSW4L          "4l"
1456b1c69cc3SEmil Constantinescu 
14576bd7aeb5SHong Zhang PETSC_EXTERN PetscErrorCode TSRosWGetType(TS, TSRosWType *);
14586bd7aeb5SHong Zhang PETSC_EXTERN PetscErrorCode TSRosWSetType(TS, TSRosWType);
1459014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS, PetscBool);
146019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
146119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal);
1462607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void);
1463560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
1464014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);
1465e27a552bSJed Brown 
1466211a84d6SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSBDFSetOrder(TS, PetscInt);
1467211a84d6SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSBDFGetOrder(TS, PetscInt *);
1468211a84d6SLisandro Dalcin 
14696bd7aeb5SHong Zhang /*J
1470195e9b02SBarry Smith   TSBasicSymplecticType - String with the name of a basic symplectic integration `TSBASICSYMPLECTIC` type
14716bd7aeb5SHong Zhang 
14726bd7aeb5SHong Zhang   Level: beginner
14736bd7aeb5SHong Zhang 
14741cc06b55SBarry Smith .seealso: [](ch_ts), `TSBasicSymplecticSetType()`, `TS`, `TSBASICSYMPLECTIC`, `TSBasicSymplecticRegister()`
14756bd7aeb5SHong Zhang J*/
1476e49d4f37SHong Zhang typedef const char *TSBasicSymplecticType;
1477e49d4f37SHong Zhang #define TSBASICSYMPLECTICSIEULER   "1"
1478e49d4f37SHong Zhang #define TSBASICSYMPLECTICVELVERLET "2"
1479e49d4f37SHong Zhang #define TSBASICSYMPLECTIC3         "3"
1480e49d4f37SHong Zhang #define TSBASICSYMPLECTIC4         "4"
1481e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticSetType(TS, TSBasicSymplecticType);
1482e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticGetType(TS, TSBasicSymplecticType *);
1483e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegister(TSBasicSymplecticType, PetscInt, PetscInt, PetscReal[], PetscReal[]);
14844bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterAll(void);
1485e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticInitializePackage(void);
1486e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticFinalizePackage(void);
1487e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterDestroy(void);
14886bd7aeb5SHong Zhang 
1489db4653c2SDaniel Finn /*J
1490195e9b02SBarry Smith   TSDISCGRAD - The Discrete Gradient integrator is a timestepper for Hamiltonian systems designed to conserve the first integral (energy),
149187497f52SBarry Smith   but also has the property for some systems of monotonicity in a functional.
1492db4653c2SDaniel Finn 
1493db4653c2SDaniel Finn   Level: beginner
1494db4653c2SDaniel Finn 
14951cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, TSDiscGradSetFormulation()`, `TSDiscGradGetFormulation()`
1496db4653c2SDaniel Finn J*/
149740be0ff1SMatthew 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 *);
14984bf303faSJacob 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 *);
1499db4653c2SDaniel Finn PETSC_EXTERN PetscErrorCode TSDiscGradIsGonzalez(TS, PetscBool *);
1500db4653c2SDaniel Finn PETSC_EXTERN PetscErrorCode TSDiscGradUseGonzalez(TS, PetscBool);
150140be0ff1SMatthew G. Knepley 
150283e2fdc7SBarry Smith /*
15030f3b3ca1SHong Zhang        PETSc interface to Sundials
150483e2fdc7SBarry Smith */
1505e808b789SPatrick Sanan #ifdef PETSC_HAVE_SUNDIALS2
15069371c9d4SSatish Balay typedef enum {
15079371c9d4SSatish Balay   SUNDIALS_ADAMS = 1,
15089371c9d4SSatish Balay   SUNDIALS_BDF   = 2
15099371c9d4SSatish Balay } TSSundialsLmmType;
15106a6fc655SJed Brown PETSC_EXTERN const char *const TSSundialsLmmTypes[];
15119371c9d4SSatish Balay typedef enum {
15129371c9d4SSatish Balay   SUNDIALS_MODIFIED_GS  = 1,
15139371c9d4SSatish Balay   SUNDIALS_CLASSICAL_GS = 2
15149371c9d4SSatish Balay } TSSundialsGramSchmidtType;
15156a6fc655SJed Brown PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
15164bf303faSJacob Faibussowitsch 
1517014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS, TSSundialsLmmType);
1518014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS, PC *);
1519014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS, PetscReal, PetscReal);
1520014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS, PetscReal);
1521014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS, PetscReal);
1522014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS, PetscInt *, PetscInt *);
1523014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS, TSSundialsGramSchmidtType);
1524014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS, PetscInt);
1525014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS, PetscReal);
1526014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS, PetscBool);
1527014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS, PetscInt);
1528c4e80e11SFlorian PETSC_EXTERN PetscErrorCode TSSundialsSetMaxord(TS, PetscInt);
1529918687b8SPatrick Sanan PETSC_EXTERN PetscErrorCode TSSundialsSetUseDense(TS, PetscBool);
15306dc17ff5SMatthew Knepley #endif
153183e2fdc7SBarry Smith 
1532014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS, PetscReal);
1533014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS, PetscReal *);
1534014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS, PetscBool *);
1535014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS, PetscBool);
15360de4c49aSJed Brown 
1537014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS, PetscReal);
1538014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS, PetscReal, PetscReal, PetscReal);
1539014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS, PetscReal *, PetscReal *, PetscReal *);
154088df8a41SLisandro Dalcin 
1541818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAlpha2SetRadius(TS, PetscReal);
1542818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAlpha2SetParams(TS, PetscReal, PetscReal, PetscReal, PetscReal);
1543818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAlpha2GetParams(TS, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
1544818efac9SLisandro Dalcin 
1545220f924aSDavid Kamensky /*S
15468434afd1SBarry Smith   TSAlpha2PredictorFn - A callback to set the predictor (i.e., the initial guess for the nonlinear solver) in
1547220f924aSDavid Kamensky   a second-order generalized-alpha time integrator.
1548220f924aSDavid Kamensky 
1549220f924aSDavid Kamensky   Calling Sequence:
1550220f924aSDavid Kamensky + ts   - the `TS` context obtained from `TSCreate()`
1551220f924aSDavid Kamensky . X0   - the previous time step's state vector
1552220f924aSDavid Kamensky . V0   - the previous time step's first derivative of the state vector
1553220f924aSDavid Kamensky . A0   - the previous time step's second derivative of the state vector
1554220f924aSDavid Kamensky . X1   - the vector into which the initial guess for the current time step will be written
1555220f924aSDavid Kamensky - ctx  - [optional] user-defined context for the predictor evaluation routine (may be `NULL`)
1556220f924aSDavid Kamensky 
1557220f924aSDavid Kamensky   Level: intermediate
1558220f924aSDavid Kamensky 
1559d1c5d1fcSBarry Smith   Note:
15608434afd1SBarry Smith   The deprecated `TSAlpha2Predictor` still works as a replacement for `TSAlpha2PredictorFn` *.
1561d1c5d1fcSBarry Smith 
1562220f924aSDavid Kamensky .seealso: [](ch_ts), `TS`, `TSAlpha2SetPredictor()`
1563220f924aSDavid Kamensky S*/
15648434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSAlpha2PredictorFn)(TS ts, Vec X0, Vec V0, Vec A0, Vec X1, void *ctx);
1565d1c5d1fcSBarry Smith 
15668434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSAlpha2PredictorFn *TSAlpha2Predictor;
1567d1c5d1fcSBarry Smith 
15688434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSAlpha2SetPredictor(TS, TSAlpha2PredictorFn *, void *ctx);
1569220f924aSDavid Kamensky 
1570014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetDM(TS, DM);
1571014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetDM(TS, DM *);
15726c699258SBarry Smith 
1573014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES, Vec, Vec, void *);
1574d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES, Vec, Mat, Mat, void *);
15750f5c6efeSJed Brown 
1576f3b1f45cSBarry Smith PETSC_EXTERN PetscErrorCode TSRHSJacobianTest(TS, PetscBool *);
1577f3b1f45cSBarry Smith PETSC_EXTERN PetscErrorCode TSRHSJacobianTestTranspose(TS, PetscBool *);
1578aad739acSMatthew G. Knepley 
15792e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSGetComputeInitialCondition(TS, PetscErrorCode (**)(TS, Vec));
15802e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSSetComputeInitialCondition(TS, PetscErrorCode (*)(TS, Vec));
15812e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSComputeInitialCondition(TS, Vec);
15822e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSGetComputeExactError(TS, PetscErrorCode (**)(TS, Vec, Vec));
1583aad739acSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSSetComputeExactError(TS, PetscErrorCode (*)(TS, Vec, Vec));
1584aad739acSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSComputeExactError(TS, Vec, Vec);
1585f2ed2dc7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscConvEstUseTS(PetscConvEst, PetscBool);
1586d60b7d5cSBarry Smith 
1587d60b7d5cSBarry Smith PETSC_EXTERN PetscErrorCode TSSetMatStructure(TS, MatStructure);
1588