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