xref: /petsc/include/petscts.h (revision aa0b44dba1dce152c77c35552895bdbefce90cc2)
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 TSPSEUDO          "pseudo"
31 #define TSCN              "cn"
32 #define TSSUNDIALS        "sundials"
33 #define TSRK              "rk"
34 #define TSPYTHON          "python"
35 #define TSTHETA           "theta"
36 #define TSALPHA           "alpha"
37 #define TSGL              "gl"
38 #define TSSSP             "ssp"
39 #define TSARKIMEX         "arkimex"
40 #define TSROSW            "rosw"
41 #define TSEIMEX           "eimex"
42 #define TSMIMEX           "mimex"
43 /*E
44     TSProblemType - Determines the type of problem this TS object is to be used to solve
45 
46    Level: beginner
47 
48 .seealso: TSCreate()
49 E*/
50 typedef enum {TS_LINEAR,TS_NONLINEAR} TSProblemType;
51 
52 /*E
53    TSEquationType - type of TS problem that is solved
54 
55    Level: beginner
56 
57    Developer Notes: this must match petsc/finclude/petscts.h
58 
59    Supported types are:
60      TS_EQ_UNSPECIFIED (default)
61      TS_EQ_EXPLICIT {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) := M(t) U_t - G(U,t) = 0
62      TS_EQ_IMPLICIT {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) = 0
63 
64 .seealso: TSGetEquationType(), TSSetEquationType()
65 E*/
66 typedef enum {
67   TS_EQ_UNSPECIFIED               = -1,
68   TS_EQ_EXPLICIT                  = 0,
69   TS_EQ_ODE_EXPLICIT              = 1,
70   TS_EQ_DAE_SEMI_EXPLICIT_INDEX1  = 100,
71   TS_EQ_DAE_SEMI_EXPLICIT_INDEX2  = 200,
72   TS_EQ_DAE_SEMI_EXPLICIT_INDEX3  = 300,
73   TS_EQ_DAE_SEMI_EXPLICIT_INDEXHI = 500,
74   TS_EQ_IMPLICIT                  = 1000,
75   TS_EQ_ODE_IMPLICIT              = 1001,
76   TS_EQ_DAE_IMPLICIT_INDEX1       = 1100,
77   TS_EQ_DAE_IMPLICIT_INDEX2       = 1200,
78   TS_EQ_DAE_IMPLICIT_INDEX3       = 1300,
79   TS_EQ_DAE_IMPLICIT_INDEXHI      = 1500
80 } TSEquationType;
81 PETSC_EXTERN const char *const*TSEquationTypes;
82 
83 /*E
84    TSConvergedReason - reason a TS method has converged or not
85 
86    Level: beginner
87 
88    Developer Notes: this must match petsc/finclude/petscts.h
89 
90    Each reason has its own manual page.
91 
92 .seealso: TSGetConvergedReason()
93 E*/
94 typedef enum {
95   TS_CONVERGED_ITERATING      = 0,
96   TS_CONVERGED_TIME           = 1,
97   TS_CONVERGED_ITS            = 2,
98   TS_CONVERGED_USER           = 3,
99   TS_CONVERGED_EVENT          = 4,
100   TS_CONVERGED_PSEUDO_FATOL   = 5,
101   TS_CONVERGED_PSEUDO_FRTOL   = 6,
102   TS_DIVERGED_NONLINEAR_SOLVE = -1,
103   TS_DIVERGED_STEP_REJECTED   = -2
104 } TSConvergedReason;
105 PETSC_EXTERN const char *const*TSConvergedReasons;
106 
107 /*MC
108    TS_CONVERGED_ITERATING - this only occurs if TSGetConvergedReason() is called during the TSSolve()
109 
110    Level: beginner
111 
112 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt()
113 M*/
114 
115 /*MC
116    TS_CONVERGED_TIME - the final time was reached
117 
118    Level: beginner
119 
120 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetDuration(), TSGetSolveTime()
121 M*/
122 
123 /*MC
124    TS_CONVERGED_ITS - the maximum number of iterations (time-steps) was reached prior to the final time
125 
126    Level: beginner
127 
128 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetDuration()
129 M*/
130 
131 /*MC
132    TS_CONVERGED_USER - user requested termination
133 
134    Level: beginner
135 
136 .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TSSetDuration()
137 M*/
138 
139 /*MC
140    TS_CONVERGED_EVENT - user requested termination on event detection
141 
142    Level: beginner
143 
144 .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TSSetDuration()
145 M*/
146 
147 /*MC
148    TS_CONVERGED_PSEUDO_FRTOL - stops when function norm decreased by a set amount, used only for TSPSEUDO
149 
150    Level: beginner
151 
152    Options Database:
153 .   -ts_pseudo_frtol <rtol>
154 
155 .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TSSetDuration(), TS_CONVERGED_PSEUDO_FATOL
156 M*/
157 
158 /*MC
159    TS_CONVERGED_PSEUDO_FATOL - stops when function norm decreases below a set amount, used only for TSPSEUDO
160 
161    Level: beginner
162 
163    Options Database:
164 .   -ts_pseudo_fatol <atol>
165 
166 .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TSSetDuration(), TS_CONVERGED_PSEUDO_FRTOL
167 M*/
168 
169 /*MC
170    TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed
171 
172    Level: beginner
173 
174    Notes: See TSSetMaxSNESFailures() for how to allow more nonlinear solver failures.
175 
176 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSGetSNES(), SNESGetConvergedReason(), TSSetMaxSNESFailures()
177 M*/
178 
179 /*MC
180    TS_DIVERGED_STEP_REJECTED - too many steps were rejected
181 
182    Level: beginner
183 
184    Notes: See TSSetMaxStepRejections() for how to allow more step rejections.
185 
186 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetMaxStepRejections()
187 M*/
188 
189 /*E
190    TSExactFinalTimeOption - option for handling of final time step
191 
192    Level: beginner
193 
194    Developer Notes: this must match petsc/finclude/petscts.h
195 
196 $  TS_EXACTFINALTIME_STEPOVER    - Don't do anything if final time is exceeded
197 $  TS_EXACTFINALTIME_INTERPOLATE - Interpolate back to final time
198 $  TS_EXACTFINALTIME_MATCHSTEP - Adapt final time step to match the final time
199 .seealso: TSGetConvergedReason(), TSSetExactFinalTime()
200 
201 E*/
202 typedef enum {TS_EXACTFINALTIME_STEPOVER=0,TS_EXACTFINALTIME_INTERPOLATE=1,TS_EXACTFINALTIME_MATCHSTEP=2} TSExactFinalTimeOption;
203 PETSC_EXTERN const char *const TSExactFinalTimeOptions[];
204 
205 
206 /* Logging support */
207 PETSC_EXTERN PetscClassId TS_CLASSID;
208 PETSC_EXTERN PetscClassId DMTS_CLASSID;
209 
210 PETSC_EXTERN PetscErrorCode TSInitializePackage(void);
211 
212 PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm,TS*);
213 PETSC_EXTERN PetscErrorCode TSClone(TS,TS*);
214 PETSC_EXTERN PetscErrorCode TSDestroy(TS*);
215 
216 PETSC_EXTERN PetscErrorCode TSSetProblemType(TS,TSProblemType);
217 PETSC_EXTERN PetscErrorCode TSGetProblemType(TS,TSProblemType*);
218 PETSC_EXTERN PetscErrorCode TSMonitor(TS,PetscInt,PetscReal,Vec);
219 PETSC_EXTERN PetscErrorCode TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**));
220 PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS);
221 
222 PETSC_EXTERN PetscErrorCode TSAdjointMonitor(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*);
223 PETSC_EXTERN PetscErrorCode TSAdjointMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *,PetscErrorCode (*)(void**));
224 PETSC_EXTERN PetscErrorCode TSAdjointMonitorCancel(TS);
225 
226 
227 PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS,const char[]);
228 PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS,const char[]);
229 PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS,const char *[]);
230 PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS);
231 PETSC_EXTERN PetscErrorCode TSSetUp(TS);
232 PETSC_EXTERN PetscErrorCode TSReset(TS);
233 
234 PETSC_EXTERN PetscErrorCode TSSetSolution(TS,Vec);
235 PETSC_EXTERN PetscErrorCode TSGetSolution(TS,Vec*);
236 
237 /*S
238      TSTrajectory - Abstract PETSc object that storing the trajectory (solution of ODE/ADE at each time step and stage)
239 
240    Level: advanced
241 
242   Concepts: ODE solvers, adjoint
243 
244 .seealso:  TSCreate(), TSSetType(), TSType, SNES, KSP, PC, TSDestroy()
245 S*/
246 typedef struct _p_TSTrajectory* TSTrajectory;
247 
248 /*J
249     TSTrajectoryType - String with the name of a PETSc TS trajectory storage method
250 
251    Level: intermediate
252 
253 .seealso: TSSetType(), TS, TSRegister(), TSTrajectoryCreate(), TSTrajectorySetType()
254 J*/
255 typedef const char* TSTrajectoryType;
256 #define TSTRAJECTORYBASIC      "basic"
257 #define TSTRAJECTORYSINGLEFILE "singlefile"
258 #define TSTRAJECTORYMEMORY     "memory"
259 
260 PETSC_EXTERN PetscFunctionList TSTrajectoryList;
261 PETSC_EXTERN PetscClassId      TSTRAJECTORY_CLASSID;
262 PETSC_EXTERN PetscBool         TSTrajectoryRegisterAllCalled;
263 
264 PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS);
265 
266 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate(MPI_Comm,TSTrajectory*);
267 PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy(TSTrajectory*);
268 PETSC_EXTERN PetscErrorCode TSTrajectorySetType(TSTrajectory,TS,const TSTrajectoryType);
269 PETSC_EXTERN PetscErrorCode TSTrajectorySet(TSTrajectory,TS,PetscInt,PetscReal,Vec);
270 PETSC_EXTERN PetscErrorCode TSTrajectoryGet(TSTrajectory,TS,PetscInt,PetscReal*);
271 PETSC_EXTERN PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory,TS);
272 PETSC_EXTERN PetscErrorCode TSTrajectoryRegisterAll(void);
273 PETSC_EXTERN PetscErrorCode TSTrajectorySetUp(TSTrajectory,TS);
274 
275 PETSC_EXTERN PetscErrorCode TSSetCostGradients(TS,PetscInt,Vec*,Vec*);
276 PETSC_EXTERN PetscErrorCode TSGetCostGradients(TS,PetscInt*,Vec**,Vec**);
277 PETSC_EXTERN PetscErrorCode TSSetCostIntegrand(TS,PetscInt,PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*),PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*),PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*),void*);
278 PETSC_EXTERN PetscErrorCode TSGetCostIntegral(TS,Vec*);
279 
280 PETSC_EXTERN PetscErrorCode TSAdjointSetRHSJacobian(TS,Mat,PetscErrorCode(*)(TS,PetscReal,Vec,Mat,void*),void*);
281 PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS);
282 PETSC_EXTERN PetscErrorCode TSAdjointSetSteps(TS,PetscInt);
283 
284 PETSC_EXTERN PetscErrorCode TSAdjointComputeRHSJacobian(TS,PetscReal,Vec,Mat);
285 PETSC_EXTERN PetscErrorCode TSAdjointStep(TS);
286 PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS);
287 PETSC_EXTERN PetscErrorCode TSAdjointComputeDRDPFunction(TS,PetscReal,Vec,Vec*);
288 PETSC_EXTERN PetscErrorCode TSAdjointComputeDRDYFunction(TS,PetscReal,Vec,Vec*);
289 PETSC_EXTERN PetscErrorCode TSAdjointComputeCostIntegrand(TS,PetscReal,Vec,Vec);
290 
291 PETSC_EXTERN PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal);
292 PETSC_EXTERN PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*);
293 PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS,TSExactFinalTimeOption);
294 
295 PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,void*);
296 
297 typedef struct _n_TSMonitorDrawCtx*  TSMonitorDrawCtx;
298 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorDrawCtx *);
299 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx*);
300 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS,PetscInt,PetscReal,Vec,void*);
301 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionPhase(TS,PetscInt,PetscReal,Vec,void*);
302 PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS,PetscInt,PetscReal,Vec,void*);
303 
304 PETSC_EXTERN PetscErrorCode TSAdjointMonitorDefault(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void *);
305 PETSC_EXTERN PetscErrorCode TSAdjointMonitorDrawSensi(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*);
306 
307 PETSC_EXTERN PetscErrorCode TSMonitorSolution(TS,PetscInt,PetscReal,Vec,void*);
308 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*);
309 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void*);
310 
311 PETSC_EXTERN PetscErrorCode TSStep(TS);
312 PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS,PetscInt,Vec,PetscBool*);
313 PETSC_EXTERN PetscErrorCode TSSolve(TS,Vec);
314 PETSC_EXTERN PetscErrorCode TSGetEquationType(TS,TSEquationType*);
315 PETSC_EXTERN PetscErrorCode TSSetEquationType(TS,TSEquationType);
316 PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*);
317 PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS,TSConvergedReason);
318 PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS,PetscReal*);
319 PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS,PetscInt*);
320 PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS,PetscInt*);
321 PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS,PetscInt*);
322 PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS,PetscInt);
323 PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS,PetscInt*);
324 PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS,PetscInt);
325 PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS,PetscBool);
326 PETSC_EXTERN PetscErrorCode TSRollBack(TS);
327 PETSC_EXTERN PetscErrorCode TSGetTotalSteps(TS,PetscInt*);
328 
329 PETSC_EXTERN PetscErrorCode TSGetStages(TS,PetscInt*,Vec**);
330 
331 PETSC_EXTERN PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal);
332 PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS,PetscReal*);
333 PETSC_EXTERN PetscErrorCode TSGetTime(TS,PetscReal*);
334 PETSC_EXTERN PetscErrorCode TSSetTime(TS,PetscReal);
335 PETSC_EXTERN PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*);
336 PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS,PetscReal);
337 PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS,PetscReal*);
338 
339 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*);
340 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat,Mat,void*);
341 PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*);
342 PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**);
343 PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*);
344 PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**);
345 PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS,PetscBool);
346 
347 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS,PetscReal,Vec,void*);
348 PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS,TSSolutionFunction,void*);
349 PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS,PetscErrorCode (*)(TS,PetscReal,Vec,void*),void*);
350 
351 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*);
352 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
353 PETSC_EXTERN PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*);
354 PETSC_EXTERN PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**);
355 PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*);
356 PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**);
357 
358 PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*);
359 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat,Mat,void*);
360 PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*);
361 PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
362 PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS,PetscReal,Vec);
363 PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS,PetscReal,Vec);
364 PETSC_EXTERN PetscErrorCode TSComputeIJacobianDefaultColor(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
365 
366 PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
367 PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS,PetscReal));
368 PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS,PetscReal,PetscInt,Vec*));
369 PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
370 PETSC_EXTERN PetscErrorCode TSPreStep(TS);
371 PETSC_EXTERN PetscErrorCode TSPreStage(TS,PetscReal);
372 PETSC_EXTERN PetscErrorCode TSPostStage(TS,PetscReal,PetscInt,Vec*);
373 PETSC_EXTERN PetscErrorCode TSPostStep(TS);
374 PETSC_EXTERN PetscErrorCode TSSetRetainStages(TS,PetscBool);
375 PETSC_EXTERN PetscErrorCode TSInterpolate(TS,PetscReal,Vec);
376 PETSC_EXTERN PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec);
377 PETSC_EXTERN PetscErrorCode TSGetTolerances(TS,PetscReal*,Vec*,PetscReal*,Vec*);
378 PETSC_EXTERN PetscErrorCode TSErrorWeightedNormInfinity(TS,Vec,Vec,PetscReal*);
379 PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm2(TS,Vec,Vec,PetscReal*);
380 PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS,Vec,Vec,NormType,PetscReal*);
381 PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal);
382 PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS,PetscReal*);
383 PETSC_EXTERN PetscErrorCode TSSetFunctionDomainError(TS, PetscErrorCode (*)(TS,PetscReal,Vec,PetscBool*));
384 PETSC_EXTERN PetscErrorCode TSFunctionDomainError(TS,PetscReal,Vec,PetscBool*);
385 
386 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*);
387 PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS,PetscReal*,void*);
388 PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *);
389 PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal);
390 
391 PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*);
392 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS,Vec,void*,PetscReal*,PetscBool *);
393 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *);
394 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal);
395 PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);
396 
397 PETSC_EXTERN PetscErrorCode TSPythonSetType(TS,const char[]);
398 
399 PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec);
400 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat,Mat);
401 PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool);
402 PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,PetscBool);
403 PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS,PetscReal,PetscReal,PetscReal*,PetscReal*);
404 
405 PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS,Vec,Vec);
406 
407 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM,TSRHSFunction,void*);
408 PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM,TSRHSFunction*,void**);
409 PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM,TSRHSJacobian,void*);
410 PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM,TSRHSJacobian*,void**);
411 PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM,TSIFunction,void*);
412 PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM,TSIFunction*,void**);
413 PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM,TSIJacobian,void*);
414 PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM,TSIJacobian*,void**);
415 PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM,TSSolutionFunction,void*);
416 PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM,TSSolutionFunction*,void**);
417 PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM,PetscErrorCode (*)(TS,PetscReal,Vec,void*),void*);
418 PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM,PetscErrorCode (**)(TS,PetscReal,Vec,void*),void**);
419 PETSC_EXTERN PetscErrorCode DMTSGetMinRadius(DM,PetscReal*);
420 PETSC_EXTERN PetscErrorCode DMTSSetMinRadius(DM,PetscReal);
421 PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **);
422 
423 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,Vec,void*),void*);
424 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*),void*);
425 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,void*),void*);
426 
427 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));
428 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));
429 
430 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*);
431 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*);
432 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*);
433 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*);
434 
435 PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*),void *);
436 PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*),void *);
437 PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*),void *);
438 PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*),void *);
439 
440 PETSC_EXTERN PetscErrorCode DMPlexTSGetGeometryFVM(DM, Vec *, Vec *, PetscReal *);
441 
442 typedef struct _n_TSMonitorLGCtx*  TSMonitorLGCtx;
443 typedef struct {
444   Vec            ray;
445   VecScatter     scatter;
446   PetscViewer    viewer;
447   TSMonitorLGCtx lgctx;
448 } TSMonitorDMDARayCtx;
449 PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void**);
450 PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS,PetscInt,PetscReal,Vec,void*);
451 PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS,PetscInt,PetscReal,Vec,void*);
452 
453 
454 /* Dynamic creation and loading functions */
455 PETSC_EXTERN PetscFunctionList TSList;
456 PETSC_EXTERN PetscErrorCode TSGetType(TS,TSType*);
457 PETSC_EXTERN PetscErrorCode TSSetType(TS,TSType);
458 PETSC_EXTERN PetscErrorCode TSRegister(const char[], PetscErrorCode (*)(TS));
459 
460 PETSC_EXTERN PetscErrorCode TSGetSNES(TS,SNES*);
461 PETSC_EXTERN PetscErrorCode TSSetSNES(TS,SNES);
462 PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*);
463 
464 PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer);
465 PETSC_EXTERN PetscErrorCode TSLoad(TS,PetscViewer);
466 PETSC_STATIC_INLINE PetscErrorCode TSViewFromOptions(TS A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
467 
468 
469 #define TS_FILE_CLASSID 1211225
470 
471 PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *);
472 PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *);
473 
474 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorLGCtx *);
475 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx*);
476 PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS,PetscInt,PetscReal,Vec,void *);
477 PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS,PetscInt,PetscReal,Vec,void *);
478 PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS,const char * const*);
479 PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS,const char *const **);
480 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx,const char * const *);
481 PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS,const char * const*);
482 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx,const char * const*);
483 PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void*);
484 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void *);
485 PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS,PetscInt,PetscReal,Vec,void *);
486 PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS,PetscInt,PetscReal,Vec,void *);
487 PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS,PetscInt,PetscReal,Vec,void *);
488 
489 typedef struct _n_TSMonitorEnvelopeCtx*  TSMonitorEnvelopeCtx;
490 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxCreate(TS,TSMonitorEnvelopeCtx*);
491 PETSC_EXTERN PetscErrorCode TSMonitorEnvelope(TS,PetscInt,PetscReal,Vec,void*);
492 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeGetBounds(TS,Vec*,Vec*);
493 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx*);
494 
495 typedef struct _n_TSMonitorSPEigCtx*  TSMonitorSPEigCtx;
496 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPEigCtx *);
497 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx*);
498 PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS,PetscInt,PetscReal,Vec,void *);
499 
500 PETSC_EXTERN PetscErrorCode TSSetEventMonitor(TS,PetscInt,PetscInt*,PetscBool*,PetscErrorCode (*)(TS,PetscReal,Vec,PetscScalar*,void*),PetscErrorCode (*)(TS,PetscInt,PetscInt[],PetscReal,Vec,PetscBool,void*),void*);
501 PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS,PetscReal,PetscReal*);
502 /*J
503    TSSSPType - string with the name of TSSSP scheme.
504 
505    Level: beginner
506 
507 .seealso: TSSSPSetType(), TS
508 J*/
509 typedef const char* TSSSPType;
510 #define TSSSPRKS2  "rks2"
511 #define TSSSPRKS3  "rks3"
512 #define TSSSPRK104 "rk104"
513 
514 PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,TSSSPType);
515 PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,TSSSPType*);
516 PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt);
517 PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*);
518 PETSC_EXTERN PetscErrorCode TSSSPFinalizePackage(void);
519 PETSC_EXTERN PetscErrorCode TSSSPInitializePackage(void);
520 
521 /*S
522    TSAdapt - Abstract object that manages time-step adaptivity
523 
524    Level: beginner
525 
526 .seealso: TS, TSAdaptCreate(), TSAdaptType
527 S*/
528 typedef struct _p_TSAdapt *TSAdapt;
529 
530 /*E
531     TSAdaptType - String with the name of TSAdapt scheme.
532 
533    Level: beginner
534 
535 .seealso: TSAdaptSetType(), TS
536 E*/
537 typedef const char *TSAdaptType;
538 #define TSADAPTBASIC "basic"
539 #define TSADAPTNONE  "none"
540 #define TSADAPTCFL   "cfl"
541 
542 PETSC_EXTERN PetscErrorCode TSGetAdapt(TS,TSAdapt*);
543 PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],PetscErrorCode (*)(TSAdapt));
544 PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void);
545 PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
546 PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*);
547 PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,TSAdaptType);
548 PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]);
549 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
550 PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool);
551 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**);
552 PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*);
553 PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscReal,Vec,PetscBool*);
554 PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer);
555 PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt,PetscViewer);
556 PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(PetscOptionItems*,TSAdapt);
557 PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt);
558 PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*);
559 PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool);
560 PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal);
561 PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt,PetscErrorCode(*)(TSAdapt,TS,PetscReal,Vec,PetscBool*));
562 
563 /*S
564    TSGLAdapt - Abstract object that manages time-step adaptivity
565 
566    Level: beginner
567 
568    Developer Notes:
569    This functionality should be replaced by the TSAdapt.
570 
571 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType
572 S*/
573 typedef struct _p_TSGLAdapt *TSGLAdapt;
574 
575 /*J
576     TSGLAdaptType - String with the name of TSGLAdapt scheme
577 
578    Level: beginner
579 
580 .seealso: TSGLAdaptSetType(), TS
581 J*/
582 typedef const char *TSGLAdaptType;
583 #define TSGLADAPT_NONE "none"
584 #define TSGLADAPT_SIZE "size"
585 #define TSGLADAPT_BOTH "both"
586 
587 PETSC_EXTERN PetscErrorCode TSGLAdaptRegister(const char[],PetscErrorCode (*)(TSGLAdapt));
588 PETSC_EXTERN PetscErrorCode TSGLAdaptInitializePackage(void);
589 PETSC_EXTERN PetscErrorCode TSGLAdaptFinalizePackage(void);
590 PETSC_EXTERN PetscErrorCode TSGLAdaptCreate(MPI_Comm,TSGLAdapt*);
591 PETSC_EXTERN PetscErrorCode TSGLAdaptSetType(TSGLAdapt,TSGLAdaptType);
592 PETSC_EXTERN PetscErrorCode TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]);
593 PETSC_EXTERN PetscErrorCode TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
594 PETSC_EXTERN PetscErrorCode TSGLAdaptView(TSGLAdapt,PetscViewer);
595 PETSC_EXTERN PetscErrorCode TSGLAdaptSetFromOptions(PetscOptionItems*,TSGLAdapt);
596 PETSC_EXTERN PetscErrorCode TSGLAdaptDestroy(TSGLAdapt*);
597 
598 /*J
599     TSGLAcceptType - String with the name of TSGLAccept scheme
600 
601    Level: beginner
602 
603 .seealso: TSGLSetAcceptType(), TS
604 J*/
605 typedef const char *TSGLAcceptType;
606 #define TSGLACCEPT_ALWAYS "always"
607 
608 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *);
609 PETSC_EXTERN PetscErrorCode TSGLAcceptRegister(const char[],TSGLAcceptFunction);
610 
611 /*J
612   TSGLType - family of time integration method within the General Linear class
613 
614   Level: beginner
615 
616 .seealso: TSGLSetType(), TSGLRegister()
617 J*/
618 typedef const char* TSGLType;
619 #define TSGL_IRKS   "irks"
620 
621 PETSC_EXTERN PetscErrorCode TSGLRegister(const char[],PetscErrorCode(*)(TS));
622 PETSC_EXTERN PetscErrorCode TSGLInitializePackage(void);
623 PETSC_EXTERN PetscErrorCode TSGLFinalizePackage(void);
624 PETSC_EXTERN PetscErrorCode TSGLSetType(TS,TSGLType);
625 PETSC_EXTERN PetscErrorCode TSGLGetAdapt(TS,TSGLAdapt*);
626 PETSC_EXTERN PetscErrorCode TSGLSetAcceptType(TS,TSGLAcceptType);
627 
628 /*J
629     TSEIMEXType - String with the name of an Extrapolated IMEX method.
630 
631    Level: beginner
632 
633 .seealso: TSEIMEXSetType(), TS, TSEIMEX, TSEIMEXRegister()
634 J*/
635 #define TSEIMEXType   char*
636 
637 PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts,PetscInt);
638 PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts,PetscInt,PetscInt);
639 PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS,PetscBool);
640 
641 /*J
642     TSRKType - String with the name of a Runge-Kutta method.
643 
644    Level: beginner
645 
646 .seealso: TSRKSetType(), TS, TSRK, TSRKRegister()
647 J*/
648 typedef const char* TSRKType;
649 #define TSRK1FE   "1fe"
650 #define TSRK2A    "2a"
651 #define TSRK3     "3"
652 #define TSRK3BS   "3bs"
653 #define TSRK4     "4"
654 #define TSRK5F    "5f"
655 #define TSRK5DP   "5dp"
656 PETSC_EXTERN PetscErrorCode TSRKGetType(TS ts,TSRKType*);
657 PETSC_EXTERN PetscErrorCode TSRKSetType(TS ts,TSRKType);
658 PETSC_EXTERN PetscErrorCode TSRKSetFullyImplicit(TS,PetscBool);
659 PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
660 PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void);
661 PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void);
662 PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void);
663 
664 /*J
665     TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX method.
666 
667    Level: beginner
668 
669 .seealso: TSARKIMEXSetType(), TS, TSARKIMEX, TSARKIMEXRegister()
670 J*/
671 typedef const char* TSARKIMEXType;
672 #define TSARKIMEX1BEE   "1bee"
673 #define TSARKIMEXA2     "a2"
674 #define TSARKIMEXL2     "l2"
675 #define TSARKIMEXARS122 "ars122"
676 #define TSARKIMEX2C     "2c"
677 #define TSARKIMEX2D     "2d"
678 #define TSARKIMEX2E     "2e"
679 #define TSARKIMEXPRSSP2 "prssp2"
680 #define TSARKIMEX3      "3"
681 #define TSARKIMEXBPR3   "bpr3"
682 #define TSARKIMEXARS443 "ars443"
683 #define TSARKIMEX4      "4"
684 #define TSARKIMEX5      "5"
685 PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType*);
686 PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType);
687 PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool);
688 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[]);
689 PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
690 PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void);
691 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);
692 
693 /*J
694     TSRosWType - String with the name of a Rosenbrock-W method.
695 
696    Level: beginner
697 
698 .seealso: TSRosWSetType(), TS, TSROSW, TSRosWRegister()
699 J*/
700 typedef const char* TSRosWType;
701 #define TSROSW2M          "2m"
702 #define TSROSW2P          "2p"
703 #define TSROSWRA3PW       "ra3pw"
704 #define TSROSWRA34PW2     "ra34pw2"
705 #define TSROSWRODAS3      "rodas3"
706 #define TSROSWSANDU3      "sandu3"
707 #define TSROSWASSP3P3S1C  "assp3p3s1c"
708 #define TSROSWLASSP3P4S2C "lassp3p4s2c"
709 #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
710 #define TSROSWARK3        "ark3"
711 #define TSROSWTHETA1      "theta1"
712 #define TSROSWTHETA2      "theta2"
713 #define TSROSWGRK4T       "grk4t"
714 #define TSROSWSHAMP4      "shamp4"
715 #define TSROSWVELDD4      "veldd4"
716 #define TSROSW4L          "4l"
717 
718 PETSC_EXTERN PetscErrorCode TSRosWGetType(TS ts,TSRosWType*);
719 PETSC_EXTERN PetscErrorCode TSRosWSetType(TS ts,TSRosWType);
720 PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool);
721 PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
722 PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
723 PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
724 PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void);
725 PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);
726 
727 /*
728        PETSc interface to Sundials
729 */
730 #ifdef PETSC_HAVE_SUNDIALS
731 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
732 PETSC_EXTERN const char *const TSSundialsLmmTypes[];
733 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
734 PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
735 PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType);
736 PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*);
737 PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal);
738 PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal);
739 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal);
740 PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
741 PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
742 PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt);
743 PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal);
744 PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool );
745 PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
746 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt);
747 #endif
748 
749 PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal);
750 PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*);
751 PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*);
752 PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool);
753 
754 PETSC_EXTERN PetscErrorCode TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*);
755 PETSC_EXTERN PetscErrorCode TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*);
756 PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal);
757 PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal);
758 PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*);
759 
760 PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM);
761 PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*);
762 
763 PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*);
764 PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat,Mat,void*);
765 
766 #endif
767