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