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