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