xref: /petsc/include/petscts.h (revision 89ae18915aa5bb3755bc46ef2d36a12ae49d3dd6)
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 
200 .seealso: TSGetConvergedReason(), TSSetExactFinalTime()
201 
202 E*/
203 typedef enum {TS_EXACTFINALTIME_UNSPECIFIED=0,TS_EXACTFINALTIME_STEPOVER=1,TS_EXACTFINALTIME_INTERPOLATE=2,TS_EXACTFINALTIME_MATCHSTEP=3} TSExactFinalTimeOption;
204 PETSC_EXTERN const char *const TSExactFinalTimeOptions[];
205 
206 
207 /* Logging support */
208 PETSC_EXTERN PetscClassId TS_CLASSID;
209 PETSC_EXTERN PetscClassId DMTS_CLASSID;
210 
211 PETSC_EXTERN PetscErrorCode TSInitializePackage(void);
212 
213 PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm,TS*);
214 PETSC_EXTERN PetscErrorCode TSClone(TS,TS*);
215 PETSC_EXTERN PetscErrorCode TSDestroy(TS*);
216 
217 PETSC_EXTERN PetscErrorCode TSSetProblemType(TS,TSProblemType);
218 PETSC_EXTERN PetscErrorCode TSGetProblemType(TS,TSProblemType*);
219 PETSC_EXTERN PetscErrorCode TSMonitor(TS,PetscInt,PetscReal,Vec);
220 PETSC_EXTERN PetscErrorCode TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**));
221 PETSC_EXTERN PetscErrorCode TSMonitorSetFromOptions(TS,const char[],const char[],const char[],PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*),PetscErrorCode (*)(TS,PetscViewerAndFormat*));
222 PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS);
223 
224 PETSC_EXTERN PetscErrorCode TSAdjointMonitor(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*);
225 PETSC_EXTERN PetscErrorCode TSAdjointMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *,PetscErrorCode (*)(void**));
226 PETSC_EXTERN PetscErrorCode TSAdjointMonitorCancel(TS);
227 PETSC_EXTERN PetscErrorCode TSAdjointMonitorSetFromOptions(TS,const char[],const char[],const char[],PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,PetscViewerAndFormat*),PetscErrorCode (*)(TS,PetscViewerAndFormat*));
228 
229 PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS,const char[]);
230 PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS,const char[]);
231 PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS,const char *[]);
232 PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS);
233 PETSC_EXTERN PetscErrorCode TSSetUp(TS);
234 PETSC_EXTERN PetscErrorCode TSReset(TS);
235 
236 PETSC_EXTERN PetscErrorCode TSSetSolution(TS,Vec);
237 PETSC_EXTERN PetscErrorCode TSGetSolution(TS,Vec*);
238 
239 /*S
240      TSTrajectory - Abstract PETSc object that storing the trajectory (solution of ODE/ADE at each time step and stage)
241 
242    Level: advanced
243 
244   Concepts: ODE solvers, adjoint
245 
246 .seealso:  TSCreate(), TSSetType(), TSType, SNES, KSP, PC, TSDestroy()
247 S*/
248 typedef struct _p_TSTrajectory* TSTrajectory;
249 
250 /*J
251     TSTrajectoryType - String with the name of a PETSc TS trajectory storage method
252 
253    Level: intermediate
254 
255 .seealso: TSSetType(), TS, TSRegister(), TSTrajectoryCreate(), TSTrajectorySetType()
256 J*/
257 typedef const char* TSTrajectoryType;
258 #define TSTRAJECTORYBASIC         "basic"
259 #define TSTRAJECTORYSINGLEFILE    "singlefile"
260 #define TSTRAJECTORYMEMORY        "memory"
261 #define TSTRAJECTORYVISUALIZATION "visualization"
262 
263 PETSC_EXTERN PetscFunctionList TSTrajectoryList;
264 PETSC_EXTERN PetscClassId      TSTRAJECTORY_CLASSID;
265 PETSC_EXTERN PetscBool         TSTrajectoryRegisterAllCalled;
266 
267 PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS);
268 
269 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate(MPI_Comm,TSTrajectory*);
270 PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy(TSTrajectory*);
271 PETSC_EXTERN PetscErrorCode TSTrajectorySetType(TSTrajectory,TS,const TSTrajectoryType);
272 PETSC_EXTERN PetscErrorCode TSTrajectorySet(TSTrajectory,TS,PetscInt,PetscReal,Vec);
273 PETSC_EXTERN PetscErrorCode TSTrajectoryGet(TSTrajectory,TS,PetscInt,PetscReal*);
274 PETSC_EXTERN PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory,TS);
275 PETSC_EXTERN PetscErrorCode TSTrajectoryRegisterAll(void);
276 PETSC_EXTERN PetscErrorCode TSTrajectorySetUp(TSTrajectory,TS);
277 
278 PETSC_EXTERN PetscErrorCode TSSetCostGradients(TS,PetscInt,Vec*,Vec*);
279 PETSC_EXTERN PetscErrorCode TSGetCostGradients(TS,PetscInt*,Vec**,Vec**);
280 PETSC_EXTERN PetscErrorCode TSSetCostIntegrand(TS,PetscInt,PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*),PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*),PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*),PetscBool,void*);
281 PETSC_EXTERN PetscErrorCode TSGetCostIntegral(TS,Vec*);
282 
283 PETSC_EXTERN PetscErrorCode TSAdjointSetRHSJacobian(TS,Mat,PetscErrorCode(*)(TS,PetscReal,Vec,Mat,void*),void*);
284 PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS);
285 PETSC_EXTERN PetscErrorCode TSAdjointSetSteps(TS,PetscInt);
286 
287 PETSC_EXTERN PetscErrorCode TSAdjointComputeRHSJacobian(TS,PetscReal,Vec,Mat);
288 PETSC_EXTERN PetscErrorCode TSAdjointStep(TS);
289 PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS);
290 PETSC_EXTERN PetscErrorCode TSAdjointComputeDRDPFunction(TS,PetscReal,Vec,Vec*);
291 PETSC_EXTERN PetscErrorCode TSAdjointComputeDRDYFunction(TS,PetscReal,Vec,Vec*);
292 PETSC_EXTERN PetscErrorCode TSAdjointComputeCostIntegrand(TS,PetscReal,Vec,Vec);
293 PETSC_EXTERN PetscErrorCode TSAdjointCostIntegral(TS);
294 PETSC_EXTERN PetscErrorCode TSForwardCostIntegral(TS);
295 
296 PETSC_EXTERN PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal);
297 PETSC_EXTERN PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*);
298 PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS,TSExactFinalTimeOption);
299 
300 PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*);
301 
302 typedef struct _n_TSMonitorDrawCtx*  TSMonitorDrawCtx;
303 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorDrawCtx *);
304 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx*);
305 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS,PetscInt,PetscReal,Vec,void*);
306 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionPhase(TS,PetscInt,PetscReal,Vec,void*);
307 PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS,PetscInt,PetscReal,Vec,void*);
308 
309 PETSC_EXTERN PetscErrorCode TSAdjointMonitorDefault(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,PetscViewerAndFormat *);
310 PETSC_EXTERN PetscErrorCode TSAdjointMonitorDrawSensi(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*);
311 
312 PETSC_EXTERN PetscErrorCode TSMonitorSolution(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*);
313 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*);
314 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void*);
315 
316 PETSC_EXTERN PetscErrorCode TSStep(TS);
317 PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS,PetscInt,Vec,PetscBool*);
318 PETSC_EXTERN PetscErrorCode TSSolve(TS,Vec);
319 PETSC_EXTERN PetscErrorCode TSGetEquationType(TS,TSEquationType*);
320 PETSC_EXTERN PetscErrorCode TSSetEquationType(TS,TSEquationType);
321 PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*);
322 PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS,TSConvergedReason);
323 PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS,PetscReal*);
324 PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS,PetscInt*);
325 PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS,PetscInt*);
326 PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS,PetscInt*);
327 PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS,PetscInt);
328 PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS,PetscInt*);
329 PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS,PetscInt);
330 PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS,PetscBool);
331 PETSC_EXTERN PetscErrorCode TSRollBack(TS);
332 PETSC_EXTERN PetscErrorCode TSGetTotalSteps(TS,PetscInt*);
333 
334 PETSC_EXTERN PetscErrorCode TSGetStages(TS,PetscInt*,Vec**);
335 
336 PETSC_EXTERN PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal);
337 PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS,PetscReal*);
338 PETSC_EXTERN PetscErrorCode TSGetTime(TS,PetscReal*);
339 PETSC_EXTERN PetscErrorCode TSSetTime(TS,PetscReal);
340 PETSC_EXTERN PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*);
341 PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS,PetscReal);
342 PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS,PetscReal*);
343 
344 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*);
345 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat,Mat,void*);
346 PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*);
347 PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**);
348 PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*);
349 PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**);
350 PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS,PetscBool);
351 
352 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS,PetscReal,Vec,void*);
353 PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS,TSSolutionFunction,void*);
354 PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS,PetscErrorCode (*)(TS,PetscReal,Vec,void*),void*);
355 
356 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*);
357 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
358 PETSC_EXTERN PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*);
359 PETSC_EXTERN PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**);
360 PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*);
361 PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**);
362 
363 PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*);
364 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat,Mat,void*);
365 PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*);
366 PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
367 PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS,PetscReal,Vec);
368 PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS,PetscReal,Vec);
369 PETSC_EXTERN PetscErrorCode TSComputeIJacobianDefaultColor(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
370 
371 PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
372 PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS,PetscReal));
373 PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS,PetscReal,PetscInt,Vec*));
374 PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
375 PETSC_EXTERN PetscErrorCode TSPreStep(TS);
376 PETSC_EXTERN PetscErrorCode TSPreStage(TS,PetscReal);
377 PETSC_EXTERN PetscErrorCode TSPostStage(TS,PetscReal,PetscInt,Vec*);
378 PETSC_EXTERN PetscErrorCode TSPostStep(TS);
379 PETSC_EXTERN PetscErrorCode TSSetRetainStages(TS,PetscBool);
380 PETSC_EXTERN PetscErrorCode TSInterpolate(TS,PetscReal,Vec);
381 PETSC_EXTERN PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec);
382 PETSC_EXTERN PetscErrorCode TSGetTolerances(TS,PetscReal*,Vec*,PetscReal*,Vec*);
383 PETSC_EXTERN PetscErrorCode TSErrorWeightedNormInfinity(TS,Vec,Vec,PetscReal*);
384 PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm2(TS,Vec,Vec,PetscReal*);
385 PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS,Vec,Vec,NormType,PetscReal*);
386 PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal);
387 PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS,PetscReal*);
388 PETSC_EXTERN PetscErrorCode TSSetFunctionDomainError(TS, PetscErrorCode (*)(TS,PetscReal,Vec,PetscBool*));
389 PETSC_EXTERN PetscErrorCode TSFunctionDomainError(TS,PetscReal,Vec,PetscBool*);
390 
391 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*);
392 PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS,PetscReal*,void*);
393 PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *);
394 PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal);
395 
396 PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*);
397 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS,Vec,void*,PetscReal*,PetscBool *);
398 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *);
399 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal);
400 PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);
401 
402 PETSC_EXTERN PetscErrorCode TSPythonSetType(TS,const char[]);
403 
404 PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec);
405 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat,Mat);
406 PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool);
407 PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,PetscBool);
408 PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS,PetscReal,PetscReal,PetscReal*,PetscReal*);
409 
410 PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS,Vec,Vec);
411 
412 PETSC_EXTERN PetscErrorCode DMTSSetBoundaryLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *);
413 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM,TSRHSFunction,void*);
414 PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM,TSRHSFunction*,void**);
415 PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM,TSRHSJacobian,void*);
416 PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM,TSRHSJacobian*,void**);
417 PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM,TSIFunction,void*);
418 PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM,TSIFunction*,void**);
419 PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM,TSIJacobian,void*);
420 PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM,TSIJacobian*,void**);
421 PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM,TSSolutionFunction,void*);
422 PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM,TSSolutionFunction*,void**);
423 PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM,PetscErrorCode (*)(TS,PetscReal,Vec,void*),void*);
424 PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM,PetscErrorCode (**)(TS,PetscReal,Vec,void*),void**);
425 PETSC_EXTERN PetscErrorCode DMTSGetMinRadius(DM,PetscReal*);
426 PETSC_EXTERN PetscErrorCode DMTSSetMinRadius(DM,PetscReal);
427 PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **);
428 
429 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,Vec,void*),void*);
430 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*),void*);
431 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,void*),void*);
432 
433 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));
434 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));
435 
436 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*);
437 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*);
438 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*);
439 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*);
440 
441 PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*),void *);
442 PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*),void *);
443 PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*),void *);
444 PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*),void *);
445 
446 PETSC_EXTERN PetscErrorCode DMPlexTSGetGeometryFVM(DM, Vec *, Vec *, PetscReal *);
447 
448 typedef struct _n_TSMonitorLGCtx*  TSMonitorLGCtx;
449 typedef struct {
450   Vec            ray;
451   VecScatter     scatter;
452   PetscViewer    viewer;
453   TSMonitorLGCtx lgctx;
454 } TSMonitorDMDARayCtx;
455 PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void**);
456 PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS,PetscInt,PetscReal,Vec,void*);
457 PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS,PetscInt,PetscReal,Vec,void*);
458 
459 
460 /* Dynamic creation and loading functions */
461 PETSC_EXTERN PetscFunctionList TSList;
462 PETSC_EXTERN PetscErrorCode TSGetType(TS,TSType*);
463 PETSC_EXTERN PetscErrorCode TSSetType(TS,TSType);
464 PETSC_EXTERN PetscErrorCode TSRegister(const char[], PetscErrorCode (*)(TS));
465 
466 PETSC_EXTERN PetscErrorCode TSGetSNES(TS,SNES*);
467 PETSC_EXTERN PetscErrorCode TSSetSNES(TS,SNES);
468 PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*);
469 
470 PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer);
471 PETSC_EXTERN PetscErrorCode TSLoad(TS,PetscViewer);
472 PETSC_STATIC_INLINE PetscErrorCode TSViewFromOptions(TS A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
473 PETSC_STATIC_INLINE PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
474 
475 #define TS_FILE_CLASSID 1211225
476 
477 PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *);
478 PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *);
479 
480 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorLGCtx *);
481 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx*);
482 PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS,PetscInt,PetscReal,Vec,void *);
483 PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS,PetscInt,PetscReal,Vec,void *);
484 PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS,const char * const*);
485 PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS,const char *const **);
486 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx,const char * const *);
487 PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS,const char * const*);
488 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx,const char * const*);
489 PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void*);
490 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void *);
491 PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS,PetscInt,PetscReal,Vec,void *);
492 PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS,PetscInt,PetscReal,Vec,void *);
493 PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS,PetscInt,PetscReal,Vec,void *);
494 
495 typedef struct _n_TSMonitorEnvelopeCtx*  TSMonitorEnvelopeCtx;
496 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxCreate(TS,TSMonitorEnvelopeCtx*);
497 PETSC_EXTERN PetscErrorCode TSMonitorEnvelope(TS,PetscInt,PetscReal,Vec,void*);
498 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeGetBounds(TS,Vec*,Vec*);
499 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx*);
500 
501 typedef struct _n_TSMonitorSPEigCtx*  TSMonitorSPEigCtx;
502 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPEigCtx *);
503 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx*);
504 PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS,PetscInt,PetscReal,Vec,void *);
505 
506 PETSC_EXTERN PetscErrorCode TSSetEventHandler(TS,PetscInt,PetscInt[],PetscBool[],PetscErrorCode (*)(TS,PetscReal,Vec,PetscScalar[],void*),PetscErrorCode (*)(TS,PetscInt,PetscInt[],PetscReal,Vec,PetscBool,void*),void*);
507 PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS,PetscReal,PetscReal[]);
508 /*J
509    TSSSPType - string with the name of TSSSP scheme.
510 
511    Level: beginner
512 
513 .seealso: TSSSPSetType(), TS
514 J*/
515 typedef const char* TSSSPType;
516 #define TSSSPRKS2  "rks2"
517 #define TSSSPRKS3  "rks3"
518 #define TSSSPRK104 "rk104"
519 
520 PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,TSSSPType);
521 PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,TSSSPType*);
522 PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt);
523 PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*);
524 PETSC_EXTERN PetscErrorCode TSSSPFinalizePackage(void);
525 PETSC_EXTERN PetscErrorCode TSSSPInitializePackage(void);
526 
527 /*S
528    TSAdapt - Abstract object that manages time-step adaptivity
529 
530    Level: beginner
531 
532 .seealso: TS, TSAdaptCreate(), TSAdaptType
533 S*/
534 typedef struct _p_TSAdapt *TSAdapt;
535 
536 /*E
537     TSAdaptType - String with the name of TSAdapt scheme.
538 
539    Level: beginner
540 
541 .seealso: TSAdaptSetType(), TS
542 E*/
543 typedef const char *TSAdaptType;
544 #define TSADAPTBASIC "basic"
545 #define TSADAPTNONE  "none"
546 #define TSADAPTCFL   "cfl"
547 
548 PETSC_EXTERN PetscErrorCode TSGetAdapt(TS,TSAdapt*);
549 PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],PetscErrorCode (*)(TSAdapt));
550 PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void);
551 PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
552 PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*);
553 PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,TSAdaptType);
554 PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]);
555 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
556 PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool);
557 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**);
558 PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*);
559 PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscReal,Vec,PetscBool*);
560 PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer);
561 PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt,PetscViewer);
562 PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(PetscOptionItems*,TSAdapt);
563 PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt);
564 PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*);
565 PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool);
566 PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal);
567 PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt,PetscErrorCode(*)(TSAdapt,TS,PetscReal,Vec,PetscBool*));
568 
569 /*S
570    TSGLAdapt - Abstract object that manages time-step adaptivity
571 
572    Level: beginner
573 
574    Developer Notes:
575    This functionality should be replaced by the TSAdapt.
576 
577 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType
578 S*/
579 typedef struct _p_TSGLAdapt *TSGLAdapt;
580 
581 /*J
582     TSGLAdaptType - String with the name of TSGLAdapt scheme
583 
584    Level: beginner
585 
586 .seealso: TSGLAdaptSetType(), TS
587 J*/
588 typedef const char *TSGLAdaptType;
589 #define TSGLADAPT_NONE "none"
590 #define TSGLADAPT_SIZE "size"
591 #define TSGLADAPT_BOTH "both"
592 
593 PETSC_EXTERN PetscErrorCode TSGLAdaptRegister(const char[],PetscErrorCode (*)(TSGLAdapt));
594 PETSC_EXTERN PetscErrorCode TSGLAdaptInitializePackage(void);
595 PETSC_EXTERN PetscErrorCode TSGLAdaptFinalizePackage(void);
596 PETSC_EXTERN PetscErrorCode TSGLAdaptCreate(MPI_Comm,TSGLAdapt*);
597 PETSC_EXTERN PetscErrorCode TSGLAdaptSetType(TSGLAdapt,TSGLAdaptType);
598 PETSC_EXTERN PetscErrorCode TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]);
599 PETSC_EXTERN PetscErrorCode TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
600 PETSC_EXTERN PetscErrorCode TSGLAdaptView(TSGLAdapt,PetscViewer);
601 PETSC_EXTERN PetscErrorCode TSGLAdaptSetFromOptions(PetscOptionItems*,TSGLAdapt);
602 PETSC_EXTERN PetscErrorCode TSGLAdaptDestroy(TSGLAdapt*);
603 
604 /*J
605     TSGLAcceptType - String with the name of TSGLAccept scheme
606 
607    Level: beginner
608 
609 .seealso: TSGLSetAcceptType(), TS
610 J*/
611 typedef const char *TSGLAcceptType;
612 #define TSGLACCEPT_ALWAYS "always"
613 
614 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *);
615 PETSC_EXTERN PetscErrorCode TSGLAcceptRegister(const char[],TSGLAcceptFunction);
616 
617 /*J
618   TSGLType - family of time integration method within the General Linear class
619 
620   Level: beginner
621 
622 .seealso: TSGLSetType(), TSGLRegister()
623 J*/
624 typedef const char* TSGLType;
625 #define TSGL_IRKS   "irks"
626 
627 PETSC_EXTERN PetscErrorCode TSGLRegister(const char[],PetscErrorCode(*)(TS));
628 PETSC_EXTERN PetscErrorCode TSGLInitializePackage(void);
629 PETSC_EXTERN PetscErrorCode TSGLFinalizePackage(void);
630 PETSC_EXTERN PetscErrorCode TSGLSetType(TS,TSGLType);
631 PETSC_EXTERN PetscErrorCode TSGLGetAdapt(TS,TSGLAdapt*);
632 PETSC_EXTERN PetscErrorCode TSGLSetAcceptType(TS,TSGLAcceptType);
633 
634 /*J
635     TSEIMEXType - String with the name of an Extrapolated IMEX method.
636 
637    Level: beginner
638 
639 .seealso: TSEIMEXSetType(), TS, TSEIMEX, TSEIMEXRegister()
640 J*/
641 #define TSEIMEXType   char*
642 
643 PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts,PetscInt);
644 PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts,PetscInt,PetscInt);
645 PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS,PetscBool);
646 
647 /*J
648     TSRKType - String with the name of a Runge-Kutta method.
649 
650    Level: beginner
651 
652 .seealso: TSRKSetType(), TS, TSRK, TSRKRegister()
653 J*/
654 typedef const char* TSRKType;
655 #define TSRK1FE   "1fe"
656 #define TSRK2A    "2a"
657 #define TSRK3     "3"
658 #define TSRK3BS   "3bs"
659 #define TSRK4     "4"
660 #define TSRK5F    "5f"
661 #define TSRK5DP   "5dp"
662 PETSC_EXTERN PetscErrorCode TSRKGetType(TS ts,TSRKType*);
663 PETSC_EXTERN PetscErrorCode TSRKSetType(TS ts,TSRKType);
664 PETSC_EXTERN PetscErrorCode TSRKSetFullyImplicit(TS,PetscBool);
665 PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
666 PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void);
667 PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void);
668 PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void);
669 
670 /*J
671     TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX method.
672 
673    Level: beginner
674 
675 .seealso: TSARKIMEXSetType(), TS, TSARKIMEX, TSARKIMEXRegister()
676 J*/
677 typedef const char* TSARKIMEXType;
678 #define TSARKIMEX1BEE   "1bee"
679 #define TSARKIMEXA2     "a2"
680 #define TSARKIMEXL2     "l2"
681 #define TSARKIMEXARS122 "ars122"
682 #define TSARKIMEX2C     "2c"
683 #define TSARKIMEX2D     "2d"
684 #define TSARKIMEX2E     "2e"
685 #define TSARKIMEXPRSSP2 "prssp2"
686 #define TSARKIMEX3      "3"
687 #define TSARKIMEXBPR3   "bpr3"
688 #define TSARKIMEXARS443 "ars443"
689 #define TSARKIMEX4      "4"
690 #define TSARKIMEX5      "5"
691 PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType*);
692 PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType);
693 PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool);
694 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[]);
695 PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
696 PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void);
697 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);
698 
699 /*J
700     TSRosWType - String with the name of a Rosenbrock-W method.
701 
702    Level: beginner
703 
704 .seealso: TSRosWSetType(), TS, TSROSW, TSRosWRegister()
705 J*/
706 typedef const char* TSRosWType;
707 #define TSROSW2M          "2m"
708 #define TSROSW2P          "2p"
709 #define TSROSWRA3PW       "ra3pw"
710 #define TSROSWRA34PW2     "ra34pw2"
711 #define TSROSWRODAS3      "rodas3"
712 #define TSROSWSANDU3      "sandu3"
713 #define TSROSWASSP3P3S1C  "assp3p3s1c"
714 #define TSROSWLASSP3P4S2C "lassp3p4s2c"
715 #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
716 #define TSROSWARK3        "ark3"
717 #define TSROSWTHETA1      "theta1"
718 #define TSROSWTHETA2      "theta2"
719 #define TSROSWGRK4T       "grk4t"
720 #define TSROSWSHAMP4      "shamp4"
721 #define TSROSWVELDD4      "veldd4"
722 #define TSROSW4L          "4l"
723 
724 PETSC_EXTERN PetscErrorCode TSRosWGetType(TS ts,TSRosWType*);
725 PETSC_EXTERN PetscErrorCode TSRosWSetType(TS ts,TSRosWType);
726 PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool);
727 PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
728 PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
729 PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
730 PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void);
731 PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);
732 
733 /*
734        PETSc interface to Sundials
735 */
736 #ifdef PETSC_HAVE_SUNDIALS
737 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
738 PETSC_EXTERN const char *const TSSundialsLmmTypes[];
739 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
740 PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
741 PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType);
742 PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*);
743 PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal);
744 PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal);
745 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal);
746 PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
747 PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
748 PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt);
749 PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal);
750 PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool );
751 PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
752 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt);
753 #endif
754 
755 PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal);
756 PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*);
757 PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*);
758 PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool);
759 
760 PETSC_EXTERN PetscErrorCode TSAlphaUseAdapt(TS,PetscBool);
761 PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal);
762 PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal);
763 PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*);
764 
765 PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM);
766 PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*);
767 
768 PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*);
769 PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat,Mat,void*);
770 
771 #endif
772