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