xref: /petsc/include/petscts.h (revision 0ed3bfb6ec21ed29e56fd524d353af30e5ebf413) !
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 TSALPHA2          "alpha2"
38 #define TSGLLE            "glle"
39 #define TSGLEE            "glee"
40 #define TSSSP             "ssp"
41 #define TSARKIMEX         "arkimex"
42 #define TSROSW            "rosw"
43 #define TSEIMEX           "eimex"
44 #define TSMIMEX           "mimex"
45 #define TSBDF             "bdf"
46 
47 /*E
48     TSProblemType - Determines the type of problem this TS object is to be used to solve
49 
50    Level: beginner
51 
52 .seealso: TSCreate()
53 E*/
54 typedef enum {TS_LINEAR,TS_NONLINEAR} TSProblemType;
55 
56 /*E
57    TSEquationType - type of TS problem that is solved
58 
59    Level: beginner
60 
61    Developer Notes: this must match petsc/finclude/petscts.h
62 
63    Supported types are:
64      TS_EQ_UNSPECIFIED (default)
65      TS_EQ_EXPLICIT {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) := M(t) U_t - G(U,t) = 0
66      TS_EQ_IMPLICIT {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) = 0
67 
68 .seealso: TSGetEquationType(), TSSetEquationType()
69 E*/
70 typedef enum {
71   TS_EQ_UNSPECIFIED               = -1,
72   TS_EQ_EXPLICIT                  = 0,
73   TS_EQ_ODE_EXPLICIT              = 1,
74   TS_EQ_DAE_SEMI_EXPLICIT_INDEX1  = 100,
75   TS_EQ_DAE_SEMI_EXPLICIT_INDEX2  = 200,
76   TS_EQ_DAE_SEMI_EXPLICIT_INDEX3  = 300,
77   TS_EQ_DAE_SEMI_EXPLICIT_INDEXHI = 500,
78   TS_EQ_IMPLICIT                  = 1000,
79   TS_EQ_ODE_IMPLICIT              = 1001,
80   TS_EQ_DAE_IMPLICIT_INDEX1       = 1100,
81   TS_EQ_DAE_IMPLICIT_INDEX2       = 1200,
82   TS_EQ_DAE_IMPLICIT_INDEX3       = 1300,
83   TS_EQ_DAE_IMPLICIT_INDEXHI      = 1500
84 } TSEquationType;
85 PETSC_EXTERN const char *const*TSEquationTypes;
86 
87 /*E
88    TSConvergedReason - reason a TS method has converged or not
89 
90    Level: beginner
91 
92    Developer Notes: this must match petsc/finclude/petscts.h
93 
94    Each reason has its own manual page.
95 
96 .seealso: TSGetConvergedReason()
97 E*/
98 typedef enum {
99   TS_CONVERGED_ITERATING      = 0,
100   TS_CONVERGED_TIME           = 1,
101   TS_CONVERGED_ITS            = 2,
102   TS_CONVERGED_USER           = 3,
103   TS_CONVERGED_EVENT          = 4,
104   TS_CONVERGED_PSEUDO_FATOL   = 5,
105   TS_CONVERGED_PSEUDO_FRTOL   = 6,
106   TS_DIVERGED_NONLINEAR_SOLVE = -1,
107   TS_DIVERGED_STEP_REJECTED   = -2
108 } TSConvergedReason;
109 PETSC_EXTERN const char *const*TSConvergedReasons;
110 
111 /*MC
112    TS_CONVERGED_ITERATING - this only occurs if TSGetConvergedReason() is called during the TSSolve()
113 
114    Level: beginner
115 
116 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt()
117 M*/
118 
119 /*MC
120    TS_CONVERGED_TIME - the final time was reached
121 
122    Level: beginner
123 
124 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetMaxTime(), TSGetMaxTime(), TSGetSolveTime()
125 M*/
126 
127 /*MC
128    TS_CONVERGED_ITS - the maximum number of iterations (time-steps) was reached prior to the final time
129 
130    Level: beginner
131 
132 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetMaxSteps(), TSGetMaxSteps()
133 M*/
134 
135 /*MC
136    TS_CONVERGED_USER - user requested termination
137 
138    Level: beginner
139 
140 .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason()
141 M*/
142 
143 /*MC
144    TS_CONVERGED_EVENT - user requested termination on event detection
145 
146    Level: beginner
147 
148 .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason()
149 M*/
150 
151 /*MC
152    TS_CONVERGED_PSEUDO_FRTOL - stops when function norm decreased by a set amount, used only for TSPSEUDO
153 
154    Level: beginner
155 
156    Options Database:
157 .   -ts_pseudo_frtol <rtol>
158 
159 .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TS_CONVERGED_PSEUDO_FATOL
160 M*/
161 
162 /*MC
163    TS_CONVERGED_PSEUDO_FATOL - stops when function norm decreases below a set amount, used only for TSPSEUDO
164 
165    Level: beginner
166 
167    Options Database:
168 .   -ts_pseudo_fatol <atol>
169 
170 .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TS_CONVERGED_PSEUDO_FRTOL
171 M*/
172 
173 /*MC
174    TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed
175 
176    Level: beginner
177 
178    Notes: See TSSetMaxSNESFailures() for how to allow more nonlinear solver failures.
179 
180 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSGetSNES(), SNESGetConvergedReason(), TSSetMaxSNESFailures()
181 M*/
182 
183 /*MC
184    TS_DIVERGED_STEP_REJECTED - too many steps were rejected
185 
186    Level: beginner
187 
188    Notes: See TSSetMaxStepRejections() for how to allow more step rejections.
189 
190 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetMaxStepRejections()
191 M*/
192 
193 /*E
194    TSExactFinalTimeOption - option for handling of final time step
195 
196    Level: beginner
197 
198    Developer Notes: this must match petsc/finclude/petscts.h
199 
200 $  TS_EXACTFINALTIME_STEPOVER    - Don't do anything if final time is exceeded
201 $  TS_EXACTFINALTIME_INTERPOLATE - Interpolate back to final time
202 $  TS_EXACTFINALTIME_MATCHSTEP - Adapt final time step to match the final time
203 
204 .seealso: TSGetConvergedReason(), TSSetExactFinalTime(), TSGetExactFinalTime()
205 
206 E*/
207 typedef enum {TS_EXACTFINALTIME_UNSPECIFIED=0,TS_EXACTFINALTIME_STEPOVER=1,TS_EXACTFINALTIME_INTERPOLATE=2,TS_EXACTFINALTIME_MATCHSTEP=3} TSExactFinalTimeOption;
208 PETSC_EXTERN const char *const TSExactFinalTimeOptions[];
209 
210 /* Logging support */
211 PETSC_EXTERN PetscClassId TS_CLASSID;
212 PETSC_EXTERN PetscClassId DMTS_CLASSID;
213 PETSC_EXTERN PetscClassId TSADAPT_CLASSID;
214 
215 PETSC_EXTERN PetscErrorCode TSInitializePackage(void);
216 PETSC_EXTERN PetscErrorCode TSFinalizePackage(void);
217 
218 PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm,TS*);
219 PETSC_EXTERN PetscErrorCode TSClone(TS,TS*);
220 PETSC_EXTERN PetscErrorCode TSDestroy(TS*);
221 
222 PETSC_EXTERN PetscErrorCode TSSetProblemType(TS,TSProblemType);
223 PETSC_EXTERN PetscErrorCode TSGetProblemType(TS,TSProblemType*);
224 PETSC_EXTERN PetscErrorCode TSMonitor(TS,PetscInt,PetscReal,Vec);
225 PETSC_EXTERN PetscErrorCode TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**));
226 PETSC_EXTERN PetscErrorCode TSMonitorSetFromOptions(TS,const char[],const char[],const char[],PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*),PetscErrorCode (*)(TS,PetscViewerAndFormat*));
227 PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS);
228 
229 PETSC_EXTERN PetscErrorCode TSAdjointMonitor(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*);
230 PETSC_EXTERN PetscErrorCode TSAdjointMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *,PetscErrorCode (*)(void**));
231 PETSC_EXTERN PetscErrorCode TSAdjointMonitorCancel(TS);
232 PETSC_EXTERN PetscErrorCode TSAdjointMonitorSetFromOptions(TS,const char[],const char[],const char[],PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,PetscViewerAndFormat*),PetscErrorCode (*)(TS,PetscViewerAndFormat*));
233 
234 PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS,const char[]);
235 PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS,const char[]);
236 PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS,const char *[]);
237 PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS);
238 PETSC_EXTERN PetscErrorCode TSSetUp(TS);
239 PETSC_EXTERN PetscErrorCode TSReset(TS);
240 
241 PETSC_EXTERN PetscErrorCode TSSetSolution(TS,Vec);
242 PETSC_EXTERN PetscErrorCode TSGetSolution(TS,Vec*);
243 
244 PETSC_EXTERN PetscErrorCode TS2SetSolution(TS,Vec,Vec);
245 PETSC_EXTERN PetscErrorCode TS2GetSolution(TS,Vec*,Vec*);
246 
247 PETSC_EXTERN PetscErrorCode TSGetSolutionComponents(TS,PetscInt*,Vec*);
248 PETSC_EXTERN PetscErrorCode TSGetAuxSolution(TS,Vec*);
249 PETSC_EXTERN PetscErrorCode TSGetTimeError(TS,PetscInt,Vec*);
250 PETSC_EXTERN PetscErrorCode TSSetTimeError(TS,Vec);
251 
252 /*S
253      TSTrajectory - Abstract PETSc object that stores the trajectory (solution of ODE/ADE at each time step)
254 
255    Level: advanced
256 
257   Concepts: ODE solvers, trajectory
258 
259 .seealso:  TSSetSaveTrajectory(), TSTrajectoryCreate(), TSTrajectorySetType(), TSTrajectoryDestroy()
260 S*/
261 typedef struct _p_TSTrajectory* TSTrajectory;
262 
263 /*J
264     TSTrajectorySetType - String with the name of a PETSc TS trajectory storage method
265 
266    Level: intermediate
267 
268 .seealso:  TSSetSaveTrajectory(), TSTrajectoryCreate(), TSTrajectoryDestroy()
269 J*/
270 typedef const char* TSTrajectoryType;
271 #define TSTRAJECTORYBASIC         "basic"
272 #define TSTRAJECTORYSINGLEFILE    "singlefile"
273 #define TSTRAJECTORYMEMORY        "memory"
274 #define TSTRAJECTORYVISUALIZATION "visualization"
275 
276 PETSC_EXTERN PetscFunctionList TSTrajectoryList;
277 PETSC_EXTERN PetscClassId      TSTRAJECTORY_CLASSID;
278 PETSC_EXTERN PetscBool         TSTrajectoryRegisterAllCalled;
279 
280 PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS);
281 
282 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate(MPI_Comm,TSTrajectory*);
283 PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy(TSTrajectory*);
284 PETSC_EXTERN PetscErrorCode TSTrajectoryView(TSTrajectory,PetscViewer);
285 PETSC_EXTERN PetscErrorCode TSTrajectorySetType(TSTrajectory,TS,const TSTrajectoryType);
286 PETSC_EXTERN PetscErrorCode TSTrajectorySet(TSTrajectory,TS,PetscInt,PetscReal,Vec);
287 PETSC_EXTERN PetscErrorCode TSTrajectoryGet(TSTrajectory,TS,PetscInt,PetscReal*);
288 PETSC_EXTERN PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory,TS);
289 PETSC_EXTERN PetscErrorCode TSTrajectoryRegister(const char[],PetscErrorCode (*)(TSTrajectory,TS));
290 PETSC_EXTERN PetscErrorCode TSTrajectoryRegisterAll(void);
291 PETSC_EXTERN PetscErrorCode TSTrajectorySetUp(TSTrajectory,TS);
292 PETSC_EXTERN PetscErrorCode TSTrajectorySetMonitor(TSTrajectory,PetscBool);
293 PETSC_EXTERN PetscErrorCode TSTrajectorySetVariableNames(TSTrajectory,const char * const*);
294 PETSC_EXTERN PetscErrorCode TSTrajectorySetTransform(TSTrajectory,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void*);
295 PETSC_EXTERN PetscErrorCode TSTrajectorySetKeepFiles(TSTrajectory,PetscBool);
296 PETSC_EXTERN PetscErrorCode TSTrajectorySetDirname(TSTrajectory,const char[]);
297 PETSC_EXTERN PetscErrorCode TSTrajectorySetFiletemplate(TSTrajectory,const char[]);
298 PETSC_EXTERN PetscErrorCode TSGetTrajectory(TS,TSTrajectory*);
299 
300 PETSC_EXTERN PetscErrorCode TSSetCostGradients(TS,PetscInt,Vec*,Vec*);
301 PETSC_EXTERN PetscErrorCode TSGetCostGradients(TS,PetscInt*,Vec**,Vec**);
302 PETSC_EXTERN PetscErrorCode TSSetCostIntegrand(TS,PetscInt,Vec,PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*),PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*),PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*),PetscBool,void*);
303 PETSC_EXTERN PetscErrorCode TSGetCostIntegral(TS,Vec*);
304 PETSC_EXTERN PetscErrorCode TSComputeCostIntegrand(TS,PetscReal,Vec,Vec);
305 
306 PETSC_EXTERN PetscErrorCode TSAdjointSetRHSJacobian(TS,Mat,PetscErrorCode(*)(TS,PetscReal,Vec,Mat,void*),void*);
307 PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS);
308 PETSC_EXTERN PetscErrorCode TSAdjointSetSteps(TS,PetscInt);
309 
310 PETSC_EXTERN PetscErrorCode TSAdjointComputeRHSJacobian(TS,PetscReal,Vec,Mat);
311 PETSC_EXTERN PetscErrorCode TSAdjointStep(TS);
312 PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS);
313 PETSC_EXTERN PetscErrorCode TSAdjointComputeDRDPFunction(TS,PetscReal,Vec,Vec*);
314 PETSC_EXTERN PetscErrorCode TSAdjointComputeDRDYFunction(TS,PetscReal,Vec,Vec*);
315 PETSC_EXTERN PetscErrorCode TSAdjointCostIntegral(TS);
316 
317 PETSC_EXTERN PetscErrorCode TSForwardSetSensitivities(TS,PetscInt,Vec*,PetscInt,Vec*);
318 PETSC_EXTERN PetscErrorCode TSForwardGetSensitivities(TS,PetscInt*,Vec**,PetscInt*,Vec**);
319 PETSC_EXTERN PetscErrorCode TSForwardSetIntegralGradients(TS,PetscInt,Vec *,Vec *);
320 PETSC_EXTERN PetscErrorCode TSForwardGetIntegralGradients(TS,PetscInt*,Vec **,Vec **);
321 PETSC_EXTERN PetscErrorCode TSForwardSetRHSJacobianP(TS,Vec*,PetscErrorCode(*)(TS,PetscReal,Vec,Vec*,void*),void*);
322 PETSC_EXTERN PetscErrorCode TSForwardComputeRHSJacobianP(TS,PetscReal,Vec,Vec*);
323 PETSC_EXTERN PetscErrorCode TSForwardSetUp(TS);
324 PETSC_EXTERN PetscErrorCode TSForwardCostIntegral(TS);
325 PETSC_EXTERN PetscErrorCode TSForwardStep(TS);
326 
327 PETSC_EXTERN PetscErrorCode TSSetMaxSteps(TS,PetscInt);
328 PETSC_EXTERN PetscErrorCode TSGetMaxSteps(TS,PetscInt*);
329 PETSC_EXTERN PetscErrorCode TSSetMaxTime(TS,PetscReal);
330 PETSC_EXTERN PetscErrorCode TSGetMaxTime(TS,PetscReal*);
331 PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS,TSExactFinalTimeOption);
332 PETSC_EXTERN PetscErrorCode TSGetExactFinalTime(TS,TSExactFinalTimeOption*);
333 
334 PETSC_EXTERN PETSC_DEPRECATED("Use TSSetTime[Step]")      PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal);
335 PETSC_EXTERN PETSC_DEPRECATED("Use TSSetMax{Steps|Time}") PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal);
336 PETSC_EXTERN PETSC_DEPRECATED("Use TSGetMax{Steps|Time}") PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*);
337 PETSC_EXTERN PETSC_DEPRECATED("Use TSGetStepNumber")      PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*);
338 PETSC_EXTERN PETSC_DEPRECATED("Use TSGetStepNumber")      PetscErrorCode TSGetTotalSteps(TS,PetscInt*);
339 
340 PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*);
341 
342 typedef struct _n_TSMonitorDrawCtx*  TSMonitorDrawCtx;
343 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorDrawCtx *);
344 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx*);
345 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS,PetscInt,PetscReal,Vec,void*);
346 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionPhase(TS,PetscInt,PetscReal,Vec,void*);
347 PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS,PetscInt,PetscReal,Vec,void*);
348 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionFunction(TS,PetscInt,PetscReal,Vec,void*);
349 
350 PETSC_EXTERN PetscErrorCode TSAdjointMonitorDefault(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,PetscViewerAndFormat *);
351 PETSC_EXTERN PetscErrorCode TSAdjointMonitorDrawSensi(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*);
352 
353 PETSC_EXTERN PetscErrorCode TSMonitorSolution(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*);
354 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*);
355 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void*);
356 
357 PETSC_EXTERN PetscErrorCode TSStep(TS);
358 PETSC_EXTERN PetscErrorCode TSEvaluateWLTE(TS,NormType,PetscInt*,PetscReal*);
359 PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS,PetscInt,Vec,PetscBool*);
360 PETSC_EXTERN PetscErrorCode TSSolve(TS,Vec);
361 PETSC_EXTERN PetscErrorCode TSGetEquationType(TS,TSEquationType*);
362 PETSC_EXTERN PetscErrorCode TSSetEquationType(TS,TSEquationType);
363 PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*);
364 PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS,TSConvergedReason);
365 PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS,PetscReal*);
366 PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS,PetscInt*);
367 PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS,PetscInt*);
368 PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS,PetscInt*);
369 PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS,PetscInt);
370 PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS,PetscInt*);
371 PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS,PetscInt);
372 PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS,PetscBool);
373 PETSC_EXTERN PetscErrorCode TSRestartStep(TS);
374 PETSC_EXTERN PetscErrorCode TSRollBack(TS);
375 
376 PETSC_EXTERN PetscErrorCode TSGetStages(TS,PetscInt*,Vec**);
377 
378 PETSC_EXTERN PetscErrorCode TSGetTime(TS,PetscReal*);
379 PETSC_EXTERN PetscErrorCode TSSetTime(TS,PetscReal);
380 PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS,PetscReal*);
381 PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS,PetscReal*);
382 PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS,PetscReal);
383 PETSC_EXTERN PetscErrorCode TSGetStepNumber(TS,PetscInt*);
384 PETSC_EXTERN PetscErrorCode TSSetStepNumber(TS,PetscInt);
385 
386 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*);
387 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat,Mat,void*);
388 PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*);
389 PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**);
390 PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*);
391 PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**);
392 PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS,PetscBool);
393 
394 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS,PetscReal,Vec,void*);
395 PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS,TSSolutionFunction,void*);
396 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSForcingFunction)(TS,PetscReal,Vec,void*);
397 PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS,TSForcingFunction,void*);
398 
399 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*);
400 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
401 PETSC_EXTERN PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*);
402 PETSC_EXTERN PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**);
403 PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*);
404 PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**);
405 
406 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Function)(TS,PetscReal,Vec,Vec,Vec,Vec,void*);
407 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Jacobian)(TS,PetscReal,Vec,Vec,Vec,PetscReal,PetscReal,Mat,Mat,void*);
408 PETSC_EXTERN PetscErrorCode TSSetI2Function(TS,Vec,TSI2Function,void*);
409 PETSC_EXTERN PetscErrorCode TSGetI2Function(TS,Vec*,TSI2Function*,void**);
410 PETSC_EXTERN PetscErrorCode TSSetI2Jacobian(TS,Mat,Mat,TSI2Jacobian,void*);
411 PETSC_EXTERN PetscErrorCode TSGetI2Jacobian(TS,Mat*,Mat*,TSI2Jacobian*,void**);
412 
413 PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*);
414 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat,Mat,void*);
415 PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*);
416 PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
417 PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS,PetscReal,Vec);
418 PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS,PetscReal,Vec);
419 PETSC_EXTERN PetscErrorCode TSComputeIJacobianDefaultColor(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
420 
421 PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
422 PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS,PetscReal));
423 PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS,PetscReal,PetscInt,Vec*));
424 PETSC_EXTERN PetscErrorCode TSSetPostEvaluate(TS, PetscErrorCode(*)(TS));
425 PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
426 PETSC_EXTERN PetscErrorCode TSPreStep(TS);
427 PETSC_EXTERN PetscErrorCode TSPreStage(TS,PetscReal);
428 PETSC_EXTERN PetscErrorCode TSPostStage(TS,PetscReal,PetscInt,Vec*);
429 PETSC_EXTERN PetscErrorCode TSPostEvaluate(TS);
430 PETSC_EXTERN PetscErrorCode TSPostStep(TS);
431 PETSC_EXTERN PetscErrorCode TSInterpolate(TS,PetscReal,Vec);
432 PETSC_EXTERN PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec);
433 PETSC_EXTERN PetscErrorCode TSGetTolerances(TS,PetscReal*,Vec*,PetscReal*,Vec*);
434 PETSC_EXTERN PetscErrorCode TSErrorWeightedNormInfinity(TS,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
435 PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm2(TS,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
436 PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS,Vec,Vec,NormType,PetscReal*,PetscReal*,PetscReal*);
437 PETSC_EXTERN PetscErrorCode TSErrorWeightedENormInfinity(TS,Vec,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
438 PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm2(TS,Vec,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
439 PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm(TS,Vec,Vec,Vec,NormType,PetscReal*,PetscReal*,PetscReal*);
440 PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal);
441 PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS,PetscReal*);
442 PETSC_EXTERN PetscErrorCode TSSetFunctionDomainError(TS, PetscErrorCode (*)(TS,PetscReal,Vec,PetscBool*));
443 PETSC_EXTERN PetscErrorCode TSFunctionDomainError(TS,PetscReal,Vec,PetscBool*);
444 
445 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*);
446 PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS,PetscReal*,void*);
447 PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *);
448 PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal);
449 PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*);
450 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS,Vec,void*,PetscReal*,PetscBool *);
451 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *);
452 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal);
453 PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);
454 
455 PETSC_EXTERN PetscErrorCode TSPythonSetType(TS,const char[]);
456 
457 PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec);
458 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat,Mat);
459 PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool);
460 PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,PetscBool);
461 PETSC_EXTERN PetscErrorCode TSComputeI2Function(TS,PetscReal,Vec,Vec,Vec,Vec);
462 PETSC_EXTERN PetscErrorCode TSComputeI2Jacobian(TS,PetscReal,Vec,Vec,Vec,PetscReal,PetscReal,Mat,Mat);
463 PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS,PetscReal,PetscReal,PetscReal*,PetscReal*);
464 
465 PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS,Vec,Vec);
466 
467 PETSC_EXTERN PetscErrorCode DMTSSetBoundaryLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *);
468 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM,TSRHSFunction,void*);
469 PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM,TSRHSFunction*,void**);
470 PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM,TSRHSJacobian,void*);
471 PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM,TSRHSJacobian*,void**);
472 PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM,TSIFunction,void*);
473 PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM,TSIFunction*,void**);
474 PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM,TSIJacobian,void*);
475 PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM,TSIJacobian*,void**);
476 PETSC_EXTERN PetscErrorCode DMTSSetI2Function(DM,TSI2Function,void*);
477 PETSC_EXTERN PetscErrorCode DMTSGetI2Function(DM,TSI2Function*,void**);
478 PETSC_EXTERN PetscErrorCode DMTSSetI2Jacobian(DM,TSI2Jacobian,void*);
479 PETSC_EXTERN PetscErrorCode DMTSGetI2Jacobian(DM,TSI2Jacobian*,void**);
480 
481 PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM,TSSolutionFunction,void*);
482 PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM,TSSolutionFunction*,void**);
483 PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM,TSForcingFunction,void*);
484 PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM,TSForcingFunction*,void**);
485 PETSC_EXTERN PetscErrorCode DMTSGetMinRadius(DM,PetscReal*);
486 PETSC_EXTERN PetscErrorCode DMTSSetMinRadius(DM,PetscReal);
487 PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **);
488 
489 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,Vec,void*),void*);
490 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*),void*);
491 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,void*),void*);
492 
493 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));
494 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));
495 
496 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*);
497 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*);
498 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*);
499 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*);
500 
501 PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*),void *);
502 PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*),void *);
503 PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*),void *);
504 PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*),void *);
505 
506 PETSC_EXTERN PetscErrorCode DMPlexTSGetGeometryFVM(DM,Vec*,Vec*,PetscReal*);
507 PETSC_EXTERN PetscErrorCode DMPlexTSGetGradientDM(DM,PetscFV,DM*);
508 
509 typedef struct _n_TSMonitorLGCtx*  TSMonitorLGCtx;
510 typedef struct {
511   Vec            ray;
512   VecScatter     scatter;
513   PetscViewer    viewer;
514   TSMonitorLGCtx lgctx;
515 } TSMonitorDMDARayCtx;
516 PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void**);
517 PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS,PetscInt,PetscReal,Vec,void*);
518 PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS,PetscInt,PetscReal,Vec,void*);
519 
520 /* Dynamic creation and loading functions */
521 PETSC_EXTERN PetscFunctionList TSList;
522 PETSC_EXTERN PetscErrorCode TSGetType(TS,TSType*);
523 PETSC_EXTERN PetscErrorCode TSSetType(TS,TSType);
524 PETSC_EXTERN PetscErrorCode TSRegister(const char[], PetscErrorCode (*)(TS));
525 
526 PETSC_EXTERN PetscErrorCode TSGetSNES(TS,SNES*);
527 PETSC_EXTERN PetscErrorCode TSSetSNES(TS,SNES);
528 PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*);
529 
530 PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer);
531 PETSC_EXTERN PetscErrorCode TSLoad(TS,PetscViewer);
532 PETSC_STATIC_INLINE PetscErrorCode TSViewFromOptions(TS A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
533 PETSC_STATIC_INLINE PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
534 
535 #define TS_FILE_CLASSID 1211225
536 
537 PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *);
538 PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *);
539 
540 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorLGCtx *);
541 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx*);
542 PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS,PetscInt,PetscReal,Vec,void *);
543 PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS,PetscInt,PetscReal,Vec,void *);
544 PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS,const char * const*);
545 PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS,const char *const **);
546 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx,const char * const *);
547 PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS,const char * const*);
548 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx,const char * const*);
549 PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void*);
550 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void *);
551 PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS,PetscInt,PetscReal,Vec,void *);
552 PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS,PetscInt,PetscReal,Vec,void *);
553 PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS,PetscInt,PetscReal,Vec,void *);
554 
555 typedef struct _n_TSMonitorEnvelopeCtx*  TSMonitorEnvelopeCtx;
556 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxCreate(TS,TSMonitorEnvelopeCtx*);
557 PETSC_EXTERN PetscErrorCode TSMonitorEnvelope(TS,PetscInt,PetscReal,Vec,void*);
558 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeGetBounds(TS,Vec*,Vec*);
559 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx*);
560 
561 typedef struct _n_TSMonitorSPEigCtx*  TSMonitorSPEigCtx;
562 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPEigCtx *);
563 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx*);
564 PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS,PetscInt,PetscReal,Vec,void *);
565 
566 PETSC_EXTERN PetscErrorCode TSSetEventHandler(TS,PetscInt,PetscInt[],PetscBool[],PetscErrorCode (*)(TS,PetscReal,Vec,PetscScalar[],void*),PetscErrorCode (*)(TS,PetscInt,PetscInt[],PetscReal,Vec,PetscBool,void*),void*);
567 PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS,PetscReal,PetscReal[]);
568 /*J
569    TSSSPType - string with the name of TSSSP scheme.
570 
571    Level: beginner
572 
573 .seealso: TSSSPSetType(), TS
574 J*/
575 typedef const char* TSSSPType;
576 #define TSSSPRKS2  "rks2"
577 #define TSSSPRKS3  "rks3"
578 #define TSSSPRK104 "rk104"
579 
580 PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,TSSSPType);
581 PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,TSSSPType*);
582 PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt);
583 PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*);
584 PETSC_EXTERN PetscErrorCode TSSSPInitializePackage(void);
585 PETSC_EXTERN PetscErrorCode TSSSPFinalizePackage(void);
586 PETSC_EXTERN PetscFunctionList TSSSPList;
587 
588 /*S
589    TSAdapt - Abstract object that manages time-step adaptivity
590 
591    Level: beginner
592 
593 .seealso: TS, TSAdaptCreate(), TSAdaptType
594 S*/
595 typedef struct _p_TSAdapt *TSAdapt;
596 
597 /*E
598     TSAdaptType - String with the name of TSAdapt scheme.
599 
600    Level: beginner
601 
602 .seealso: TSAdaptSetType(), TS
603 E*/
604 typedef const char *TSAdaptType;
605 #define TSADAPTNONE  "none"
606 #define TSADAPTBASIC "basic"
607 #define TSADAPTDSP   "dsp"
608 #define TSADAPTCFL   "cfl"
609 #define TSADAPTGLEE  "glee"
610 
611 PETSC_EXTERN PetscErrorCode TSGetAdapt(TS,TSAdapt*);
612 PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],PetscErrorCode (*)(TSAdapt));
613 PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void);
614 PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
615 PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*);
616 PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,TSAdaptType);
617 PETSC_EXTERN PetscErrorCode TSAdaptGetType(TSAdapt,TSAdaptType*);
618 PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]);
619 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
620 PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool);
621 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**);
622 PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*);
623 PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscReal,Vec,PetscBool*);
624 PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer);
625 PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt,PetscViewer);
626 PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(PetscOptionItems*,TSAdapt);
627 PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt);
628 PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*);
629 PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool);
630 PETSC_EXTERN PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt,PetscBool);
631 PETSC_EXTERN PetscErrorCode TSAdaptSetSafety(TSAdapt,PetscReal,PetscReal);
632 PETSC_EXTERN PetscErrorCode TSAdaptGetSafety(TSAdapt,PetscReal*,PetscReal*);
633 PETSC_EXTERN PetscErrorCode TSAdaptSetClip(TSAdapt,PetscReal,PetscReal);
634 PETSC_EXTERN PetscErrorCode TSAdaptGetClip(TSAdapt,PetscReal*,PetscReal*);
635 PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal);
636 PETSC_EXTERN PetscErrorCode TSAdaptGetStepLimits(TSAdapt,PetscReal*,PetscReal*);
637 PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt,PetscErrorCode(*)(TSAdapt,TS,PetscReal,Vec,PetscBool*));
638 PETSC_EXTERN PetscErrorCode TSAdaptDSPSetFilter(TSAdapt,const char *);
639 PETSC_EXTERN PetscErrorCode TSAdaptDSPSetPID(TSAdapt,PetscReal,PetscReal,PetscReal);
640 
641 /*S
642    TSGLLEAdapt - Abstract object that manages time-step adaptivity
643 
644    Level: beginner
645 
646    Developer Notes:
647    This functionality should be replaced by the TSAdapt.
648 
649 .seealso: TSGLLE, TSGLLEAdaptCreate(), TSGLLEAdaptType
650 S*/
651 typedef struct _p_TSGLLEAdapt *TSGLLEAdapt;
652 
653 /*J
654     TSGLLEAdaptType - String with the name of TSGLLEAdapt scheme
655 
656    Level: beginner
657 
658 .seealso: TSGLLEAdaptSetType(), TS
659 J*/
660 typedef const char *TSGLLEAdaptType;
661 #define TSGLLEADAPT_NONE "none"
662 #define TSGLLEADAPT_SIZE "size"
663 #define TSGLLEADAPT_BOTH "both"
664 
665 PETSC_EXTERN PetscErrorCode TSGLLEAdaptRegister(const char[],PetscErrorCode (*)(TSGLLEAdapt));
666 PETSC_EXTERN PetscErrorCode TSGLLEAdaptInitializePackage(void);
667 PETSC_EXTERN PetscErrorCode TSGLLEAdaptFinalizePackage(void);
668 PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate(MPI_Comm,TSGLLEAdapt*);
669 PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt,TSGLLEAdaptType);
670 PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt,const char[]);
671 PETSC_EXTERN PetscErrorCode TSGLLEAdaptChoose(TSGLLEAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
672 PETSC_EXTERN PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt,PetscViewer);
673 PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetFromOptions(PetscOptionItems*,TSGLLEAdapt);
674 PETSC_EXTERN PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt*);
675 
676 /*J
677     TSGLLEAcceptType - String with the name of TSGLLEAccept scheme
678 
679    Level: beginner
680 
681 .seealso: TSGLLESetAcceptType(), TS
682 J*/
683 typedef const char *TSGLLEAcceptType;
684 #define TSGLLEACCEPT_ALWAYS "always"
685 
686 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLLEAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *);
687 PETSC_EXTERN PetscErrorCode TSGLLEAcceptRegister(const char[],TSGLLEAcceptFunction);
688 
689 /*J
690   TSGLLEType - family of time integration method within the General Linear class
691 
692   Level: beginner
693 
694 .seealso: TSGLLESetType(), TSGLLERegister()
695 J*/
696 typedef const char* TSGLLEType;
697 #define TSGLLE_IRKS   "irks"
698 
699 PETSC_EXTERN PetscErrorCode TSGLLERegister(const char[],PetscErrorCode(*)(TS));
700 PETSC_EXTERN PetscErrorCode TSGLLEInitializePackage(void);
701 PETSC_EXTERN PetscErrorCode TSGLLEFinalizePackage(void);
702 PETSC_EXTERN PetscErrorCode TSGLLESetType(TS,TSGLLEType);
703 PETSC_EXTERN PetscErrorCode TSGLLEGetAdapt(TS,TSGLLEAdapt*);
704 PETSC_EXTERN PetscErrorCode TSGLLESetAcceptType(TS,TSGLLEAcceptType);
705 
706 /*J
707     TSEIMEXType - String with the name of an Extrapolated IMEX method.
708 
709    Level: beginner
710 
711 .seealso: TSEIMEXSetType(), TS, TSEIMEX, TSEIMEXRegister()
712 J*/
713 #define TSEIMEXType   char*
714 
715 PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts,PetscInt);
716 PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts,PetscInt,PetscInt);
717 PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS,PetscBool);
718 
719 /*J
720     TSRKType - String with the name of a Runge-Kutta method.
721 
722    Level: beginner
723 
724 .seealso: TSRKSetType(), TS, TSRK, TSRKRegister()
725 J*/
726 typedef const char* TSRKType;
727 #define TSRK1FE   "1fe"
728 #define TSRK2A    "2a"
729 #define TSRK3     "3"
730 #define TSRK3BS   "3bs"
731 #define TSRK4     "4"
732 #define TSRK5F    "5f"
733 #define TSRK5DP   "5dp"
734 #define TSRK5BS   "5bs"
735 
736 PETSC_EXTERN PetscErrorCode TSRKGetType(TS ts,TSRKType*);
737 PETSC_EXTERN PetscErrorCode TSRKSetType(TS ts,TSRKType);
738 PETSC_EXTERN PetscErrorCode TSRKSetFullyImplicit(TS,PetscBool);
739 PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
740 PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void);
741 PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void);
742 PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void);
743 
744 /*J
745     TSGLEEType - String with the name of a General Linear with Error Estimation method.
746 
747    Level: beginner
748 
749 .seealso: TSGLEESetType(), TS, TSGLEE, TSGLEERegister()
750 J*/
751 typedef const char* TSGLEEType;
752 #define TSGLEEi1      "BE1"
753 #define TSGLEE23      "23"
754 #define TSGLEE24      "24"
755 #define TSGLEE25I     "25i"
756 #define TSGLEE35      "35"
757 #define TSGLEEEXRK2A  "exrk2a"
758 #define TSGLEERK32G1  "rk32g1"
759 #define TSGLEERK285EX "rk285ex"
760 /*J
761     TSGLEEMode - String with the mode of error estimation for a General Linear with Error Estimation method.
762 
763    Level: beginner
764 
765 .seealso: TSGLEESetMode(), TS, TSGLEE, TSGLEERegister()
766 J*/
767 PETSC_EXTERN PetscErrorCode TSGLEEGetType(TS ts,TSGLEEType*);
768 PETSC_EXTERN PetscErrorCode TSGLEESetType(TS ts,TSGLEEType);
769 PETSC_EXTERN PetscErrorCode TSGLEERegister(TSGLEEType,PetscInt,PetscInt,PetscInt,PetscReal,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
770 PETSC_EXTERN PetscErrorCode TSGLEEFinalizePackage(void);
771 PETSC_EXTERN PetscErrorCode TSGLEEInitializePackage(void);
772 PETSC_EXTERN PetscErrorCode TSGLEERegisterDestroy(void);
773 
774 /*J
775     TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX method.
776 
777    Level: beginner
778 
779 .seealso: TSARKIMEXSetType(), TS, TSARKIMEX, TSARKIMEXRegister()
780 J*/
781 typedef const char* TSARKIMEXType;
782 #define TSARKIMEX1BEE   "1bee"
783 #define TSARKIMEXA2     "a2"
784 #define TSARKIMEXL2     "l2"
785 #define TSARKIMEXARS122 "ars122"
786 #define TSARKIMEX2C     "2c"
787 #define TSARKIMEX2D     "2d"
788 #define TSARKIMEX2E     "2e"
789 #define TSARKIMEXPRSSP2 "prssp2"
790 #define TSARKIMEX3      "3"
791 #define TSARKIMEXBPR3   "bpr3"
792 #define TSARKIMEXARS443 "ars443"
793 #define TSARKIMEX4      "4"
794 #define TSARKIMEX5      "5"
795 PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType*);
796 PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType);
797 PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool);
798 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[]);
799 PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void);
800 PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
801 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);
802 
803 /*J
804     TSRosWType - String with the name of a Rosenbrock-W method.
805 
806    Level: beginner
807 
808 .seealso: TSRosWSetType(), TS, TSROSW, TSRosWRegister()
809 J*/
810 typedef const char* TSRosWType;
811 #define TSROSW2M          "2m"
812 #define TSROSW2P          "2p"
813 #define TSROSWRA3PW       "ra3pw"
814 #define TSROSWRA34PW2     "ra34pw2"
815 #define TSROSWRODAS3      "rodas3"
816 #define TSROSWSANDU3      "sandu3"
817 #define TSROSWASSP3P3S1C  "assp3p3s1c"
818 #define TSROSWLASSP3P4S2C "lassp3p4s2c"
819 #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
820 #define TSROSWARK3        "ark3"
821 #define TSROSWTHETA1      "theta1"
822 #define TSROSWTHETA2      "theta2"
823 #define TSROSWGRK4T       "grk4t"
824 #define TSROSWSHAMP4      "shamp4"
825 #define TSROSWVELDD4      "veldd4"
826 #define TSROSW4L          "4l"
827 
828 PETSC_EXTERN PetscErrorCode TSRosWGetType(TS ts,TSRosWType*);
829 PETSC_EXTERN PetscErrorCode TSRosWSetType(TS ts,TSRosWType);
830 PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool);
831 PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
832 PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
833 PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void);
834 PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
835 PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);
836 
837 PETSC_EXTERN PetscErrorCode TSBDFSetOrder(TS,PetscInt);
838 PETSC_EXTERN PetscErrorCode TSBDFGetOrder(TS,PetscInt*);
839 
840 /*
841        PETSc interface to Sundials
842 */
843 #ifdef PETSC_HAVE_SUNDIALS
844 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
845 PETSC_EXTERN const char *const TSSundialsLmmTypes[];
846 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
847 PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
848 PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType);
849 PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*);
850 PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal);
851 PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal);
852 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal);
853 PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
854 PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
855 PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt);
856 PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal);
857 PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool );
858 PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
859 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt);
860 #endif
861 
862 PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal);
863 PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*);
864 PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*);
865 PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool);
866 
867 PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal);
868 PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal);
869 PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*);
870 
871 PETSC_EXTERN PetscErrorCode TSAlpha2SetRadius(TS,PetscReal);
872 PETSC_EXTERN PetscErrorCode TSAlpha2SetParams(TS,PetscReal,PetscReal,PetscReal,PetscReal);
873 PETSC_EXTERN PetscErrorCode TSAlpha2GetParams(TS,PetscReal*,PetscReal*,PetscReal*,PetscReal*);
874 
875 PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM);
876 PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*);
877 
878 PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*);
879 PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat,Mat,void*);
880 
881 #endif
882