xref: /petsc/include/petscts.h (revision 60e16b1b373f0d7eb06432ae904f4dc6ff882c7e)
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 
11ac09b921SBarry Smith /* SUBMANSEC = TS */
12ac09b921SBarry Smith 
13435da068SBarry Smith /*S
14435da068SBarry Smith      TS - Abstract PETSc object that manages all time-steppers (ODE integrators)
15435da068SBarry Smith 
16435da068SBarry Smith    Level: beginner
17435da068SBarry Smith 
18db781477SPatrick Sanan .seealso: `TSCreate()`, `TSSetType()`, `TSType`, `SNES`, `KSP`, `PC`, `TSDestroy()`
19435da068SBarry Smith S*/
20f09e8eb9SSatish Balay typedef struct _p_TS *TS;
21435da068SBarry Smith 
2276bdecfbSBarry Smith /*J
2387497f52SBarry Smith     TSType - String with the name of a PETSc `TS` method.
24435da068SBarry Smith 
25435da068SBarry Smith    Level: beginner
26435da068SBarry Smith 
27db781477SPatrick Sanan .seealso: `TSSetType()`, `TS`, `TSRegister()`
2876bdecfbSBarry Smith J*/
2919fd82e9SBarry Smith typedef const char *TSType;
309596e0b4SJed Brown #define TSEULER           "euler"
319596e0b4SJed Brown #define TSBEULER          "beuler"
32e49d4f37SHong Zhang #define TSBASICSYMPLECTIC "basicsymplectic"
339596e0b4SJed Brown #define TSPSEUDO          "pseudo"
344d91e141SJed Brown #define TSCN              "cn"
359596e0b4SJed Brown #define TSSUNDIALS        "sundials"
364d91e141SJed Brown #define TSRK              "rk"
379596e0b4SJed Brown #define TSPYTHON          "python"
389596e0b4SJed Brown #define TSTHETA           "theta"
3988df8a41SLisandro Dalcin #define TSALPHA           "alpha"
40818efac9SLisandro Dalcin #define TSALPHA2          "alpha2"
4126d28e4eSEmil Constantinescu #define TSGLLE            "glle"
42b6a60446SDebojyoti Ghosh #define TSGLEE            "glee"
43ef7bb5aaSJed Brown #define TSSSP             "ssp"
448a381b04SJed Brown #define TSARKIMEX         "arkimex"
45e27a552bSJed Brown #define TSROSW            "rosw"
467d5697caSHong Zhang #define TSEIMEX           "eimex"
47abadfa0fSMatthew G. Knepley #define TSMIMEX           "mimex"
48211a84d6SLisandro Dalcin #define TSBDF             "bdf"
49d249a78fSBarry Smith #define TSRADAU5          "radau5"
504b84eec9SHong Zhang #define TSMPRK            "mprk"
5140be0ff1SMatthew G. Knepley #define TSDISCGRAD        "discgrad"
52d2567f34SHong Zhang #define TSIRK             "irk"
5340be0ff1SMatthew G. Knepley 
54435da068SBarry Smith /*E
5587497f52SBarry Smith     TSProblemType - Determines the type of problem this `TS` object is to be used to solve
56435da068SBarry Smith 
57435da068SBarry Smith    Level: beginner
58435da068SBarry Smith 
59db781477SPatrick Sanan .seealso: `TSCreate()`
60435da068SBarry Smith E*/
619371c9d4SSatish Balay typedef enum {
629371c9d4SSatish Balay   TS_LINEAR,
639371c9d4SSatish Balay   TS_NONLINEAR
649371c9d4SSatish Balay } TSProblemType;
65818ad0c1SBarry Smith 
66b93fea0eSJed Brown /*E
6787497f52SBarry Smith    TSEquationType - type of `TS` problem that is solved
68e817cc15SEmil Constantinescu 
69e817cc15SEmil Constantinescu    Level: beginner
70e817cc15SEmil Constantinescu 
7195452b02SPatrick Sanan    Developer Notes:
7295452b02SPatrick Sanan     this must match petsc/finclude/petscts.h
73e817cc15SEmil Constantinescu 
74e817cc15SEmil Constantinescu    Supported types are:
7587497f52SBarry Smith      `TS_EQ_UNSPECIFIED` (default)
7687497f52SBarry 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
7787497f52SBarry Smith      `TS_EQ_IMPLICIT` {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) = 0
78e817cc15SEmil Constantinescu 
79db781477SPatrick Sanan .seealso: `TSGetEquationType()`, `TSSetEquationType()`
80e817cc15SEmil Constantinescu E*/
81e817cc15SEmil Constantinescu typedef enum {
82e817cc15SEmil Constantinescu   TS_EQ_UNSPECIFIED               = -1,
83e817cc15SEmil Constantinescu   TS_EQ_EXPLICIT                  = 0,
84e817cc15SEmil Constantinescu   TS_EQ_ODE_EXPLICIT              = 1,
85e817cc15SEmil Constantinescu   TS_EQ_DAE_SEMI_EXPLICIT_INDEX1  = 100,
86e817cc15SEmil Constantinescu   TS_EQ_DAE_SEMI_EXPLICIT_INDEX2  = 200,
87e817cc15SEmil Constantinescu   TS_EQ_DAE_SEMI_EXPLICIT_INDEX3  = 300,
88e817cc15SEmil Constantinescu   TS_EQ_DAE_SEMI_EXPLICIT_INDEXHI = 500,
89e817cc15SEmil Constantinescu   TS_EQ_IMPLICIT                  = 1000,
90e817cc15SEmil Constantinescu   TS_EQ_ODE_IMPLICIT              = 1001,
91e817cc15SEmil Constantinescu   TS_EQ_DAE_IMPLICIT_INDEX1       = 1100,
92e817cc15SEmil Constantinescu   TS_EQ_DAE_IMPLICIT_INDEX2       = 1200,
93e817cc15SEmil Constantinescu   TS_EQ_DAE_IMPLICIT_INDEX3       = 1300,
9419436ca2SJed Brown   TS_EQ_DAE_IMPLICIT_INDEXHI      = 1500
95e817cc15SEmil Constantinescu } TSEquationType;
96e817cc15SEmil Constantinescu PETSC_EXTERN const char *const *TSEquationTypes;
97e817cc15SEmil Constantinescu 
98e817cc15SEmil Constantinescu /*E
9987497f52SBarry Smith    TSConvergedReason - reason a `TS` method has converged or not
100b93fea0eSJed Brown 
101b93fea0eSJed Brown    Level: beginner
102b93fea0eSJed Brown 
10395452b02SPatrick Sanan    Developer Notes:
10495452b02SPatrick Sanan     this must match petsc/finclude/petscts.h
105b93fea0eSJed Brown 
106b93fea0eSJed Brown    Each reason has its own manual page.
107b93fea0eSJed Brown 
108db781477SPatrick Sanan .seealso: `TSGetConvergedReason()`
109b93fea0eSJed Brown E*/
110193ac0bcSJed Brown typedef enum {
111193ac0bcSJed Brown   TS_CONVERGED_ITERATING          = 0,
112193ac0bcSJed Brown   TS_CONVERGED_TIME               = 1,
113193ac0bcSJed Brown   TS_CONVERGED_ITS                = 2,
114d6ad946cSShri Abhyankar   TS_CONVERGED_USER               = 3,
1152b7db910SShri Abhyankar   TS_CONVERGED_EVENT              = 4,
1163118ae5eSBarry Smith   TS_CONVERGED_PSEUDO_FATOL       = 5,
1173118ae5eSBarry Smith   TS_CONVERGED_PSEUDO_FRTOL       = 6,
118193ac0bcSJed Brown   TS_DIVERGED_NONLINEAR_SOLVE     = -1,
119b5b4017aSHong Zhang   TS_DIVERGED_STEP_REJECTED       = -2,
12005c9950eSHong Zhang   TSFORWARD_DIVERGED_LINEAR_SOLVE = -3,
12105c9950eSHong Zhang   TSADJOINT_DIVERGED_LINEAR_SOLVE = -4
122193ac0bcSJed Brown } TSConvergedReason;
123014dd563SJed Brown PETSC_EXTERN const char *const *TSConvergedReasons;
1241957e957SBarry Smith 
125b93fea0eSJed Brown /*MC
12687497f52SBarry Smith    TS_CONVERGED_ITERATING - this only occurs if `TSGetConvergedReason()` is called during the `TSSolve()`
127b93fea0eSJed Brown 
128b93fea0eSJed Brown    Level: beginner
129b93fea0eSJed Brown 
130db781477SPatrick Sanan .seealso: `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`
131b93fea0eSJed Brown M*/
132b93fea0eSJed Brown 
133b93fea0eSJed Brown /*MC
134b93fea0eSJed Brown    TS_CONVERGED_TIME - the final time was reached
135b93fea0eSJed Brown 
136b93fea0eSJed Brown    Level: beginner
137b93fea0eSJed Brown 
138db781477SPatrick Sanan .seealso: `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxTime()`, `TSGetMaxTime()`, `TSGetSolveTime()`
139b93fea0eSJed Brown M*/
140b93fea0eSJed Brown 
141b93fea0eSJed Brown /*MC
1421957e957SBarry Smith    TS_CONVERGED_ITS - the maximum number of iterations (time-steps) was reached prior to the final time
143b93fea0eSJed Brown 
144b93fea0eSJed Brown    Level: beginner
145b93fea0eSJed Brown 
146db781477SPatrick Sanan .seealso: `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxSteps()`, `TSGetMaxSteps()`
147b93fea0eSJed Brown M*/
1481957e957SBarry Smith 
149d6ad946cSShri Abhyankar /*MC
150d6ad946cSShri Abhyankar    TS_CONVERGED_USER - user requested termination
151d6ad946cSShri Abhyankar 
152d6ad946cSShri Abhyankar    Level: beginner
153d6ad946cSShri Abhyankar 
154db781477SPatrick Sanan .seealso: `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`
155b93fea0eSJed Brown M*/
156b93fea0eSJed Brown 
157b93fea0eSJed Brown /*MC
158d76fc68cSShri Abhyankar    TS_CONVERGED_EVENT - user requested termination on event detection
159d76fc68cSShri Abhyankar 
160d76fc68cSShri Abhyankar    Level: beginner
161d76fc68cSShri Abhyankar 
162db781477SPatrick Sanan .seealso: `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`
163d76fc68cSShri Abhyankar M*/
164d76fc68cSShri Abhyankar 
165d76fc68cSShri Abhyankar /*MC
16687497f52SBarry Smith    TS_CONVERGED_PSEUDO_FRTOL - stops when function norm decreased by a set amount, used only for `TSPSEUDO`
1673118ae5eSBarry Smith 
1683118ae5eSBarry Smith    Level: beginner
1693118ae5eSBarry Smith 
1703c7db156SBarry Smith    Options Database Key:
171d5b43468SJose E. Roman .   -ts_pseudo_frtol <rtol> - use specified rtol
1723118ae5eSBarry Smith 
173db781477SPatrick Sanan .seealso: `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`, `TS_CONVERGED_PSEUDO_FATOL`
1743118ae5eSBarry Smith M*/
1753118ae5eSBarry Smith 
1763118ae5eSBarry Smith /*MC
17787497f52SBarry Smith    TS_CONVERGED_PSEUDO_FATOL - stops when function norm decreases below a set amount, used only for `TSPSEUDO`
1783118ae5eSBarry Smith 
1793118ae5eSBarry Smith    Level: beginner
1803118ae5eSBarry Smith 
1813c7db156SBarry Smith    Options Database Key:
18267b8a455SSatish Balay .   -ts_pseudo_fatol <atol> - use specified atol
1833118ae5eSBarry Smith 
184db781477SPatrick Sanan .seealso: `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`, `TS_CONVERGED_PSEUDO_FRTOL`
1853118ae5eSBarry Smith M*/
1863118ae5eSBarry Smith 
1873118ae5eSBarry Smith /*MC
188b93fea0eSJed Brown    TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed
189b93fea0eSJed Brown 
190b93fea0eSJed Brown    Level: beginner
191b93fea0eSJed Brown 
19295452b02SPatrick Sanan    Notes:
19395452b02SPatrick Sanan     See TSSetMaxSNESFailures() for how to allow more nonlinear solver failures.
1941957e957SBarry Smith 
195db781477SPatrick Sanan .seealso: `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSGetSNES()`, `SNESGetConvergedReason()`, `TSSetMaxSNESFailures()`
196b93fea0eSJed Brown M*/
197b93fea0eSJed Brown 
198b93fea0eSJed Brown /*MC
199b93fea0eSJed Brown    TS_DIVERGED_STEP_REJECTED - too many steps were rejected
200b93fea0eSJed Brown 
201b93fea0eSJed Brown    Level: beginner
202b93fea0eSJed Brown 
20395452b02SPatrick Sanan    Notes:
20495452b02SPatrick Sanan     See TSSetMaxStepRejections() for how to allow more step rejections.
2051957e957SBarry Smith 
206db781477SPatrick Sanan .seealso: `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxStepRejections()`
207b93fea0eSJed Brown M*/
208b93fea0eSJed Brown 
209d6ad946cSShri Abhyankar /*E
210d6ad946cSShri Abhyankar    TSExactFinalTimeOption - option for handling of final time step
211d6ad946cSShri Abhyankar 
212d6ad946cSShri Abhyankar    Level: beginner
213d6ad946cSShri Abhyankar 
21495452b02SPatrick Sanan    Developer Notes:
21595452b02SPatrick Sanan     this must match petsc/finclude/petscts.h
216d6ad946cSShri Abhyankar 
21787497f52SBarry Smith $  `TS_EXACTFINALTIME_STEPOVER`    - Don't do anything if final time is exceeded
21887497f52SBarry Smith $  `TS_EXACTFINALTIME_INTERPOLATE` - Interpolate back to final time
21987497f52SBarry Smith $  `TS_EXACTFINALTIME_MATCHSTEP` - Adapt final time step to match the final time
220feed9e9dSBarry Smith 
221db781477SPatrick Sanan .seealso: `TSGetConvergedReason()`, `TSSetExactFinalTime()`, `TSGetExactFinalTime()`
2224ceb6c04SBarry Smith 
223d6ad946cSShri Abhyankar E*/
2249371c9d4SSatish Balay typedef enum {
2259371c9d4SSatish Balay   TS_EXACTFINALTIME_UNSPECIFIED = 0,
2269371c9d4SSatish Balay   TS_EXACTFINALTIME_STEPOVER    = 1,
2279371c9d4SSatish Balay   TS_EXACTFINALTIME_INTERPOLATE = 2,
2289371c9d4SSatish Balay   TS_EXACTFINALTIME_MATCHSTEP   = 3
2299371c9d4SSatish Balay } TSExactFinalTimeOption;
230d6ad946cSShri Abhyankar PETSC_EXTERN const char *const TSExactFinalTimeOptions[];
231d6ad946cSShri Abhyankar 
232000e7ae3SMatthew Knepley /* Logging support */
233014dd563SJed Brown PETSC_EXTERN PetscClassId TS_CLASSID;
234d74926cbSBarry Smith PETSC_EXTERN PetscClassId DMTS_CLASSID;
2351b9b13dfSLisandro Dalcin PETSC_EXTERN PetscClassId TSADAPT_CLASSID;
236000e7ae3SMatthew Knepley 
237607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSInitializePackage(void);
238560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSFinalizePackage(void);
2398ba1e511SMatthew Knepley 
240014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm, TS *);
241baa10174SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSClone(TS, TS *);
242014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSDestroy(TS *);
243818ad0c1SBarry Smith 
244014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetProblemType(TS, TSProblemType);
245014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetProblemType(TS, TSProblemType *);
246014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSMonitor(TS, PetscInt, PetscReal, Vec);
247014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSMonitorSet(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, void *), void *, PetscErrorCode (*)(void **));
248721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorSetFromOptions(TS, const char[], const char[], const char[], PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *), PetscErrorCode (*)(TS, PetscViewerAndFormat *));
249014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS);
250818ad0c1SBarry Smith 
251014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS, const char[]);
252014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS, const char[]);
253014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS, const char *[]);
254014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS);
255014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetUp(TS);
256014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSReset(TS);
257818ad0c1SBarry Smith 
258014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetSolution(TS, Vec);
259014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetSolution(TS, Vec *);
260818ad0c1SBarry Smith 
261efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TS2SetSolution(TS, Vec, Vec);
262efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TS2GetSolution(TS, Vec *, Vec *);
26357df6a1bSDebojyoti Ghosh 
264642f3702SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetSolutionComponents(TS, PetscInt *, Vec *);
265642f3702SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetAuxSolution(TS, Vec *);
2660a01e1b2SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetTimeError(TS, PetscInt, Vec *);
26757df6a1bSDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSSetTimeError(TS, Vec);
2686e2e232bSDebojyoti Ghosh 
269a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSSetRHSJacobianP(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), void *);
270cd4cee2dSHong Zhang PETSC_EXTERN PetscErrorCode TSGetRHSJacobianP(TS, Mat *, PetscErrorCode (**)(TS, PetscReal, Vec, Mat, void *), void **);
271a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianP(TS, PetscReal, Vec, Mat);
27233f72c5dSHong Zhang PETSC_EXTERN PetscErrorCode TSSetIJacobianP(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Vec, PetscReal, Mat, void *), void *);
27333f72c5dSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIJacobianP(TS, PetscReal, Vec, Vec, PetscReal, Mat, PetscBool);
27425ef9dfeSBarry Smith PETSC_EXTERN                PETSC_DEPRECATED_FUNCTION("Use TSGetQuadratureTS then TSComputeRHSJacobianP") PetscErrorCode TSComputeDRDPFunction(TS, PetscReal, Vec, Vec *);
27525ef9dfeSBarry Smith PETSC_EXTERN                PETSC_DEPRECATED_FUNCTION("Use TSGetQuadratureTS then TSComputeRHSJacobian") PetscErrorCode TSComputeDRDUFunction(TS, PetscReal, Vec, Vec *);
2769371c9d4SSatish 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 *);
277b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUU(TS, PetscReal, Vec, Vec *, Vec, Vec *);
278b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUP(TS, PetscReal, Vec, Vec *, Vec, Vec *);
279b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPU(TS, PetscReal, Vec, Vec *, Vec, Vec *);
280b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPP(TS, PetscReal, Vec, Vec *, Vec, Vec *);
2819371c9d4SSatish 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 *);
282b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS, PetscReal, Vec, Vec *, Vec, Vec *);
283b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS, PetscReal, Vec, Vec *, Vec, Vec *);
284b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS, PetscReal, Vec, Vec *, Vec, Vec *);
285b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS, PetscReal, Vec, Vec *, Vec, Vec *);
286b5b4017aSHong Zhang PETSC_EXTERN PetscErrorCode TSSetCostHessianProducts(TS, PetscInt, Vec *, Vec *, Vec);
287b5b4017aSHong Zhang PETSC_EXTERN PetscErrorCode TSGetCostHessianProducts(TS, PetscInt *, Vec **, Vec **, Vec *);
288ba0675f6SHong Zhang PETSC_EXTERN PetscErrorCode TSComputeSNESJacobian(TS, Vec, Mat, Mat);
289a05bf03eSHong Zhang 
290bc952696SBarry Smith /*S
291e91eccc2SStefano Zampini      TSTrajectory - Abstract PETSc object that stores the trajectory (solution of ODE/DAE at each time step)
292bc952696SBarry Smith 
293bc952696SBarry Smith    Level: advanced
294bc952696SBarry Smith 
295db781477SPatrick Sanan .seealso: `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectorySetType()`, `TSTrajectoryDestroy()`, `TSTrajectoryReset()`
296bc952696SBarry Smith S*/
297bc952696SBarry Smith typedef struct _p_TSTrajectory *TSTrajectory;
298bc952696SBarry Smith 
299bc952696SBarry Smith /*J
30087497f52SBarry Smith     TSTrajectorySetType - String with the name of a PETSc `TS` trajectory storage method
301bc952696SBarry Smith 
302bc952696SBarry Smith    Level: intermediate
303bc952696SBarry Smith 
304db781477SPatrick Sanan .seealso: `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectoryDestroy()`
305bc952696SBarry Smith J*/
306bc952696SBarry Smith typedef const char *TSTrajectoryType;
307bc952696SBarry Smith #define TSTRAJECTORYBASIC         "basic"
3081c8c567eSBarry Smith #define TSTRAJECTORYSINGLEFILE    "singlefile"
3099a53571cSHong Zhang #define TSTRAJECTORYMEMORY        "memory"
3102b043167SHong Zhang #define TSTRAJECTORYVISUALIZATION "visualization"
311bc952696SBarry Smith 
312bc952696SBarry Smith PETSC_EXTERN PetscFunctionList TSTrajectoryList;
313bc952696SBarry Smith PETSC_EXTERN PetscClassId      TSTRAJECTORY_CLASSID;
314bc952696SBarry Smith PETSC_EXTERN PetscBool         TSTrajectoryRegisterAllCalled;
315bc952696SBarry Smith 
316bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS);
3172d29f1f2SStefano Zampini PETSC_EXTERN PetscErrorCode TSResetTrajectory(TS);
31867a3cfb0SHong Zhang PETSC_EXTERN PetscErrorCode TSRemoveTrajectory(TS);
319bc952696SBarry Smith 
320bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectoryCreate(MPI_Comm, TSTrajectory *);
3219a992471SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryReset(TSTrajectory);
322bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy(TSTrajectory *);
323560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSTrajectoryView(TSTrajectory, PetscViewer);
324fd9d3c67SJed Brown PETSC_EXTERN PetscErrorCode TSTrajectorySetType(TSTrajectory, TS, TSTrajectoryType);
325881c1a9bSHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryGetType(TSTrajectory, TS, TSTrajectoryType *);
326bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySet(TSTrajectory, TS, PetscInt, PetscReal, Vec);
327c679fc15SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryGet(TSTrajectory, TS, PetscInt, PetscReal *);
328fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetVecs(TSTrajectory, TS, PetscInt, PetscReal *, Vec, Vec);
329fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetUpdatedHistoryVecs(TSTrajectory, TS, PetscReal, Vec *, Vec *);
330fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetNumSteps(TSTrajectory, PetscInt *);
331fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryRestoreUpdatedHistoryVecs(TSTrajectory, Vec *, Vec *);
332972caf09SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory, TS);
333560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSTrajectoryRegister(const char[], PetscErrorCode (*)(TSTrajectory, TS));
334bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectoryRegisterAll(void);
33568bece0bSHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetUp(TSTrajectory, TS);
3369ffb3502SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetUseHistory(TSTrajectory, PetscBool);
3372bee684fSHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetMonitor(TSTrajectory, PetscBool);
33878fbdcc8SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySetVariableNames(TSTrajectory, const char *const *);
33908347785SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySetTransform(TSTrajectory, PetscErrorCode (*)(void *, Vec, Vec *), PetscErrorCode (*)(void *), void *);
340fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectorySetSolutionOnly(TSTrajectory, PetscBool);
341fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetSolutionOnly(TSTrajectory, PetscBool *);
34264fc91eeSBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySetKeepFiles(TSTrajectory, PetscBool);
34364e38db7SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetDirname(TSTrajectory, const char[]);
34464e38db7SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetFiletemplate(TSTrajectory, const char[]);
34578fbdcc8SBarry Smith PETSC_EXTERN PetscErrorCode TSGetTrajectory(TS, TSTrajectory *);
346bc952696SBarry Smith 
347dfb21088SHong Zhang PETSC_EXTERN PetscErrorCode TSSetCostGradients(TS, PetscInt, Vec *, Vec *);
348dfb21088SHong Zhang PETSC_EXTERN PetscErrorCode TSGetCostGradients(TS, PetscInt *, Vec **, Vec **);
34925ef9dfeSBarry Smith PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSCreateQuadratureTS() then set up the sub-TS (since version 3.12)") 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 *);
350dfb21088SHong Zhang PETSC_EXTERN PetscErrorCode TSGetCostIntegral(TS, Vec *);
351022c081aSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeCostIntegrand(TS, PetscReal, Vec, Vec);
352cd4cee2dSHong Zhang PETSC_EXTERN PetscErrorCode TSCreateQuadratureTS(TS, PetscBool, TS *);
353cd4cee2dSHong Zhang PETSC_EXTERN PetscErrorCode TSGetQuadratureTS(TS, PetscBool *, TS *);
354d4aa0a58SBarry Smith 
355dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSetFromOptions(TS, PetscOptionItems *);
356a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitor(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *);
357a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitorSet(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *), void *, PetscErrorCode (*)(void **));
358a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitorCancel(TS);
359a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitorSetFromOptions(TS, const char[], const char[], const char[], PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, PetscViewerAndFormat *), PetscErrorCode (*)(TS, PetscViewerAndFormat *));
360a05bf03eSHong Zhang 
36125ef9dfeSBarry Smith PETSC_EXTERN                PETSC_DEPRECATED_FUNCTION("Use TSSetRHSJacobianP()") PetscErrorCode TSAdjointSetRHSJacobian(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), void *);
36225ef9dfeSBarry Smith PETSC_EXTERN                PETSC_DEPRECATED_FUNCTION("Use TSComputeRHSJacobianP()") PetscErrorCode TSAdjointComputeRHSJacobian(TS, PetscReal, Vec, Mat);
36325ef9dfeSBarry Smith PETSC_EXTERN                PETSC_DEPRECATED_FUNCTION("Use TSGetQuadratureTS then TSComputeRHSJacobianP") PetscErrorCode TSAdjointComputeDRDPFunction(TS, PetscReal, Vec, Vec *);
36425ef9dfeSBarry Smith PETSC_EXTERN                PETSC_DEPRECATED_FUNCTION("Use TSGetQuadratureTS then TSComputeRHSJacobian") PetscErrorCode TSAdjointComputeDRDYFunction(TS, PetscReal, Vec, Vec *);
3652c39e106SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS);
36673b18844SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSetSteps(TS, PetscInt);
367d4aa0a58SBarry Smith 
36873b18844SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointStep(TS);
369f6a906c0SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS);
370ece46509SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointReset(TS);
371999a3721SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointCostIntegral(TS);
372ecf68647SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointSetForward(TS, Mat);
373ecf68647SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointResetForward(TS);
374715f1b00SHong Zhang 
37513b1b0ffSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardSetSensitivities(TS, PetscInt, Mat);
37613b1b0ffSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardGetSensitivities(TS, PetscInt *, Mat *);
37725ef9dfeSBarry Smith PETSC_EXTERN                PETSC_DEPRECATED_FUNCTION("Use TSCreateQuadratureTS() and TSForwardSetSensitivities() (since version 3.12)") PetscErrorCode TSForwardSetIntegralGradients(TS, PetscInt, Vec *);
37825ef9dfeSBarry Smith PETSC_EXTERN                PETSC_DEPRECATED_FUNCTION("Use TSForwardGetSensitivities()") PetscErrorCode TSForwardGetIntegralGradients(TS, PetscInt *, Vec **);
379715f1b00SHong Zhang PETSC_EXTERN PetscErrorCode TSForwardSetUp(TS);
3807adebddeSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardReset(TS);
381999a3721SHong Zhang PETSC_EXTERN PetscErrorCode TSForwardCostIntegral(TS);
382715f1b00SHong Zhang PETSC_EXTERN PetscErrorCode TSForwardStep(TS);
383b5b4017aSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardSetInitialSensitivities(TS, Mat);
384cc4f23bcSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardGetStages(TS, PetscInt *, Mat *[]);
385c235aa8dSHong Zhang 
386618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetMaxSteps(TS, PetscInt);
387618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetMaxSteps(TS, PetscInt *);
388618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetMaxTime(TS, PetscReal);
389618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetMaxTime(TS, PetscReal *);
39049354f04SShri Abhyankar PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS, TSExactFinalTimeOption);
391f6953c82SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetExactFinalTime(TS, TSExactFinalTimeOption *);
3924a658b32SHong Zhang PETSC_EXTERN PetscErrorCode TSSetTimeSpan(TS, PetscInt, PetscReal *);
3934a658b32SHong Zhang PETSC_EXTERN PetscErrorCode TSGetTimeSpan(TS, PetscInt *, const PetscReal **);
3944a658b32SHong Zhang PETSC_EXTERN PetscErrorCode TSGetTimeSpanSolutions(TS, PetscInt *, Vec **);
395818ad0c1SBarry Smith 
39625ef9dfeSBarry Smith PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSSetTime[Step]() (since version 3.8)") PetscErrorCode TSSetInitialTimeStep(TS, PetscReal, PetscReal);
39725ef9dfeSBarry Smith PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSSetMax{Steps|Time}() (since version 3.8)") PetscErrorCode TSSetDuration(TS, PetscInt, PetscReal);
39825ef9dfeSBarry Smith PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetMax{Steps|Time}() (since version 3.8)") PetscErrorCode TSGetDuration(TS, PetscInt *, PetscReal *);
39925ef9dfeSBarry Smith PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetStepNumber() (since version 3.8)") PetscErrorCode TSGetTimeStepNumber(TS, PetscInt *);
40025ef9dfeSBarry Smith PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetStepNumber() (since version 3.8)") PetscErrorCode TSGetTotalSteps(TS, PetscInt *);
40119eac22cSLisandro Dalcin 
402721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
403cc9c3a59SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorExtreme(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
40483a4ac43SBarry Smith 
40583a4ac43SBarry Smith typedef struct _n_TSMonitorDrawCtx *TSMonitorDrawCtx;
40683a4ac43SBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorDrawCtx *);
40783a4ac43SBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *);
40883a4ac43SBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawSolution(TS, PetscInt, PetscReal, Vec, void *);
4092d5ee99bSBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawSolutionPhase(TS, PetscInt, PetscReal, Vec, void *);
41083a4ac43SBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawError(TS, PetscInt, PetscReal, Vec, void *);
4110ed3bfb6SBarry Smith PETSC_EXTERN PetscErrorCode         TSMonitorDrawSolutionFunction(TS, PetscInt, PetscReal, Vec, void *);
41283a4ac43SBarry Smith 
413721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointMonitorDefault(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, PetscViewerAndFormat *);
4149110b2e7SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitorDrawSensi(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *);
4159110b2e7SHong Zhang 
416721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorSolution(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
417014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS, PetscInt, PetscReal, Vec, void *);
418014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void *);
419fb1732b5SBarry Smith 
420014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSStep(TS);
4217cbde773SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSEvaluateWLTE(TS, NormType, PetscInt *, PetscReal *);
422014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS, PetscInt, Vec, PetscBool *);
423cc708dedSBarry Smith PETSC_EXTERN PetscErrorCode TSSolve(TS, Vec);
424e817cc15SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetEquationType(TS, TSEquationType *);
425e817cc15SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSSetEquationType(TS, TSEquationType);
426014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS, TSConvergedReason *);
427d6ad946cSShri Abhyankar PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS, TSConvergedReason);
428487e0bb9SJed Brown PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS, PetscReal *);
4295ef26d82SJed Brown PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS, PetscInt *);
4305ef26d82SJed Brown PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS, PetscInt *);
431cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS, PetscInt *);
432cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS, PetscInt);
433cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS, PetscInt *);
434cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS, PetscInt);
435cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS, PetscBool);
436dcb233daSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRestartStep(TS);
43724655328SShri PETSC_EXTERN PetscErrorCode TSRollBack(TS);
438d2daff3dSHong Zhang 
439cc4f23bcSHong Zhang PETSC_EXTERN PetscErrorCode TSGetStages(TS, PetscInt *, Vec *[]);
440818ad0c1SBarry Smith 
441014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetTime(TS, PetscReal *);
442014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetTime(TS, PetscReal);
443e5e524a1SHong Zhang PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS, PetscReal *);
44480275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS, PetscReal *);
44580275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS, PetscReal);
44680275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetStepNumber(TS, PetscInt *);
44780275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetStepNumber(TS, PetscInt);
448818ad0c1SBarry Smith 
449638cfed1SJed Brown PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS, PetscReal, Vec, Vec, void *);
450d1e9a80fSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS, PetscReal, Vec, Mat, Mat, void *);
451cd4cee2dSHong Zhang PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobianP)(TS, PetscReal, Vec, Mat, void *);
452014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS, Vec, TSRHSFunction, void *);
453014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS, Vec *, TSRHSFunction *, void **);
454014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS, Mat, Mat, TSRHSJacobian, void *);
455014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS, Mat *, Mat *, TSRHSJacobian *, void **);
456e1244c69SJed Brown PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS, PetscBool);
457818ad0c1SBarry Smith 
458ef20d060SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS, PetscReal, Vec, void *);
459ef20d060SBarry Smith PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS, TSSolutionFunction, void *);
460d56366bfSLisandro Dalcin PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSForcingFunction)(TS, PetscReal, Vec, void *);
461d56366bfSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS, TSForcingFunction, void *);
462ef20d060SBarry Smith 
463638cfed1SJed Brown PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS, PetscReal, Vec, Vec, Vec, void *);
464d1e9a80fSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
465014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetIFunction(TS, Vec, TSIFunction, void *);
466014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetIFunction(TS, Vec *, TSIFunction *, void **);
467014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS, Mat, Mat, TSIJacobian, void *);
468014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS, Mat *, Mat *, TSIJacobian *, void **);
469316643e7SJed Brown 
470efe9872eSLisandro Dalcin PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Function)(TS, PetscReal, Vec, Vec, Vec, Vec, void *);
471efe9872eSLisandro Dalcin PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Jacobian)(TS, PetscReal, Vec, Vec, Vec, PetscReal, PetscReal, Mat, Mat, void *);
472efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetI2Function(TS, Vec, TSI2Function, void *);
473efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetI2Function(TS, Vec *, TSI2Function *, void **);
474efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetI2Jacobian(TS, Mat, Mat, TSI2Jacobian, void *);
475efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetI2Jacobian(TS, Mat *, Mat *, TSI2Jacobian *, void **);
476efe9872eSLisandro Dalcin 
477545aaa6fSHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitSetIS(TS, const char[], IS);
478545aaa6fSHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetIS(TS, const char[], IS *);
479545aaa6fSHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitSetRHSFunction(TS, const char[], Vec, TSRHSFunction, void *);
4801d06f6b3SHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTS(TS, const char[], TS *);
4811d06f6b3SHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTSs(TS, PetscInt *, TS *[]);
4820fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSSetUseSplitRHSFunction(TS, PetscBool);
4830fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSGetUseSplitRHSFunction(TS, PetscBool *);
484d9194312SHong Zhang 
485014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS, PetscReal, Vec, Vec, void *);
486d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS, PetscReal, Vec, Mat, Mat, void *);
487014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS, PetscReal, Vec, Vec, Vec, void *);
488d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
4891c8a10a0SJed Brown PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS, PetscReal, Vec);
4909b7cd975SBarry Smith PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS, PetscReal, Vec);
491847ff0e1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSComputeIJacobianDefaultColor(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
492d7cfae9bSHong Zhang PETSC_EXTERN PetscErrorCode TSPruneIJacobianColor(TS, Mat, Mat);
493e34be4c2SBarry Smith 
494014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
495b8123daeSJed Brown PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS, PetscReal));
4969be3e283SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS, PetscReal, PetscInt, Vec *));
497c688d042SShri Abhyankar PETSC_EXTERN PetscErrorCode TSSetPostEvaluate(TS, PetscErrorCode (*)(TS));
498014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
499014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPreStep(TS);
500b8123daeSJed Brown PETSC_EXTERN PetscErrorCode TSPreStage(TS, PetscReal);
5019be3e283SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSPostStage(TS, PetscReal, PetscInt, Vec *);
502c688d042SShri Abhyankar PETSC_EXTERN PetscErrorCode TSPostEvaluate(TS);
503014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPostStep(TS);
504014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSInterpolate(TS, PetscReal, Vec);
505014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetTolerances(TS, PetscReal, Vec, PetscReal, Vec);
506c5033834SJed Brown PETSC_EXTERN PetscErrorCode TSGetTolerances(TS, PetscReal *, Vec *, PetscReal *, Vec *);
5077453f775SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSErrorWeightedNormInfinity(TS, Vec, Vec, PetscReal *, PetscReal *, PetscReal *);
5087453f775SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm2(TS, Vec, Vec, PetscReal *, PetscReal *, PetscReal *);
5097453f775SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *);
5108a175baeSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSErrorWeightedENormInfinity(TS, Vec, Vec, Vec, PetscReal *, PetscReal *, PetscReal *);
5118a175baeSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm2(TS, Vec, Vec, Vec, PetscReal *, PetscReal *, PetscReal *);
5128a175baeSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm(TS, Vec, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *);
513014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS, PetscReal);
514014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS, PetscReal *);
515d183316bSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSSetFunctionDomainError(TS, PetscErrorCode (*)(TS, PetscReal, Vec, PetscBool *));
516d183316bSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSFunctionDomainError(TS, PetscReal, Vec, PetscBool *);
517000e7ae3SMatthew Knepley 
518014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS, PetscErrorCode (*)(TS, PetscReal *, void *), void *);
5198d359177SBarry Smith PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS, PetscReal *, void *);
520014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS, PetscReal *);
521014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS, PetscReal);
522014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS, PetscErrorCode (*)(TS, Vec, void *, PetscReal *, PetscBool *), void *);
5238d359177SBarry Smith PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS, Vec, void *, PetscReal *, PetscBool *);
524014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS, Vec, PetscReal *, PetscBool *);
525014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS, PetscReal);
526014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);
52721c89e3eSBarry Smith 
528014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPythonSetType(TS, const char[]);
529ebead697SStefano Zampini PETSC_EXTERN PetscErrorCode TSPythonGetType(TS, const char *[]);
5301d6018f0SLisandro Dalcin 
531014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS, PetscReal, Vec, Vec);
532d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS, PetscReal, Vec, Mat, Mat);
533014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS, PetscReal, Vec, Vec, Vec, PetscBool);
534d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscBool);
535efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSComputeI2Function(TS, PetscReal, Vec, Vec, Vec, Vec);
536efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSComputeI2Jacobian(TS, PetscReal, Vec, Vec, Vec, PetscReal, PetscReal, Mat, Mat);
537f9c1d6abSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS, PetscReal, PetscReal, PetscReal *, PetscReal *);
538818ad0c1SBarry Smith 
539014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS, Vec, Vec);
540d6ebe24aSShri Abhyankar 
541967bb25aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetBoundaryLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *);
54224989b8cSPeter Brune PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM, TSRHSFunction, void *);
54324989b8cSPeter Brune PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM, TSRHSFunction *, void **);
544800f99ffSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionContextDestroy(DM, PetscErrorCode (*)(void *));
54524989b8cSPeter Brune PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM, TSRHSJacobian, void *);
54624989b8cSPeter Brune PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM, TSRHSJacobian *, void **);
547800f99ffSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobianContextDestroy(DM, PetscErrorCode (*)(void *));
54824989b8cSPeter Brune PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM, TSIFunction, void *);
54924989b8cSPeter Brune PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM, TSIFunction *, void **);
550800f99ffSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMTSSetIFunctionContextDestroy(DM, PetscErrorCode (*)(void *));
55124989b8cSPeter Brune PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM, TSIJacobian, void *);
55224989b8cSPeter Brune PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM, TSIJacobian *, void **);
553800f99ffSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMTSSetIJacobianContextDestroy(DM, PetscErrorCode (*)(void *));
554efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode DMTSSetI2Function(DM, TSI2Function, void *);
555efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode DMTSGetI2Function(DM, TSI2Function *, void **);
556800f99ffSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMTSSetI2FunctionContextDestroy(DM, PetscErrorCode (*)(void *));
557efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode DMTSSetI2Jacobian(DM, TSI2Jacobian, void *);
558efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode DMTSGetI2Jacobian(DM, TSI2Jacobian *, void **);
559800f99ffSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMTSSetI2JacobianContextDestroy(DM, PetscErrorCode (*)(void *));
560efe9872eSLisandro Dalcin 
561e3c11fc1SJed Brown PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSTransientVariable)(TS, Vec, Vec, void *);
562438f35afSJed Brown PETSC_EXTERN PetscErrorCode TSSetTransientVariable(TS, TSTransientVariable, void *);
563e3c11fc1SJed Brown PETSC_EXTERN PetscErrorCode DMTSSetTransientVariable(DM, TSTransientVariable, void *);
564e3c11fc1SJed Brown PETSC_EXTERN PetscErrorCode DMTSGetTransientVariable(DM, TSTransientVariable *, void *);
565e3c11fc1SJed Brown PETSC_EXTERN PetscErrorCode TSComputeTransientVariable(TS, Vec, Vec);
566e3c11fc1SJed Brown PETSC_EXTERN PetscErrorCode TSHasTransientVariable(TS, PetscBool *);
567e3c11fc1SJed Brown 
568ef20d060SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM, TSSolutionFunction, void *);
569ef20d060SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM, TSSolutionFunction *, void **);
570d56366bfSLisandro Dalcin PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM, TSForcingFunction, void *);
571d56366bfSLisandro Dalcin PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM, TSForcingFunction *, void **);
572e17bd7bbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCheckResidual(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscReal *);
573e17bd7bbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCheckJacobian(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscBool *, PetscReal *);
5747f96f943SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec);
5759b7cd975SBarry Smith 
5760c9409e4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSGetIFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, Vec, void *), void **);
577d433e6cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, Vec, void *), void *);
5780c9409e4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSGetIJacobianLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void **);
579d433e6cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void *);
5800c9409e4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSGetRHSFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, void *), void **);
581d433e6cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *);
582b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrix(DM);
583b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrixLumped(DM);
584b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSDestroyRHSMassMatrix(DM);
5856c6b9e74SPeter Brune 
5866c6b9e74SPeter Brune PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM, PetscErrorCode (*)(void *, PetscViewer), PetscErrorCode (*)(void **, PetscViewer));
5876c6b9e74SPeter Brune PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM, PetscErrorCode (*)(void *, PetscViewer), PetscErrorCode (*)(void **, PetscViewer));
5886c6b9e74SPeter Brune 
5896c6b9e74SPeter Brune PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo *, PetscReal, void *, void *, void *);
590d1e9a80fSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo *, PetscReal, void *, Mat, Mat, void *);
5916c6b9e74SPeter Brune PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo *, PetscReal, void *, void *, void *, void *);
592d1e9a80fSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo *, PetscReal, void *, void *, PetscReal, Mat, Mat, void *);
5936c6b9e74SPeter Brune 
5946c6b9e74SPeter Brune PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM, InsertMode, PetscErrorCode (*)(DMDALocalInfo *, PetscReal, void *, void *, void *), void *);
595d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM, PetscErrorCode (*)(DMDALocalInfo *, PetscReal, void *, Mat, Mat, void *), void *);
5966c6b9e74SPeter Brune PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM, InsertMode, PetscErrorCode (*)(DMDALocalInfo *, PetscReal, void *, void *, void *, void *), void *);
597d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM, PetscErrorCode (*)(DMDALocalInfo *, PetscReal, void *, void *, PetscReal, Mat, Mat, void *), void *);
5986c6b9e74SPeter Brune 
599b6aca0f9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTSGetGeometryFVM(DM, Vec *, Vec *, PetscReal *);
6006dbbd306SMatthew G. Knepley 
601be1b0d75SMatthew G. Knepley typedef struct _n_TSMonitorLGCtx *TSMonitorLGCtx;
602d1212d36SBarry Smith typedef struct {
603d1212d36SBarry Smith   Vec            ray;
604d1212d36SBarry Smith   VecScatter     scatter;
605d1212d36SBarry Smith   PetscViewer    viewer;
60651b4a12fSMatthew G. Knepley   TSMonitorLGCtx lgctx;
607d1212d36SBarry Smith } TSMonitorDMDARayCtx;
608d1212d36SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void **);
609d1212d36SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS, PetscInt, PetscReal, Vec, void *);
61051b4a12fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS, PetscInt, PetscReal, Vec, void *);
611d1212d36SBarry Smith 
612bdad233fSMatthew Knepley /* Dynamic creation and loading functions */
613140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList TSList;
61419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode    TSGetType(TS, TSType *);
61519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode    TSSetType(TS, TSType);
616bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode    TSRegister(const char[], PetscErrorCode (*)(TS));
61730de9b25SBarry Smith 
618014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetSNES(TS, SNES *);
619deb2cd25SJed Brown PETSC_EXTERN PetscErrorCode TSSetSNES(TS, SNES);
620014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetKSP(TS, KSP *);
621818ad0c1SBarry Smith 
622014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSView(TS, PetscViewer);
62355849f57SBarry Smith PETSC_EXTERN PetscErrorCode TSLoad(TS, PetscViewer);
624fe2efc57SMark PETSC_EXTERN PetscErrorCode TSViewFromOptions(TS, PetscObject, const char[]);
625fe2efc57SMark PETSC_EXTERN PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory, PetscObject, const char[]);
62655849f57SBarry Smith 
62755849f57SBarry Smith #define TS_FILE_CLASSID 1211225
62821c89e3eSBarry Smith 
629014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS, void *);
630014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS, void *);
63121c89e3eSBarry Smith 
632a9f9c1f6SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtx *);
633a9f9c1f6SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *);
6344f09c107SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS, PetscInt, PetscReal, Vec, void *);
6354f09c107SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS, PetscInt, PetscReal, Vec, void *);
63631152f8aSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS, const char *const *);
637b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS, const char *const **);
638e673d494SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx, const char *const *);
639a66092f1SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS, const char *const *);
640a66092f1SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx, const char *const *);
6417684fa3eSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS, PetscErrorCode (*)(void *, Vec, Vec *), PetscErrorCode (*)(void *), void *);
6427684fa3eSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx, PetscErrorCode (*)(void *, Vec, Vec *), PetscErrorCode (*)(void *), void *);
6434f09c107SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS, PetscInt, PetscReal, Vec, void *);
644201da799SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS, PetscInt, PetscReal, Vec, void *);
645201da799SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS, PetscInt, PetscReal, Vec, void *);
646edbaebb3SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorError(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
647d0c080abSJoseph Pusztay PETSC_EXTERN PetscErrorCode TSDMSwarmMonitorMoments(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *);
648ef20d060SBarry Smith 
649e669de00SBarry Smith struct _n_TSMonitorLGCtxNetwork {
650e669de00SBarry Smith   PetscInt     nlg;
651e669de00SBarry Smith   PetscDrawLG *lg;
652e669de00SBarry Smith   PetscBool    semilogy;
653e669de00SBarry Smith   PetscInt     howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
654e669de00SBarry Smith };
655e669de00SBarry Smith typedef struct _n_TSMonitorLGCtxNetwork *TSMonitorLGCtxNetwork;
656e669de00SBarry Smith PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkDestroy(TSMonitorLGCtxNetwork *);
657e669de00SBarry Smith PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkCreate(TS, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtxNetwork *);
658e669de00SBarry Smith PETSC_EXTERN PetscErrorCode              TSMonitorLGCtxNetworkSolution(TS, PetscInt, PetscReal, Vec, void *);
659e669de00SBarry Smith 
660b3d3934dSBarry Smith typedef struct _n_TSMonitorEnvelopeCtx *TSMonitorEnvelopeCtx;
661b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeCtxCreate(TS, TSMonitorEnvelopeCtx *);
662b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode             TSMonitorEnvelope(TS, PetscInt, PetscReal, Vec, void *);
663b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeGetBounds(TS, Vec *, Vec *);
664b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode             TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *);
665b3d3934dSBarry Smith 
6668189c53fSBarry Smith typedef struct _n_TSMonitorSPEigCtx *TSMonitorSPEigCtx;
6678189c53fSBarry Smith PETSC_EXTERN PetscErrorCode          TSMonitorSPEigCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorSPEigCtx *);
6688189c53fSBarry Smith PETSC_EXTERN PetscErrorCode          TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *);
6698189c53fSBarry Smith PETSC_EXTERN PetscErrorCode          TSMonitorSPEig(TS, PetscInt, PetscReal, Vec, void *);
6703914022bSBarry Smith 
6711b575b74SJoseph Pusztay typedef struct _n_TSMonitorSPCtx *TSMonitorSPCtx;
672*60e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorSPCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscBool, PetscBool, TSMonitorSPCtx *);
6731b575b74SJoseph Pusztay PETSC_EXTERN PetscErrorCode       TSMonitorSPCtxDestroy(TSMonitorSPCtx *);
6740ec8ee2bSJoseph Pusztay PETSC_EXTERN PetscErrorCode       TSMonitorSPSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
6751b575b74SJoseph Pusztay 
676*60e16b1bSMatthew G. Knepley typedef struct _n_TSMonitorHGCtx *TSMonitorHGCtx;
677*60e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorHGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscInt, PetscBool, TSMonitorHGCtx *);
678*60e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
679*60e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorHGCtxDestroy(TSMonitorHGCtx *);
680*60e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode       TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *);
681*60e16b1bSMatthew G. Knepley 
6822589fa24SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetEventHandler(TS, PetscInt, PetscInt[], PetscBool[], PetscErrorCode (*)(TS, PetscReal, Vec, PetscScalar[], void *), PetscErrorCode (*)(TS, PetscInt, PetscInt[], PetscReal, Vec, PetscBool, void *), void *);
6830e37a75fSKarl Rupp PETSC_EXTERN PetscErrorCode TSSetPostEventIntervalStep(TS, PetscReal);
6846427ac75SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS, PetscReal, PetscReal[]);
6851ea83e56SMiguel PETSC_EXTERN PetscErrorCode TSGetNumEvents(TS, PetscInt *);
6861ea83e56SMiguel 
68776bdecfbSBarry Smith /*J
68887497f52SBarry Smith    TSSSPType - string with the name of `TSSSP` scheme.
689815f1ad5SJed Brown 
690815f1ad5SJed Brown    Level: beginner
691815f1ad5SJed Brown 
692db781477SPatrick Sanan .seealso: `TSSSPSetType()`, `TS`
69376bdecfbSBarry Smith J*/
69419fd82e9SBarry Smith typedef const char *TSSSPType;
695815f1ad5SJed Brown #define TSSSPRKS2  "rks2"
696815f1ad5SJed Brown #define TSSSPRKS3  "rks3"
697815f1ad5SJed Brown #define TSSSPRK104 "rk104"
698815f1ad5SJed Brown 
69919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode    TSSSPSetType(TS, TSSSPType);
70019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode    TSSSPGetType(TS, TSSSPType *);
701014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSSPSetNumStages(TS, PetscInt);
702014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSSPGetNumStages(TS, PetscInt *);
703787849ffSJed Brown PETSC_EXTERN PetscErrorCode    TSSSPInitializePackage(void);
704560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode    TSSSPFinalizePackage(void);
705013797aaSBarry Smith PETSC_EXTERN PetscFunctionList TSSSPList;
706815f1ad5SJed Brown 
707ed657a08SJed Brown /*S
70884df9cb4SJed Brown    TSAdapt - Abstract object that manages time-step adaptivity
70984df9cb4SJed Brown 
71084df9cb4SJed Brown    Level: beginner
71184df9cb4SJed Brown 
712db781477SPatrick Sanan .seealso: `TS`, `TSAdaptCreate()`, `TSAdaptType`
71384df9cb4SJed Brown S*/
71484df9cb4SJed Brown typedef struct _p_TSAdapt *TSAdapt;
71584df9cb4SJed Brown 
716f0fc11ceSJed Brown /*J
71787497f52SBarry Smith     TSAdaptType - String with the name of `TSAdapt` scheme.
71884df9cb4SJed Brown 
71984df9cb4SJed Brown    Level: beginner
72084df9cb4SJed Brown 
721db781477SPatrick Sanan .seealso: `TSAdaptSetType()`, `TS`
722f0fc11ceSJed Brown J*/
723f8963224SJed Brown typedef const char *TSAdaptType;
724b0f836d7SJed Brown #define TSADAPTNONE    "none"
7251917a363SLisandro Dalcin #define TSADAPTBASIC   "basic"
726cb7ab186SLisandro Dalcin #define TSADAPTDSP     "dsp"
7278d59e960SJed Brown #define TSADAPTCFL     "cfl"
7281917a363SLisandro Dalcin #define TSADAPTGLEE    "glee"
729d949e4a4SStefano Zampini #define TSADAPTHISTORY "history"
73084df9cb4SJed Brown 
731552698daSJed Brown PETSC_EXTERN PetscErrorCode TSGetAdapt(TS, TSAdapt *);
732bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[], PetscErrorCode (*)(TSAdapt));
733607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void);
734014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
735014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm, TSAdapt *);
73619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt, TSAdaptType);
737d0288e4fSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetType(TSAdapt, TSAdaptType *);
738014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt, const char[]);
739014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
740014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt, const char[], PetscInt, PetscInt, PetscReal, PetscReal, PetscBool);
741014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt, PetscInt *, const PetscInt **, const PetscInt **, const PetscReal **, const PetscReal **);
742014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt, TS, PetscReal, PetscInt *, PetscReal *, PetscBool *);
743b295832fSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt, TS, PetscReal, Vec, PetscBool *);
744014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt, PetscViewer);
745ad6bc421SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt, PetscViewer);
746dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(TSAdapt, PetscOptionItems *);
747099af0a3SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt);
748014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt *);
749014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt, PetscBool);
750bf997491SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt, PetscBool);
7511917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptSetSafety(TSAdapt, PetscReal, PetscReal);
7521917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetSafety(TSAdapt, PetscReal *, PetscReal *);
75376cddca1SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt, PetscReal);
75476cddca1SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt, PetscReal *);
7551917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptSetClip(TSAdapt, PetscReal, PetscReal);
7561917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetClip(TSAdapt, PetscReal *, PetscReal *);
75762c23b28SMark PETSC_EXTERN PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt, PetscReal);
75862c23b28SMark PETSC_EXTERN PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt, PetscReal *);
759014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt, PetscReal, PetscReal);
7601917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetStepLimits(TSAdapt, PetscReal *, PetscReal *);
761b295832fSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt, PetscErrorCode (*)(TSAdapt, TS, PetscReal, Vec, PetscBool *));
762d949e4a4SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptHistorySetHistory(TSAdapt, PetscInt n, PetscReal hist[], PetscBool);
76375017583SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptHistorySetTrajectory(TSAdapt, TSTrajectory, PetscBool);
76475017583SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptHistoryGetStep(TSAdapt, PetscInt, PetscReal *, PetscReal *);
765de50f1caSBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt, PetscInt);
766cb7ab186SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptDSPSetFilter(TSAdapt, const char *);
767cb7ab186SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptDSPSetPID(TSAdapt, PetscReal, PetscReal, PetscReal);
76884df9cb4SJed Brown 
76984df9cb4SJed Brown /*S
77087497f52SBarry Smith    TSGLLEAdapt - Abstract object that manages time-step adaptivity for `TSGLLE`
771ed657a08SJed Brown 
772ed657a08SJed Brown    Level: beginner
773ed657a08SJed Brown 
77484df9cb4SJed Brown    Developer Notes:
77587497f52SBarry Smith    This functionality should be replaced by the `TSAdapt`.
77684df9cb4SJed Brown 
777db781477SPatrick Sanan .seealso: `TSGLLE`, `TSGLLEAdaptCreate()`, `TSGLLEAdaptType`
778ed657a08SJed Brown S*/
77926d28e4eSEmil Constantinescu typedef struct _p_TSGLLEAdapt *TSGLLEAdapt;
780ed657a08SJed Brown 
78176bdecfbSBarry Smith /*J
78287497f52SBarry Smith     TSGLLEAdaptType - String with the name of `TSGLLEAdapt` scheme
783ed657a08SJed Brown 
784ed657a08SJed Brown    Level: beginner
785ed657a08SJed Brown 
786db781477SPatrick Sanan .seealso: `TSGLLEAdaptSetType()`, `TS`
78776bdecfbSBarry Smith J*/
78826d28e4eSEmil Constantinescu typedef const char *TSGLLEAdaptType;
78926d28e4eSEmil Constantinescu #define TSGLLEADAPT_NONE "none"
79026d28e4eSEmil Constantinescu #define TSGLLEADAPT_SIZE "size"
79126d28e4eSEmil Constantinescu #define TSGLLEADAPT_BOTH "both"
792ed657a08SJed Brown 
79326d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptRegister(const char[], PetscErrorCode (*)(TSGLLEAdapt));
79426d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptInitializePackage(void);
79526d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptFinalizePackage(void);
79626d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate(MPI_Comm, TSGLLEAdapt *);
79726d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt, TSGLLEAdaptType);
79826d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt, const char[]);
79926d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptChoose(TSGLLEAdapt, PetscInt, const PetscInt[], const PetscReal[], const PetscReal[], PetscInt, PetscReal, PetscReal, PetscInt *, PetscReal *, PetscBool *);
80026d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt, PetscViewer);
801dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetFromOptions(TSGLLEAdapt, PetscOptionItems *);
80226d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt *);
803ed657a08SJed Brown 
80476bdecfbSBarry Smith /*J
80587497f52SBarry Smith     TSGLLEAcceptType - String with the name of `TSGLLEAccept` scheme
806ed657a08SJed Brown 
807ed657a08SJed Brown    Level: beginner
808ed657a08SJed Brown 
809db781477SPatrick Sanan .seealso: `TSGLLESetAcceptType()`, `TS`
81076bdecfbSBarry Smith J*/
81126d28e4eSEmil Constantinescu typedef const char *TSGLLEAcceptType;
81226d28e4eSEmil Constantinescu #define TSGLLEACCEPT_ALWAYS "always"
813ed657a08SJed Brown 
81426d28e4eSEmil Constantinescu PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLLEAcceptFunction)(TS, PetscReal, PetscReal, const PetscReal[], PetscBool *);
81526d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAcceptRegister(const char[], TSGLLEAcceptFunction);
816ed657a08SJed Brown 
81776bdecfbSBarry Smith /*J
81826d28e4eSEmil Constantinescu   TSGLLEType - family of time integration method within the General Linear class
81918b56cb9SJed Brown 
82018b56cb9SJed Brown   Level: beginner
82118b56cb9SJed Brown 
822db781477SPatrick Sanan .seealso: `TSGLLESetType()`, `TSGLLERegister()`
82376bdecfbSBarry Smith J*/
82426d28e4eSEmil Constantinescu typedef const char *TSGLLEType;
82526d28e4eSEmil Constantinescu #define TSGLLE_IRKS "irks"
82618b56cb9SJed Brown 
82726d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLERegister(const char[], PetscErrorCode (*)(TS));
82826d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEInitializePackage(void);
82926d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEFinalizePackage(void);
83026d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLESetType(TS, TSGLLEType);
83126d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEGetAdapt(TS, TSGLLEAdapt *);
83226d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLESetAcceptType(TS, TSGLLEAcceptType);
83318b56cb9SJed Brown 
834020d8f30SJed Brown /*J
8357d5697caSHong Zhang     TSEIMEXType - String with the name of an Extrapolated IMEX method.
8367d5697caSHong Zhang 
8377d5697caSHong Zhang    Level: beginner
8387d5697caSHong Zhang 
839db781477SPatrick Sanan .seealso: `TSEIMEXSetType()`, `TS`, `TSEIMEX`, `TSEIMEXRegister()`
8407d5697caSHong Zhang J*/
8417d5697caSHong Zhang #define TSEIMEXType char *
8427d5697caSHong Zhang 
8437d5697caSHong Zhang PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts, PetscInt);
8447d5697caSHong Zhang PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts, PetscInt, PetscInt);
8457d5697caSHong Zhang PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS, PetscBool);
8467d5697caSHong Zhang 
8477d5697caSHong Zhang /*J
848f68a32c8SEmil Constantinescu     TSRKType - String with the name of a Runge-Kutta method.
849f68a32c8SEmil Constantinescu 
850f68a32c8SEmil Constantinescu    Level: beginner
851f68a32c8SEmil Constantinescu 
852db781477SPatrick Sanan .seealso: `TSRKSetType()`, `TS`, `TSRK`, `TSRKRegister()`
853f68a32c8SEmil Constantinescu J*/
854f68a32c8SEmil Constantinescu typedef const char *TSRKType;
855f68a32c8SEmil Constantinescu #define TSRK1FE "1fe"
856fdefd5e5SDebojyoti Ghosh #define TSRK2A  "2a"
857e7685601SHong Zhang #define TSRK2B  "2b"
858f68a32c8SEmil Constantinescu #define TSRK3   "3"
859fdefd5e5SDebojyoti Ghosh #define TSRK3BS "3bs"
860f68a32c8SEmil Constantinescu #define TSRK4   "4"
861fdefd5e5SDebojyoti Ghosh #define TSRK5F  "5f"
862fdefd5e5SDebojyoti Ghosh #define TSRK5DP "5dp"
86305e23783SLisandro Dalcin #define TSRK5BS "5bs"
86437acaa02SHendrik Ranocha #define TSRK6VR "6vr"
86537acaa02SHendrik Ranocha #define TSRK7VR "7vr"
86637acaa02SHendrik Ranocha #define TSRK8VR "8vr"
86705e23783SLisandro Dalcin 
8682ea87ba9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRKGetOrder(TS, PetscInt *);
8690f7a28cdSStefano Zampini PETSC_EXTERN PetscErrorCode TSRKGetType(TS, TSRKType *);
8700f7a28cdSStefano Zampini PETSC_EXTERN PetscErrorCode TSRKSetType(TS, TSRKType);
8710f7a28cdSStefano Zampini PETSC_EXTERN PetscErrorCode TSRKGetTableau(TS, PetscInt *, const PetscReal **, const PetscReal **, const PetscReal **, const PetscReal **, PetscInt *, const PetscReal **, PetscBool *);
8720fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSRKSetMultirate(TS, PetscBool);
8730fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSRKGetMultirate(TS, PetscBool *);
874f68a32c8SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
875f68a32c8SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void);
876560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void);
877f68a32c8SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void);
878f68a32c8SEmil Constantinescu 
879f68a32c8SEmil Constantinescu /*J
8804b84eec9SHong Zhang    TSMPRKType - String with the name of a Partitioned Runge-Kutta method
8819f537349Sluzhanghpp 
8829f537349Sluzhanghpp    Level: beginner
8839f537349Sluzhanghpp 
884db781477SPatrick Sanan .seealso: `TSMPRKSetType()`, `TS`, `TSMPRK`, `TSMPRKRegister()`
8859f537349Sluzhanghpp J*/
8864b84eec9SHong Zhang typedef const char *TSMPRKType;
88719c2959aSHong Zhang #define TSMPRK2A22 "2a22"
88819c2959aSHong Zhang #define TSMPRK2A23 "2a23"
88919c2959aSHong Zhang #define TSMPRK2A32 "2a32"
89019c2959aSHong Zhang #define TSMPRK2A33 "2a33"
8914b84eec9SHong Zhang #define TSMPRKP2   "p2"
8924b84eec9SHong Zhang #define TSMPRKP3   "p3"
8939f537349Sluzhanghpp 
8944b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKGetType(TS ts, TSMPRKType *);
8954b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKSetType(TS ts, TSMPRKType);
89679bc8a77SHong 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[]);
8974b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKInitializePackage(void);
8984b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKFinalizePackage(void);
8994b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKRegisterDestroy(void);
9009f537349Sluzhanghpp 
9019f537349Sluzhanghpp /*J
902d2567f34SHong Zhang     TSIRKType - String with the name of an implicit Runge-Kutta method.
903d2567f34SHong Zhang 
904d2567f34SHong Zhang    Level: beginner
905d2567f34SHong Zhang 
906db781477SPatrick Sanan .seealso: `TSIRKSetType()`, `TS`, `TSIRK`, `TSIRKRegister()`
907d2567f34SHong Zhang J*/
908d2567f34SHong Zhang typedef const char *TSIRKType;
909d2567f34SHong Zhang #define TSIRKGAUSS "gauss"
910d2567f34SHong Zhang 
911d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKGetOrder(TS, PetscInt *);
912d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKGetType(TS, TSIRKType *);
913d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKSetType(TS, TSIRKType);
914d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKGetNumStages(TS, PetscInt *);
915d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKSetNumStages(TS, PetscInt);
916d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKRegister(const char[], PetscErrorCode (*function)(TS));
917d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKTableauCreate(TS, PetscInt, const PetscReal *, const PetscReal *, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, const PetscScalar *);
918d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKInitializePackage(void);
919d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKFinalizePackage(void);
920d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKRegisterDestroy(void);
921d2567f34SHong Zhang 
922d2567f34SHong Zhang /*J
923b6a60446SDebojyoti Ghosh     TSGLEEType - String with the name of a General Linear with Error Estimation method.
924b6a60446SDebojyoti Ghosh 
925b6a60446SDebojyoti Ghosh    Level: beginner
926b6a60446SDebojyoti Ghosh 
927db781477SPatrick Sanan .seealso: `TSGLEESetType()`, `TS`, `TSGLEE`, `TSGLEERegister()`
928b6a60446SDebojyoti Ghosh J*/
929b6a60446SDebojyoti Ghosh typedef const char *TSGLEEType;
930ab8c5912SEmil Constantinescu #define TSGLEEi1      "BE1"
931b6a60446SDebojyoti Ghosh #define TSGLEE23      "23"
932b6a60446SDebojyoti Ghosh #define TSGLEE24      "24"
933b6a60446SDebojyoti Ghosh #define TSGLEE25I     "25i"
934b6a60446SDebojyoti Ghosh #define TSGLEE35      "35"
935b6a60446SDebojyoti Ghosh #define TSGLEEEXRK2A  "exrk2a"
936b6a60446SDebojyoti Ghosh #define TSGLEERK32G1  "rk32g1"
937b6a60446SDebojyoti Ghosh #define TSGLEERK285EX "rk285ex"
938b6a60446SDebojyoti Ghosh /*J
939b6a60446SDebojyoti Ghosh     TSGLEEMode - String with the mode of error estimation for a General Linear with Error Estimation method.
940b6a60446SDebojyoti Ghosh 
941b6a60446SDebojyoti Ghosh    Level: beginner
942b6a60446SDebojyoti Ghosh 
943db781477SPatrick Sanan .seealso: `TSGLEESetMode()`, `TS`, `TSGLEE`, `TSGLEERegister()`
944b6a60446SDebojyoti Ghosh J*/
945b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEEGetType(TS ts, TSGLEEType *);
946b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEESetType(TS ts, TSGLEEType);
94757df6a1bSDebojyoti 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[]);
948b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEEFinalizePackage(void);
949b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEEInitializePackage(void);
950b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEERegisterDestroy(void);
951b6a60446SDebojyoti Ghosh 
952b6a60446SDebojyoti Ghosh /*J
953020d8f30SJed Brown     TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX method.
954020d8f30SJed Brown 
955020d8f30SJed Brown    Level: beginner
956020d8f30SJed Brown 
957db781477SPatrick Sanan .seealso: `TSARKIMEXSetType()`, `TS`, `TSARKIMEX`, `TSARKIMEXRegister()`
958020d8f30SJed Brown J*/
95919fd82e9SBarry Smith typedef const char *TSARKIMEXType;
960e817cc15SEmil Constantinescu #define TSARKIMEX1BEE   "1bee"
9611f80e275SEmil Constantinescu #define TSARKIMEXA2     "a2"
9621f80e275SEmil Constantinescu #define TSARKIMEXL2     "l2"
9631f80e275SEmil Constantinescu #define TSARKIMEXARS122 "ars122"
9641f80e275SEmil Constantinescu #define TSARKIMEX2C     "2c"
9658a381b04SJed Brown #define TSARKIMEX2D     "2d"
966a3a57f36SJed Brown #define TSARKIMEX2E     "2e"
9676cf0794eSJed Brown #define TSARKIMEXPRSSP2 "prssp2"
9688a381b04SJed Brown #define TSARKIMEX3      "3"
9696cf0794eSJed Brown #define TSARKIMEXBPR3   "bpr3"
9706cf0794eSJed Brown #define TSARKIMEXARS443 "ars443"
9718a381b04SJed Brown #define TSARKIMEX4      "4"
9728a381b04SJed Brown #define TSARKIMEX5      "5"
97319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts, TSARKIMEXType *);
97419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts, TSARKIMEXType);
975014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS, PetscBool);
9763a28c0e5SStefano Zampini PETSC_EXTERN PetscErrorCode TSARKIMEXGetFullyImplicit(TS, PetscBool *);
97719fd82e9SBarry 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[]);
978607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void);
979560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
980014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);
9818a381b04SJed Brown 
982020d8f30SJed Brown /*J
983020d8f30SJed Brown     TSRosWType - String with the name of a Rosenbrock-W method.
984020d8f30SJed Brown 
985020d8f30SJed Brown    Level: beginner
986020d8f30SJed Brown 
987db781477SPatrick Sanan .seealso: `TSRosWSetType()`, `TS`, `TSROSW`, `TSRosWRegister()`
988020d8f30SJed Brown J*/
98919fd82e9SBarry Smith typedef const char *TSRosWType;
99061692a83SJed Brown #define TSROSW2M          "2m"
99161692a83SJed Brown #define TSROSW2P          "2p"
992fe7e6d57SJed Brown #define TSROSWRA3PW       "ra3pw"
993fe7e6d57SJed Brown #define TSROSWRA34PW2     "ra34pw2"
994ef3c5b88SJed Brown #define TSROSWRODAS3      "rodas3"
995ef3c5b88SJed Brown #define TSROSWSANDU3      "sandu3"
996961f28d0SJed Brown #define TSROSWASSP3P3S1C  "assp3p3s1c"
997961f28d0SJed Brown #define TSROSWLASSP3P4S2C "lassp3p4s2c"
99843b21953SEmil Constantinescu #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
999753f8adbSEmil Constantinescu #define TSROSWARK3        "ark3"
10003606a31eSEmil Constantinescu #define TSROSWTHETA1      "theta1"
10013606a31eSEmil Constantinescu #define TSROSWTHETA2      "theta2"
100242faf41dSJed Brown #define TSROSWGRK4T       "grk4t"
100342faf41dSJed Brown #define TSROSWSHAMP4      "shamp4"
100442faf41dSJed Brown #define TSROSWVELDD4      "veldd4"
100542faf41dSJed Brown #define TSROSW4L          "4l"
1006b1c69cc3SEmil Constantinescu 
10076bd7aeb5SHong Zhang PETSC_EXTERN PetscErrorCode TSRosWGetType(TS, TSRosWType *);
10086bd7aeb5SHong Zhang PETSC_EXTERN PetscErrorCode TSRosWSetType(TS, TSRosWType);
1009014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS, PetscBool);
101019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]);
101119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal);
1012607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void);
1013560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
1014014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);
1015e27a552bSJed Brown 
1016211a84d6SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSBDFSetOrder(TS, PetscInt);
1017211a84d6SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSBDFGetOrder(TS, PetscInt *);
1018211a84d6SLisandro Dalcin 
10196bd7aeb5SHong Zhang /*J
1020e49d4f37SHong Zhang   TSBasicSymplecticType - String with the name of a basic symplectic integration method.
10216bd7aeb5SHong Zhang 
10226bd7aeb5SHong Zhang   Level: beginner
10236bd7aeb5SHong Zhang 
1024db781477SPatrick Sanan   .seealso: `TSBasicSymplecticSetType()`, `TS`, `TSBASICSYMPLECTIC`, `TSBasicSymplecticRegister()`
10256bd7aeb5SHong Zhang J*/
1026e49d4f37SHong Zhang typedef const char *TSBasicSymplecticType;
1027e49d4f37SHong Zhang #define TSBASICSYMPLECTICSIEULER   "1"
1028e49d4f37SHong Zhang #define TSBASICSYMPLECTICVELVERLET "2"
1029e49d4f37SHong Zhang #define TSBASICSYMPLECTIC3         "3"
1030e49d4f37SHong Zhang #define TSBASICSYMPLECTIC4         "4"
1031e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticSetType(TS, TSBasicSymplecticType);
1032e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticGetType(TS, TSBasicSymplecticType *);
1033e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegister(TSBasicSymplecticType, PetscInt, PetscInt, PetscReal[], PetscReal[]);
1034e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticInitializePackage(void);
1035e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticFinalizePackage(void);
1036e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterDestroy(void);
10376bd7aeb5SHong Zhang 
1038db4653c2SDaniel Finn /*J
103987497f52SBarry Smith   TSDiscreteGradient - The Discrete Gradient integrator is a timestepper for Hamiltonian systems designed to conserve the first integral (energy),
104087497f52SBarry Smith   but also has the property for some systems of monotonicity in a functional.
1041db4653c2SDaniel Finn 
1042db4653c2SDaniel Finn   Level: beginner
1043db4653c2SDaniel Finn 
1044db781477SPatrick Sanan   .seealso: `TS`, `TSDISCGRAD`, `TSDiscGradSetFormulation()`, `TSDiscGradGetFormulation`
1045db4653c2SDaniel Finn J*/
104640be0ff1SMatthew 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 *);
1047db4653c2SDaniel Finn PETSC_EXTERN PetscErrorCode TSDiscGradIsGonzalez(TS, PetscBool *);
1048db4653c2SDaniel Finn PETSC_EXTERN PetscErrorCode TSDiscGradUseGonzalez(TS, PetscBool);
104940be0ff1SMatthew G. Knepley 
105083e2fdc7SBarry Smith /*
10510f3b3ca1SHong Zhang        PETSc interface to Sundials
105283e2fdc7SBarry Smith */
1053e808b789SPatrick Sanan #ifdef PETSC_HAVE_SUNDIALS2
10549371c9d4SSatish Balay typedef enum {
10559371c9d4SSatish Balay   SUNDIALS_ADAMS = 1,
10569371c9d4SSatish Balay   SUNDIALS_BDF   = 2
10579371c9d4SSatish Balay } TSSundialsLmmType;
10586a6fc655SJed Brown PETSC_EXTERN const char *const TSSundialsLmmTypes[];
10599371c9d4SSatish Balay typedef enum {
10609371c9d4SSatish Balay   SUNDIALS_MODIFIED_GS  = 1,
10619371c9d4SSatish Balay   SUNDIALS_CLASSICAL_GS = 2
10629371c9d4SSatish Balay } TSSundialsGramSchmidtType;
10636a6fc655SJed Brown PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
1064014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSundialsSetType(TS, TSSundialsLmmType);
1065014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSundialsGetPC(TS, PC *);
1066014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSundialsSetTolerance(TS, PetscReal, PetscReal);
1067014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSundialsSetMinTimeStep(TS, PetscReal);
1068014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSundialsSetMaxTimeStep(TS, PetscReal);
1069014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSundialsGetIterations(TS, PetscInt *, PetscInt *);
1070014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSundialsSetGramSchmidtType(TS, TSSundialsGramSchmidtType);
1071014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSundialsSetGMRESRestart(TS, PetscInt);
1072014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSundialsSetLinearTolerance(TS, PetscReal);
1073014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSundialsMonitorInternalSteps(TS, PetscBool);
1074014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSundialsGetParameters(TS, PetscInt *, long *[], double *[]);
1075014dd563SJed Brown PETSC_EXTERN PetscErrorCode    TSSundialsSetMaxl(TS, PetscInt);
1076c4e80e11SFlorian PETSC_EXTERN PetscErrorCode    TSSundialsSetMaxord(TS, PetscInt);
1077918687b8SPatrick Sanan PETSC_EXTERN PetscErrorCode    TSSundialsSetUseDense(TS, PetscBool);
10786dc17ff5SMatthew Knepley #endif
107983e2fdc7SBarry Smith 
1080014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS, PetscReal);
1081014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS, PetscReal *);
1082014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS, PetscBool *);
1083014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS, PetscBool);
10840de4c49aSJed Brown 
1085014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS, PetscReal);
1086014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS, PetscReal, PetscReal, PetscReal);
1087014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS, PetscReal *, PetscReal *, PetscReal *);
108888df8a41SLisandro Dalcin 
1089818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAlpha2SetRadius(TS, PetscReal);
1090818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAlpha2SetParams(TS, PetscReal, PetscReal, PetscReal, PetscReal);
1091818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAlpha2GetParams(TS, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
1092818efac9SLisandro Dalcin 
1093014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetDM(TS, DM);
1094014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetDM(TS, DM *);
10956c699258SBarry Smith 
1096014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES, Vec, Vec, void *);
1097d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES, Vec, Mat, Mat, void *);
10980f5c6efeSJed Brown 
1099f3b1f45cSBarry Smith PETSC_EXTERN PetscErrorCode TSRHSJacobianTest(TS, PetscBool *);
1100f3b1f45cSBarry Smith PETSC_EXTERN PetscErrorCode TSRHSJacobianTestTranspose(TS, PetscBool *);
1101aad739acSMatthew G. Knepley 
11022e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSGetComputeInitialCondition(TS, PetscErrorCode (**)(TS, Vec));
11032e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSSetComputeInitialCondition(TS, PetscErrorCode (*)(TS, Vec));
11042e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSComputeInitialCondition(TS, Vec);
11052e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSGetComputeExactError(TS, PetscErrorCode (**)(TS, Vec, Vec));
1106aad739acSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSSetComputeExactError(TS, PetscErrorCode (*)(TS, Vec, Vec));
1107aad739acSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSComputeExactError(TS, Vec, Vec);
1108f2ed2dc7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscConvEstUseTS(PetscConvEst, PetscBool);
1109d60b7d5cSBarry Smith 
1110d60b7d5cSBarry Smith PETSC_EXTERN PetscErrorCode TSSetMatStructure(TS, MatStructure);
1111818ad0c1SBarry Smith #endif
1112