xref: /petsc/include/petscts.h (revision e468013211bb8e2e81ded0f6b4ec8e550e081d7e)
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 TSSetSensitivity(TS,Vec*,PetscInt);
207 PETSC_EXTERN PetscErrorCode TSGetSensitivity(TS,Vec**,PetscInt*);
208 PETSC_EXTERN PetscErrorCode TSSetSensitivityP(TS,Vec*,PetscInt);
209 PETSC_EXTERN PetscErrorCode TSGetSensitivityP(TS,Vec**,PetscInt*);
210 PETSC_EXTERN PetscErrorCode TSSetRHSJacobianP(TS,Mat,PetscErrorCode(*)(TS,PetscReal,Vec,Mat,void*),void*);
211 PETSC_EXTERN PetscErrorCode TSRHSJacobianP(TS,PetscReal,Vec,Mat);
212 PETSC_EXTERN PetscErrorCode TSSetCostIntegrand(TS,PetscInt,Vec,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,void*),void*);
213 PETSC_EXTERN PetscErrorCode TSComputeCostIntegrand(TS,PetscReal,Vec,Vec,void*);
214 
215 PETSC_EXTERN PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal);
216 PETSC_EXTERN PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*);
217 PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS,TSExactFinalTimeOption);
218 
219 PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,void*);
220 
221 typedef struct _n_TSMonitorDrawCtx*  TSMonitorDrawCtx;
222 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorDrawCtx *);
223 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx*);
224 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS,PetscInt,PetscReal,Vec,void*);
225 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionPhase(TS,PetscInt,PetscReal,Vec,void*);
226 PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS,PetscInt,PetscReal,Vec,void*);
227 
228 PETSC_EXTERN PetscErrorCode TSMonitorSolutionBinary(TS,PetscInt,PetscReal,Vec,void*);
229 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*);
230 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void*);
231 
232 PETSC_EXTERN PetscErrorCode TSStep(TS);
233 PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS,PetscInt,Vec,PetscBool*);
234 PETSC_EXTERN PetscErrorCode TSSolve(TS,Vec);
235 PETSC_EXTERN PetscErrorCode TSGetEquationType(TS,TSEquationType*);
236 PETSC_EXTERN PetscErrorCode TSSetEquationType(TS,TSEquationType);
237 PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*);
238 PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS,TSConvergedReason);
239 PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS,PetscReal*);
240 PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS,PetscInt*);
241 PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS,PetscInt*);
242 PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS,PetscInt*);
243 PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS,PetscInt);
244 PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS,PetscInt*);
245 PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS,PetscInt);
246 PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS,PetscBool);
247 PETSC_EXTERN PetscErrorCode TSRollBack(TS);
248 
249 PETSC_EXTERN PetscErrorCode TSGetStages(TS,PetscInt*,Vec**);
250 PETSC_EXTERN PetscErrorCode TSSetCheckpoint(TS,PetscBool);
251 PETSC_EXTERN PetscErrorCode TSSetReverseMode(TS,PetscBool);
252 
253 PETSC_EXTERN PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal);
254 PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS,PetscReal*);
255 PETSC_EXTERN PetscErrorCode TSGetTime(TS,PetscReal*);
256 PETSC_EXTERN PetscErrorCode TSSetTime(TS,PetscReal);
257 PETSC_EXTERN PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*);
258 PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS,PetscReal);
259 PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS,PetscReal*);
260 
261 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*);
262 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat,Mat,void*);
263 PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*);
264 PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**);
265 PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*);
266 PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**);
267 PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS,PetscBool);
268 
269 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS,PetscReal,Vec,void*);
270 PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS,TSSolutionFunction,void*);
271 PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS,PetscErrorCode (*)(TS,PetscReal,Vec,void*),void*);
272 
273 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*);
274 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
275 PETSC_EXTERN PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*);
276 PETSC_EXTERN PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**);
277 PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*);
278 PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**);
279 
280 PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*);
281 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat,Mat,void*);
282 PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*);
283 PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
284 PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS,PetscReal,Vec);
285 PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS,PetscReal,Vec);
286 
287 PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
288 PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS,PetscReal));
289 PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS,PetscReal,PetscInt,Vec*));
290 PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
291 PETSC_EXTERN PetscErrorCode TSPreStep(TS);
292 PETSC_EXTERN PetscErrorCode TSPreStage(TS,PetscReal);
293 PETSC_EXTERN PetscErrorCode TSPostStage(TS,PetscReal,PetscInt,Vec*);
294 PETSC_EXTERN PetscErrorCode TSPostStep(TS);
295 PETSC_EXTERN PetscErrorCode TSSetRetainStages(TS,PetscBool);
296 PETSC_EXTERN PetscErrorCode TSInterpolate(TS,PetscReal,Vec);
297 PETSC_EXTERN PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec);
298 PETSC_EXTERN PetscErrorCode TSGetTolerances(TS,PetscReal*,Vec*,PetscReal*,Vec*);
299 PETSC_EXTERN PetscErrorCode TSErrorNormWRMS(TS,Vec,PetscReal*);
300 PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal);
301 PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS,PetscReal*);
302 
303 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*);
304 PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS,PetscReal*,void*);
305 PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *);
306 PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal);
307 
308 PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*);
309 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS,Vec,void*,PetscReal*,PetscBool *);
310 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *);
311 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal);
312 PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);
313 
314 PETSC_EXTERN PetscErrorCode TSPythonSetType(TS,const char[]);
315 
316 PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec);
317 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat,Mat);
318 PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool);
319 PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,PetscBool);
320 PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS,PetscReal,PetscReal,PetscReal*,PetscReal*);
321 
322 PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS,Vec,Vec);
323 
324 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM,TSRHSFunction,void*);
325 PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM,TSRHSFunction*,void**);
326 PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM,TSRHSJacobian,void*);
327 PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM,TSRHSJacobian*,void**);
328 PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM,TSIFunction,void*);
329 PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM,TSIFunction*,void**);
330 PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM,TSIJacobian,void*);
331 PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM,TSIJacobian*,void**);
332 PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM,TSSolutionFunction,void*);
333 PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM,TSSolutionFunction*,void**);
334 PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM,PetscErrorCode (*)(TS,PetscReal,Vec,void*),void*);
335 PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM,PetscErrorCode (**)(TS,PetscReal,Vec,void*),void**);
336 PETSC_EXTERN PetscErrorCode DMTSGetMinRadius(DM,PetscReal*);
337 PETSC_EXTERN PetscErrorCode DMTSSetMinRadius(DM,PetscReal);
338 PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec, void (**)(const PetscReal[], PetscScalar *, void *), void **);
339 
340 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,Vec,void*),void*);
341 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*),void*);
342 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,void*),void*);
343 
344 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));
345 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));
346 
347 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*);
348 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*);
349 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*);
350 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*);
351 
352 PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*),void *);
353 PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*),void *);
354 PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*),void *);
355 PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*),void *);
356 
357 PETSC_EXTERN PetscErrorCode DMPlexTSGetGeometryFVM(DM, Vec *, Vec *, PetscReal *);
358 
359 typedef struct _n_TSMonitorLGCtx*  TSMonitorLGCtx;
360 typedef struct {
361   Vec            ray;
362   VecScatter     scatter;
363   PetscViewer    viewer;
364   TSMonitorLGCtx lgctx;
365 } TSMonitorDMDARayCtx;
366 PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void**);
367 PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS,PetscInt,PetscReal,Vec,void*);
368 PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS,PetscInt,PetscReal,Vec,void*);
369 
370 
371 /* Dynamic creation and loading functions */
372 PETSC_EXTERN PetscFunctionList TSList;
373 PETSC_EXTERN PetscBool         TSRegisterAllCalled;
374 PETSC_EXTERN PetscErrorCode TSGetType(TS,TSType*);
375 PETSC_EXTERN PetscErrorCode TSSetType(TS,TSType);
376 PETSC_EXTERN PetscErrorCode TSRegister(const char[], PetscErrorCode (*)(TS));
377 PETSC_EXTERN PetscErrorCode TSRegisterAll(void);
378 
379 PETSC_EXTERN PetscErrorCode TSGetSNES(TS,SNES*);
380 PETSC_EXTERN PetscErrorCode TSSetSNES(TS,SNES);
381 PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*);
382 
383 PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer);
384 PETSC_EXTERN PetscErrorCode TSLoad(TS,PetscViewer);
385 PETSC_STATIC_INLINE PetscErrorCode TSViewFromOptions(TS A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);}
386 
387 
388 #define TS_FILE_CLASSID 1211225
389 
390 PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *);
391 PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *);
392 
393 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorLGCtx *);
394 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx*);
395 PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS,PetscInt,PetscReal,Vec,void *);
396 PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS,PetscInt,PetscReal,Vec,void *);
397 PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS,PetscInt,PetscReal,Vec,void *);
398 PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS,PetscInt,PetscReal,Vec,void *);
399 PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS,PetscInt,PetscReal,Vec,void *);
400 
401 typedef struct _n_TSMonitorSPEigCtx*  TSMonitorSPEigCtx;
402 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPEigCtx *);
403 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx*);
404 PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS,PetscInt,PetscReal,Vec,void *);
405 
406 PETSC_EXTERN PetscErrorCode TSSetEventMonitor(TS,PetscInt,PetscInt*,PetscBool*,PetscErrorCode (*)(TS,PetscReal,Vec,PetscScalar*,void*),PetscErrorCode (*)(TS,PetscInt,PetscInt[],PetscReal,Vec,void*),void*);
407 /*J
408    TSSSPType - string with the name of TSSSP scheme.
409 
410    Level: beginner
411 
412 .seealso: TSSSPSetType(), TS
413 J*/
414 typedef const char* TSSSPType;
415 #define TSSSPRKS2  "rks2"
416 #define TSSSPRKS3  "rks3"
417 #define TSSSPRK104 "rk104"
418 
419 PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,TSSSPType);
420 PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,TSSSPType*);
421 PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt);
422 PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*);
423 PETSC_EXTERN PetscErrorCode TSSSPFinalizePackage(void);
424 PETSC_EXTERN PetscErrorCode TSSSPInitializePackage(void);
425 
426 /*S
427    TSAdapt - Abstract object that manages time-step adaptivity
428 
429    Level: beginner
430 
431 .seealso: TS, TSAdaptCreate(), TSAdaptType
432 S*/
433 typedef struct _p_TSAdapt *TSAdapt;
434 
435 /*E
436     TSAdaptType - String with the name of TSAdapt scheme.
437 
438    Level: beginner
439 
440 .seealso: TSAdaptSetType(), TS
441 E*/
442 typedef const char *TSAdaptType;
443 #define TSADAPTBASIC "basic"
444 #define TSADAPTNONE  "none"
445 #define TSADAPTCFL   "cfl"
446 
447 PETSC_EXTERN PetscErrorCode TSGetAdapt(TS,TSAdapt*);
448 PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],PetscErrorCode (*)(TSAdapt));
449 PETSC_EXTERN PetscErrorCode TSAdaptRegisterAll(void);
450 PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void);
451 PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
452 PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*);
453 PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,TSAdaptType);
454 PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]);
455 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
456 PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool);
457 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**);
458 PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*);
459 PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscBool*);
460 PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer);
461 PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt,PetscViewer);
462 PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(PetscOptions*,TSAdapt);
463 PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*);
464 PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool);
465 PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal);
466 PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt,PetscErrorCode(*)(TSAdapt,TS,PetscBool*));
467 
468 /*S
469    TSGLAdapt - Abstract object that manages time-step adaptivity
470 
471    Level: beginner
472 
473    Developer Notes:
474    This functionality should be replaced by the TSAdapt.
475 
476 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType
477 S*/
478 typedef struct _p_TSGLAdapt *TSGLAdapt;
479 
480 /*J
481     TSGLAdaptType - String with the name of TSGLAdapt scheme
482 
483    Level: beginner
484 
485 .seealso: TSGLAdaptSetType(), TS
486 J*/
487 typedef const char *TSGLAdaptType;
488 #define TSGLADAPT_NONE "none"
489 #define TSGLADAPT_SIZE "size"
490 #define TSGLADAPT_BOTH "both"
491 
492 PETSC_EXTERN PetscErrorCode TSGLAdaptRegister(const char[],PetscErrorCode (*)(TSGLAdapt));
493 PETSC_EXTERN PetscErrorCode TSGLAdaptRegisterAll(void);
494 PETSC_EXTERN PetscErrorCode TSGLAdaptInitializePackage(void);
495 PETSC_EXTERN PetscErrorCode TSGLAdaptFinalizePackage(void);
496 PETSC_EXTERN PetscErrorCode TSGLAdaptCreate(MPI_Comm,TSGLAdapt*);
497 PETSC_EXTERN PetscErrorCode TSGLAdaptSetType(TSGLAdapt,TSGLAdaptType);
498 PETSC_EXTERN PetscErrorCode TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]);
499 PETSC_EXTERN PetscErrorCode TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
500 PETSC_EXTERN PetscErrorCode TSGLAdaptView(TSGLAdapt,PetscViewer);
501 PETSC_EXTERN PetscErrorCode TSGLAdaptSetFromOptions(PetscOptions*,TSGLAdapt);
502 PETSC_EXTERN PetscErrorCode TSGLAdaptDestroy(TSGLAdapt*);
503 
504 /*J
505     TSGLAcceptType - String with the name of TSGLAccept scheme
506 
507    Level: beginner
508 
509 .seealso: TSGLSetAcceptType(), TS
510 J*/
511 typedef const char *TSGLAcceptType;
512 #define TSGLACCEPT_ALWAYS "always"
513 
514 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *);
515 PETSC_EXTERN PetscErrorCode TSGLAcceptRegister(const char[],TSGLAcceptFunction);
516 
517 /*J
518   TSGLType - family of time integration method within the General Linear class
519 
520   Level: beginner
521 
522 .seealso: TSGLSetType(), TSGLRegister()
523 J*/
524 typedef const char* TSGLType;
525 #define TSGL_IRKS   "irks"
526 
527 PETSC_EXTERN PetscErrorCode TSGLRegister(const char[],PetscErrorCode(*)(TS));
528 PETSC_EXTERN PetscErrorCode TSGLRegisterAll(void);
529 PETSC_EXTERN PetscErrorCode TSGLInitializePackage(void);
530 PETSC_EXTERN PetscErrorCode TSGLFinalizePackage(void);
531 PETSC_EXTERN PetscErrorCode TSGLSetType(TS,TSGLType);
532 PETSC_EXTERN PetscErrorCode TSGLGetAdapt(TS,TSGLAdapt*);
533 PETSC_EXTERN PetscErrorCode TSGLSetAcceptType(TS,TSGLAcceptType);
534 
535 /*J
536     TSEIMEXType - String with the name of an Extrapolated IMEX method.
537 
538    Level: beginner
539 
540 .seealso: TSEIMEXSetType(), TS, TSEIMEX, TSEIMEXRegister()
541 J*/
542 #define TSEIMEXType   char*
543 
544 PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts,PetscInt);
545 PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts,PetscInt,PetscInt);
546 PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS,PetscBool);
547 
548 /*J
549     TSRKType - String with the name of a Runge-Kutta method.
550 
551    Level: beginner
552 
553 .seealso: TSRKSetType(), TS, TSRK, TSRKRegister()
554 J*/
555 typedef const char* TSRKType;
556 #define TSRK1FE   "1fe"
557 #define TSRK2A    "2a"
558 #define TSRK3     "3"
559 #define TSRK3BS   "3bs"
560 #define TSRK4     "4"
561 #define TSRK5F    "5f"
562 #define TSRK5DP   "5dp"
563 PETSC_EXTERN PetscErrorCode TSRKGetType(TS ts,TSRKType*);
564 PETSC_EXTERN PetscErrorCode TSRKSetType(TS ts,TSRKType);
565 PETSC_EXTERN PetscErrorCode TSRKSetFullyImplicit(TS,PetscBool);
566 PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
567 PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void);
568 PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void);
569 PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void);
570 PETSC_EXTERN PetscErrorCode TSRKRegisterAll(void);
571 
572 /*J
573     TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX method.
574 
575    Level: beginner
576 
577 .seealso: TSARKIMEXSetType(), TS, TSARKIMEX, TSARKIMEXRegister()
578 J*/
579 typedef const char* TSARKIMEXType;
580 #define TSARKIMEX1BEE   "1bee"
581 #define TSARKIMEXA2     "a2"
582 #define TSARKIMEXL2     "l2"
583 #define TSARKIMEXARS122 "ars122"
584 #define TSARKIMEX2C     "2c"
585 #define TSARKIMEX2D     "2d"
586 #define TSARKIMEX2E     "2e"
587 #define TSARKIMEXPRSSP2 "prssp2"
588 #define TSARKIMEX3      "3"
589 #define TSARKIMEXBPR3   "bpr3"
590 #define TSARKIMEXARS443 "ars443"
591 #define TSARKIMEX4      "4"
592 #define TSARKIMEX5      "5"
593 PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType*);
594 PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType);
595 PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool);
596 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[]);
597 PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
598 PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void);
599 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);
600 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterAll(void);
601 
602 /*J
603     TSRosWType - String with the name of a Rosenbrock-W method.
604 
605    Level: beginner
606 
607 .seealso: TSRosWSetType(), TS, TSROSW, TSRosWRegister()
608 J*/
609 typedef const char* TSRosWType;
610 #define TSROSW2M          "2m"
611 #define TSROSW2P          "2p"
612 #define TSROSWRA3PW       "ra3pw"
613 #define TSROSWRA34PW2     "ra34pw2"
614 #define TSROSWRODAS3      "rodas3"
615 #define TSROSWSANDU3      "sandu3"
616 #define TSROSWASSP3P3S1C  "assp3p3s1c"
617 #define TSROSWLASSP3P4S2C "lassp3p4s2c"
618 #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
619 #define TSROSWARK3        "ark3"
620 #define TSROSWTHETA1      "theta1"
621 #define TSROSWTHETA2      "theta2"
622 #define TSROSWGRK4T       "grk4t"
623 #define TSROSWSHAMP4      "shamp4"
624 #define TSROSWVELDD4      "veldd4"
625 #define TSROSW4L          "4l"
626 
627 PETSC_EXTERN PetscErrorCode TSRosWGetType(TS ts,TSRosWType*);
628 PETSC_EXTERN PetscErrorCode TSRosWSetType(TS ts,TSRosWType);
629 PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool);
630 PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
631 PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
632 PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
633 PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void);
634 PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);
635 PETSC_EXTERN PetscErrorCode TSRosWRegisterAll(void);
636 
637 /*
638        PETSc interface to Sundials
639 */
640 #ifdef PETSC_HAVE_SUNDIALS
641 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
642 PETSC_EXTERN const char *const TSSundialsLmmTypes[];
643 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
644 PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
645 PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType);
646 PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*);
647 PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal);
648 PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal);
649 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal);
650 PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
651 PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
652 PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt);
653 PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal);
654 PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool );
655 PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
656 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt);
657 #endif
658 
659 PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal);
660 PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*);
661 PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*);
662 PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool);
663 
664 PETSC_EXTERN PetscErrorCode TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*);
665 PETSC_EXTERN PetscErrorCode TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*);
666 PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal);
667 PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal);
668 PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*);
669 
670 PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM);
671 PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*);
672 
673 PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*);
674 PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat,Mat,void*);
675 
676 #endif
677