xref: /petsc/include/petscts.h (revision 76be6f4ff3bd4e251c19fc00ebbebfd58b6e7589)
1 /*
2    User interface for the timestepping package. This package
3    is for use in solving time-dependent PDEs.
4 */
5 #if !defined(PETSCTS_H)
6 #define PETSCTS_H
7 
8 #include <petscsnes.h>
9 #include <petscconvest.h>
10 
11 /* SUBMANSEC = TS */
12 
13 /*S
14      TS - Abstract PETSc object that manages all time-steppers (ODE integrators)
15 
16    Level: beginner
17 
18 .seealso: `TSCreate()`, `TSSetType()`, `TSType`, `SNES`, `KSP`, `PC`, `TSDestroy()`
19 S*/
20 typedef struct _p_TS* TS;
21 
22 /*J
23     TSType - String with the name of a PETSc TS method.
24 
25    Level: beginner
26 
27 .seealso: `TSSetType()`, `TS`, `TSRegister()`
28 J*/
29 typedef const char* TSType;
30 #define TSEULER           "euler"
31 #define TSBEULER          "beuler"
32 #define TSBASICSYMPLECTIC "basicsymplectic"
33 #define TSPSEUDO          "pseudo"
34 #define TSCN              "cn"
35 #define TSSUNDIALS        "sundials"
36 #define TSRK              "rk"
37 #define TSPYTHON          "python"
38 #define TSTHETA           "theta"
39 #define TSALPHA           "alpha"
40 #define TSALPHA2          "alpha2"
41 #define TSGLLE            "glle"
42 #define TSGLEE            "glee"
43 #define TSSSP             "ssp"
44 #define TSARKIMEX         "arkimex"
45 #define TSROSW            "rosw"
46 #define TSEIMEX           "eimex"
47 #define TSMIMEX           "mimex"
48 #define TSBDF             "bdf"
49 #define TSRADAU5          "radau5"
50 #define TSMPRK            "mprk"
51 #define TSDISCGRAD        "discgrad"
52 #define TSIRK             "irk"
53 
54 /*E
55     TSProblemType - Determines the type of problem this TS object is to be used to solve
56 
57    Level: beginner
58 
59 .seealso: `TSCreate()`
60 E*/
61 typedef enum {TS_LINEAR,TS_NONLINEAR} TSProblemType;
62 
63 /*E
64    TSEquationType - type of TS problem that is solved
65 
66    Level: beginner
67 
68    Developer Notes:
69     this must match petsc/finclude/petscts.h
70 
71    Supported types are:
72      TS_EQ_UNSPECIFIED (default)
73      TS_EQ_EXPLICIT {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) := M(t) U_t - G(U,t) = 0
74      TS_EQ_IMPLICIT {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) = 0
75 
76 .seealso: `TSGetEquationType()`, `TSSetEquationType()`
77 E*/
78 typedef enum {
79   TS_EQ_UNSPECIFIED               = -1,
80   TS_EQ_EXPLICIT                  = 0,
81   TS_EQ_ODE_EXPLICIT              = 1,
82   TS_EQ_DAE_SEMI_EXPLICIT_INDEX1  = 100,
83   TS_EQ_DAE_SEMI_EXPLICIT_INDEX2  = 200,
84   TS_EQ_DAE_SEMI_EXPLICIT_INDEX3  = 300,
85   TS_EQ_DAE_SEMI_EXPLICIT_INDEXHI = 500,
86   TS_EQ_IMPLICIT                  = 1000,
87   TS_EQ_ODE_IMPLICIT              = 1001,
88   TS_EQ_DAE_IMPLICIT_INDEX1       = 1100,
89   TS_EQ_DAE_IMPLICIT_INDEX2       = 1200,
90   TS_EQ_DAE_IMPLICIT_INDEX3       = 1300,
91   TS_EQ_DAE_IMPLICIT_INDEXHI      = 1500
92 } TSEquationType;
93 PETSC_EXTERN const char *const*TSEquationTypes;
94 
95 /*E
96    TSConvergedReason - reason a TS method has converged or not
97 
98    Level: beginner
99 
100    Developer Notes:
101     this must match petsc/finclude/petscts.h
102 
103    Each reason has its own manual page.
104 
105 .seealso: `TSGetConvergedReason()`
106 E*/
107 typedef enum {
108   TS_CONVERGED_ITERATING          = 0,
109   TS_CONVERGED_TIME               = 1,
110   TS_CONVERGED_ITS                = 2,
111   TS_CONVERGED_USER               = 3,
112   TS_CONVERGED_EVENT              = 4,
113   TS_CONVERGED_PSEUDO_FATOL       = 5,
114   TS_CONVERGED_PSEUDO_FRTOL       = 6,
115   TS_DIVERGED_NONLINEAR_SOLVE     = -1,
116   TS_DIVERGED_STEP_REJECTED       = -2,
117   TSFORWARD_DIVERGED_LINEAR_SOLVE = -3,
118   TSADJOINT_DIVERGED_LINEAR_SOLVE = -4
119 } TSConvergedReason;
120 PETSC_EXTERN const char *const*TSConvergedReasons;
121 
122 /*MC
123    TS_CONVERGED_ITERATING - this only occurs if TSGetConvergedReason() is called during the TSSolve()
124 
125    Level: beginner
126 
127 .seealso: `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`
128 M*/
129 
130 /*MC
131    TS_CONVERGED_TIME - the final time was reached
132 
133    Level: beginner
134 
135 .seealso: `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxTime()`, `TSGetMaxTime()`, `TSGetSolveTime()`
136 M*/
137 
138 /*MC
139    TS_CONVERGED_ITS - the maximum number of iterations (time-steps) was reached prior to the final time
140 
141    Level: beginner
142 
143 .seealso: `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxSteps()`, `TSGetMaxSteps()`
144 M*/
145 
146 /*MC
147    TS_CONVERGED_USER - user requested termination
148 
149    Level: beginner
150 
151 .seealso: `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`
152 M*/
153 
154 /*MC
155    TS_CONVERGED_EVENT - user requested termination on event detection
156 
157    Level: beginner
158 
159 .seealso: `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`
160 M*/
161 
162 /*MC
163    TS_CONVERGED_PSEUDO_FRTOL - stops when function norm decreased by a set amount, used only for TSPSEUDO
164 
165    Level: beginner
166 
167    Options Database:
168 .   -ts_pseudo_frtol <rtol> - use sepcified rtol
169 
170 .seealso: `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`, `TS_CONVERGED_PSEUDO_FATOL`
171 M*/
172 
173 /*MC
174    TS_CONVERGED_PSEUDO_FATOL - stops when function norm decreases below a set amount, used only for TSPSEUDO
175 
176    Level: beginner
177 
178    Options Database:
179 .   -ts_pseudo_fatol <atol> - use specified atol
180 
181 .seealso: `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`, `TS_CONVERGED_PSEUDO_FRTOL`
182 M*/
183 
184 /*MC
185    TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed
186 
187    Level: beginner
188 
189    Notes:
190     See TSSetMaxSNESFailures() for how to allow more nonlinear solver failures.
191 
192 .seealso: `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSGetSNES()`, `SNESGetConvergedReason()`, `TSSetMaxSNESFailures()`
193 M*/
194 
195 /*MC
196    TS_DIVERGED_STEP_REJECTED - too many steps were rejected
197 
198    Level: beginner
199 
200    Notes:
201     See TSSetMaxStepRejections() for how to allow more step rejections.
202 
203 .seealso: `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxStepRejections()`
204 M*/
205 
206 /*E
207    TSExactFinalTimeOption - option for handling of final time step
208 
209    Level: beginner
210 
211    Developer Notes:
212     this must match petsc/finclude/petscts.h
213 
214 $  TS_EXACTFINALTIME_STEPOVER    - Don't do anything if final time is exceeded
215 $  TS_EXACTFINALTIME_INTERPOLATE - Interpolate back to final time
216 $  TS_EXACTFINALTIME_MATCHSTEP - Adapt final time step to match the final time
217 
218 .seealso: `TSGetConvergedReason()`, `TSSetExactFinalTime()`, `TSGetExactFinalTime()`
219 
220 E*/
221 typedef enum {TS_EXACTFINALTIME_UNSPECIFIED=0,TS_EXACTFINALTIME_STEPOVER=1,TS_EXACTFINALTIME_INTERPOLATE=2,TS_EXACTFINALTIME_MATCHSTEP=3} TSExactFinalTimeOption;
222 PETSC_EXTERN const char *const TSExactFinalTimeOptions[];
223 
224 /* Logging support */
225 PETSC_EXTERN PetscClassId TS_CLASSID;
226 PETSC_EXTERN PetscClassId DMTS_CLASSID;
227 PETSC_EXTERN PetscClassId TSADAPT_CLASSID;
228 
229 PETSC_EXTERN PetscErrorCode TSInitializePackage(void);
230 PETSC_EXTERN PetscErrorCode TSFinalizePackage(void);
231 
232 PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm,TS*);
233 PETSC_EXTERN PetscErrorCode TSClone(TS,TS*);
234 PETSC_EXTERN PetscErrorCode TSDestroy(TS*);
235 
236 PETSC_EXTERN PetscErrorCode TSSetProblemType(TS,TSProblemType);
237 PETSC_EXTERN PetscErrorCode TSGetProblemType(TS,TSProblemType*);
238 PETSC_EXTERN PetscErrorCode TSMonitor(TS,PetscInt,PetscReal,Vec);
239 PETSC_EXTERN PetscErrorCode TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**));
240 PETSC_EXTERN PetscErrorCode TSMonitorSetFromOptions(TS,const char[],const char[],const char[],PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*),PetscErrorCode (*)(TS,PetscViewerAndFormat*));
241 PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS);
242 
243 PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS,const char[]);
244 PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS,const char[]);
245 PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS,const char *[]);
246 PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS);
247 PETSC_EXTERN PetscErrorCode TSSetUp(TS);
248 PETSC_EXTERN PetscErrorCode TSReset(TS);
249 
250 PETSC_EXTERN PetscErrorCode TSSetSolution(TS,Vec);
251 PETSC_EXTERN PetscErrorCode TSGetSolution(TS,Vec*);
252 
253 PETSC_EXTERN PetscErrorCode TS2SetSolution(TS,Vec,Vec);
254 PETSC_EXTERN PetscErrorCode TS2GetSolution(TS,Vec*,Vec*);
255 
256 PETSC_EXTERN PetscErrorCode TSGetSolutionComponents(TS,PetscInt*,Vec*);
257 PETSC_EXTERN PetscErrorCode TSGetAuxSolution(TS,Vec*);
258 PETSC_EXTERN PetscErrorCode TSGetTimeError(TS,PetscInt,Vec*);
259 PETSC_EXTERN PetscErrorCode TSSetTimeError(TS,Vec);
260 
261 PETSC_EXTERN PetscErrorCode TSSetRHSJacobianP(TS,Mat,PetscErrorCode(*)(TS,PetscReal,Vec,Mat,void*),void*);
262 PETSC_EXTERN PetscErrorCode TSGetRHSJacobianP(TS,Mat*,PetscErrorCode(**)(TS,PetscReal,Vec,Mat,void*),void**);
263 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianP(TS,PetscReal,Vec,Mat);
264 PETSC_EXTERN PetscErrorCode TSSetIJacobianP(TS,Mat,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal,Mat,void*),void*);
265 PETSC_EXTERN PetscErrorCode TSComputeIJacobianP(TS,PetscReal,Vec,Vec,PetscReal,Mat,PetscBool);
266 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetQuadratureTS then TSComputeRHSJacobianP") PetscErrorCode TSComputeDRDPFunction(TS,PetscReal,Vec,Vec*);
267 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetQuadratureTS then TSComputeRHSJacobian") PetscErrorCode TSComputeDRDUFunction(TS,PetscReal,Vec,Vec*);
268 PETSC_EXTERN PetscErrorCode TSSetIHessianProduct(TS,Vec*,PetscErrorCode(*)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
269                                                     Vec*,PetscErrorCode(*)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
270                                                     Vec*,PetscErrorCode(*)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
271                                                     Vec*,PetscErrorCode(*)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
272                                                     void*);
273 PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUU(TS,PetscReal,Vec,Vec*,Vec,Vec*);
274 PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUP(TS,PetscReal,Vec,Vec*,Vec,Vec*);
275 PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPU(TS,PetscReal,Vec,Vec*,Vec,Vec*);
276 PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPP(TS,PetscReal,Vec,Vec*,Vec,Vec*);
277 PETSC_EXTERN PetscErrorCode TSSetRHSHessianProduct(TS,Vec*,PetscErrorCode(*)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
278                                                     Vec*,PetscErrorCode(*)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
279                                                     Vec*,PetscErrorCode(*)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
280                                                     Vec*,PetscErrorCode(*)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
281                                                     void*);
282 PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS,PetscReal,Vec,Vec*,Vec,Vec*);
283 PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS,PetscReal,Vec,Vec*,Vec,Vec*);
284 PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS,PetscReal,Vec,Vec*,Vec,Vec*);
285 PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS,PetscReal,Vec,Vec*,Vec,Vec*);
286 PETSC_EXTERN PetscErrorCode TSSetCostHessianProducts(TS,PetscInt,Vec*,Vec*,Vec);
287 PETSC_EXTERN PetscErrorCode TSGetCostHessianProducts(TS,PetscInt*,Vec**,Vec**,Vec*);
288 PETSC_EXTERN PetscErrorCode TSComputeSNESJacobian(TS,Vec,Mat,Mat);
289 
290 /*S
291      TSTrajectory - Abstract PETSc object that stores the trajectory (solution of ODE/DAE at each time step)
292 
293    Level: advanced
294 
295 .seealso: `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectorySetType()`, `TSTrajectoryDestroy()`, `TSTrajectoryReset()`
296 S*/
297 typedef struct _p_TSTrajectory* TSTrajectory;
298 
299 /*J
300     TSTrajectorySetType - String with the name of a PETSc TS trajectory storage method
301 
302    Level: intermediate
303 
304 .seealso: `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectoryDestroy()`
305 J*/
306 typedef const char* TSTrajectoryType;
307 #define TSTRAJECTORYBASIC         "basic"
308 #define TSTRAJECTORYSINGLEFILE    "singlefile"
309 #define TSTRAJECTORYMEMORY        "memory"
310 #define TSTRAJECTORYVISUALIZATION "visualization"
311 
312 PETSC_EXTERN PetscFunctionList TSTrajectoryList;
313 PETSC_EXTERN PetscClassId      TSTRAJECTORY_CLASSID;
314 PETSC_EXTERN PetscBool         TSTrajectoryRegisterAllCalled;
315 
316 PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS);
317 PETSC_EXTERN PetscErrorCode TSResetTrajectory(TS);
318 PETSC_EXTERN PetscErrorCode TSRemoveTrajectory(TS);
319 
320 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate(MPI_Comm,TSTrajectory*);
321 PETSC_EXTERN PetscErrorCode TSTrajectoryReset(TSTrajectory);
322 PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy(TSTrajectory*);
323 PETSC_EXTERN PetscErrorCode TSTrajectoryView(TSTrajectory,PetscViewer);
324 PETSC_EXTERN PetscErrorCode TSTrajectorySetType(TSTrajectory,TS,TSTrajectoryType);
325 PETSC_EXTERN PetscErrorCode TSTrajectoryGetType(TSTrajectory,TS,TSTrajectoryType*);
326 PETSC_EXTERN PetscErrorCode TSTrajectorySet(TSTrajectory,TS,PetscInt,PetscReal,Vec);
327 PETSC_EXTERN PetscErrorCode TSTrajectoryGet(TSTrajectory,TS,PetscInt,PetscReal*);
328 PETSC_EXTERN PetscErrorCode TSTrajectoryGetVecs(TSTrajectory,TS,PetscInt,PetscReal*,Vec,Vec);
329 PETSC_EXTERN PetscErrorCode TSTrajectoryGetUpdatedHistoryVecs(TSTrajectory,TS,PetscReal,Vec*,Vec*);
330 PETSC_EXTERN PetscErrorCode TSTrajectoryGetNumSteps(TSTrajectory,PetscInt*);
331 PETSC_EXTERN PetscErrorCode TSTrajectoryRestoreUpdatedHistoryVecs(TSTrajectory,Vec*,Vec*);
332 PETSC_EXTERN PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory,TS);
333 PETSC_EXTERN PetscErrorCode TSTrajectoryRegister(const char[],PetscErrorCode (*)(TSTrajectory,TS));
334 PETSC_EXTERN PetscErrorCode TSTrajectoryRegisterAll(void);
335 PETSC_EXTERN PetscErrorCode TSTrajectorySetUp(TSTrajectory,TS);
336 PETSC_EXTERN PetscErrorCode TSTrajectorySetUseHistory(TSTrajectory,PetscBool);
337 PETSC_EXTERN PetscErrorCode TSTrajectorySetMonitor(TSTrajectory,PetscBool);
338 PETSC_EXTERN PetscErrorCode TSTrajectorySetVariableNames(TSTrajectory,const char * const*);
339 PETSC_EXTERN PetscErrorCode TSTrajectorySetTransform(TSTrajectory,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void*);
340 PETSC_EXTERN PetscErrorCode TSTrajectorySetSolutionOnly(TSTrajectory,PetscBool);
341 PETSC_EXTERN PetscErrorCode TSTrajectoryGetSolutionOnly(TSTrajectory,PetscBool*);
342 PETSC_EXTERN PetscErrorCode TSTrajectorySetKeepFiles(TSTrajectory,PetscBool);
343 PETSC_EXTERN PetscErrorCode TSTrajectorySetDirname(TSTrajectory,const char[]);
344 PETSC_EXTERN PetscErrorCode TSTrajectorySetFiletemplate(TSTrajectory,const char[]);
345 PETSC_EXTERN PetscErrorCode TSGetTrajectory(TS,TSTrajectory*);
346 
347 PETSC_EXTERN PetscErrorCode TSSetCostGradients(TS,PetscInt,Vec*,Vec*);
348 PETSC_EXTERN PetscErrorCode TSGetCostGradients(TS,PetscInt*,Vec**,Vec**);
349 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*);
350 PETSC_EXTERN PetscErrorCode TSGetCostIntegral(TS,Vec*);
351 PETSC_EXTERN PetscErrorCode TSComputeCostIntegrand(TS,PetscReal,Vec,Vec);
352 PETSC_EXTERN PetscErrorCode TSCreateQuadratureTS(TS,PetscBool,TS*);
353 PETSC_EXTERN PetscErrorCode TSGetQuadratureTS(TS,PetscBool*,TS*);
354 
355 PETSC_EXTERN PetscErrorCode TSAdjointSetFromOptions(PetscOptionItems*,TS);
356 PETSC_EXTERN PetscErrorCode TSAdjointMonitor(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*);
357 PETSC_EXTERN PetscErrorCode TSAdjointMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *,PetscErrorCode (*)(void**));
358 PETSC_EXTERN PetscErrorCode TSAdjointMonitorCancel(TS);
359 PETSC_EXTERN PetscErrorCode TSAdjointMonitorSetFromOptions(TS,const char[],const char[],const char[],PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,PetscViewerAndFormat*),PetscErrorCode (*)(TS,PetscViewerAndFormat*));
360 
361 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSSetRHSJacobianP()")     PetscErrorCode TSAdjointSetRHSJacobian(TS,Mat,PetscErrorCode(*)(TS,PetscReal,Vec,Mat,void*),void*);
362 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSComputeRHSJacobianP()") PetscErrorCode TSAdjointComputeRHSJacobian(TS,PetscReal,Vec,Mat);
363 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetQuadratureTS then TSComputeRHSJacobianP") PetscErrorCode TSAdjointComputeDRDPFunction(TS,PetscReal,Vec,Vec*);
364 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetQuadratureTS then TSComputeRHSJacobian") PetscErrorCode TSAdjointComputeDRDYFunction(TS,PetscReal,Vec,Vec*);
365 PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS);
366 PETSC_EXTERN PetscErrorCode TSAdjointSetSteps(TS,PetscInt);
367 
368 PETSC_EXTERN PetscErrorCode TSAdjointStep(TS);
369 PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS);
370 PETSC_EXTERN PetscErrorCode TSAdjointReset(TS);
371 PETSC_EXTERN PetscErrorCode TSAdjointCostIntegral(TS);
372 PETSC_EXTERN PetscErrorCode TSAdjointSetForward(TS,Mat);
373 PETSC_EXTERN PetscErrorCode TSAdjointResetForward(TS);
374 
375 PETSC_EXTERN PetscErrorCode TSForwardSetSensitivities(TS,PetscInt,Mat);
376 PETSC_EXTERN PetscErrorCode TSForwardGetSensitivities(TS,PetscInt*,Mat*);
377 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSCreateQuadratureTS() and TSForwardSetSensitivities() (since version 3.12)") PetscErrorCode TSForwardSetIntegralGradients(TS,PetscInt,Vec *);
378 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSForwardGetSensitivities()")                          PetscErrorCode TSForwardGetIntegralGradients(TS,PetscInt*,Vec **);
379 PETSC_EXTERN PetscErrorCode TSForwardSetUp(TS);
380 PETSC_EXTERN PetscErrorCode TSForwardReset(TS);
381 PETSC_EXTERN PetscErrorCode TSForwardCostIntegral(TS);
382 PETSC_EXTERN PetscErrorCode TSForwardStep(TS);
383 PETSC_EXTERN PetscErrorCode TSForwardSetInitialSensitivities(TS,Mat);
384 PETSC_EXTERN PetscErrorCode TSForwardGetStages(TS,PetscInt*,Mat*[]);
385 
386 PETSC_EXTERN PetscErrorCode TSSetMaxSteps(TS,PetscInt);
387 PETSC_EXTERN PetscErrorCode TSGetMaxSteps(TS,PetscInt*);
388 PETSC_EXTERN PetscErrorCode TSSetMaxTime(TS,PetscReal);
389 PETSC_EXTERN PetscErrorCode TSGetMaxTime(TS,PetscReal*);
390 PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS,TSExactFinalTimeOption);
391 PETSC_EXTERN PetscErrorCode TSGetExactFinalTime(TS,TSExactFinalTimeOption*);
392 PETSC_EXTERN PetscErrorCode TSSetTimeSpan(TS,PetscInt,PetscReal*);
393 PETSC_EXTERN PetscErrorCode TSGetTimeSpan(TS,PetscInt*,const PetscReal**);
394 PETSC_EXTERN PetscErrorCode TSGetTimeSpanSolutions(TS,PetscInt*,Vec**);
395 
396 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSSetTime[Step]() (since version 3.8)")      PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal);
397 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSSetMax{Steps|Time}() (since version 3.8)") PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal);
398 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetMax{Steps|Time}() (since version 3.8)") PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*);
399 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetStepNumber() (since version 3.8)")      PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*);
400 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use TSGetStepNumber() (since version 3.8)")      PetscErrorCode TSGetTotalSteps(TS,PetscInt*);
401 
402 PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*);
403 PETSC_EXTERN PetscErrorCode TSMonitorExtreme(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*);
404 
405 typedef struct _n_TSMonitorDrawCtx*  TSMonitorDrawCtx;
406 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorDrawCtx *);
407 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx*);
408 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS,PetscInt,PetscReal,Vec,void*);
409 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionPhase(TS,PetscInt,PetscReal,Vec,void*);
410 PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS,PetscInt,PetscReal,Vec,void*);
411 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionFunction(TS,PetscInt,PetscReal,Vec,void*);
412 
413 PETSC_EXTERN PetscErrorCode TSAdjointMonitorDefault(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,PetscViewerAndFormat *);
414 PETSC_EXTERN PetscErrorCode TSAdjointMonitorDrawSensi(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*);
415 
416 PETSC_EXTERN PetscErrorCode TSMonitorSolution(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*);
417 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*);
418 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void*);
419 
420 PETSC_EXTERN PetscErrorCode TSStep(TS);
421 PETSC_EXTERN PetscErrorCode TSEvaluateWLTE(TS,NormType,PetscInt*,PetscReal*);
422 PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS,PetscInt,Vec,PetscBool*);
423 PETSC_EXTERN PetscErrorCode TSSolve(TS,Vec);
424 PETSC_EXTERN PetscErrorCode TSGetEquationType(TS,TSEquationType*);
425 PETSC_EXTERN PetscErrorCode TSSetEquationType(TS,TSEquationType);
426 PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*);
427 PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS,TSConvergedReason);
428 PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS,PetscReal*);
429 PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS,PetscInt*);
430 PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS,PetscInt*);
431 PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS,PetscInt*);
432 PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS,PetscInt);
433 PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS,PetscInt*);
434 PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS,PetscInt);
435 PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS,PetscBool);
436 PETSC_EXTERN PetscErrorCode TSRestartStep(TS);
437 PETSC_EXTERN PetscErrorCode TSRollBack(TS);
438 
439 PETSC_EXTERN PetscErrorCode TSGetStages(TS,PetscInt*,Vec*[]);
440 
441 PETSC_EXTERN PetscErrorCode TSGetTime(TS,PetscReal*);
442 PETSC_EXTERN PetscErrorCode TSSetTime(TS,PetscReal);
443 PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS,PetscReal*);
444 PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS,PetscReal*);
445 PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS,PetscReal);
446 PETSC_EXTERN PetscErrorCode TSGetStepNumber(TS,PetscInt*);
447 PETSC_EXTERN PetscErrorCode TSSetStepNumber(TS,PetscInt);
448 
449 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*);
450 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat,Mat,void*);
451 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobianP)(TS,PetscReal,Vec,Mat,void*);
452 PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*);
453 PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**);
454 PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*);
455 PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**);
456 PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS,PetscBool);
457 
458 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS,PetscReal,Vec,void*);
459 PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS,TSSolutionFunction,void*);
460 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSForcingFunction)(TS,PetscReal,Vec,void*);
461 PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS,TSForcingFunction,void*);
462 
463 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*);
464 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
465 PETSC_EXTERN PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*);
466 PETSC_EXTERN PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**);
467 PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*);
468 PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**);
469 
470 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Function)(TS,PetscReal,Vec,Vec,Vec,Vec,void*);
471 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Jacobian)(TS,PetscReal,Vec,Vec,Vec,PetscReal,PetscReal,Mat,Mat,void*);
472 PETSC_EXTERN PetscErrorCode TSSetI2Function(TS,Vec,TSI2Function,void*);
473 PETSC_EXTERN PetscErrorCode TSGetI2Function(TS,Vec*,TSI2Function*,void**);
474 PETSC_EXTERN PetscErrorCode TSSetI2Jacobian(TS,Mat,Mat,TSI2Jacobian,void*);
475 PETSC_EXTERN PetscErrorCode TSGetI2Jacobian(TS,Mat*,Mat*,TSI2Jacobian*,void**);
476 
477 PETSC_EXTERN PetscErrorCode TSRHSSplitSetIS(TS,const char[],IS);
478 PETSC_EXTERN PetscErrorCode TSRHSSplitGetIS(TS,const char[],IS*);
479 PETSC_EXTERN PetscErrorCode TSRHSSplitSetRHSFunction(TS,const char[],Vec,TSRHSFunction,void*);
480 PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTS(TS,const char[],TS*);
481 PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTSs(TS,PetscInt*,TS*[]);
482 PETSC_EXTERN PetscErrorCode TSSetUseSplitRHSFunction(TS, PetscBool);
483 PETSC_EXTERN PetscErrorCode TSGetUseSplitRHSFunction(TS, PetscBool*);
484 
485 PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*);
486 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat,Mat,void*);
487 PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*);
488 PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
489 PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS,PetscReal,Vec);
490 PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS,PetscReal,Vec);
491 PETSC_EXTERN PetscErrorCode TSComputeIJacobianDefaultColor(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
492 
493 PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
494 PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS,PetscReal));
495 PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS,PetscReal,PetscInt,Vec*));
496 PETSC_EXTERN PetscErrorCode TSSetPostEvaluate(TS, PetscErrorCode(*)(TS));
497 PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
498 PETSC_EXTERN PetscErrorCode TSPreStep(TS);
499 PETSC_EXTERN PetscErrorCode TSPreStage(TS,PetscReal);
500 PETSC_EXTERN PetscErrorCode TSPostStage(TS,PetscReal,PetscInt,Vec*);
501 PETSC_EXTERN PetscErrorCode TSPostEvaluate(TS);
502 PETSC_EXTERN PetscErrorCode TSPostStep(TS);
503 PETSC_EXTERN PetscErrorCode TSInterpolate(TS,PetscReal,Vec);
504 PETSC_EXTERN PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec);
505 PETSC_EXTERN PetscErrorCode TSGetTolerances(TS,PetscReal*,Vec*,PetscReal*,Vec*);
506 PETSC_EXTERN PetscErrorCode TSErrorWeightedNormInfinity(TS,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
507 PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm2(TS,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
508 PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS,Vec,Vec,NormType,PetscReal*,PetscReal*,PetscReal*);
509 PETSC_EXTERN PetscErrorCode TSErrorWeightedENormInfinity(TS,Vec,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
510 PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm2(TS,Vec,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
511 PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm(TS,Vec,Vec,Vec,NormType,PetscReal*,PetscReal*,PetscReal*);
512 PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal);
513 PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS,PetscReal*);
514 PETSC_EXTERN PetscErrorCode TSSetFunctionDomainError(TS, PetscErrorCode (*)(TS,PetscReal,Vec,PetscBool*));
515 PETSC_EXTERN PetscErrorCode TSFunctionDomainError(TS,PetscReal,Vec,PetscBool*);
516 
517 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*);
518 PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS,PetscReal*,void*);
519 PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *);
520 PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal);
521 PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*);
522 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS,Vec,void*,PetscReal*,PetscBool *);
523 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *);
524 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal);
525 PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);
526 
527 PETSC_EXTERN PetscErrorCode TSPythonSetType(TS,const char[]);
528 
529 PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec);
530 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat,Mat);
531 PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool);
532 PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,PetscBool);
533 PETSC_EXTERN PetscErrorCode TSComputeI2Function(TS,PetscReal,Vec,Vec,Vec,Vec);
534 PETSC_EXTERN PetscErrorCode TSComputeI2Jacobian(TS,PetscReal,Vec,Vec,Vec,PetscReal,PetscReal,Mat,Mat);
535 PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS,PetscReal,PetscReal,PetscReal*,PetscReal*);
536 
537 PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS,Vec,Vec);
538 
539 PETSC_EXTERN PetscErrorCode DMTSSetBoundaryLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *);
540 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM,TSRHSFunction,void*);
541 PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM,TSRHSFunction*,void**);
542 PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM,TSRHSJacobian,void*);
543 PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM,TSRHSJacobian*,void**);
544 PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM,TSIFunction,void*);
545 PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM,TSIFunction*,void**);
546 PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM,TSIJacobian,void*);
547 PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM,TSIJacobian*,void**);
548 PETSC_EXTERN PetscErrorCode DMTSSetI2Function(DM,TSI2Function,void*);
549 PETSC_EXTERN PetscErrorCode DMTSGetI2Function(DM,TSI2Function*,void**);
550 PETSC_EXTERN PetscErrorCode DMTSSetI2Jacobian(DM,TSI2Jacobian,void*);
551 PETSC_EXTERN PetscErrorCode DMTSGetI2Jacobian(DM,TSI2Jacobian*,void**);
552 
553 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSTransientVariable)(TS,Vec,Vec,void*);
554 PETSC_EXTERN PetscErrorCode TSSetTransientVariable(TS,TSTransientVariable,void*);
555 PETSC_EXTERN PetscErrorCode DMTSSetTransientVariable(DM,TSTransientVariable,void*);
556 PETSC_EXTERN PetscErrorCode DMTSGetTransientVariable(DM,TSTransientVariable*,void*);
557 PETSC_EXTERN PetscErrorCode TSComputeTransientVariable(TS,Vec,Vec);
558 PETSC_EXTERN PetscErrorCode TSHasTransientVariable(TS,PetscBool*);
559 
560 PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM,TSSolutionFunction,void*);
561 PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM,TSSolutionFunction*,void**);
562 PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM,TSForcingFunction,void*);
563 PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM,TSForcingFunction*,void**);
564 PETSC_EXTERN PetscErrorCode DMTSCheckResidual(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscReal *);
565 PETSC_EXTERN PetscErrorCode DMTSCheckJacobian(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscBool *, PetscReal *);
566 PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec);
567 
568 PETSC_EXTERN PetscErrorCode DMTSGetIFunctionLocal(DM,PetscErrorCode (**)(DM,PetscReal,Vec,Vec,Vec,void*),void**);
569 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,Vec,void*),void*);
570 PETSC_EXTERN PetscErrorCode DMTSGetIJacobianLocal(DM,PetscErrorCode (**)(DM,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*),void**);
571 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*),void*);
572 PETSC_EXTERN PetscErrorCode DMTSGetRHSFunctionLocal(DM,PetscErrorCode (**)(DM,PetscReal,Vec,Vec,void*),void**);
573 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,void*),void*);
574 PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrix(DM);
575 PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrixLumped(DM);
576 PETSC_EXTERN PetscErrorCode DMTSDestroyRHSMassMatrix(DM);
577 
578 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));
579 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));
580 
581 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*);
582 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*);
583 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*);
584 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*);
585 
586 PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*),void *);
587 PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*),void *);
588 PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*),void *);
589 PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*),void *);
590 
591 PETSC_EXTERN PetscErrorCode DMPlexTSGetGeometryFVM(DM,Vec*,Vec*,PetscReal*);
592 
593 typedef struct _n_TSMonitorLGCtx*  TSMonitorLGCtx;
594 typedef struct {
595   Vec            ray;
596   VecScatter     scatter;
597   PetscViewer    viewer;
598   TSMonitorLGCtx lgctx;
599 } TSMonitorDMDARayCtx;
600 PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void**);
601 PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS,PetscInt,PetscReal,Vec,void*);
602 PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS,PetscInt,PetscReal,Vec,void*);
603 
604 /* Dynamic creation and loading functions */
605 PETSC_EXTERN PetscFunctionList TSList;
606 PETSC_EXTERN PetscErrorCode TSGetType(TS,TSType*);
607 PETSC_EXTERN PetscErrorCode TSSetType(TS,TSType);
608 PETSC_EXTERN PetscErrorCode TSRegister(const char[], PetscErrorCode (*)(TS));
609 
610 PETSC_EXTERN PetscErrorCode TSGetSNES(TS,SNES*);
611 PETSC_EXTERN PetscErrorCode TSSetSNES(TS,SNES);
612 PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*);
613 
614 PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer);
615 PETSC_EXTERN PetscErrorCode TSLoad(TS,PetscViewer);
616 PETSC_EXTERN PetscErrorCode TSViewFromOptions(TS,PetscObject,const char[]);
617 PETSC_EXTERN PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory,PetscObject,const char[]);
618 
619 #define TS_FILE_CLASSID 1211225
620 
621 PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *);
622 PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *);
623 
624 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorLGCtx *);
625 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx*);
626 PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS,PetscInt,PetscReal,Vec,void *);
627 PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS,PetscInt,PetscReal,Vec,void *);
628 PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS,const char * const*);
629 PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS,const char *const **);
630 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx,const char * const *);
631 PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS,const char * const*);
632 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx,const char * const*);
633 PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void*);
634 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void *);
635 PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS,PetscInt,PetscReal,Vec,void *);
636 PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS,PetscInt,PetscReal,Vec,void *);
637 PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS,PetscInt,PetscReal,Vec,void *);
638 PETSC_EXTERN PetscErrorCode TSMonitorError(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat *);
639 PETSC_EXTERN PetscErrorCode TSDMSwarmMonitorMoments(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat *);
640 
641 struct _n_TSMonitorLGCtxNetwork {
642   PetscInt       nlg;
643   PetscDrawLG    *lg;
644   PetscBool      semilogy;
645   PetscInt       howoften;  /* when > 0 uses step % howoften, when negative only final solution plotted */
646 };
647 typedef struct _n_TSMonitorLGCtxNetwork*  TSMonitorLGCtxNetwork;
648 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxNetworkDestroy(TSMonitorLGCtxNetwork*);
649 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxNetworkCreate(TS,const char[],const char[],int,int,int,int,PetscInt,TSMonitorLGCtxNetwork*);
650 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxNetworkSolution(TS,PetscInt,PetscReal,Vec,void*);
651 
652 typedef struct _n_TSMonitorEnvelopeCtx*  TSMonitorEnvelopeCtx;
653 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxCreate(TS,TSMonitorEnvelopeCtx*);
654 PETSC_EXTERN PetscErrorCode TSMonitorEnvelope(TS,PetscInt,PetscReal,Vec,void*);
655 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeGetBounds(TS,Vec*,Vec*);
656 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx*);
657 
658 typedef struct _n_TSMonitorSPEigCtx*  TSMonitorSPEigCtx;
659 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPEigCtx *);
660 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx*);
661 PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS,PetscInt,PetscReal,Vec,void *);
662 
663 typedef struct _n_TSMonitorSPCtx* TSMonitorSPCtx;
664 PETSC_EXTERN PetscErrorCode TSMonitorSPCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,PetscInt,PetscBool,TSMonitorSPCtx*);
665 PETSC_EXTERN PetscErrorCode TSMonitorSPCtxDestroy(TSMonitorSPCtx*);
666 PETSC_EXTERN PetscErrorCode TSMonitorSPSwarmSolution(TS,PetscInt,PetscReal,Vec,void*);
667 
668 PETSC_EXTERN PetscErrorCode TSSetEventHandler(TS,PetscInt,PetscInt[],PetscBool[],PetscErrorCode (*)(TS,PetscReal,Vec,PetscScalar[],void*),PetscErrorCode (*)(TS,PetscInt,PetscInt[],PetscReal,Vec,PetscBool,void*),void*);
669 PETSC_EXTERN PetscErrorCode TSSetPostEventIntervalStep(TS,PetscReal);
670 PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS,PetscReal,PetscReal[]);
671 PETSC_EXTERN PetscErrorCode TSGetNumEvents(TS,PetscInt*);
672 
673 /*J
674    TSSSPType - string with the name of TSSSP scheme.
675 
676    Level: beginner
677 
678 .seealso: `TSSSPSetType()`, `TS`
679 J*/
680 typedef const char* TSSSPType;
681 #define TSSSPRKS2  "rks2"
682 #define TSSSPRKS3  "rks3"
683 #define TSSSPRK104 "rk104"
684 
685 PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,TSSSPType);
686 PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,TSSSPType*);
687 PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt);
688 PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*);
689 PETSC_EXTERN PetscErrorCode TSSSPInitializePackage(void);
690 PETSC_EXTERN PetscErrorCode TSSSPFinalizePackage(void);
691 PETSC_EXTERN PetscFunctionList TSSSPList;
692 
693 /*S
694    TSAdapt - Abstract object that manages time-step adaptivity
695 
696    Level: beginner
697 
698 .seealso: `TS`, `TSAdaptCreate()`, `TSAdaptType`
699 S*/
700 typedef struct _p_TSAdapt *TSAdapt;
701 
702 /*J
703     TSAdaptType - String with the name of TSAdapt scheme.
704 
705    Level: beginner
706 
707 .seealso: `TSAdaptSetType()`, `TS`
708 J*/
709 typedef const char *TSAdaptType;
710 #define TSADAPTNONE    "none"
711 #define TSADAPTBASIC   "basic"
712 #define TSADAPTDSP     "dsp"
713 #define TSADAPTCFL     "cfl"
714 #define TSADAPTGLEE    "glee"
715 #define TSADAPTHISTORY "history"
716 
717 PETSC_EXTERN PetscErrorCode TSGetAdapt(TS,TSAdapt*);
718 PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],PetscErrorCode (*)(TSAdapt));
719 PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void);
720 PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
721 PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*);
722 PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,TSAdaptType);
723 PETSC_EXTERN PetscErrorCode TSAdaptGetType(TSAdapt,TSAdaptType*);
724 PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]);
725 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
726 PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool);
727 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**);
728 PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*);
729 PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscReal,Vec,PetscBool*);
730 PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer);
731 PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt,PetscViewer);
732 PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(PetscOptionItems*,TSAdapt);
733 PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt);
734 PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*);
735 PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool);
736 PETSC_EXTERN PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt,PetscBool);
737 PETSC_EXTERN PetscErrorCode TSAdaptSetSafety(TSAdapt,PetscReal,PetscReal);
738 PETSC_EXTERN PetscErrorCode TSAdaptGetSafety(TSAdapt,PetscReal*,PetscReal*);
739 PETSC_EXTERN PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt,PetscReal);
740 PETSC_EXTERN PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt,PetscReal*);
741 PETSC_EXTERN PetscErrorCode TSAdaptSetClip(TSAdapt,PetscReal,PetscReal);
742 PETSC_EXTERN PetscErrorCode TSAdaptGetClip(TSAdapt,PetscReal*,PetscReal*);
743 PETSC_EXTERN PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt,PetscReal);
744 PETSC_EXTERN PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt,PetscReal*);
745 PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal);
746 PETSC_EXTERN PetscErrorCode TSAdaptGetStepLimits(TSAdapt,PetscReal*,PetscReal*);
747 PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt,PetscErrorCode(*)(TSAdapt,TS,PetscReal,Vec,PetscBool*));
748 PETSC_EXTERN PetscErrorCode TSAdaptHistorySetHistory(TSAdapt,PetscInt n,PetscReal hist[],PetscBool);
749 PETSC_EXTERN PetscErrorCode TSAdaptHistorySetTrajectory(TSAdapt,TSTrajectory,PetscBool);
750 PETSC_EXTERN PetscErrorCode TSAdaptHistoryGetStep(TSAdapt,PetscInt,PetscReal*,PetscReal*);
751 PETSC_EXTERN PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt,PetscInt);
752 PETSC_EXTERN PetscErrorCode TSAdaptDSPSetFilter(TSAdapt,const char *);
753 PETSC_EXTERN PetscErrorCode TSAdaptDSPSetPID(TSAdapt,PetscReal,PetscReal,PetscReal);
754 
755 /*S
756    TSGLLEAdapt - Abstract object that manages time-step adaptivity
757 
758    Level: beginner
759 
760    Developer Notes:
761    This functionality should be replaced by the TSAdapt.
762 
763 .seealso: `TSGLLE`, `TSGLLEAdaptCreate()`, `TSGLLEAdaptType`
764 S*/
765 typedef struct _p_TSGLLEAdapt *TSGLLEAdapt;
766 
767 /*J
768     TSGLLEAdaptType - String with the name of TSGLLEAdapt scheme
769 
770    Level: beginner
771 
772 .seealso: `TSGLLEAdaptSetType()`, `TS`
773 J*/
774 typedef const char *TSGLLEAdaptType;
775 #define TSGLLEADAPT_NONE "none"
776 #define TSGLLEADAPT_SIZE "size"
777 #define TSGLLEADAPT_BOTH "both"
778 
779 PETSC_EXTERN PetscErrorCode TSGLLEAdaptRegister(const char[],PetscErrorCode (*)(TSGLLEAdapt));
780 PETSC_EXTERN PetscErrorCode TSGLLEAdaptInitializePackage(void);
781 PETSC_EXTERN PetscErrorCode TSGLLEAdaptFinalizePackage(void);
782 PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate(MPI_Comm,TSGLLEAdapt*);
783 PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt,TSGLLEAdaptType);
784 PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt,const char[]);
785 PETSC_EXTERN PetscErrorCode TSGLLEAdaptChoose(TSGLLEAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
786 PETSC_EXTERN PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt,PetscViewer);
787 PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetFromOptions(PetscOptionItems*,TSGLLEAdapt);
788 PETSC_EXTERN PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt*);
789 
790 /*J
791     TSGLLEAcceptType - String with the name of TSGLLEAccept scheme
792 
793    Level: beginner
794 
795 .seealso: `TSGLLESetAcceptType()`, `TS`
796 J*/
797 typedef const char *TSGLLEAcceptType;
798 #define TSGLLEACCEPT_ALWAYS "always"
799 
800 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLLEAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *);
801 PETSC_EXTERN PetscErrorCode TSGLLEAcceptRegister(const char[],TSGLLEAcceptFunction);
802 
803 /*J
804   TSGLLEType - family of time integration method within the General Linear class
805 
806   Level: beginner
807 
808 .seealso: `TSGLLESetType()`, `TSGLLERegister()`
809 J*/
810 typedef const char* TSGLLEType;
811 #define TSGLLE_IRKS   "irks"
812 
813 PETSC_EXTERN PetscErrorCode TSGLLERegister(const char[],PetscErrorCode(*)(TS));
814 PETSC_EXTERN PetscErrorCode TSGLLEInitializePackage(void);
815 PETSC_EXTERN PetscErrorCode TSGLLEFinalizePackage(void);
816 PETSC_EXTERN PetscErrorCode TSGLLESetType(TS,TSGLLEType);
817 PETSC_EXTERN PetscErrorCode TSGLLEGetAdapt(TS,TSGLLEAdapt*);
818 PETSC_EXTERN PetscErrorCode TSGLLESetAcceptType(TS,TSGLLEAcceptType);
819 
820 /*J
821     TSEIMEXType - String with the name of an Extrapolated IMEX method.
822 
823    Level: beginner
824 
825 .seealso: `TSEIMEXSetType()`, `TS`, `TSEIMEX`, `TSEIMEXRegister()`
826 J*/
827 #define TSEIMEXType   char*
828 
829 PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts,PetscInt);
830 PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts,PetscInt,PetscInt);
831 PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS,PetscBool);
832 
833 /*J
834     TSRKType - String with the name of a Runge-Kutta method.
835 
836    Level: beginner
837 
838 .seealso: `TSRKSetType()`, `TS`, `TSRK`, `TSRKRegister()`
839 J*/
840 typedef const char* TSRKType;
841 #define TSRK1FE   "1fe"
842 #define TSRK2A    "2a"
843 #define TSRK2B    "2b"
844 #define TSRK3     "3"
845 #define TSRK3BS   "3bs"
846 #define TSRK4     "4"
847 #define TSRK5F    "5f"
848 #define TSRK5DP   "5dp"
849 #define TSRK5BS   "5bs"
850 #define TSRK6VR   "6vr"
851 #define TSRK7VR   "7vr"
852 #define TSRK8VR   "8vr"
853 
854 PETSC_EXTERN PetscErrorCode TSRKGetOrder(TS,PetscInt*);
855 PETSC_EXTERN PetscErrorCode TSRKGetType(TS,TSRKType*);
856 PETSC_EXTERN PetscErrorCode TSRKSetType(TS,TSRKType);
857 PETSC_EXTERN PetscErrorCode TSRKGetTableau(TS,PetscInt*,const PetscReal**,const PetscReal**,const PetscReal**,const PetscReal**,PetscInt*,const PetscReal**,PetscBool*);
858 PETSC_EXTERN PetscErrorCode TSRKSetMultirate(TS,PetscBool);
859 PETSC_EXTERN PetscErrorCode TSRKGetMultirate(TS,PetscBool*);
860 PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
861 PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void);
862 PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void);
863 PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void);
864 
865 /*J
866    TSMPRKType - String with the name of a Partitioned Runge-Kutta method
867 
868    Level: beginner
869 
870 .seealso: `TSMPRKSetType()`, `TS`, `TSMPRK`, `TSMPRKRegister()`
871 J*/
872 typedef const char* TSMPRKType;
873 #define TSMPRK2A22  "2a22"
874 #define TSMPRK2A23  "2a23"
875 #define TSMPRK2A32  "2a32"
876 #define TSMPRK2A33  "2a33"
877 #define TSMPRKP2    "p2"
878 #define TSMPRKP3    "p3"
879 
880 PETSC_EXTERN PetscErrorCode TSMPRKGetType(TS ts,TSMPRKType*);
881 PETSC_EXTERN PetscErrorCode TSMPRKSetType(TS ts,TSMPRKType);
882 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[]);
883 PETSC_EXTERN PetscErrorCode TSMPRKInitializePackage(void);
884 PETSC_EXTERN PetscErrorCode TSMPRKFinalizePackage(void);
885 PETSC_EXTERN PetscErrorCode TSMPRKRegisterDestroy(void);
886 
887 /*J
888     TSIRKType - String with the name of an implicit Runge-Kutta method.
889 
890    Level: beginner
891 
892 .seealso: `TSIRKSetType()`, `TS`, `TSIRK`, `TSIRKRegister()`
893 J*/
894 typedef const char* TSIRKType;
895 #define TSIRKGAUSS   "gauss"
896 
897 PETSC_EXTERN PetscErrorCode TSIRKGetOrder(TS,PetscInt*);
898 PETSC_EXTERN PetscErrorCode TSIRKGetType(TS,TSIRKType*);
899 PETSC_EXTERN PetscErrorCode TSIRKSetType(TS,TSIRKType);
900 PETSC_EXTERN PetscErrorCode TSIRKGetNumStages(TS,PetscInt*);
901 PETSC_EXTERN PetscErrorCode TSIRKSetNumStages(TS,PetscInt);
902 PETSC_EXTERN PetscErrorCode TSIRKRegister(const char[],PetscErrorCode (*function)(TS));
903 PETSC_EXTERN PetscErrorCode TSIRKTableauCreate(TS,PetscInt,const PetscReal*,const PetscReal*,const PetscReal*,const PetscReal*,const PetscScalar*,const PetscScalar*,const PetscScalar*);
904 PETSC_EXTERN PetscErrorCode TSIRKInitializePackage(void);
905 PETSC_EXTERN PetscErrorCode TSIRKFinalizePackage(void);
906 PETSC_EXTERN PetscErrorCode TSIRKRegisterDestroy(void);
907 
908 /*J
909     TSGLEEType - String with the name of a General Linear with Error Estimation method.
910 
911    Level: beginner
912 
913 .seealso: `TSGLEESetType()`, `TS`, `TSGLEE`, `TSGLEERegister()`
914 J*/
915 typedef const char* TSGLEEType;
916 #define TSGLEEi1      "BE1"
917 #define TSGLEE23      "23"
918 #define TSGLEE24      "24"
919 #define TSGLEE25I     "25i"
920 #define TSGLEE35      "35"
921 #define TSGLEEEXRK2A  "exrk2a"
922 #define TSGLEERK32G1  "rk32g1"
923 #define TSGLEERK285EX "rk285ex"
924 /*J
925     TSGLEEMode - String with the mode of error estimation for a General Linear with Error Estimation method.
926 
927    Level: beginner
928 
929 .seealso: `TSGLEESetMode()`, `TS`, `TSGLEE`, `TSGLEERegister()`
930 J*/
931 PETSC_EXTERN PetscErrorCode TSGLEEGetType(TS ts,TSGLEEType*);
932 PETSC_EXTERN PetscErrorCode TSGLEESetType(TS ts,TSGLEEType);
933 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[]);
934 PETSC_EXTERN PetscErrorCode TSGLEEFinalizePackage(void);
935 PETSC_EXTERN PetscErrorCode TSGLEEInitializePackage(void);
936 PETSC_EXTERN PetscErrorCode TSGLEERegisterDestroy(void);
937 
938 /*J
939     TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX method.
940 
941    Level: beginner
942 
943 .seealso: `TSARKIMEXSetType()`, `TS`, `TSARKIMEX`, `TSARKIMEXRegister()`
944 J*/
945 typedef const char* TSARKIMEXType;
946 #define TSARKIMEX1BEE   "1bee"
947 #define TSARKIMEXA2     "a2"
948 #define TSARKIMEXL2     "l2"
949 #define TSARKIMEXARS122 "ars122"
950 #define TSARKIMEX2C     "2c"
951 #define TSARKIMEX2D     "2d"
952 #define TSARKIMEX2E     "2e"
953 #define TSARKIMEXPRSSP2 "prssp2"
954 #define TSARKIMEX3      "3"
955 #define TSARKIMEXBPR3   "bpr3"
956 #define TSARKIMEXARS443 "ars443"
957 #define TSARKIMEX4      "4"
958 #define TSARKIMEX5      "5"
959 PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType*);
960 PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType);
961 PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool);
962 PETSC_EXTERN PetscErrorCode TSARKIMEXGetFullyImplicit(TS,PetscBool*);
963 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[]);
964 PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void);
965 PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
966 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);
967 
968 /*J
969     TSRosWType - String with the name of a Rosenbrock-W method.
970 
971    Level: beginner
972 
973 .seealso: `TSRosWSetType()`, `TS`, `TSROSW`, `TSRosWRegister()`
974 J*/
975 typedef const char* TSRosWType;
976 #define TSROSW2M          "2m"
977 #define TSROSW2P          "2p"
978 #define TSROSWRA3PW       "ra3pw"
979 #define TSROSWRA34PW2     "ra34pw2"
980 #define TSROSWRODAS3      "rodas3"
981 #define TSROSWSANDU3      "sandu3"
982 #define TSROSWASSP3P3S1C  "assp3p3s1c"
983 #define TSROSWLASSP3P4S2C "lassp3p4s2c"
984 #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
985 #define TSROSWARK3        "ark3"
986 #define TSROSWTHETA1      "theta1"
987 #define TSROSWTHETA2      "theta2"
988 #define TSROSWGRK4T       "grk4t"
989 #define TSROSWSHAMP4      "shamp4"
990 #define TSROSWVELDD4      "veldd4"
991 #define TSROSW4L          "4l"
992 
993 PETSC_EXTERN PetscErrorCode TSRosWGetType(TS,TSRosWType*);
994 PETSC_EXTERN PetscErrorCode TSRosWSetType(TS,TSRosWType);
995 PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool);
996 PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
997 PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
998 PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void);
999 PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
1000 PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);
1001 
1002 PETSC_EXTERN PetscErrorCode TSBDFSetOrder(TS,PetscInt);
1003 PETSC_EXTERN PetscErrorCode TSBDFGetOrder(TS,PetscInt*);
1004 
1005 /*J
1006   TSBasicSymplecticType - String with the name of a basic symplectic integration method.
1007 
1008   Level: beginner
1009 
1010   .seealso: `TSBasicSymplecticSetType()`, `TS`, `TSBASICSYMPLECTIC`, `TSBasicSymplecticRegister()`
1011 J*/
1012 typedef const char* TSBasicSymplecticType;
1013 #define TSBASICSYMPLECTICSIEULER   "1"
1014 #define TSBASICSYMPLECTICVELVERLET "2"
1015 #define TSBASICSYMPLECTIC3         "3"
1016 #define TSBASICSYMPLECTIC4         "4"
1017 PETSC_EXTERN PetscErrorCode TSBasicSymplecticSetType(TS,TSBasicSymplecticType);
1018 PETSC_EXTERN PetscErrorCode TSBasicSymplecticGetType(TS,TSBasicSymplecticType*);
1019 PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegister(TSBasicSymplecticType,PetscInt,PetscInt,PetscReal[],PetscReal[]);
1020 PETSC_EXTERN PetscErrorCode TSBasicSymplecticInitializePackage(void);
1021 PETSC_EXTERN PetscErrorCode TSBasicSymplecticFinalizePackage(void);
1022 PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterDestroy(void);
1023 
1024 /*J
1025   TSDiscreteGradient - The Discrete Gradient integrator is a timestepper for hamiltonian systems designed to conserve the first integral (energy), but also has the property for some systems of monotonicity in a functional.
1026 
1027   Level: beginner
1028 
1029   .seealso: `TS`, `TSDISCGRAD`, `TSDiscGradSetFormulation()`, `TSDiscGradGetFormulation`
1030 J*/
1031 PETSC_EXTERN PetscErrorCode TSDiscGradSetFormulation(TS, PetscErrorCode(*)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode(*)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode(*)(TS, PetscReal, Vec, Vec, void *), void *);
1032 PETSC_EXTERN PetscErrorCode TSDiscGradIsGonzalez(TS,PetscBool*);
1033 PETSC_EXTERN PetscErrorCode TSDiscGradUseGonzalez(TS,PetscBool);
1034 
1035 /*
1036        PETSc interface to Sundials
1037 */
1038 #ifdef PETSC_HAVE_SUNDIALS2
1039 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
1040 PETSC_EXTERN const char *const TSSundialsLmmTypes[];
1041 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
1042 PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
1043 PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType);
1044 PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*);
1045 PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal);
1046 PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal);
1047 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal);
1048 PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
1049 PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
1050 PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt);
1051 PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal);
1052 PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool);
1053 PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
1054 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt);
1055 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxord(TS,PetscInt);
1056 PETSC_EXTERN PetscErrorCode TSSundialsSetUseDense(TS,PetscBool);
1057 #endif
1058 
1059 PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal);
1060 PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*);
1061 PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*);
1062 PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool);
1063 
1064 PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal);
1065 PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal);
1066 PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*);
1067 
1068 PETSC_EXTERN PetscErrorCode TSAlpha2SetRadius(TS,PetscReal);
1069 PETSC_EXTERN PetscErrorCode TSAlpha2SetParams(TS,PetscReal,PetscReal,PetscReal,PetscReal);
1070 PETSC_EXTERN PetscErrorCode TSAlpha2GetParams(TS,PetscReal*,PetscReal*,PetscReal*,PetscReal*);
1071 
1072 PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM);
1073 PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*);
1074 
1075 PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*);
1076 PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat,Mat,void*);
1077 
1078 PETSC_EXTERN PetscErrorCode TSRHSJacobianTest(TS,PetscBool*);
1079 PETSC_EXTERN PetscErrorCode TSRHSJacobianTestTranspose(TS,PetscBool*);
1080 
1081 PETSC_EXTERN PetscErrorCode TSGetComputeInitialCondition(TS, PetscErrorCode (**)(TS, Vec));
1082 PETSC_EXTERN PetscErrorCode TSSetComputeInitialCondition(TS, PetscErrorCode (*)(TS, Vec));
1083 PETSC_EXTERN PetscErrorCode TSComputeInitialCondition(TS, Vec);
1084 PETSC_EXTERN PetscErrorCode TSGetComputeExactError(TS, PetscErrorCode (**)(TS, Vec, Vec));
1085 PETSC_EXTERN PetscErrorCode TSSetComputeExactError(TS, PetscErrorCode (*)(TS, Vec, Vec));
1086 PETSC_EXTERN PetscErrorCode TSComputeExactError(TS, Vec, Vec);
1087 PETSC_EXTERN PetscErrorCode PetscConvEstUseTS(PetscConvEst, PetscBool);
1088 
1089 PETSC_EXTERN PetscErrorCode TSSetMatStructure(TS,MatStructure);
1090 #endif
1091