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