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