xref: /petsc/include/petscts.h (revision d27334e213e4e789cf4a49a6816302757d1e7cf0)
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 */
56524c165SJacob Faibussowitsch #ifndef PETSCTS_H
626bd1501SBarry Smith #define PETSCTS_H
7ac09b921SBarry Smith 
82c8e378dSBarry Smith #include <petscsnes.h>
9900f6b5bSMatthew G. Knepley #include <petscconvest.h>
10818ad0c1SBarry Smith 
1114d0ab18SJacob Faibussowitsch /*I <petscts.h> I*/
1214d0ab18SJacob Faibussowitsch 
13ac09b921SBarry Smith /* SUBMANSEC = TS */
14ac09b921SBarry Smith 
15435da068SBarry Smith /*S
16435da068SBarry Smith      TS - Abstract PETSc object that manages all time-steppers (ODE integrators)
17435da068SBarry Smith 
18435da068SBarry Smith    Level: beginner
19435da068SBarry Smith 
201cc06b55SBarry Smith .seealso: [](integrator_table), [](ch_ts), `TSCreate()`, `TSSetType()`, `TSType`, `SNES`, `KSP`, `PC`, `TSDestroy()`
21435da068SBarry Smith S*/
22f09e8eb9SSatish Balay typedef struct _p_TS *TS;
23435da068SBarry Smith 
2476bdecfbSBarry Smith /*J
2587497f52SBarry Smith     TSType - String with the name of a PETSc `TS` method.
26435da068SBarry Smith 
27435da068SBarry Smith    Level: beginner
28435da068SBarry Smith 
291cc06b55SBarry Smith .seealso: [](integrator_table), [](ch_ts), `TSSetType()`, `TS`, `TSRegister()`
3076bdecfbSBarry Smith J*/
3119fd82e9SBarry Smith typedef const char *TSType;
329596e0b4SJed Brown #define TSEULER           "euler"
339596e0b4SJed Brown #define TSBEULER          "beuler"
34e49d4f37SHong Zhang #define TSBASICSYMPLECTIC "basicsymplectic"
359596e0b4SJed Brown #define TSPSEUDO          "pseudo"
364d91e141SJed Brown #define TSCN              "cn"
379596e0b4SJed Brown #define TSSUNDIALS        "sundials"
384d91e141SJed Brown #define TSRK              "rk"
399596e0b4SJed Brown #define TSPYTHON          "python"
409596e0b4SJed Brown #define TSTHETA           "theta"
4188df8a41SLisandro Dalcin #define TSALPHA           "alpha"
42818efac9SLisandro Dalcin #define TSALPHA2          "alpha2"
4326d28e4eSEmil Constantinescu #define TSGLLE            "glle"
44b6a60446SDebojyoti Ghosh #define TSGLEE            "glee"
45ef7bb5aaSJed Brown #define TSSSP             "ssp"
468a381b04SJed Brown #define TSARKIMEX         "arkimex"
47e27a552bSJed Brown #define TSROSW            "rosw"
487d5697caSHong Zhang #define TSEIMEX           "eimex"
49abadfa0fSMatthew G. Knepley #define TSMIMEX           "mimex"
50211a84d6SLisandro Dalcin #define TSBDF             "bdf"
51d249a78fSBarry Smith #define TSRADAU5          "radau5"
524b84eec9SHong Zhang #define TSMPRK            "mprk"
5340be0ff1SMatthew G. Knepley #define TSDISCGRAD        "discgrad"
54d2567f34SHong Zhang #define TSIRK             "irk"
55*d27334e2SStefano Zampini #define TSDIRK            "dirk"
5640be0ff1SMatthew G. Knepley 
57435da068SBarry Smith /*E
5887497f52SBarry Smith     TSProblemType - Determines the type of problem this `TS` object is to be used to solve
59435da068SBarry Smith 
60195e9b02SBarry Smith    Values:
61195e9b02SBarry Smith  + `TS_LINEAR` - a linear ODE or DAE
62195e9b02SBarry Smith  - `TS_NONLINEAR` - a nonlinear ODE or DAE
63195e9b02SBarry Smith 
64435da068SBarry Smith    Level: beginner
65435da068SBarry Smith 
661cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSCreate()`
67435da068SBarry Smith E*/
689371c9d4SSatish Balay typedef enum {
699371c9d4SSatish Balay   TS_LINEAR,
709371c9d4SSatish Balay   TS_NONLINEAR
719371c9d4SSatish Balay } TSProblemType;
72818ad0c1SBarry Smith 
73b93fea0eSJed Brown /*E
7487497f52SBarry Smith    TSEquationType - type of `TS` problem that is solved
75e817cc15SEmil Constantinescu 
76e817cc15SEmil Constantinescu    Level: beginner
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 
831cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGetEquationType()`, `TSSetEquationType()`
84e817cc15SEmil Constantinescu E*/
85e817cc15SEmil Constantinescu typedef enum {
86e817cc15SEmil Constantinescu   TS_EQ_UNSPECIFIED               = -1,
87e817cc15SEmil Constantinescu   TS_EQ_EXPLICIT                  = 0,
88e817cc15SEmil Constantinescu   TS_EQ_ODE_EXPLICIT              = 1,
89e817cc15SEmil Constantinescu   TS_EQ_DAE_SEMI_EXPLICIT_INDEX1  = 100,
90e817cc15SEmil Constantinescu   TS_EQ_DAE_SEMI_EXPLICIT_INDEX2  = 200,
91e817cc15SEmil Constantinescu   TS_EQ_DAE_SEMI_EXPLICIT_INDEX3  = 300,
92e817cc15SEmil Constantinescu   TS_EQ_DAE_SEMI_EXPLICIT_INDEXHI = 500,
93e817cc15SEmil Constantinescu   TS_EQ_IMPLICIT                  = 1000,
94e817cc15SEmil Constantinescu   TS_EQ_ODE_IMPLICIT              = 1001,
95e817cc15SEmil Constantinescu   TS_EQ_DAE_IMPLICIT_INDEX1       = 1100,
96e817cc15SEmil Constantinescu   TS_EQ_DAE_IMPLICIT_INDEX2       = 1200,
97e817cc15SEmil Constantinescu   TS_EQ_DAE_IMPLICIT_INDEX3       = 1300,
9819436ca2SJed Brown   TS_EQ_DAE_IMPLICIT_INDEXHI      = 1500
99e817cc15SEmil Constantinescu } TSEquationType;
100e817cc15SEmil Constantinescu PETSC_EXTERN const char *const *TSEquationTypes;
101e817cc15SEmil Constantinescu 
102e817cc15SEmil Constantinescu /*E
10316a05f60SBarry Smith    TSConvergedReason - reason a `TS` method has converged (integrated to the requested time) or not
104b93fea0eSJed Brown 
105195e9b02SBarry Smith    Values:
106195e9b02SBarry Smith +  `TS_CONVERGED_ITERATING`          - this only occurs if `TSGetConvergedReason()` is called during the `TSSolve()`
107195e9b02SBarry Smith .  `TS_CONVERGED_TIME`               - the final time was reached
108195e9b02SBarry Smith .  `TS_CONVERGED_ITS`                - the maximum number of iterations (time-steps) was reached prior to the final time
109195e9b02SBarry Smith .  `TS_CONVERGED_USER`               - user requested termination
110195e9b02SBarry Smith .  `TS_CONVERGED_EVENT`              - user requested termination on event detection
111195e9b02SBarry Smith .  `TS_CONVERGED_PSEUDO_FATOL`       - stops when function norm decreased by a set amount, used only for `TSPSEUDO`
112195e9b02SBarry Smith .  `TS_CONVERGED_PSEUDO_FRTOL`       - stops when function norm decreases below a set amount, used only for `TSPSEUDO`
113195e9b02SBarry Smith .  `TS_DIVERGED_NONLINEAR_SOLVE`     - too many nonlinear solve failures have occurred
114195e9b02SBarry Smith .  `TS_DIVERGED_STEP_REJECTED`       - too many steps were rejected
115195e9b02SBarry Smith .  `TSFORWARD_DIVERGED_LINEAR_SOLVE` - tangent linear solve failed
116195e9b02SBarry Smith -  `TSADJOINT_DIVERGED_LINEAR_SOLVE` - transposed linear solve failed
117195e9b02SBarry Smith 
118b93fea0eSJed Brown    Level: beginner
119b93fea0eSJed Brown 
1201cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGetConvergedReason()`
121b93fea0eSJed Brown E*/
122193ac0bcSJed Brown typedef enum {
123193ac0bcSJed Brown   TS_CONVERGED_ITERATING          = 0,
124193ac0bcSJed Brown   TS_CONVERGED_TIME               = 1,
125193ac0bcSJed Brown   TS_CONVERGED_ITS                = 2,
126d6ad946cSShri Abhyankar   TS_CONVERGED_USER               = 3,
1272b7db910SShri Abhyankar   TS_CONVERGED_EVENT              = 4,
1283118ae5eSBarry Smith   TS_CONVERGED_PSEUDO_FATOL       = 5,
1293118ae5eSBarry Smith   TS_CONVERGED_PSEUDO_FRTOL       = 6,
130193ac0bcSJed Brown   TS_DIVERGED_NONLINEAR_SOLVE     = -1,
131b5b4017aSHong Zhang   TS_DIVERGED_STEP_REJECTED       = -2,
13205c9950eSHong Zhang   TSFORWARD_DIVERGED_LINEAR_SOLVE = -3,
13305c9950eSHong Zhang   TSADJOINT_DIVERGED_LINEAR_SOLVE = -4
134193ac0bcSJed Brown } TSConvergedReason;
135014dd563SJed Brown PETSC_EXTERN const char *const *TSConvergedReasons;
1361957e957SBarry Smith 
137b93fea0eSJed Brown /*MC
13887497f52SBarry Smith    TS_CONVERGED_ITERATING - this only occurs if `TSGetConvergedReason()` is called during the `TSSolve()`
139b93fea0eSJed Brown 
140b93fea0eSJed Brown    Level: beginner
141b93fea0eSJed Brown 
1421cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`
143b93fea0eSJed Brown M*/
144b93fea0eSJed Brown 
145b93fea0eSJed Brown /*MC
146b93fea0eSJed Brown    TS_CONVERGED_TIME - the final time was reached
147b93fea0eSJed Brown 
148b93fea0eSJed Brown    Level: beginner
149b93fea0eSJed Brown 
1501cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxTime()`, `TSGetMaxTime()`, `TSGetSolveTime()`
151b93fea0eSJed Brown M*/
152b93fea0eSJed Brown 
153b93fea0eSJed Brown /*MC
1541957e957SBarry Smith    TS_CONVERGED_ITS - the maximum number of iterations (time-steps) was reached prior to the final time
155b93fea0eSJed Brown 
156b93fea0eSJed Brown    Level: beginner
157b93fea0eSJed Brown 
1581cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxSteps()`, `TSGetMaxSteps()`
159b93fea0eSJed Brown M*/
1601957e957SBarry Smith 
161d6ad946cSShri Abhyankar /*MC
162d6ad946cSShri Abhyankar    TS_CONVERGED_USER - user requested termination
163d6ad946cSShri Abhyankar 
164d6ad946cSShri Abhyankar    Level: beginner
165d6ad946cSShri Abhyankar 
1661cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`
167b93fea0eSJed Brown M*/
168b93fea0eSJed Brown 
169b93fea0eSJed Brown /*MC
170d76fc68cSShri Abhyankar    TS_CONVERGED_EVENT - user requested termination on event detection
171d76fc68cSShri Abhyankar 
172d76fc68cSShri Abhyankar    Level: beginner
173d76fc68cSShri Abhyankar 
1741cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`
175d76fc68cSShri Abhyankar M*/
176d76fc68cSShri Abhyankar 
177d76fc68cSShri Abhyankar /*MC
17887497f52SBarry Smith    TS_CONVERGED_PSEUDO_FRTOL - stops when function norm decreased by a set amount, used only for `TSPSEUDO`
1793118ae5eSBarry Smith 
1803c7db156SBarry Smith    Options Database Key:
181d5b43468SJose E. Roman .   -ts_pseudo_frtol <rtol> - use specified rtol
1823118ae5eSBarry Smith 
183195e9b02SBarry Smith    Level: beginner
184195e9b02SBarry Smith 
1851cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`, `TS_CONVERGED_PSEUDO_FATOL`
1863118ae5eSBarry Smith M*/
1873118ae5eSBarry Smith 
1883118ae5eSBarry Smith /*MC
18987497f52SBarry Smith    TS_CONVERGED_PSEUDO_FATOL - stops when function norm decreases below a set amount, used only for `TSPSEUDO`
1903118ae5eSBarry Smith 
1913c7db156SBarry Smith    Options Database Key:
19267b8a455SSatish Balay .   -ts_pseudo_fatol <atol> - use specified atol
1933118ae5eSBarry Smith 
194195e9b02SBarry Smith    Level: beginner
195195e9b02SBarry Smith 
1961cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`, `TS_CONVERGED_PSEUDO_FRTOL`
1973118ae5eSBarry Smith M*/
1983118ae5eSBarry Smith 
1993118ae5eSBarry Smith /*MC
200b93fea0eSJed Brown    TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed
201b93fea0eSJed Brown 
202b93fea0eSJed Brown    Level: beginner
203b93fea0eSJed Brown 
204195e9b02SBarry Smith    Note:
205195e9b02SBarry Smith     See `TSSetMaxSNESFailures()` for how to allow more nonlinear solver failures.
2061957e957SBarry Smith 
2071cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSGetSNES()`, `SNESGetConvergedReason()`, `TSSetMaxSNESFailures()`
208b93fea0eSJed Brown M*/
209b93fea0eSJed Brown 
210b93fea0eSJed Brown /*MC
211b93fea0eSJed Brown    TS_DIVERGED_STEP_REJECTED - too many steps were rejected
212b93fea0eSJed Brown 
213b93fea0eSJed Brown    Level: beginner
214b93fea0eSJed Brown 
21595452b02SPatrick Sanan    Notes:
216195e9b02SBarry Smith     See `TSSetMaxStepRejections()` for how to allow more step rejections.
2171957e957SBarry Smith 
2181cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxStepRejections()`
219b93fea0eSJed Brown M*/
220b93fea0eSJed Brown 
221d6ad946cSShri Abhyankar /*E
222d6ad946cSShri Abhyankar    TSExactFinalTimeOption - option for handling of final time step
223d6ad946cSShri Abhyankar 
224195e9b02SBarry Smith    Values:
22516a05f60SBarry Smith +  `TS_EXACTFINALTIME_STEPOVER`    - Don't do anything if requested final time is exceeded
226195e9b02SBarry Smith .  `TS_EXACTFINALTIME_INTERPOLATE` - Interpolate back to final time
22716a05f60SBarry Smith -  `TS_EXACTFINALTIME_MATCHSTEP` - Adapt final time step to match the final time requested
228195e9b02SBarry Smith 
229d6ad946cSShri Abhyankar    Level: beginner
230d6ad946cSShri Abhyankar 
2311cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGetConvergedReason()`, `TSSetExactFinalTime()`, `TSGetExactFinalTime()`
232d6ad946cSShri Abhyankar E*/
2339371c9d4SSatish Balay typedef enum {
2349371c9d4SSatish Balay   TS_EXACTFINALTIME_UNSPECIFIED = 0,
2359371c9d4SSatish Balay   TS_EXACTFINALTIME_STEPOVER    = 1,
2369371c9d4SSatish Balay   TS_EXACTFINALTIME_INTERPOLATE = 2,
2379371c9d4SSatish Balay   TS_EXACTFINALTIME_MATCHSTEP   = 3
2389371c9d4SSatish Balay } TSExactFinalTimeOption;
239d6ad946cSShri Abhyankar PETSC_EXTERN const char *const TSExactFinalTimeOptions[];
240d6ad946cSShri Abhyankar 
241000e7ae3SMatthew Knepley /* Logging support */
242014dd563SJed Brown PETSC_EXTERN PetscClassId TS_CLASSID;
243d74926cbSBarry Smith PETSC_EXTERN PetscClassId DMTS_CLASSID;
2441b9b13dfSLisandro Dalcin PETSC_EXTERN PetscClassId TSADAPT_CLASSID;
245000e7ae3SMatthew Knepley 
246607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSInitializePackage(void);
247560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSFinalizePackage(void);
2488ba1e511SMatthew Knepley 
249014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm, TS *);
250baa10174SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSClone(TS, TS *);
251014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSDestroy(TS *);
252818ad0c1SBarry Smith 
253014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetProblemType(TS, TSProblemType);
254014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetProblemType(TS, TSProblemType *);
255014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSMonitor(TS, PetscInt, PetscReal, Vec);
256014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSMonitorSet(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, void *), void *, PetscErrorCode (*)(void **));
257721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorSetFromOptions(TS, const char[], const char[], const char[], PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *), PetscErrorCode (*)(TS, PetscViewerAndFormat *));
258014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS);
259818ad0c1SBarry Smith 
260014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS, const char[]);
261014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS, const char[]);
262014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS, const char *[]);
263014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS);
264014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetUp(TS);
265014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSReset(TS);
266818ad0c1SBarry Smith 
267014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetSolution(TS, Vec);
268014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetSolution(TS, Vec *);
269818ad0c1SBarry Smith 
270efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TS2SetSolution(TS, Vec, Vec);
271efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TS2GetSolution(TS, Vec *, Vec *);
27257df6a1bSDebojyoti Ghosh 
273642f3702SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetSolutionComponents(TS, PetscInt *, Vec *);
274642f3702SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetAuxSolution(TS, Vec *);
2750a01e1b2SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetTimeError(TS, PetscInt, Vec *);
27657df6a1bSDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSSetTimeError(TS, Vec);
2776e2e232bSDebojyoti Ghosh 
278a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSSetRHSJacobianP(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), void *);
279cd4cee2dSHong Zhang PETSC_EXTERN PetscErrorCode TSGetRHSJacobianP(TS, Mat *, PetscErrorCode (**)(TS, PetscReal, Vec, Mat, void *), void **);
280a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianP(TS, PetscReal, Vec, Mat);
28133f72c5dSHong Zhang PETSC_EXTERN PetscErrorCode TSSetIJacobianP(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Vec, PetscReal, Mat, void *), void *);
28233f72c5dSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIJacobianP(TS, PetscReal, Vec, Vec, PetscReal, Mat, PetscBool);
283edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS() then TSComputeRHSJacobianP()", ) PetscErrorCode TSComputeDRDPFunction(TS, PetscReal, Vec, Vec *);
284edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS() then TSComputeRHSJacobian()", ) PetscErrorCode TSComputeDRDUFunction(TS, PetscReal, Vec, Vec *);
2859371c9d4SSatish 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 *);
286b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUU(TS, PetscReal, Vec, Vec *, Vec, Vec *);
287b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUP(TS, PetscReal, Vec, Vec *, Vec, Vec *);
288b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPU(TS, PetscReal, Vec, Vec *, Vec, Vec *);
289b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPP(TS, PetscReal, Vec, Vec *, Vec, Vec *);
2909371c9d4SSatish 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 *);
291b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS, PetscReal, Vec, Vec *, Vec, Vec *);
292b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS, PetscReal, Vec, Vec *, Vec, Vec *);
293b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS, PetscReal, Vec, Vec *, Vec, Vec *);
294b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS, PetscReal, Vec, Vec *, Vec, Vec *);
295b5b4017aSHong Zhang PETSC_EXTERN PetscErrorCode TSSetCostHessianProducts(TS, PetscInt, Vec *, Vec *, Vec);
296b5b4017aSHong Zhang PETSC_EXTERN PetscErrorCode TSGetCostHessianProducts(TS, PetscInt *, Vec **, Vec **, Vec *);
297ba0675f6SHong Zhang PETSC_EXTERN PetscErrorCode TSComputeSNESJacobian(TS, Vec, Mat, Mat);
298a05bf03eSHong Zhang 
299bc952696SBarry Smith /*S
300e91eccc2SStefano Zampini      TSTrajectory - Abstract PETSc object that stores the trajectory (solution of ODE/DAE at each time step)
301bc952696SBarry Smith 
302bc952696SBarry Smith    Level: advanced
303bc952696SBarry Smith 
3041cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectorySetType()`, `TSTrajectoryDestroy()`, `TSTrajectoryReset()`
305bc952696SBarry Smith S*/
306bc952696SBarry Smith typedef struct _p_TSTrajectory *TSTrajectory;
307bc952696SBarry Smith 
308bc952696SBarry Smith /*J
309195e9b02SBarry Smith     TSTrajectoryType - String with the name of a PETSc `TS` trajectory storage method
310bc952696SBarry Smith 
311bc952696SBarry Smith    Level: intermediate
312bc952696SBarry Smith 
3131cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectoryDestroy()`
314bc952696SBarry Smith J*/
315bc952696SBarry Smith typedef const char *TSTrajectoryType;
316bc952696SBarry Smith #define TSTRAJECTORYBASIC         "basic"
3171c8c567eSBarry Smith #define TSTRAJECTORYSINGLEFILE    "singlefile"
3189a53571cSHong Zhang #define TSTRAJECTORYMEMORY        "memory"
3192b043167SHong Zhang #define TSTRAJECTORYVISUALIZATION "visualization"
320bc952696SBarry Smith 
321bc952696SBarry Smith PETSC_EXTERN PetscFunctionList TSTrajectoryList;
322bc952696SBarry Smith PETSC_EXTERN PetscClassId      TSTRAJECTORY_CLASSID;
323bc952696SBarry Smith PETSC_EXTERN PetscBool         TSTrajectoryRegisterAllCalled;
324bc952696SBarry Smith 
325bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS);
3262d29f1f2SStefano Zampini PETSC_EXTERN PetscErrorCode TSResetTrajectory(TS);
32767a3cfb0SHong Zhang PETSC_EXTERN PetscErrorCode TSRemoveTrajectory(TS);
328bc952696SBarry Smith 
329bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectoryCreate(MPI_Comm, TSTrajectory *);
3309a992471SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryReset(TSTrajectory);
331bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy(TSTrajectory *);
332560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSTrajectoryView(TSTrajectory, PetscViewer);
333fd9d3c67SJed Brown PETSC_EXTERN PetscErrorCode TSTrajectorySetType(TSTrajectory, TS, TSTrajectoryType);
334881c1a9bSHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryGetType(TSTrajectory, TS, TSTrajectoryType *);
335bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySet(TSTrajectory, TS, PetscInt, PetscReal, Vec);
336c679fc15SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryGet(TSTrajectory, TS, PetscInt, PetscReal *);
337fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetVecs(TSTrajectory, TS, PetscInt, PetscReal *, Vec, Vec);
338fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetUpdatedHistoryVecs(TSTrajectory, TS, PetscReal, Vec *, Vec *);
339fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetNumSteps(TSTrajectory, PetscInt *);
340fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryRestoreUpdatedHistoryVecs(TSTrajectory, Vec *, Vec *);
341972caf09SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory, TS);
342560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSTrajectoryRegister(const char[], PetscErrorCode (*)(TSTrajectory, TS));
343bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectoryRegisterAll(void);
34468bece0bSHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetUp(TSTrajectory, TS);
3459ffb3502SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetUseHistory(TSTrajectory, PetscBool);
3462bee684fSHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetMonitor(TSTrajectory, PetscBool);
34778fbdcc8SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySetVariableNames(TSTrajectory, const char *const *);
34808347785SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySetTransform(TSTrajectory, PetscErrorCode (*)(void *, Vec, Vec *), PetscErrorCode (*)(void *), void *);
349fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectorySetSolutionOnly(TSTrajectory, PetscBool);
350fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetSolutionOnly(TSTrajectory, PetscBool *);
35164fc91eeSBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySetKeepFiles(TSTrajectory, PetscBool);
35264e38db7SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetDirname(TSTrajectory, const char[]);
35364e38db7SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetFiletemplate(TSTrajectory, const char[]);
35478fbdcc8SBarry Smith PETSC_EXTERN PetscErrorCode TSGetTrajectory(TS, TSTrajectory *);
355bc952696SBarry Smith 
356dfb21088SHong Zhang PETSC_EXTERN PetscErrorCode TSSetCostGradients(TS, PetscInt, Vec *, Vec *);
357dfb21088SHong Zhang PETSC_EXTERN PetscErrorCode TSGetCostGradients(TS, PetscInt *, Vec **, Vec **);
358edd03b47SJacob 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 *);
359dfb21088SHong Zhang PETSC_EXTERN PetscErrorCode TSGetCostIntegral(TS, Vec *);
360022c081aSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeCostIntegrand(TS, PetscReal, Vec, Vec);
361cd4cee2dSHong Zhang PETSC_EXTERN PetscErrorCode TSCreateQuadratureTS(TS, PetscBool, TS *);
362cd4cee2dSHong Zhang PETSC_EXTERN PetscErrorCode TSGetQuadratureTS(TS, PetscBool *, TS *);
363d4aa0a58SBarry Smith 
364dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSetFromOptions(TS, PetscOptionItems *);
365a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitor(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *);
366a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitorSet(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *), void *, PetscErrorCode (*)(void **));
367a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitorCancel(TS);
368a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitorSetFromOptions(TS, const char[], const char[], const char[], PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, PetscViewerAndFormat *), PetscErrorCode (*)(TS, PetscViewerAndFormat *));
369a05bf03eSHong Zhang 
370edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSSetRHSJacobianP()", ) PetscErrorCode TSAdjointSetRHSJacobian(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), void *);
371edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSComputeRHSJacobianP()", ) PetscErrorCode TSAdjointComputeRHSJacobian(TS, PetscReal, Vec, Mat);
372edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS()", ) PetscErrorCode TSAdjointComputeDRDPFunction(TS, PetscReal, Vec, Vec *);
373edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS()", ) PetscErrorCode TSAdjointComputeDRDYFunction(TS, PetscReal, Vec, Vec *);
3742c39e106SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS);
37573b18844SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSetSteps(TS, PetscInt);
376d4aa0a58SBarry Smith 
37773b18844SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointStep(TS);
378f6a906c0SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS);
379ece46509SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointReset(TS);
380999a3721SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointCostIntegral(TS);
381ecf68647SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointSetForward(TS, Mat);
382ecf68647SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointResetForward(TS);
383715f1b00SHong Zhang 
38413b1b0ffSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardSetSensitivities(TS, PetscInt, Mat);
38513b1b0ffSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardGetSensitivities(TS, PetscInt *, Mat *);
386edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 12, 0, "TSCreateQuadratureTS()", ) PetscErrorCode TSForwardSetIntegralGradients(TS, PetscInt, Vec *);
387edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 12, 0, "TSForwardGetSensitivities()", ) PetscErrorCode TSForwardGetIntegralGradients(TS, PetscInt *, Vec **);
388715f1b00SHong Zhang PETSC_EXTERN PetscErrorCode TSForwardSetUp(TS);
3897adebddeSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardReset(TS);
390999a3721SHong Zhang PETSC_EXTERN PetscErrorCode TSForwardCostIntegral(TS);
391715f1b00SHong Zhang PETSC_EXTERN PetscErrorCode TSForwardStep(TS);
392b5b4017aSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardSetInitialSensitivities(TS, Mat);
393cc4f23bcSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardGetStages(TS, PetscInt *, Mat *[]);
394c235aa8dSHong Zhang 
395618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetMaxSteps(TS, PetscInt);
396618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetMaxSteps(TS, PetscInt *);
397618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetMaxTime(TS, PetscReal);
398618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetMaxTime(TS, PetscReal *);
39949354f04SShri Abhyankar PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS, TSExactFinalTimeOption);
400f6953c82SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetExactFinalTime(TS, TSExactFinalTimeOption *);
4014a658b32SHong Zhang PETSC_EXTERN PetscErrorCode TSSetTimeSpan(TS, PetscInt, PetscReal *);
4024a658b32SHong Zhang PETSC_EXTERN PetscErrorCode TSGetTimeSpan(TS, PetscInt *, const PetscReal **);
4034a658b32SHong Zhang PETSC_EXTERN PetscErrorCode TSGetTimeSpanSolutions(TS, PetscInt *, Vec **);
404818ad0c1SBarry Smith 
405edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSSetTime()", ) PetscErrorCode TSSetInitialTimeStep(TS, PetscReal, PetscReal);
406edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSSetMax()", ) PetscErrorCode TSSetDuration(TS, PetscInt, PetscReal);
407edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetMax()", ) PetscErrorCode TSGetDuration(TS, PetscInt *, PetscReal *);
408edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetStepNumber()", ) PetscErrorCode TSGetTimeStepNumber(TS, PetscInt *);
409edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetStepNumber()", ) PetscErrorCode TSGetTotalSteps(TS, PetscInt *);
41019eac22cSLisandro Dalcin 
411721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
412cc9c3a59SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorExtreme(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
41383a4ac43SBarry Smith 
41483a4ac43SBarry Smith typedef struct _n_TSMonitorDrawCtx *TSMonitorDrawCtx;
41583a4ac43SBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorDrawCtx *);
41683a4ac43SBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *);
41783a4ac43SBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawSolution(TS, PetscInt, PetscReal, Vec, void *);
4182d5ee99bSBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawSolutionPhase(TS, PetscInt, PetscReal, Vec, void *);
41983a4ac43SBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawError(TS, PetscInt, PetscReal, Vec, void *);
4200ed3bfb6SBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawSolutionFunction(TS, PetscInt, PetscReal, Vec, void *);
42183a4ac43SBarry Smith 
422721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointMonitorDefault(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, PetscViewerAndFormat *);
4239110b2e7SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitorDrawSensi(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *);
4249110b2e7SHong Zhang 
425721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorSolution(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
426014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS, PetscInt, PetscReal, Vec, void *);
427014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void *);
428fb1732b5SBarry Smith 
429014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSStep(TS);
4307cbde773SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSEvaluateWLTE(TS, NormType, PetscInt *, PetscReal *);
431014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS, PetscInt, Vec, PetscBool *);
432cc708dedSBarry Smith PETSC_EXTERN PetscErrorCode TSSolve(TS, Vec);
433e817cc15SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetEquationType(TS, TSEquationType *);
434e817cc15SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSSetEquationType(TS, TSEquationType);
435014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS, TSConvergedReason *);
436d6ad946cSShri Abhyankar PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS, TSConvergedReason);
437487e0bb9SJed Brown PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS, PetscReal *);
4385ef26d82SJed Brown PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS, PetscInt *);
4395ef26d82SJed Brown PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS, PetscInt *);
440cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS, PetscInt *);
441cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS, PetscInt);
442cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS, PetscInt *);
443cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS, PetscInt);
444cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS, PetscBool);
445dcb233daSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRestartStep(TS);
44624655328SShri PETSC_EXTERN PetscErrorCode TSRollBack(TS);
447d2daff3dSHong Zhang 
448cc4f23bcSHong Zhang PETSC_EXTERN PetscErrorCode TSGetStages(TS, PetscInt *, Vec *[]);
449818ad0c1SBarry Smith 
450014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetTime(TS, PetscReal *);
451014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetTime(TS, PetscReal);
452e5e524a1SHong Zhang PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS, PetscReal *);
45380275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS, PetscReal *);
45480275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS, PetscReal);
45580275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetStepNumber(TS, PetscInt *);
45680275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetStepNumber(TS, PetscInt);
457818ad0c1SBarry Smith 
45814d0ab18SJacob Faibussowitsch /*S
45914d0ab18SJacob Faibussowitsch   TSRHSFunction - An `TS` right-hand-side evaluation function
46014d0ab18SJacob Faibussowitsch 
46114d0ab18SJacob Faibussowitsch   Calling Sequence:
46214d0ab18SJacob Faibussowitsch + ts  - timestep context
46314d0ab18SJacob Faibussowitsch . t   - current time
46414d0ab18SJacob Faibussowitsch . u   - input vector
46514d0ab18SJacob Faibussowitsch . F   - function vector
46614d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context
46714d0ab18SJacob Faibussowitsch 
46814d0ab18SJacob Faibussowitsch   Level: beginner
46914d0ab18SJacob Faibussowitsch 
47014d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetRHSFunction()`, `DMTSSetRHSFunction()`, `TSIFunction()`,
47114d0ab18SJacob Faibussowitsch `TSIJacobian()`, `TSRHSJacobian()`
47214d0ab18SJacob Faibussowitsch S*/
47314d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS ts, PetscReal t, Vec u, Vec F, void *ctx);
47414d0ab18SJacob Faibussowitsch 
47514d0ab18SJacob Faibussowitsch /*S
47614d0ab18SJacob Faibussowitsch   TSRHSJacobian - A `TS` right-hand-side Jacobian evaluation function
47714d0ab18SJacob Faibussowitsch 
47814d0ab18SJacob Faibussowitsch   Calling Sequence:
47914d0ab18SJacob Faibussowitsch + ts   - the `TS` context obtained from `TSCreate()`
48014d0ab18SJacob Faibussowitsch . t    - current time
48114d0ab18SJacob Faibussowitsch . u    - input vector
48214d0ab18SJacob Faibussowitsch . Amat - (approximate) Jacobian matrix
48314d0ab18SJacob Faibussowitsch . Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat)
48414d0ab18SJacob Faibussowitsch - ctx  - [optional] user-defined context for matrix evaluation routine
48514d0ab18SJacob Faibussowitsch 
48614d0ab18SJacob Faibussowitsch   Level: beginner
48714d0ab18SJacob Faibussowitsch 
48814d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetRHSJacobian()`, `DMTSSetRHSJacobian()`, `TSRHSFunction()`,
48914d0ab18SJacob Faibussowitsch `TSIFunction()`, `TSIJacobian()`
49014d0ab18SJacob Faibussowitsch S*/
49114d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS ts, PetscReal t, Vec u, Mat Amat, Mat Pmat, void *ctx);
49214d0ab18SJacob Faibussowitsch 
49314d0ab18SJacob Faibussowitsch /*S
49414d0ab18SJacob Faibussowitsch   TSRHSJacobianP - A function that computes the Jacobian of G w.r.t. the parameters P where U_t
49514d0ab18SJacob Faibussowitsch   = G(U,P,t), as well as the location to store the matrix.
49614d0ab18SJacob Faibussowitsch 
49714d0ab18SJacob Faibussowitsch   Calling Sequence:
49814d0ab18SJacob Faibussowitsch + ts  - the `TS` context
49914d0ab18SJacob Faibussowitsch . t   - current timestep
50014d0ab18SJacob Faibussowitsch . U   - input vector (current ODE solution)
50114d0ab18SJacob Faibussowitsch . A   - output matrix
50214d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context
50314d0ab18SJacob Faibussowitsch 
50414d0ab18SJacob Faibussowitsch   Level: beginner
50514d0ab18SJacob Faibussowitsch 
50614d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetRHSJacobianP()`, `TSGetRHSJacobianP()`
50714d0ab18SJacob Faibussowitsch S*/
50814d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobianP)(TS ts, PetscReal t, Vec U, Mat A, void *ctx);
50914d0ab18SJacob Faibussowitsch 
510014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS, Vec, TSRHSFunction, void *);
511014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS, Vec *, TSRHSFunction *, void **);
512014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS, Mat, Mat, TSRHSJacobian, void *);
513014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS, Mat *, Mat *, TSRHSJacobian *, void **);
514e1244c69SJed Brown PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS, PetscBool);
515818ad0c1SBarry Smith 
51614d0ab18SJacob Faibussowitsch /*S
51714d0ab18SJacob Faibussowitsch   TSSolutionFunction - A `TS` solution evaluation function
51814d0ab18SJacob Faibussowitsch 
51914d0ab18SJacob Faibussowitsch   Calling Sequence:
52014d0ab18SJacob Faibussowitsch + ts  - timestep context
52114d0ab18SJacob Faibussowitsch . t   - current time
52214d0ab18SJacob Faibussowitsch . u   - output vector
52314d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context
52414d0ab18SJacob Faibussowitsch 
52514d0ab18SJacob Faibussowitsch   Level: advanced
52614d0ab18SJacob Faibussowitsch 
52714d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetSolutionFunction()`, `DMTSSetSolutionFunction()`
52814d0ab18SJacob Faibussowitsch S*/
52914d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS ts, PetscReal t, Vec u, void *ctx);
53014d0ab18SJacob Faibussowitsch 
531ef20d060SBarry Smith PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS, TSSolutionFunction, void *);
53214d0ab18SJacob Faibussowitsch 
53314d0ab18SJacob Faibussowitsch /*S
53414d0ab18SJacob Faibussowitsch   TSForcingFunction - A `TS` forcing function evaluation function
53514d0ab18SJacob Faibussowitsch 
53614d0ab18SJacob Faibussowitsch   Calling Sequence:
53714d0ab18SJacob Faibussowitsch + ts  - timestep context
53814d0ab18SJacob Faibussowitsch . t   - current time
53914d0ab18SJacob Faibussowitsch . f   - output vector
54014d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context
54114d0ab18SJacob Faibussowitsch 
54214d0ab18SJacob Faibussowitsch   Level: advanced
54314d0ab18SJacob Faibussowitsch 
54414d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetForcingFunction()`, `DMTSSetForcingFunction()`
54514d0ab18SJacob Faibussowitsch S*/
54614d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSForcingFunction)(TS ts, PetscReal t, Vec f, void *ctx);
54714d0ab18SJacob Faibussowitsch 
548d56366bfSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS, TSForcingFunction, void *);
549ef20d060SBarry Smith 
55014d0ab18SJacob Faibussowitsch /*S
55114d0ab18SJacob Faibussowitsch   TSIFunction - A `TS` implicit function evaluation function
55214d0ab18SJacob Faibussowitsch 
55314d0ab18SJacob Faibussowitsch   Calling Sequence:
55414d0ab18SJacob Faibussowitsch + ts  - the `TS` context obtained from `TSCreate()`
55514d0ab18SJacob Faibussowitsch . t   - time at step/stage being solved
55614d0ab18SJacob Faibussowitsch . U   - state vector
55714d0ab18SJacob Faibussowitsch . U_t - time derivative of state vector
55814d0ab18SJacob Faibussowitsch . F   - function vector
55914d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined context for function
56014d0ab18SJacob Faibussowitsch 
56114d0ab18SJacob Faibussowitsch   Level: beginner
56214d0ab18SJacob Faibussowitsch 
56314d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetIFunction()`, `DMTSSetIFunction()`, `TSIJacobian()`, `TSRHSFunction()`, `TSRHSJacobian()`
56414d0ab18SJacob Faibussowitsch S*/
56514d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS ts, PetscReal t, Vec U, Vec U_tt, Vec F, void *ctx);
56614d0ab18SJacob Faibussowitsch 
56714d0ab18SJacob Faibussowitsch /*S
56814d0ab18SJacob Faibussowitsch   TSIJacobian - A `TS` Jacobian evaluation function
56914d0ab18SJacob Faibussowitsch 
57014d0ab18SJacob Faibussowitsch   Calling Sequence:
57114d0ab18SJacob Faibussowitsch + ts   - the `TS` context obtained from `TSCreate()`
57214d0ab18SJacob Faibussowitsch . t    - time at step/stage being solved
57314d0ab18SJacob Faibussowitsch . U    - state vector
57414d0ab18SJacob Faibussowitsch . U_t  - time derivative of state vector
57514d0ab18SJacob Faibussowitsch . a    - shift
57614d0ab18SJacob Faibussowitsch . Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
57714d0ab18SJacob Faibussowitsch . Pmat - matrix used for constructing preconditioner, usually the same as `Amat`
57814d0ab18SJacob Faibussowitsch - ctx  - [optional] user-defined context for Jacobian evaluation routine
57914d0ab18SJacob Faibussowitsch 
58014d0ab18SJacob Faibussowitsch   Level: beginner
58114d0ab18SJacob Faibussowitsch 
58214d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TSSetIJacobian()`, `DMTSSetIJacobian()`, `TSIFunction()`, `TSRHSFunction()`, `TSRHSJacobian()`
58314d0ab18SJacob Faibussowitsch S*/
58414d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS ts, PetscReal t, Vec U, Vec U_t, PetscReal a, Mat Amat, Mat Pmat, void *ctx);
58514d0ab18SJacob Faibussowitsch 
586014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetIFunction(TS, Vec, TSIFunction, void *);
587014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetIFunction(TS, Vec *, TSIFunction *, void **);
588014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS, Mat, Mat, TSIJacobian, void *);
589014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS, Mat *, Mat *, TSIJacobian *, void **);
590316643e7SJed Brown 
59114d0ab18SJacob Faibussowitsch /*S
59214d0ab18SJacob Faibussowitsch   TSI2Function - A `TS` implicit function evaluation function for 2nd order systems
59314d0ab18SJacob Faibussowitsch 
59414d0ab18SJacob Faibussowitsch   Calling Sequence:
59514d0ab18SJacob Faibussowitsch + ts   - the `TS` context obtained from `TSCreate()`
59614d0ab18SJacob Faibussowitsch . t    - time at step/stage being solved
59714d0ab18SJacob Faibussowitsch . U    - state vector
59814d0ab18SJacob Faibussowitsch . U_t  - time derivative of state vector
59914d0ab18SJacob Faibussowitsch . U_tt - second time derivative of state vector
60014d0ab18SJacob Faibussowitsch . F    - function vector
60114d0ab18SJacob Faibussowitsch - ctx  - [optional] user-defined context for matrix evaluation routine (may be `NULL`)
60214d0ab18SJacob Faibussowitsch 
60314d0ab18SJacob Faibussowitsch   Level: advanced
60414d0ab18SJacob Faibussowitsch 
60514d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetI2Function()`, `DMTSSetI2Function()`
60614d0ab18SJacob Faibussowitsch S*/
60714d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Function)(TS ts, PetscReal t, Vec U, Vec U_t, Vec U_tt, Vec F, void *ctx);
60814d0ab18SJacob Faibussowitsch 
60914d0ab18SJacob Faibussowitsch /*S
61014d0ab18SJacob Faibussowitsch   TSI2Jacobian - A `TS` implicit Jacobian evaluation function for 2nd order systems
61114d0ab18SJacob Faibussowitsch 
61214d0ab18SJacob Faibussowitsch   Calling Sequence:
61314d0ab18SJacob Faibussowitsch + ts   - the `TS` context obtained from `TSCreate()`
61414d0ab18SJacob Faibussowitsch . t    - time at step/stage being solved
61514d0ab18SJacob Faibussowitsch . U    - state vector
61614d0ab18SJacob Faibussowitsch . U_t  - time derivative of state vector
61714d0ab18SJacob Faibussowitsch . U_tt - second time derivative of state vector
61814d0ab18SJacob Faibussowitsch . v    - shift for U_t
61914d0ab18SJacob Faibussowitsch . a    - shift for U_tt
62014d0ab18SJacob 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
62114d0ab18SJacob Faibussowitsch . jac  - matrix from which to construct the preconditioner, may be same as `J`
62214d0ab18SJacob Faibussowitsch - ctx  - [optional] user-defined context for matrix evaluation routine
62314d0ab18SJacob Faibussowitsch 
62414d0ab18SJacob Faibussowitsch   Level: advanced
62514d0ab18SJacob Faibussowitsch 
62614d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetI2Jacobian()`, `DMTSSetI2Jacobian()`, `TSIFunction()`, `TSIJacobian()`, `TSRHSFunction()`, `TSRHSJacobian()`
62714d0ab18SJacob Faibussowitsch S*/
62814d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Jacobian)(TS ts, PetscReal t, Vec U, Vec U_t, Vec U_tt, PetscReal v, PetscReal a, Mat J, Mat Jac, void *ctx);
62914d0ab18SJacob Faibussowitsch 
630efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetI2Function(TS, Vec, TSI2Function, void *);
631efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetI2Function(TS, Vec *, TSI2Function *, void **);
632efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetI2Jacobian(TS, Mat, Mat, TSI2Jacobian, void *);
633efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetI2Jacobian(TS, Mat *, Mat *, TSI2Jacobian *, void **);
634efe9872eSLisandro Dalcin 
635545aaa6fSHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitSetIS(TS, const char[], IS);
636545aaa6fSHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetIS(TS, const char[], IS *);
637545aaa6fSHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitSetRHSFunction(TS, const char[], Vec, TSRHSFunction, void *);
6381d06f6b3SHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTS(TS, const char[], TS *);
6391d06f6b3SHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTSs(TS, PetscInt *, TS *[]);
6400fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSSetUseSplitRHSFunction(TS, PetscBool);
6410fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSGetUseSplitRHSFunction(TS, PetscBool *);
642d9194312SHong Zhang 
643014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS, PetscReal, Vec, Vec, void *);
644d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS, PetscReal, Vec, Mat, Mat, void *);
645014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS, PetscReal, Vec, Vec, Vec, void *);
646d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
6471c8a10a0SJed Brown PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS, PetscReal, Vec);
6489b7cd975SBarry Smith PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS, PetscReal, Vec);
649847ff0e1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSComputeIJacobianDefaultColor(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
650d7cfae9bSHong Zhang PETSC_EXTERN PetscErrorCode TSPruneIJacobianColor(TS, Mat, Mat);
651e34be4c2SBarry Smith 
652014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
653b8123daeSJed Brown PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS, PetscReal));
6549be3e283SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS, PetscReal, PetscInt, Vec *));
655c688d042SShri Abhyankar PETSC_EXTERN PetscErrorCode TSSetPostEvaluate(TS, PetscErrorCode (*)(TS));
656014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
6576bd3a4fdSStefano Zampini PETSC_EXTERN PetscErrorCode TSSetResize(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscBool *, void *), PetscErrorCode (*)(TS, PetscInt, Vec[], Vec[], void *), void *);
658014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPreStep(TS);
659b8123daeSJed Brown PETSC_EXTERN PetscErrorCode TSPreStage(TS, PetscReal);
6609be3e283SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSPostStage(TS, PetscReal, PetscInt, Vec *);
661c688d042SShri Abhyankar PETSC_EXTERN PetscErrorCode TSPostEvaluate(TS);
662014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPostStep(TS);
6636bd3a4fdSStefano Zampini PETSC_EXTERN PetscErrorCode TSResize(TS);
6646bd3a4fdSStefano Zampini PETSC_EXTERN PetscErrorCode TSResizeRetrieveVec(TS, const char *, Vec *);
6656bd3a4fdSStefano Zampini PETSC_EXTERN PetscErrorCode TSResizeRegisterVec(TS, const char *, Vec);
6666bd3a4fdSStefano Zampini 
667014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSInterpolate(TS, PetscReal, Vec);
668014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetTolerances(TS, PetscReal, Vec, PetscReal, Vec);
669c5033834SJed Brown PETSC_EXTERN PetscErrorCode TSGetTolerances(TS, PetscReal *, Vec *, PetscReal *, Vec *);
6707453f775SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *);
6718a175baeSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm(TS, Vec, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *);
672014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS, PetscReal);
673014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS, PetscReal *);
674d183316bSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSSetFunctionDomainError(TS, PetscErrorCode (*)(TS, PetscReal, Vec, PetscBool *));
675d183316bSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSFunctionDomainError(TS, PetscReal, Vec, PetscBool *);
676000e7ae3SMatthew Knepley 
677014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS, PetscErrorCode (*)(TS, PetscReal *, void *), void *);
6788d359177SBarry Smith PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS, PetscReal *, void *);
679014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS, PetscReal *);
680014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS, PetscReal);
681014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS, PetscErrorCode (*)(TS, Vec, void *, PetscReal *, PetscBool *), void *);
6828d359177SBarry Smith PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS, Vec, void *, PetscReal *, PetscBool *);
683014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS, Vec, PetscReal *, PetscBool *);
684014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS, PetscReal);
685014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);
68621c89e3eSBarry Smith 
687014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPythonSetType(TS, const char[]);
688ebead697SStefano Zampini PETSC_EXTERN PetscErrorCode TSPythonGetType(TS, const char *[]);
6891d6018f0SLisandro Dalcin 
690014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS, PetscReal, Vec, Vec);
691d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS, PetscReal, Vec, Mat, Mat);
692014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS, PetscReal, Vec, Vec, Vec, PetscBool);
693d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscBool);
694efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSComputeI2Function(TS, PetscReal, Vec, Vec, Vec, Vec);
695efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSComputeI2Jacobian(TS, PetscReal, Vec, Vec, Vec, PetscReal, PetscReal, Mat, Mat);
696f9c1d6abSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS, PetscReal, PetscReal, PetscReal *, PetscReal *);
697818ad0c1SBarry Smith 
698014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS, Vec, Vec);
699d6ebe24aSShri Abhyankar 
700967bb25aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetBoundaryLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *);
70124989b8cSPeter Brune PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM, TSRHSFunction, void *);
70224989b8cSPeter Brune PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM, TSRHSFunction *, void **);
703800f99ffSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionContextDestroy(DM, PetscErrorCode (*)(void *));
70424989b8cSPeter Brune PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM, TSRHSJacobian, void *);
70524989b8cSPeter Brune PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM, TSRHSJacobian *, void **);
706800f99ffSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobianContextDestroy(DM, PetscErrorCode (*)(void *));
70724989b8cSPeter Brune PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM, TSIFunction, void *);
70824989b8cSPeter Brune PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM, TSIFunction *, void **);
709800f99ffSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMTSSetIFunctionContextDestroy(DM, PetscErrorCode (*)(void *));
71024989b8cSPeter Brune PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM, TSIJacobian, void *);
71124989b8cSPeter Brune PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM, TSIJacobian *, void **);
712800f99ffSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMTSSetIJacobianContextDestroy(DM, PetscErrorCode (*)(void *));
713efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode DMTSSetI2Function(DM, TSI2Function, void *);
714efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode DMTSGetI2Function(DM, TSI2Function *, void **);
715800f99ffSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMTSSetI2FunctionContextDestroy(DM, PetscErrorCode (*)(void *));
716efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode DMTSSetI2Jacobian(DM, TSI2Jacobian, void *);
717efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode DMTSGetI2Jacobian(DM, TSI2Jacobian *, void **);
718800f99ffSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMTSSetI2JacobianContextDestroy(DM, PetscErrorCode (*)(void *));
719efe9872eSLisandro Dalcin 
72014d0ab18SJacob Faibussowitsch /*S
72114d0ab18SJacob Faibussowitsch   TSTransientVariable - A function to transform from state to transient variables
72214d0ab18SJacob Faibussowitsch 
72314d0ab18SJacob Faibussowitsch   Calling Sequence:
72414d0ab18SJacob Faibussowitsch + ts  - timestep context
72514d0ab18SJacob Faibussowitsch . p   - input vector (primitive form)
72614d0ab18SJacob Faibussowitsch . c   - output vector, transient variables (conservative form)
72714d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context
72814d0ab18SJacob Faibussowitsch 
72914d0ab18SJacob Faibussowitsch   Level: advanced
73014d0ab18SJacob Faibussowitsch 
73114d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetTransientVariable()`, `DMTSSetTransientVariable()`
73214d0ab18SJacob Faibussowitsch S*/
73314d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSTransientVariable)(TS ts, Vec p, Vec c, void *ctx);
73414d0ab18SJacob Faibussowitsch 
735438f35afSJed Brown PETSC_EXTERN PetscErrorCode TSSetTransientVariable(TS, TSTransientVariable, void *);
736e3c11fc1SJed Brown PETSC_EXTERN PetscErrorCode DMTSSetTransientVariable(DM, TSTransientVariable, void *);
737e3c11fc1SJed Brown PETSC_EXTERN PetscErrorCode DMTSGetTransientVariable(DM, TSTransientVariable *, void *);
738e3c11fc1SJed Brown PETSC_EXTERN PetscErrorCode TSComputeTransientVariable(TS, Vec, Vec);
739e3c11fc1SJed Brown PETSC_EXTERN PetscErrorCode TSHasTransientVariable(TS, PetscBool *);
740e3c11fc1SJed Brown 
741ef20d060SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM, TSSolutionFunction, void *);
742ef20d060SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM, TSSolutionFunction *, void **);
743d56366bfSLisandro Dalcin PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM, TSForcingFunction, void *);
744d56366bfSLisandro Dalcin PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM, TSForcingFunction *, void **);
745e17bd7bbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCheckResidual(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscReal *);
746e17bd7bbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCheckJacobian(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscBool *, PetscReal *);
7477f96f943SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec);
7489b7cd975SBarry Smith 
7490c9409e4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSGetIFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, Vec, void *), void **);
750d433e6cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, Vec, void *), void *);
7510c9409e4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSGetIJacobianLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void **);
752d433e6cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void *);
7530c9409e4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSGetRHSFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, void *), void **);
754d433e6cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *);
755b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrix(DM);
756b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrixLumped(DM);
757b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSDestroyRHSMassMatrix(DM);
7586c6b9e74SPeter Brune 
7596c6b9e74SPeter Brune PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM, PetscErrorCode (*)(void *, PetscViewer), PetscErrorCode (*)(void **, PetscViewer));
7606c6b9e74SPeter Brune PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM, PetscErrorCode (*)(void *, PetscViewer), PetscErrorCode (*)(void **, PetscViewer));
7616c6b9e74SPeter Brune 
76214d0ab18SJacob Faibussowitsch /*S
76314d0ab18SJacob Faibussowitsch   DMDATSRHSFunctionLocal - A local `TS` right hand side residual evaluation function for use with `DMDA`
7646c6b9e74SPeter Brune 
76514d0ab18SJacob Faibussowitsch   Calling Sequence:
76614d0ab18SJacob Faibussowitsch + info - defines the subdomain to evaluate the residual on
76714d0ab18SJacob Faibussowitsch . t    - time at which to evaluate residual
76814d0ab18SJacob Faibussowitsch . x    - array of local state information
76914d0ab18SJacob Faibussowitsch . f    - output array of local residual information
77014d0ab18SJacob Faibussowitsch - ctx  - optional user context
77114d0ab18SJacob Faibussowitsch 
77214d0ab18SJacob Faibussowitsch   Level: beginner
77314d0ab18SJacob Faibussowitsch 
77414d0ab18SJacob Faibussowitsch .seealso: `DMDA`, `DMDATSSetRHSFunctionLocal()`, `TSRHSFunction()`, `DMDATSRHSJacobianLocal()`, `DMDATSIJacobianLocal()`, `DMDATSIFunctionLocal()`
77514d0ab18SJacob Faibussowitsch S*/
77614d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo *info, PetscReal t, void *x, void *f, void *ctx);
77714d0ab18SJacob Faibussowitsch 
77814d0ab18SJacob Faibussowitsch /*S
77914d0ab18SJacob Faibussowitsch   DMDATSRHSJacobianLocal - A local residual evaluation function for use with `DMDA`
78014d0ab18SJacob Faibussowitsch 
78114d0ab18SJacob Faibussowitsch   Calling Sequence:
78214d0ab18SJacob Faibussowitsch + info - defines the subdomain to evaluate the residual on
78314d0ab18SJacob Faibussowitsch . t    - time at which to evaluate residual
78414d0ab18SJacob Faibussowitsch . x    - array of local state information
78514d0ab18SJacob Faibussowitsch . J    - Jacobian matrix
78614d0ab18SJacob Faibussowitsch . B    - matrix from which to construct the preconditioner; often same as `J`
78714d0ab18SJacob Faibussowitsch - ctx  - optional context
78814d0ab18SJacob Faibussowitsch 
78914d0ab18SJacob Faibussowitsch   Level: beginner
79014d0ab18SJacob Faibussowitsch 
79114d0ab18SJacob Faibussowitsch .seealso: `DMDA`, `DMDATSSetRHSJacobianLocal()`, `TSRHSJacobian()`, `DMDATSRHSFunctionLocal()`, `DMDATSIJacobianLocal()`, `DMDATSIFunctionLocal()`
79214d0ab18SJacob Faibussowitsch S*/
79314d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo *info, PetscReal t, void *x, Mat J, Mat B, void *ctx);
79414d0ab18SJacob Faibussowitsch 
79514d0ab18SJacob Faibussowitsch /*S
79614d0ab18SJacob Faibussowitsch   DMDATSIFunctionLocal - A local residual evaluation function for use with `DMDA`
79714d0ab18SJacob Faibussowitsch 
79814d0ab18SJacob Faibussowitsch   Calling Sequence:
79914d0ab18SJacob Faibussowitsch + info  - defines the subdomain to evaluate the residual on
80014d0ab18SJacob Faibussowitsch . t     - time at which to evaluate residual
80114d0ab18SJacob Faibussowitsch . x     - array of local state information
80214d0ab18SJacob Faibussowitsch . xdot  - array of local time derivative information
80314d0ab18SJacob Faibussowitsch . imode - output array of local function evaluation information
80414d0ab18SJacob Faibussowitsch - ctx   - optional context
80514d0ab18SJacob Faibussowitsch 
80614d0ab18SJacob Faibussowitsch   Level: beginner
80714d0ab18SJacob Faibussowitsch 
80814d0ab18SJacob Faibussowitsch .seealso: `DMDA`, `DMDATSSetIFunctionLocal()`, `DMDATSIJacobianLocal()`, `TSIFunction()`
80914d0ab18SJacob Faibussowitsch S*/
81014d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo *info, PetscReal t, void *x, void *xdot, void *imode, void *ctx);
81114d0ab18SJacob Faibussowitsch 
81214d0ab18SJacob Faibussowitsch /*S
81314d0ab18SJacob Faibussowitsch   DMDATSIJacobianLocal - A local residual evaluation function for use with `DMDA`
81414d0ab18SJacob Faibussowitsch 
81514d0ab18SJacob Faibussowitsch   Calling Sequence:
81614d0ab18SJacob Faibussowitsch + info  - defines the subdomain to evaluate the residual on
81714d0ab18SJacob Faibussowitsch . t     - time at which to evaluate the jacobian
81814d0ab18SJacob Faibussowitsch . x     - array of local state information
81914d0ab18SJacob Faibussowitsch . xdot  - time derivative at this state
82014d0ab18SJacob Faibussowitsch . shift - see `TSSetIJacobian()` for the meaning of this parameter
82114d0ab18SJacob Faibussowitsch . J     - Jacobian matrix
82214d0ab18SJacob Faibussowitsch . B     - matrix from which to construct the preconditioner; often same as `J`
82314d0ab18SJacob Faibussowitsch - ctx   - optional context
82414d0ab18SJacob Faibussowitsch 
82514d0ab18SJacob Faibussowitsch   Level: beginner
82614d0ab18SJacob Faibussowitsch 
82714d0ab18SJacob Faibussowitsch .seealso: `DMDA` `DMDATSSetIJacobianLocal()`, `TSIJacobian()`, `DMDATSIFunctionLocal()`, `DMDATSRHSFunctionLocal()`,  `DMDATSRHSJacob ianlocal()`
82814d0ab18SJacob Faibussowitsch S*/
82914d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo *info, PetscReal t, void *x, void *xdot, PetscReal shift, Mat J, Mat B, void *ctx);
83014d0ab18SJacob Faibussowitsch 
83114d0ab18SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM, InsertMode, DMDATSRHSFunctionLocal, void *);
83214d0ab18SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM, DMDATSRHSJacobianLocal, void *);
83314d0ab18SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM, InsertMode, DMDATSIFunctionLocal, void *);
83414d0ab18SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM, DMDATSIJacobianLocal, void *);
8356c6b9e74SPeter Brune 
836be1b0d75SMatthew G. Knepley typedef struct _n_TSMonitorLGCtx *TSMonitorLGCtx;
837d1212d36SBarry Smith typedef struct {
838d1212d36SBarry Smith   Vec            ray;
839d1212d36SBarry Smith   VecScatter     scatter;
840d1212d36SBarry Smith   PetscViewer    viewer;
84151b4a12fSMatthew G. Knepley   TSMonitorLGCtx lgctx;
842d1212d36SBarry Smith } TSMonitorDMDARayCtx;
843d1212d36SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void **);
844d1212d36SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS, PetscInt, PetscReal, Vec, void *);
84551b4a12fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS, PetscInt, PetscReal, Vec, void *);
846d1212d36SBarry Smith 
847bdad233fSMatthew Knepley /* Dynamic creation and loading functions */
848140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList TSList;
84919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode    TSGetType(TS, TSType *);
85019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode    TSSetType(TS, TSType);
851bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode    TSRegister(const char[], PetscErrorCode (*)(TS));
85230de9b25SBarry Smith 
853014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetSNES(TS, SNES *);
854deb2cd25SJed Brown PETSC_EXTERN PetscErrorCode TSSetSNES(TS, SNES);
855014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetKSP(TS, KSP *);
856818ad0c1SBarry Smith 
857014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSView(TS, PetscViewer);
85855849f57SBarry Smith PETSC_EXTERN PetscErrorCode TSLoad(TS, PetscViewer);
859fe2efc57SMark PETSC_EXTERN PetscErrorCode TSViewFromOptions(TS, PetscObject, const char[]);
860fe2efc57SMark PETSC_EXTERN PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory, PetscObject, const char[]);
86155849f57SBarry Smith 
86255849f57SBarry Smith #define TS_FILE_CLASSID 1211225
86321c89e3eSBarry Smith 
864014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS, void *);
865014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS, void *);
86621c89e3eSBarry Smith 
867a9f9c1f6SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtx *);
868a9f9c1f6SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *);
8694f09c107SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS, PetscInt, PetscReal, Vec, void *);
8704f09c107SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS, PetscInt, PetscReal, Vec, void *);
87131152f8aSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS, const char *const *);
872b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS, const char *const **);
873e673d494SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx, const char *const *);
874a66092f1SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS, const char *const *);
875a66092f1SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx, const char *const *);
8767684fa3eSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS, PetscErrorCode (*)(void *, Vec, Vec *), PetscErrorCode (*)(void *), void *);
8777684fa3eSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx, PetscErrorCode (*)(void *, Vec, Vec *), PetscErrorCode (*)(void *), void *);
8784f09c107SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS, PetscInt, PetscReal, Vec, void *);
879201da799SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS, PetscInt, PetscReal, Vec, void *);
880201da799SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS, PetscInt, PetscReal, Vec, void *);
881edbaebb3SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorError(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
882d0c080abSJoseph Pusztay PETSC_EXTERN PetscErrorCode TSDMSwarmMonitorMoments(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
883ef20d060SBarry Smith 
884e669de00SBarry Smith struct _n_TSMonitorLGCtxNetwork {
885e669de00SBarry Smith   PetscInt     nlg;
886e669de00SBarry Smith   PetscDrawLG *lg;
887e669de00SBarry Smith   PetscBool    semilogy;
888e669de00SBarry Smith   PetscInt     howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
889e669de00SBarry Smith };
890e669de00SBarry Smith typedef struct _n_TSMonitorLGCtxNetwork *TSMonitorLGCtxNetwork;
891e669de00SBarry Smith PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkDestroy(TSMonitorLGCtxNetwork *);
892e669de00SBarry Smith PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkCreate(TS, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtxNetwork *);
893e669de00SBarry Smith PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkSolution(TS, PetscInt, PetscReal, Vec, void *);
894e669de00SBarry Smith 
895b3d3934dSBarry Smith typedef struct _n_TSMonitorEnvelopeCtx *TSMonitorEnvelopeCtx;
896b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeCtxCreate(TS, TSMonitorEnvelopeCtx *);
897b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode             TSMonitorEnvelope(TS, PetscInt, PetscReal, Vec, void *);
898b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeGetBounds(TS, Vec *, Vec *);
899b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *);
900b3d3934dSBarry Smith 
9018189c53fSBarry Smith typedef struct _n_TSMonitorSPEigCtx *TSMonitorSPEigCtx;
9028189c53fSBarry Smith PETSC_EXTERN PetscErrorCode          TSMonitorSPEigCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorSPEigCtx *);
9038189c53fSBarry Smith PETSC_EXTERN PetscErrorCode          TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *);
9048189c53fSBarry Smith PETSC_EXTERN PetscErrorCode          TSMonitorSPEig(TS, PetscInt, PetscReal, Vec, void *);
9053914022bSBarry Smith 
9061b575b74SJoseph Pusztay typedef struct _n_TSMonitorSPCtx *TSMonitorSPCtx;
90760e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorSPCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscBool, PetscBool, TSMonitorSPCtx *);
9081b575b74SJoseph Pusztay PETSC_EXTERN PetscErrorCode       TSMonitorSPCtxDestroy(TSMonitorSPCtx *);
9090ec8ee2bSJoseph Pusztay PETSC_EXTERN PetscErrorCode       TSMonitorSPSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
9101b575b74SJoseph Pusztay 
91160e16b1bSMatthew G. Knepley typedef struct _n_TSMonitorHGCtx *TSMonitorHGCtx;
91260e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorHGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscInt, PetscBool, TSMonitorHGCtx *);
91360e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
91460e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorHGCtxDestroy(TSMonitorHGCtx *);
91560e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
91660e16b1bSMatthew G. Knepley 
9172589fa24SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetEventHandler(TS, PetscInt, PetscInt[], PetscBool[], PetscErrorCode (*)(TS, PetscReal, Vec, PetscScalar[], void *), PetscErrorCode (*)(TS, PetscInt, PetscInt[], PetscReal, Vec, PetscBool, void *), void *);
9180e37a75fSKarl Rupp PETSC_EXTERN PetscErrorCode TSSetPostEventIntervalStep(TS, PetscReal);
9196427ac75SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS, PetscReal, PetscReal[]);
9201ea83e56SMiguel PETSC_EXTERN PetscErrorCode TSGetNumEvents(TS, PetscInt *);
9211ea83e56SMiguel 
92276bdecfbSBarry Smith /*J
923195e9b02SBarry Smith    TSSSPType - string with the name of a `TSSSP` scheme.
924815f1ad5SJed Brown 
925815f1ad5SJed Brown    Level: beginner
926815f1ad5SJed Brown 
9271cc06b55SBarry Smith .seealso: [](ch_ts), `TSSSPSetType()`, `TS`, `TSSSP`
92876bdecfbSBarry Smith J*/
92919fd82e9SBarry Smith typedef const char *TSSSPType;
930815f1ad5SJed Brown #define TSSSPRKS2  "rks2"
931815f1ad5SJed Brown #define TSSSPRKS3  "rks3"
932815f1ad5SJed Brown #define TSSSPRK104 "rk104"
933815f1ad5SJed Brown 
93419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode    TSSSPSetType(TS, TSSSPType);
93519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode    TSSSPGetType(TS, TSSSPType *);
936014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSSPSetNumStages(TS, PetscInt);
937014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSSPGetNumStages(TS, PetscInt *);
938787849ffSJed Brown PETSC_EXTERN PetscErrorCode    TSSSPInitializePackage(void);
939560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode    TSSSPFinalizePackage(void);
940013797aaSBarry Smith PETSC_EXTERN PetscFunctionList TSSSPList;
941815f1ad5SJed Brown 
942ed657a08SJed Brown /*S
94384df9cb4SJed Brown    TSAdapt - Abstract object that manages time-step adaptivity
94484df9cb4SJed Brown 
94584df9cb4SJed Brown    Level: beginner
94684df9cb4SJed Brown 
9471cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSAdaptCreate()`, `TSAdaptType`
94884df9cb4SJed Brown S*/
94984df9cb4SJed Brown typedef struct _p_TSAdapt *TSAdapt;
95084df9cb4SJed Brown 
951f0fc11ceSJed Brown /*J
95287497f52SBarry Smith     TSAdaptType - String with the name of `TSAdapt` scheme.
95384df9cb4SJed Brown 
95484df9cb4SJed Brown    Level: beginner
95584df9cb4SJed Brown 
9561cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdaptSetType()`, `TS`, `TSAdapt`
957f0fc11ceSJed Brown J*/
958f8963224SJed Brown typedef const char *TSAdaptType;
959b0f836d7SJed Brown #define TSADAPTNONE    "none"
9601917a363SLisandro Dalcin #define TSADAPTBASIC   "basic"
961cb7ab186SLisandro Dalcin #define TSADAPTDSP     "dsp"
9628d59e960SJed Brown #define TSADAPTCFL     "cfl"
9631917a363SLisandro Dalcin #define TSADAPTGLEE    "glee"
964d949e4a4SStefano Zampini #define TSADAPTHISTORY "history"
96584df9cb4SJed Brown 
966552698daSJed Brown PETSC_EXTERN PetscErrorCode TSGetAdapt(TS, TSAdapt *);
967bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[], PetscErrorCode (*)(TSAdapt));
968607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void);
969014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
970014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm, TSAdapt *);
97119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt, TSAdaptType);
972d0288e4fSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetType(TSAdapt, TSAdaptType *);
973014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt, const char[]);
974014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
975014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt, const char[], PetscInt, PetscInt, PetscReal, PetscReal, PetscBool);
976014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt, PetscInt *, const PetscInt **, const PetscInt **, const PetscReal **, const PetscReal **);
977014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt, TS, PetscReal, PetscInt *, PetscReal *, PetscBool *);
978b295832fSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt, TS, PetscReal, Vec, PetscBool *);
979014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt, PetscViewer);
980ad6bc421SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt, PetscViewer);
981dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(TSAdapt, PetscOptionItems *);
982099af0a3SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt);
983014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt *);
984014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt, PetscBool);
985bf997491SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt, PetscBool);
9861917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptSetSafety(TSAdapt, PetscReal, PetscReal);
9871917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetSafety(TSAdapt, PetscReal *, PetscReal *);
98876cddca1SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt, PetscReal);
98976cddca1SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt, PetscReal *);
9901917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptSetClip(TSAdapt, PetscReal, PetscReal);
9911917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetClip(TSAdapt, PetscReal *, PetscReal *);
99262c23b28SMark PETSC_EXTERN PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt, PetscReal);
99362c23b28SMark PETSC_EXTERN PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt, PetscReal *);
994014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt, PetscReal, PetscReal);
9951917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetStepLimits(TSAdapt, PetscReal *, PetscReal *);
996b295832fSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt, PetscErrorCode (*)(TSAdapt, TS, PetscReal, Vec, PetscBool *));
997d949e4a4SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptHistorySetHistory(TSAdapt, PetscInt n, PetscReal hist[], PetscBool);
99875017583SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptHistorySetTrajectory(TSAdapt, TSTrajectory, PetscBool);
99975017583SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptHistoryGetStep(TSAdapt, PetscInt, PetscReal *, PetscReal *);
1000de50f1caSBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt, PetscInt);
1001cb7ab186SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptDSPSetFilter(TSAdapt, const char *);
1002cb7ab186SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptDSPSetPID(TSAdapt, PetscReal, PetscReal, PetscReal);
100384df9cb4SJed Brown 
100484df9cb4SJed Brown /*S
100587497f52SBarry Smith    TSGLLEAdapt - Abstract object that manages time-step adaptivity for `TSGLLE`
1006ed657a08SJed Brown 
1007ed657a08SJed Brown    Level: beginner
1008ed657a08SJed Brown 
1009195e9b02SBarry Smith    Developer Note:
101087497f52SBarry Smith    This functionality should be replaced by the `TSAdapt`.
101184df9cb4SJed Brown 
10121cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLEAdaptCreate()`, `TSGLLEAdaptType`
1013ed657a08SJed Brown S*/
101426d28e4eSEmil Constantinescu typedef struct _p_TSGLLEAdapt *TSGLLEAdapt;
1015ed657a08SJed Brown 
101676bdecfbSBarry Smith /*J
101787497f52SBarry Smith     TSGLLEAdaptType - String with the name of `TSGLLEAdapt` scheme
1018ed657a08SJed Brown 
1019ed657a08SJed Brown    Level: beginner
1020ed657a08SJed Brown 
1021195e9b02SBarry Smith    Developer Note:
1022195e9b02SBarry Smith    This functionality should be replaced by the `TSAdaptType`.
1023195e9b02SBarry Smith 
10241cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLLEAdaptSetType()`, `TS`
102576bdecfbSBarry Smith J*/
102626d28e4eSEmil Constantinescu typedef const char *TSGLLEAdaptType;
102726d28e4eSEmil Constantinescu #define TSGLLEADAPT_NONE "none"
102826d28e4eSEmil Constantinescu #define TSGLLEADAPT_SIZE "size"
102926d28e4eSEmil Constantinescu #define TSGLLEADAPT_BOTH "both"
1030ed657a08SJed Brown 
103126d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptRegister(const char[], PetscErrorCode (*)(TSGLLEAdapt));
103226d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptInitializePackage(void);
103326d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptFinalizePackage(void);
103426d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate(MPI_Comm, TSGLLEAdapt *);
103526d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt, TSGLLEAdaptType);
103626d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt, const char[]);
103726d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptChoose(TSGLLEAdapt, PetscInt, const PetscInt[], const PetscReal[], const PetscReal[], PetscInt, PetscReal, PetscReal, PetscInt *, PetscReal *, PetscBool *);
103826d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt, PetscViewer);
1039dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetFromOptions(TSGLLEAdapt, PetscOptionItems *);
104026d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt *);
1041ed657a08SJed Brown 
104276bdecfbSBarry Smith /*J
104387497f52SBarry Smith     TSGLLEAcceptType - String with the name of `TSGLLEAccept` scheme
1044ed657a08SJed Brown 
1045ed657a08SJed Brown    Level: beginner
1046ed657a08SJed Brown 
10471cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLLESetAcceptType()`, `TS`, `TSGLLEAccept`
104876bdecfbSBarry Smith J*/
104926d28e4eSEmil Constantinescu typedef const char *TSGLLEAcceptType;
105026d28e4eSEmil Constantinescu #define TSGLLEACCEPT_ALWAYS "always"
1051ed657a08SJed Brown 
105226d28e4eSEmil Constantinescu PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLLEAcceptFunction)(TS, PetscReal, PetscReal, const PetscReal[], PetscBool *);
105326d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAcceptRegister(const char[], TSGLLEAcceptFunction);
1054ed657a08SJed Brown 
105576bdecfbSBarry Smith /*J
1056195e9b02SBarry Smith   TSGLLEType - string with the name of a General Linear `TSGLLE` type
105718b56cb9SJed Brown 
105818b56cb9SJed Brown   Level: beginner
105918b56cb9SJed Brown 
10601cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLESetType()`, `TSGLLERegister()`, `TSGLLEAccept`
106176bdecfbSBarry Smith J*/
106226d28e4eSEmil Constantinescu typedef const char *TSGLLEType;
106326d28e4eSEmil Constantinescu #define TSGLLE_IRKS "irks"
106418b56cb9SJed Brown 
106526d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLERegister(const char[], PetscErrorCode (*)(TS));
106626d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEInitializePackage(void);
106726d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEFinalizePackage(void);
106826d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLESetType(TS, TSGLLEType);
106926d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEGetAdapt(TS, TSGLLEAdapt *);
107026d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLESetAcceptType(TS, TSGLLEAcceptType);
107118b56cb9SJed Brown 
1072020d8f30SJed Brown /*J
1073195e9b02SBarry Smith     TSEIMEXType - String with the name of an Extrapolated IMEX `TSEIMEX` type
10747d5697caSHong Zhang 
10757d5697caSHong Zhang    Level: beginner
10767d5697caSHong Zhang 
10771cc06b55SBarry Smith .seealso: [](ch_ts), `TSEIMEXSetType()`, `TS`, `TSEIMEX`, `TSEIMEXRegister()`
10787d5697caSHong Zhang J*/
10797d5697caSHong Zhang #define TSEIMEXType char *
10807d5697caSHong Zhang 
10817d5697caSHong Zhang PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts, PetscInt);
10827d5697caSHong Zhang PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts, PetscInt, PetscInt);
10837d5697caSHong Zhang PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS, PetscBool);
10847d5697caSHong Zhang 
10857d5697caSHong Zhang /*J
1086195e9b02SBarry Smith     TSRKType - String with the name of a Runge-Kutta `TSRK` type
1087f68a32c8SEmil Constantinescu 
1088f68a32c8SEmil Constantinescu    Level: beginner
1089f68a32c8SEmil Constantinescu 
10901cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSRKSetType()`, `TS`, `TSRK`, `TSRKRegister()`
1091f68a32c8SEmil Constantinescu J*/
1092f68a32c8SEmil Constantinescu typedef const char *TSRKType;
1093f68a32c8SEmil Constantinescu #define TSRK1FE "1fe"
1094fdefd5e5SDebojyoti Ghosh #define TSRK2A  "2a"
1095e7685601SHong Zhang #define TSRK2B  "2b"
1096f68a32c8SEmil Constantinescu #define TSRK3   "3"
1097fdefd5e5SDebojyoti Ghosh #define TSRK3BS "3bs"
1098f68a32c8SEmil Constantinescu #define TSRK4   "4"
1099fdefd5e5SDebojyoti Ghosh #define TSRK5F  "5f"
1100fdefd5e5SDebojyoti Ghosh #define TSRK5DP "5dp"
110105e23783SLisandro Dalcin #define TSRK5BS "5bs"
110237acaa02SHendrik Ranocha #define TSRK6VR "6vr"
110337acaa02SHendrik Ranocha #define TSRK7VR "7vr"
110437acaa02SHendrik Ranocha #define TSRK8VR "8vr"
110505e23783SLisandro Dalcin 
11062ea87ba9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRKGetOrder(TS, PetscInt *);
11070f7a28cdSStefano Zampini PETSC_EXTERN PetscErrorCode TSRKGetType(TS, TSRKType *);
11080f7a28cdSStefano Zampini PETSC_EXTERN PetscErrorCode TSRKSetType(TS, TSRKType);
11090f7a28cdSStefano Zampini PETSC_EXTERN PetscErrorCode TSRKGetTableau(TS, PetscInt *, const PetscReal **, const PetscReal **, const PetscReal **, const PetscReal **, PetscInt *, const PetscReal **, PetscBool *);
11100fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSRKSetMultirate(TS, PetscBool);
11110fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSRKGetMultirate(TS, PetscBool *);
1112f68a32c8SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
1113f68a32c8SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void);
1114560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void);
1115f68a32c8SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void);
1116f68a32c8SEmil Constantinescu 
1117f68a32c8SEmil Constantinescu /*J
1118195e9b02SBarry Smith    TSMPRKType - String with the name of a Partitioned Runge-Kutta `TSMPRK` type
11199f537349Sluzhanghpp 
11209f537349Sluzhanghpp    Level: beginner
11219f537349Sluzhanghpp 
11221cc06b55SBarry Smith .seealso: [](ch_ts), `TSMPRKSetType()`, `TS`, `TSMPRK`, `TSMPRKRegister()`
11239f537349Sluzhanghpp J*/
11244b84eec9SHong Zhang typedef const char *TSMPRKType;
112519c2959aSHong Zhang #define TSMPRK2A22 "2a22"
112619c2959aSHong Zhang #define TSMPRK2A23 "2a23"
112719c2959aSHong Zhang #define TSMPRK2A32 "2a32"
112819c2959aSHong Zhang #define TSMPRK2A33 "2a33"
11294b84eec9SHong Zhang #define TSMPRKP2   "p2"
11304b84eec9SHong Zhang #define TSMPRKP3   "p3"
11319f537349Sluzhanghpp 
11324b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKGetType(TS ts, TSMPRKType *);
11334b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKSetType(TS ts, TSMPRKType);
113479bc8a77SHong 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[]);
11354b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKInitializePackage(void);
11364b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKFinalizePackage(void);
11374b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKRegisterDestroy(void);
11389f537349Sluzhanghpp 
11399f537349Sluzhanghpp /*J
1140195e9b02SBarry Smith     TSIRKType - String with the name of an implicit Runge-Kutta `TSIRK` type
1141d2567f34SHong Zhang 
1142d2567f34SHong Zhang    Level: beginner
1143d2567f34SHong Zhang 
11441cc06b55SBarry Smith .seealso: [](ch_ts), `TSIRKSetType()`, `TS`, `TSIRK`, `TSIRKRegister()`
1145d2567f34SHong Zhang J*/
1146d2567f34SHong Zhang typedef const char *TSIRKType;
1147d2567f34SHong Zhang #define TSIRKGAUSS "gauss"
1148d2567f34SHong Zhang 
1149d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKGetType(TS, TSIRKType *);
1150d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKSetType(TS, TSIRKType);
1151d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKGetNumStages(TS, PetscInt *);
1152d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKSetNumStages(TS, PetscInt);
1153d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKRegister(const char[], PetscErrorCode (*function)(TS));
1154d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKTableauCreate(TS, PetscInt, const PetscReal *, const PetscReal *, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, const PetscScalar *);
1155d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKInitializePackage(void);
1156d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKFinalizePackage(void);
1157d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKRegisterDestroy(void);
1158d2567f34SHong Zhang 
1159d2567f34SHong Zhang /*J
1160195e9b02SBarry Smith     TSGLEEType - String with the name of a General Linear with Error Estimation `TSGLEE` type
1161b6a60446SDebojyoti Ghosh 
1162b6a60446SDebojyoti Ghosh    Level: beginner
1163b6a60446SDebojyoti Ghosh 
11641cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEESetType()`, `TS`, `TSGLEE`, `TSGLEERegister()`
1165b6a60446SDebojyoti Ghosh J*/
1166b6a60446SDebojyoti Ghosh typedef const char *TSGLEEType;
1167ab8c5912SEmil Constantinescu #define TSGLEEi1      "BE1"
1168b6a60446SDebojyoti Ghosh #define TSGLEE23      "23"
1169b6a60446SDebojyoti Ghosh #define TSGLEE24      "24"
1170b6a60446SDebojyoti Ghosh #define TSGLEE25I     "25i"
1171b6a60446SDebojyoti Ghosh #define TSGLEE35      "35"
1172b6a60446SDebojyoti Ghosh #define TSGLEEEXRK2A  "exrk2a"
1173b6a60446SDebojyoti Ghosh #define TSGLEERK32G1  "rk32g1"
1174b6a60446SDebojyoti Ghosh #define TSGLEERK285EX "rk285ex"
117516a05f60SBarry Smith 
1176b6a60446SDebojyoti Ghosh /*J
1177195e9b02SBarry Smith     TSGLEEMode - String with the mode of error estimation for a General Linear with Error Estimation `TSGLEE` type
1178b6a60446SDebojyoti Ghosh 
1179b6a60446SDebojyoti Ghosh    Level: beginner
1180b6a60446SDebojyoti Ghosh 
11811cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEESetMode()`, `TS`, `TSGLEE`, `TSGLEERegister()`
1182b6a60446SDebojyoti Ghosh J*/
1183b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEEGetType(TS ts, TSGLEEType *);
1184b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEESetType(TS ts, TSGLEEType);
118557df6a1bSDebojyoti 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[]);
1186b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEEFinalizePackage(void);
1187b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEEInitializePackage(void);
1188b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEERegisterDestroy(void);
1189b6a60446SDebojyoti Ghosh 
1190b6a60446SDebojyoti Ghosh /*J
1191195e9b02SBarry Smith     TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX `TSARKIMEX` type
1192020d8f30SJed Brown 
1193020d8f30SJed Brown    Level: beginner
1194020d8f30SJed Brown 
11951cc06b55SBarry Smith .seealso: [](ch_ts), `TSARKIMEXSetType()`, `TS`, `TSARKIMEX`, `TSARKIMEXRegister()`
1196020d8f30SJed Brown J*/
119719fd82e9SBarry Smith typedef const char *TSARKIMEXType;
1198e817cc15SEmil Constantinescu #define TSARKIMEX1BEE   "1bee"
11991f80e275SEmil Constantinescu #define TSARKIMEXA2     "a2"
12001f80e275SEmil Constantinescu #define TSARKIMEXL2     "l2"
12011f80e275SEmil Constantinescu #define TSARKIMEXARS122 "ars122"
12021f80e275SEmil Constantinescu #define TSARKIMEX2C     "2c"
12038a381b04SJed Brown #define TSARKIMEX2D     "2d"
1204a3a57f36SJed Brown #define TSARKIMEX2E     "2e"
12056cf0794eSJed Brown #define TSARKIMEXPRSSP2 "prssp2"
12068a381b04SJed Brown #define TSARKIMEX3      "3"
12076cf0794eSJed Brown #define TSARKIMEXBPR3   "bpr3"
12086cf0794eSJed Brown #define TSARKIMEXARS443 "ars443"
12098a381b04SJed Brown #define TSARKIMEX4      "4"
12108a381b04SJed Brown #define TSARKIMEX5      "5"
121119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts, TSARKIMEXType *);
121219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts, TSARKIMEXType);
1213014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS, PetscBool);
12143a28c0e5SStefano Zampini PETSC_EXTERN PetscErrorCode TSARKIMEXGetFullyImplicit(TS, PetscBool *);
121519fd82e9SBarry 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[]);
1216607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void);
1217560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
1218014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);
12198a381b04SJed Brown 
1220020d8f30SJed Brown /*J
1221*d27334e2SStefano Zampini     TSDIRKType - String with the name of a Diagonally Implicit Runge-Kutta `TSDIRK` type
1222*d27334e2SStefano Zampini 
1223*d27334e2SStefano Zampini    Level: beginner
1224*d27334e2SStefano Zampini 
1225*d27334e2SStefano Zampini .seealso: [](ch_ts), `TSDIRKSetType()`, `TS`, `TSDIRK`, `TSDIRKRegister()`
1226*d27334e2SStefano Zampini J*/
1227*d27334e2SStefano Zampini typedef const char *TSDIRKType;
1228*d27334e2SStefano Zampini #define TSDIRKS212      "s212"
1229*d27334e2SStefano Zampini #define TSDIRKES212     "es212"
1230*d27334e2SStefano Zampini #define TSDIRKES213SAL  "es213sal"
1231*d27334e2SStefano Zampini #define TSDIRKES324SAL  "es324sal"
1232*d27334e2SStefano Zampini #define TSDIRKES325SAL  "es325sal"
1233*d27334e2SStefano Zampini #define TSDIRK657A      "657a"
1234*d27334e2SStefano Zampini #define TSDIRKES648SA   "es648sa"
1235*d27334e2SStefano Zampini #define TSDIRK658A      "658a"
1236*d27334e2SStefano Zampini #define TSDIRKS659A     "s659a"
1237*d27334e2SStefano Zampini #define TSDIRK7510SAL   "7510sal"
1238*d27334e2SStefano Zampini #define TSDIRKES7510SA  "es7510sa"
1239*d27334e2SStefano Zampini #define TSDIRK759A      "759a"
1240*d27334e2SStefano Zampini #define TSDIRKS7511SAL  "s7511sal"
1241*d27334e2SStefano Zampini #define TSDIRK8614A     "8614a"
1242*d27334e2SStefano Zampini #define TSDIRK8616SAL   "8616sal"
1243*d27334e2SStefano Zampini #define TSDIRKES8516SAL "es8516sal"
1244*d27334e2SStefano Zampini PETSC_EXTERN PetscErrorCode TSDIRKGetType(TS ts, TSDIRKType *);
1245*d27334e2SStefano Zampini PETSC_EXTERN PetscErrorCode TSDIRKSetType(TS ts, TSDIRKType);
1246*d27334e2SStefano Zampini PETSC_EXTERN PetscErrorCode TSDIRKRegister(TSDIRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
1247*d27334e2SStefano Zampini 
1248*d27334e2SStefano Zampini /*J
1249195e9b02SBarry Smith     TSRosWType - String with the name of a Rosenbrock-W `TSROSW` type
1250020d8f30SJed Brown 
1251020d8f30SJed Brown    Level: beginner
1252020d8f30SJed Brown 
12531cc06b55SBarry Smith .seealso: [](ch_ts), `TSRosWSetType()`, `TS`, `TSROSW`, `TSRosWRegister()`
1254020d8f30SJed Brown J*/
125519fd82e9SBarry Smith typedef const char *TSRosWType;
125661692a83SJed Brown #define TSROSW2M          "2m"
125761692a83SJed Brown #define TSROSW2P          "2p"
1258fe7e6d57SJed Brown #define TSROSWRA3PW       "ra3pw"
1259fe7e6d57SJed Brown #define TSROSWRA34PW2     "ra34pw2"
1260ef3c5b88SJed Brown #define TSROSWRODAS3      "rodas3"
1261ef3c5b88SJed Brown #define TSROSWSANDU3      "sandu3"
1262961f28d0SJed Brown #define TSROSWASSP3P3S1C  "assp3p3s1c"
1263961f28d0SJed Brown #define TSROSWLASSP3P4S2C "lassp3p4s2c"
126443b21953SEmil Constantinescu #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
1265753f8adbSEmil Constantinescu #define TSROSWARK3        "ark3"
12663606a31eSEmil Constantinescu #define TSROSWTHETA1      "theta1"
12673606a31eSEmil Constantinescu #define TSROSWTHETA2      "theta2"
126842faf41dSJed Brown #define TSROSWGRK4T       "grk4t"
126942faf41dSJed Brown #define TSROSWSHAMP4      "shamp4"
127042faf41dSJed Brown #define TSROSWVELDD4      "veldd4"
127142faf41dSJed Brown #define TSROSW4L          "4l"
1272b1c69cc3SEmil Constantinescu 
12736bd7aeb5SHong Zhang PETSC_EXTERN PetscErrorCode TSRosWGetType(TS, TSRosWType *);
12746bd7aeb5SHong Zhang PETSC_EXTERN PetscErrorCode TSRosWSetType(TS, TSRosWType);
1275014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS, PetscBool);
127619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
127719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal);
1278607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void);
1279560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
1280014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);
1281e27a552bSJed Brown 
1282211a84d6SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSBDFSetOrder(TS, PetscInt);
1283211a84d6SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSBDFGetOrder(TS, PetscInt *);
1284211a84d6SLisandro Dalcin 
12856bd7aeb5SHong Zhang /*J
1286195e9b02SBarry Smith   TSBasicSymplecticType - String with the name of a basic symplectic integration `TSBASICSYMPLECTIC` type
12876bd7aeb5SHong Zhang 
12886bd7aeb5SHong Zhang   Level: beginner
12896bd7aeb5SHong Zhang 
12901cc06b55SBarry Smith .seealso: [](ch_ts), `TSBasicSymplecticSetType()`, `TS`, `TSBASICSYMPLECTIC`, `TSBasicSymplecticRegister()`
12916bd7aeb5SHong Zhang J*/
1292e49d4f37SHong Zhang typedef const char *TSBasicSymplecticType;
1293e49d4f37SHong Zhang #define TSBASICSYMPLECTICSIEULER   "1"
1294e49d4f37SHong Zhang #define TSBASICSYMPLECTICVELVERLET "2"
1295e49d4f37SHong Zhang #define TSBASICSYMPLECTIC3         "3"
1296e49d4f37SHong Zhang #define TSBASICSYMPLECTIC4         "4"
1297e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticSetType(TS, TSBasicSymplecticType);
1298e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticGetType(TS, TSBasicSymplecticType *);
1299e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegister(TSBasicSymplecticType, PetscInt, PetscInt, PetscReal[], PetscReal[]);
1300e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticInitializePackage(void);
1301e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticFinalizePackage(void);
1302e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterDestroy(void);
13036bd7aeb5SHong Zhang 
1304db4653c2SDaniel Finn /*J
1305195e9b02SBarry Smith   TSDISCGRAD - The Discrete Gradient integrator is a timestepper for Hamiltonian systems designed to conserve the first integral (energy),
130687497f52SBarry Smith   but also has the property for some systems of monotonicity in a functional.
1307db4653c2SDaniel Finn 
1308db4653c2SDaniel Finn   Level: beginner
1309db4653c2SDaniel Finn 
13101cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, TSDiscGradSetFormulation()`, `TSDiscGradGetFormulation()`
1311db4653c2SDaniel Finn J*/
131240be0ff1SMatthew 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 *);
1313db4653c2SDaniel Finn PETSC_EXTERN PetscErrorCode TSDiscGradIsGonzalez(TS, PetscBool *);
1314db4653c2SDaniel Finn PETSC_EXTERN PetscErrorCode TSDiscGradUseGonzalez(TS, PetscBool);
131540be0ff1SMatthew G. Knepley 
131683e2fdc7SBarry Smith /*
13170f3b3ca1SHong Zhang        PETSc interface to Sundials
131883e2fdc7SBarry Smith */
1319e808b789SPatrick Sanan #ifdef PETSC_HAVE_SUNDIALS2
13209371c9d4SSatish Balay typedef enum {
13219371c9d4SSatish Balay   SUNDIALS_ADAMS = 1,
13229371c9d4SSatish Balay   SUNDIALS_BDF   = 2
13239371c9d4SSatish Balay } TSSundialsLmmType;
13246a6fc655SJed Brown PETSC_EXTERN const char *const TSSundialsLmmTypes[];
13259371c9d4SSatish Balay typedef enum {
13269371c9d4SSatish Balay   SUNDIALS_MODIFIED_GS  = 1,
13279371c9d4SSatish Balay   SUNDIALS_CLASSICAL_GS = 2
13289371c9d4SSatish Balay } TSSundialsGramSchmidtType;
13296a6fc655SJed Brown PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
1330014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSundialsSetType(TS, TSSundialsLmmType);
1331014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSundialsGetPC(TS, PC *);
1332014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSundialsSetTolerance(TS, PetscReal, PetscReal);
1333014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSundialsSetMinTimeStep(TS, PetscReal);
1334014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSundialsSetMaxTimeStep(TS, PetscReal);
1335014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSundialsGetIterations(TS, PetscInt *, PetscInt *);
1336014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSundialsSetGramSchmidtType(TS, TSSundialsGramSchmidtType);
1337014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSundialsSetGMRESRestart(TS, PetscInt);
1338014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSundialsSetLinearTolerance(TS, PetscReal);
1339014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSundialsMonitorInternalSteps(TS, PetscBool);
1340014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSundialsSetMaxl(TS, PetscInt);
1341c4e80e11SFlorian PETSC_EXTERN PetscErrorCode    TSSundialsSetMaxord(TS, PetscInt);
1342918687b8SPatrick Sanan PETSC_EXTERN PetscErrorCode    TSSundialsSetUseDense(TS, PetscBool);
13436dc17ff5SMatthew Knepley #endif
134483e2fdc7SBarry Smith 
1345014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS, PetscReal);
1346014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS, PetscReal *);
1347014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS, PetscBool *);
1348014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS, PetscBool);
13490de4c49aSJed Brown 
1350014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS, PetscReal);
1351014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS, PetscReal, PetscReal, PetscReal);
1352014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS, PetscReal *, PetscReal *, PetscReal *);
135388df8a41SLisandro Dalcin 
1354818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAlpha2SetRadius(TS, PetscReal);
1355818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAlpha2SetParams(TS, PetscReal, PetscReal, PetscReal, PetscReal);
1356818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAlpha2GetParams(TS, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
1357818efac9SLisandro Dalcin 
1358014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetDM(TS, DM);
1359014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetDM(TS, DM *);
13606c699258SBarry Smith 
1361014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES, Vec, Vec, void *);
1362d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES, Vec, Mat, Mat, void *);
13630f5c6efeSJed Brown 
1364f3b1f45cSBarry Smith PETSC_EXTERN PetscErrorCode TSRHSJacobianTest(TS, PetscBool *);
1365f3b1f45cSBarry Smith PETSC_EXTERN PetscErrorCode TSRHSJacobianTestTranspose(TS, PetscBool *);
1366aad739acSMatthew G. Knepley 
13672e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSGetComputeInitialCondition(TS, PetscErrorCode (**)(TS, Vec));
13682e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSSetComputeInitialCondition(TS, PetscErrorCode (*)(TS, Vec));
13692e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSComputeInitialCondition(TS, Vec);
13702e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSGetComputeExactError(TS, PetscErrorCode (**)(TS, Vec, Vec));
1371aad739acSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSSetComputeExactError(TS, PetscErrorCode (*)(TS, Vec, Vec));
1372aad739acSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSComputeExactError(TS, Vec, Vec);
1373f2ed2dc7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscConvEstUseTS(PetscConvEst, PetscBool);
1374d60b7d5cSBarry Smith 
1375d60b7d5cSBarry Smith PETSC_EXTERN PetscErrorCode TSSetMatStructure(TS, MatStructure);
1376818ad0c1SBarry Smith #endif
1377