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