xref: /petsc/include/petscts.h (revision e5e524a124c0c207e0d57149ef18b3d5980bbea4)
1 /*
2    User interface for the timestepping package. This package
3    is for use in solving time-dependent PDEs.
4 */
5 #if !defined(__PETSCTS_H)
6 #define __PETSCTS_H
7 #include <petscsnes.h>
8 
9 /*S
10      TS - Abstract PETSc object that manages all time-steppers (ODE integrators)
11 
12    Level: beginner
13 
14   Concepts: ODE solvers
15 
16 .seealso:  TSCreate(), TSSetType(), TSType, SNES, KSP, PC, TSDestroy()
17 S*/
18 typedef struct _p_TS* TS;
19 
20 /*J
21     TSType - String with the name of a PETSc TS method.
22 
23    Level: beginner
24 
25 .seealso: TSSetType(), TS, TSRegister()
26 J*/
27 typedef const char* TSType;
28 #define TSEULER           "euler"
29 #define TSBEULER          "beuler"
30 #define TSPSEUDO          "pseudo"
31 #define TSCN              "cn"
32 #define TSSUNDIALS        "sundials"
33 #define TSRK              "rk"
34 #define TSPYTHON          "python"
35 #define TSTHETA           "theta"
36 #define TSALPHA           "alpha"
37 #define TSGL              "gl"
38 #define TSSSP             "ssp"
39 #define TSARKIMEX         "arkimex"
40 #define TSROSW            "rosw"
41 #define TSEIMEX           "eimex"
42 /*E
43     TSProblemType - Determines the type of problem this TS object is to be used to solve
44 
45    Level: beginner
46 
47 .seealso: TSCreate()
48 E*/
49 typedef enum {TS_LINEAR,TS_NONLINEAR} TSProblemType;
50 
51 /*E
52    TSEquationType - type of TS problem that is solved
53 
54    Level: beginner
55 
56    Developer Notes: this must match petsc-finclude/petscts.h
57 
58    Supported types are:
59      TS_EQ_UNSPECIFIED (default)
60      TS_EQ_EXPLICIT {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) := M(t) U_t - G(U,t) = 0
61      TS_EQ_IMPLICIT {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) = 0
62 
63 .seealso: TSGetEquationType(), TSSetEquationType()
64 E*/
65 typedef enum {
66   TS_EQ_UNSPECIFIED               = -1,
67   TS_EQ_EXPLICIT                  = 0,
68   TS_EQ_ODE_EXPLICIT              = 1,
69   TS_EQ_DAE_SEMI_EXPLICIT_INDEX1  = 100,
70   TS_EQ_DAE_SEMI_EXPLICIT_INDEX2  = 200,
71   TS_EQ_DAE_SEMI_EXPLICIT_INDEX3  = 300,
72   TS_EQ_DAE_SEMI_EXPLICIT_INDEXHI = 500,
73   TS_EQ_IMPLICIT                  = 1000,
74   TS_EQ_ODE_IMPLICIT              = 1001,
75   TS_EQ_DAE_IMPLICIT_INDEX1       = 1100,
76   TS_EQ_DAE_IMPLICIT_INDEX2       = 1200,
77   TS_EQ_DAE_IMPLICIT_INDEX3       = 1300,
78   TS_EQ_DAE_IMPLICIT_INDEXHI      = 1500
79 } TSEquationType;
80 PETSC_EXTERN const char *const*TSEquationTypes;
81 
82 /*E
83    TSConvergedReason - reason a TS method has converged or not
84 
85    Level: beginner
86 
87    Developer Notes: this must match petsc-finclude/petscts.h
88 
89    Each reason has its own manual page.
90 
91 .seealso: TSGetConvergedReason()
92 E*/
93 typedef enum {
94   TS_CONVERGED_ITERATING      = 0,
95   TS_CONVERGED_TIME           = 1,
96   TS_CONVERGED_ITS            = 2,
97   TS_CONVERGED_USER           = 3,
98   TS_CONVERGED_EVENT          = 4,
99   TS_DIVERGED_NONLINEAR_SOLVE = -1,
100   TS_DIVERGED_STEP_REJECTED   = -2
101 } TSConvergedReason;
102 PETSC_EXTERN const char *const*TSConvergedReasons;
103 
104 /*MC
105    TS_CONVERGED_ITERATING - this only occurs if TSGetConvergedReason() is called during the TSSolve()
106 
107    Level: beginner
108 
109 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt()
110 M*/
111 
112 /*MC
113    TS_CONVERGED_TIME - the final time was reached
114 
115    Level: beginner
116 
117 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetDuration(), TSGetSolveTime()
118 M*/
119 
120 /*MC
121    TS_CONVERGED_ITS - the maximum number of iterations (time-steps) was reached prior to the final time
122 
123    Level: beginner
124 
125 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetDuration()
126 M*/
127 
128 /*MC
129    TS_CONVERGED_USER - user requested termination
130 
131    Level: beginner
132 
133 .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TSSetDuration()
134 M*/
135 
136 /*MC
137    TS_CONVERGED_EVENT - user requested termination on event detection
138 
139    Level: beginner
140 
141 .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TSSetDuration()
142 M*/
143 
144 /*MC
145    TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed
146 
147    Level: beginner
148 
149    Notes: See TSSetMaxSNESFailures() for how to allow more nonlinear solver failures.
150 
151 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSGetSNES(), SNESGetConvergedReason(), TSSetMaxSNESFailures()
152 M*/
153 
154 /*MC
155    TS_DIVERGED_STEP_REJECTED - too many steps were rejected
156 
157    Level: beginner
158 
159    Notes: See TSSetMaxStepRejections() for how to allow more step rejections.
160 
161 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetMaxStepRejections()
162 M*/
163 
164 /*E
165    TSExactFinalTimeOption - option for handling of final time step
166 
167    Level: beginner
168 
169    Developer Notes: this must match petsc-finclude/petscts.h
170 
171 $  TS_EXACTFINALTIME_STEPOVER    - Don't do anything if final time is exceeded
172 $  TS_EXACTFINALTIME_INTERPOLATE - Interpolate back to final time
173 $  TS_EXACTFINALTIME_MATCHSTEP - Adapt final time step to match the final time
174 .seealso: TSGetConvergedReason(), TSSetExactFinalTime()
175 
176 E*/
177 typedef enum {TS_EXACTFINALTIME_STEPOVER=0,TS_EXACTFINALTIME_INTERPOLATE=1,TS_EXACTFINALTIME_MATCHSTEP=2} TSExactFinalTimeOption;
178 PETSC_EXTERN const char *const TSExactFinalTimeOptions[];
179 
180 
181 /* Logging support */
182 PETSC_EXTERN PetscClassId TS_CLASSID;
183 PETSC_EXTERN PetscClassId DMTS_CLASSID;
184 
185 PETSC_EXTERN PetscErrorCode TSInitializePackage(void);
186 
187 PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm,TS*);
188 PETSC_EXTERN PetscErrorCode TSDestroy(TS*);
189 
190 PETSC_EXTERN PetscErrorCode TSSetProblemType(TS,TSProblemType);
191 PETSC_EXTERN PetscErrorCode TSGetProblemType(TS,TSProblemType*);
192 PETSC_EXTERN PetscErrorCode TSMonitor(TS,PetscInt,PetscReal,Vec);
193 PETSC_EXTERN PetscErrorCode TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**));
194 PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS);
195 
196 PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS,const char[]);
197 PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS,const char[]);
198 PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS,const char *[]);
199 PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS);
200 PETSC_EXTERN PetscErrorCode TSSetUp(TS);
201 PETSC_EXTERN PetscErrorCode TSReset(TS);
202 
203 PETSC_EXTERN PetscErrorCode TSSetSolution(TS,Vec);
204 PETSC_EXTERN PetscErrorCode TSGetSolution(TS,Vec*);
205 
206 PETSC_EXTERN PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal);
207 PETSC_EXTERN PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*);
208 PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS,TSExactFinalTimeOption);
209 
210 PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,void*);
211 
212 typedef struct _n_TSMonitorDrawCtx*  TSMonitorDrawCtx;
213 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorDrawCtx *);
214 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx*);
215 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS,PetscInt,PetscReal,Vec,void*);
216 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionPhase(TS,PetscInt,PetscReal,Vec,void*);
217 PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS,PetscInt,PetscReal,Vec,void*);
218 
219 PETSC_EXTERN PetscErrorCode TSMonitorSolutionBinary(TS,PetscInt,PetscReal,Vec,void*);
220 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*);
221 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void*);
222 
223 PETSC_EXTERN PetscErrorCode TSStep(TS);
224 PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS,PetscInt,Vec,PetscBool*);
225 PETSC_EXTERN PetscErrorCode TSSolve(TS,Vec);
226 PETSC_EXTERN PetscErrorCode TSSolveADJ(TS,Vec);
227 PETSC_EXTERN PetscErrorCode TSGetEquationType(TS,TSEquationType*);
228 PETSC_EXTERN PetscErrorCode TSSetEquationType(TS,TSEquationType);
229 PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*);
230 PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS,TSConvergedReason);
231 PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS,PetscReal*);
232 PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS,PetscInt*);
233 PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS,PetscInt*);
234 PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS,PetscInt*);
235 PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS,PetscInt);
236 PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS,PetscInt*);
237 PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS,PetscInt);
238 PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS,PetscBool);
239 PETSC_EXTERN PetscErrorCode TSRollBack(TS);
240 PETSC_EXTERN PetscErrorCode TSGetStages(TS,PetscInt*,Vec**);
241 
242 PETSC_EXTERN PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal);
243 PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS,PetscReal*);
244 PETSC_EXTERN PetscErrorCode TSGetTime(TS,PetscReal*);
245 PETSC_EXTERN PetscErrorCode TSSetTime(TS,PetscReal);
246 PETSC_EXTERN PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*);
247 PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS,PetscReal);
248 PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS,PetscReal*);
249 
250 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*);
251 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat,Mat,void*);
252 PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*);
253 PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**);
254 PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*);
255 PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**);
256 PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS,PetscBool);
257 
258 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS,PetscReal,Vec,void*);
259 PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS,TSSolutionFunction,void*);
260 PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS,PetscErrorCode (*)(TS,PetscReal,Vec,void*),void*);
261 
262 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*);
263 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
264 PETSC_EXTERN PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*);
265 PETSC_EXTERN PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**);
266 PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*);
267 PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**);
268 
269 PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*);
270 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat,Mat,void*);
271 PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*);
272 PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
273 PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS,PetscReal,Vec);
274 PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS,PetscReal,Vec);
275 
276 PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
277 PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS,PetscReal));
278 PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS,PetscReal,PetscInt,Vec*));
279 PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
280 PETSC_EXTERN PetscErrorCode TSPreStep(TS);
281 PETSC_EXTERN PetscErrorCode TSPreStage(TS,PetscReal);
282 PETSC_EXTERN PetscErrorCode TSPostStage(TS,PetscReal,PetscInt,Vec*);
283 PETSC_EXTERN PetscErrorCode TSPostStep(TS);
284 PETSC_EXTERN PetscErrorCode TSSetRetainStages(TS,PetscBool);
285 PETSC_EXTERN PetscErrorCode TSInterpolate(TS,PetscReal,Vec);
286 PETSC_EXTERN PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec);
287 PETSC_EXTERN PetscErrorCode TSGetTolerances(TS,PetscReal*,Vec*,PetscReal*,Vec*);
288 PETSC_EXTERN PetscErrorCode TSErrorNormWRMS(TS,Vec,PetscReal*);
289 PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal);
290 PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS,PetscReal*);
291 
292 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*);
293 PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS,PetscReal*,void*);
294 PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *);
295 PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal);
296 
297 PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*);
298 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS,Vec,void*,PetscReal*,PetscBool *);
299 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *);
300 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal);
301 PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);
302 
303 PETSC_EXTERN PetscErrorCode TSPythonSetType(TS,const char[]);
304 
305 PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec);
306 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat,Mat);
307 PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool);
308 PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,PetscBool);
309 PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS,PetscReal,PetscReal,PetscReal*,PetscReal*);
310 
311 PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS,Vec,Vec);
312 
313 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM,TSRHSFunction,void*);
314 PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM,TSRHSFunction*,void**);
315 PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM,TSRHSJacobian,void*);
316 PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM,TSRHSJacobian*,void**);
317 PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM,TSIFunction,void*);
318 PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM,TSIFunction*,void**);
319 PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM,TSIJacobian,void*);
320 PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM,TSIJacobian*,void**);
321 PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM,TSSolutionFunction,void*);
322 PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM,TSSolutionFunction*,void**);
323 PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM,PetscErrorCode (*)(TS,PetscReal,Vec,void*),void*);
324 PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM,PetscErrorCode (**)(TS,PetscReal,Vec,void*),void**);
325 PETSC_EXTERN PetscErrorCode DMTSGetMinRadius(DM,PetscReal*);
326 PETSC_EXTERN PetscErrorCode DMTSSetMinRadius(DM,PetscReal);
327 PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec, void (**)(const PetscReal[], PetscScalar *, void *), void **);
328 
329 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,Vec,void*),void*);
330 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*),void*);
331 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,void*),void*);
332 
333 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));
334 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));
335 
336 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*);
337 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*);
338 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*);
339 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*);
340 
341 PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*),void *);
342 PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*),void *);
343 PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*),void *);
344 PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*),void *);
345 
346 PETSC_EXTERN PetscErrorCode DMPlexTSGetGeometryFVM(DM, Vec *, Vec *, PetscReal *);
347 
348 typedef struct _n_TSMonitorLGCtx*  TSMonitorLGCtx;
349 typedef struct {
350   Vec            ray;
351   VecScatter     scatter;
352   PetscViewer    viewer;
353   TSMonitorLGCtx lgctx;
354 } TSMonitorDMDARayCtx;
355 PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void**);
356 PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS,PetscInt,PetscReal,Vec,void*);
357 PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS,PetscInt,PetscReal,Vec,void*);
358 
359 
360 /* Dynamic creation and loading functions */
361 PETSC_EXTERN PetscFunctionList TSList;
362 PETSC_EXTERN PetscBool         TSRegisterAllCalled;
363 PETSC_EXTERN PetscErrorCode TSGetType(TS,TSType*);
364 PETSC_EXTERN PetscErrorCode TSSetType(TS,TSType);
365 PETSC_EXTERN PetscErrorCode TSRegister(const char[], PetscErrorCode (*)(TS));
366 PETSC_EXTERN PetscErrorCode TSRegisterAll(void);
367 
368 PETSC_EXTERN PetscErrorCode TSGetSNES(TS,SNES*);
369 PETSC_EXTERN PetscErrorCode TSSetSNES(TS,SNES);
370 PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*);
371 
372 PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer);
373 PETSC_EXTERN PetscErrorCode TSLoad(TS,PetscViewer);
374 PETSC_STATIC_INLINE PetscErrorCode TSViewFromOptions(TS A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);}
375 
376 
377 #define TS_FILE_CLASSID 1211225
378 
379 PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *);
380 PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *);
381 
382 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorLGCtx *);
383 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx*);
384 PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS,PetscInt,PetscReal,Vec,void *);
385 PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS,PetscInt,PetscReal,Vec,void *);
386 PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS,PetscInt,PetscReal,Vec,void *);
387 PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS,PetscInt,PetscReal,Vec,void *);
388 PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS,PetscInt,PetscReal,Vec,void *);
389 
390 typedef struct _n_TSMonitorSPEigCtx*  TSMonitorSPEigCtx;
391 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPEigCtx *);
392 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx*);
393 PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS,PetscInt,PetscReal,Vec,void *);
394 
395 PETSC_EXTERN PetscErrorCode TSSetEventMonitor(TS,PetscInt,PetscInt*,PetscBool*,PetscErrorCode (*)(TS,PetscReal,Vec,PetscScalar*,void*),PetscErrorCode (*)(TS,PetscInt,PetscInt[],PetscReal,Vec,void*),void*);
396 /*J
397    TSSSPType - string with the name of TSSSP scheme.
398 
399    Level: beginner
400 
401 .seealso: TSSSPSetType(), TS
402 J*/
403 typedef const char* TSSSPType;
404 #define TSSSPRKS2  "rks2"
405 #define TSSSPRKS3  "rks3"
406 #define TSSSPRK104 "rk104"
407 
408 PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,TSSSPType);
409 PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,TSSSPType*);
410 PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt);
411 PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*);
412 PETSC_EXTERN PetscErrorCode TSSSPFinalizePackage(void);
413 PETSC_EXTERN PetscErrorCode TSSSPInitializePackage(void);
414 
415 /*S
416    TSAdapt - Abstract object that manages time-step adaptivity
417 
418    Level: beginner
419 
420 .seealso: TS, TSAdaptCreate(), TSAdaptType
421 S*/
422 typedef struct _p_TSAdapt *TSAdapt;
423 
424 /*E
425     TSAdaptType - String with the name of TSAdapt scheme.
426 
427    Level: beginner
428 
429 .seealso: TSAdaptSetType(), TS
430 E*/
431 typedef const char *TSAdaptType;
432 #define TSADAPTBASIC "basic"
433 #define TSADAPTNONE  "none"
434 #define TSADAPTCFL   "cfl"
435 
436 PETSC_EXTERN PetscErrorCode TSGetAdapt(TS,TSAdapt*);
437 PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],PetscErrorCode (*)(TSAdapt));
438 PETSC_EXTERN PetscErrorCode TSAdaptRegisterAll(void);
439 PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void);
440 PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
441 PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*);
442 PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,TSAdaptType);
443 PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]);
444 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
445 PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool);
446 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**);
447 PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*);
448 PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscBool*);
449 PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer);
450 PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt,PetscViewer);
451 PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(PetscOptions*,TSAdapt);
452 PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*);
453 PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool);
454 PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal);
455 PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt,PetscErrorCode(*)(TSAdapt,TS,PetscBool*));
456 
457 /*S
458    TSGLAdapt - Abstract object that manages time-step adaptivity
459 
460    Level: beginner
461 
462    Developer Notes:
463    This functionality should be replaced by the TSAdapt.
464 
465 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType
466 S*/
467 typedef struct _p_TSGLAdapt *TSGLAdapt;
468 
469 /*J
470     TSGLAdaptType - String with the name of TSGLAdapt scheme
471 
472    Level: beginner
473 
474 .seealso: TSGLAdaptSetType(), TS
475 J*/
476 typedef const char *TSGLAdaptType;
477 #define TSGLADAPT_NONE "none"
478 #define TSGLADAPT_SIZE "size"
479 #define TSGLADAPT_BOTH "both"
480 
481 PETSC_EXTERN PetscErrorCode TSGLAdaptRegister(const char[],PetscErrorCode (*)(TSGLAdapt));
482 PETSC_EXTERN PetscErrorCode TSGLAdaptRegisterAll(void);
483 PETSC_EXTERN PetscErrorCode TSGLAdaptInitializePackage(void);
484 PETSC_EXTERN PetscErrorCode TSGLAdaptFinalizePackage(void);
485 PETSC_EXTERN PetscErrorCode TSGLAdaptCreate(MPI_Comm,TSGLAdapt*);
486 PETSC_EXTERN PetscErrorCode TSGLAdaptSetType(TSGLAdapt,TSGLAdaptType);
487 PETSC_EXTERN PetscErrorCode TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]);
488 PETSC_EXTERN PetscErrorCode TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
489 PETSC_EXTERN PetscErrorCode TSGLAdaptView(TSGLAdapt,PetscViewer);
490 PETSC_EXTERN PetscErrorCode TSGLAdaptSetFromOptions(PetscOptions*,TSGLAdapt);
491 PETSC_EXTERN PetscErrorCode TSGLAdaptDestroy(TSGLAdapt*);
492 
493 /*J
494     TSGLAcceptType - String with the name of TSGLAccept scheme
495 
496    Level: beginner
497 
498 .seealso: TSGLSetAcceptType(), TS
499 J*/
500 typedef const char *TSGLAcceptType;
501 #define TSGLACCEPT_ALWAYS "always"
502 
503 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *);
504 PETSC_EXTERN PetscErrorCode TSGLAcceptRegister(const char[],TSGLAcceptFunction);
505 
506 /*J
507   TSGLType - family of time integration method within the General Linear class
508 
509   Level: beginner
510 
511 .seealso: TSGLSetType(), TSGLRegister()
512 J*/
513 typedef const char* TSGLType;
514 #define TSGL_IRKS   "irks"
515 
516 PETSC_EXTERN PetscErrorCode TSGLRegister(const char[],PetscErrorCode(*)(TS));
517 PETSC_EXTERN PetscErrorCode TSGLRegisterAll(void);
518 PETSC_EXTERN PetscErrorCode TSGLInitializePackage(void);
519 PETSC_EXTERN PetscErrorCode TSGLFinalizePackage(void);
520 PETSC_EXTERN PetscErrorCode TSGLSetType(TS,TSGLType);
521 PETSC_EXTERN PetscErrorCode TSGLGetAdapt(TS,TSGLAdapt*);
522 PETSC_EXTERN PetscErrorCode TSGLSetAcceptType(TS,TSGLAcceptType);
523 
524 /*J
525     TSEIMEXType - String with the name of an Extrapolated IMEX method.
526 
527    Level: beginner
528 
529 .seealso: TSEIMEXSetType(), TS, TSEIMEX, TSEIMEXRegister()
530 J*/
531 #define TSEIMEXType   char*
532 
533 PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts,PetscInt);
534 PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts,PetscInt,PetscInt);
535 PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS,PetscBool);
536 
537 /*J
538     TSRKType - String with the name of a Runge-Kutta method.
539 
540    Level: beginner
541 
542 .seealso: TSRKSetType(), TS, TSRK, TSRKRegister()
543 J*/
544 typedef const char* TSRKType;
545 #define TSRK1FE   "1fe"
546 #define TSRK2A    "2a"
547 #define TSRK3     "3"
548 #define TSRK3BS   "3bs"
549 #define TSRK4     "4"
550 #define TSRK5F    "5f"
551 #define TSRK5DP   "5dp"
552 PETSC_EXTERN PetscErrorCode TSRKGetType(TS ts,TSRKType*);
553 PETSC_EXTERN PetscErrorCode TSRKSetType(TS ts,TSRKType);
554 PETSC_EXTERN PetscErrorCode TSRKSetFullyImplicit(TS,PetscBool);
555 PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
556 PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void);
557 PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void);
558 PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void);
559 PETSC_EXTERN PetscErrorCode TSRKRegisterAll(void);
560 
561 /*J
562     TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX method.
563 
564    Level: beginner
565 
566 .seealso: TSARKIMEXSetType(), TS, TSARKIMEX, TSARKIMEXRegister()
567 J*/
568 typedef const char* TSARKIMEXType;
569 #define TSARKIMEX1BEE   "1bee"
570 #define TSARKIMEXA2     "a2"
571 #define TSARKIMEXL2     "l2"
572 #define TSARKIMEXARS122 "ars122"
573 #define TSARKIMEX2C     "2c"
574 #define TSARKIMEX2D     "2d"
575 #define TSARKIMEX2E     "2e"
576 #define TSARKIMEXPRSSP2 "prssp2"
577 #define TSARKIMEX3      "3"
578 #define TSARKIMEXBPR3   "bpr3"
579 #define TSARKIMEXARS443 "ars443"
580 #define TSARKIMEX4      "4"
581 #define TSARKIMEX5      "5"
582 PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType*);
583 PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType);
584 PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool);
585 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[]);
586 PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
587 PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void);
588 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);
589 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterAll(void);
590 
591 /*J
592     TSRosWType - String with the name of a Rosenbrock-W method.
593 
594    Level: beginner
595 
596 .seealso: TSRosWSetType(), TS, TSROSW, TSRosWRegister()
597 J*/
598 typedef const char* TSRosWType;
599 #define TSROSW2M          "2m"
600 #define TSROSW2P          "2p"
601 #define TSROSWRA3PW       "ra3pw"
602 #define TSROSWRA34PW2     "ra34pw2"
603 #define TSROSWRODAS3      "rodas3"
604 #define TSROSWSANDU3      "sandu3"
605 #define TSROSWASSP3P3S1C  "assp3p3s1c"
606 #define TSROSWLASSP3P4S2C "lassp3p4s2c"
607 #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
608 #define TSROSWARK3        "ark3"
609 #define TSROSWTHETA1      "theta1"
610 #define TSROSWTHETA2      "theta2"
611 #define TSROSWGRK4T       "grk4t"
612 #define TSROSWSHAMP4      "shamp4"
613 #define TSROSWVELDD4      "veldd4"
614 #define TSROSW4L          "4l"
615 
616 PETSC_EXTERN PetscErrorCode TSRosWGetType(TS ts,TSRosWType*);
617 PETSC_EXTERN PetscErrorCode TSRosWSetType(TS ts,TSRosWType);
618 PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool);
619 PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
620 PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
621 PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
622 PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void);
623 PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);
624 PETSC_EXTERN PetscErrorCode TSRosWRegisterAll(void);
625 
626 /*
627        PETSc interface to Sundials
628 */
629 #ifdef PETSC_HAVE_SUNDIALS
630 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
631 PETSC_EXTERN const char *const TSSundialsLmmTypes[];
632 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
633 PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
634 PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType);
635 PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*);
636 PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal);
637 PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal);
638 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal);
639 PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
640 PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
641 PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt);
642 PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal);
643 PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool );
644 PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
645 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt);
646 #endif
647 
648 PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal);
649 PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*);
650 PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*);
651 PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool);
652 
653 PETSC_EXTERN PetscErrorCode TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*);
654 PETSC_EXTERN PetscErrorCode TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*);
655 PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal);
656 PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal);
657 PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*);
658 
659 PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM);
660 PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*);
661 
662 PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*);
663 PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat,Mat,void*);
664 
665 #endif
666