xref: /petsc/include/petscts.h (revision bfef2c865f0b1fc45dc9e1c1c1a6a07f48bb2825)
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 PETSC_EXTERN_CXX_BEGIN
9 
10 /*S
11      TS - Abstract PETSc object that manages all time-steppers (ODE integrators)
12 
13    Level: beginner
14 
15   Concepts: ODE solvers
16 
17 .seealso:  TSCreate(), TSSetType(), TSType, SNES, KSP, PC
18 S*/
19 typedef struct _p_TS* TS;
20 
21 /*E
22     TSType - String with the name of a PETSc TS method or the creation function
23        with an optional dynamic library name, for example
24        http://www.mcs.anl.gov/petsc/lib.a:mytscreate()
25 
26    Level: beginner
27 
28 .seealso: TSSetType(), TS
29 E*/
30 #define TSType char*
31 #define TSEULER           "euler"
32 #define TSBEULER          "beuler"
33 #define TSPSEUDO          "pseudo"
34 #define TSCN              "cn"
35 #define TSSUNDIALS        "sundials"
36 #define TSRK              "rk"
37 #define TSPYTHON          "python"
38 #define TSTHETA           "theta"
39 #define TSALPHA           "alpha"
40 #define TSGL              "gl"
41 #define TSSSP             "ssp"
42 
43 /*E
44     TSProblemType - Determines the type of problem this TS object is to be used to solve
45 
46    Level: beginner
47 
48 .seealso: TSCreate()
49 E*/
50 typedef enum {TS_LINEAR,TS_NONLINEAR} TSProblemType;
51 
52 /* Logging support */
53 extern PetscClassId  TS_CLASSID;
54 
55 extern PetscErrorCode   TSInitializePackage(const char[]);
56 
57 extern PetscErrorCode   TSCreate(MPI_Comm,TS*);
58 extern PetscErrorCode   TSDestroy(TS);
59 
60 extern PetscErrorCode   TSSetProblemType(TS,TSProblemType);
61 extern PetscErrorCode   TSGetProblemType(TS,TSProblemType*);
62 extern PetscErrorCode   TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void*));
63 extern PetscErrorCode   TSMonitorCancel(TS);
64 
65 extern PetscErrorCode   TSSetOptionsPrefix(TS,const char[]);
66 extern PetscErrorCode   TSAppendOptionsPrefix(TS,const char[]);
67 extern PetscErrorCode   TSGetOptionsPrefix(TS,const char *[]);
68 extern PetscErrorCode   TSSetFromOptions(TS);
69 extern PetscErrorCode   TSSetUp(TS);
70 
71 extern PetscErrorCode   TSSetSolution(TS,Vec);
72 extern PetscErrorCode   TSGetSolution(TS,Vec*);
73 
74 extern PetscErrorCode   TSSetDuration(TS,PetscInt,PetscReal);
75 extern PetscErrorCode   TSGetDuration(TS,PetscInt*,PetscReal*);
76 
77 extern PetscErrorCode   TSMonitorDefault(TS,PetscInt,PetscReal,Vec,void*);
78 extern PetscErrorCode   TSMonitorSolution(TS,PetscInt,PetscReal,Vec,void*);
79 extern PetscErrorCode   TSStep(TS,PetscInt *,PetscReal*);
80 extern PetscErrorCode   TSSolve(TS,Vec);
81 
82 
83 extern PetscErrorCode   TSSetInitialTimeStep(TS,PetscReal,PetscReal);
84 extern PetscErrorCode   TSGetTimeStep(TS,PetscReal*);
85 extern PetscErrorCode   TSGetTime(TS,PetscReal*);
86 extern PetscErrorCode   TSSetTime(TS,PetscReal);
87 extern PetscErrorCode   TSGetTimeStepNumber(TS,PetscInt*);
88 extern PetscErrorCode   TSSetTimeStep(TS,PetscReal);
89 
90 extern PetscErrorCode   TSSetRHSFunction(TS,PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*),void*);
91 extern PetscErrorCode   TSSetMatrices(TS,Mat,PetscErrorCode (*)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),Mat,PetscErrorCode (*)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),MatStructure,void*);
92 extern PetscErrorCode   TSGetMatrices(TS,Mat*,Mat*,void**);
93 extern PetscErrorCode   TSSetRHSJacobian(TS,Mat,Mat,PetscErrorCode (*)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*),void*);
94 extern PetscErrorCode   TSGetRHSJacobian(TS,Mat*,Mat*,void**);
95 
96 typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*);
97 typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*);
98 extern PetscErrorCode   TSSetIFunction(TS,TSIFunction,void*);
99 extern PetscErrorCode   TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*);
100 extern PetscErrorCode   TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**);
101 
102 extern PetscErrorCode   TSDefaultComputeJacobianColor(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
103 extern PetscErrorCode   TSDefaultComputeJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
104 
105 extern PetscErrorCode   TSSetPreStep(TS, PetscErrorCode (*)(TS));
106 extern PetscErrorCode   TSSetPostStep(TS, PetscErrorCode (*)(TS));
107 extern PetscErrorCode   TSPreStep(TS);
108 extern PetscErrorCode   TSPostStep(TS);
109 
110 extern PetscErrorCode   TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*);
111 extern PetscErrorCode   TSPseudoDefaultTimeStep(TS,PetscReal*,void*);
112 extern PetscErrorCode   TSPseudoComputeTimeStep(TS,PetscReal *);
113 
114 extern PetscErrorCode   TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*);
115 extern PetscErrorCode   TSPseudoDefaultVerifyTimeStep(TS,Vec,void*,PetscReal*,PetscBool *);
116 extern PetscErrorCode   TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *);
117 extern PetscErrorCode   TSPseudoSetTimeStepIncrement(TS,PetscReal);
118 extern PetscErrorCode   TSPseudoIncrementDtFromInitialDt(TS);
119 
120 extern PetscErrorCode   TSPythonSetType(TS,const char[]);
121 
122 extern PetscErrorCode   TSComputeRHSFunction(TS,PetscReal,Vec,Vec);
123 extern PetscErrorCode   TSComputeRHSJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*);
124 extern PetscErrorCode   TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec);
125 extern PetscErrorCode   TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*);
126 
127 /* Dynamic creation and loading functions */
128 extern PetscFList TSList;
129 extern PetscBool  TSRegisterAllCalled;
130 extern PetscErrorCode   TSGetType(TS,const TSType*);
131 extern PetscErrorCode   TSSetType(TS,const TSType);
132 extern PetscErrorCode   TSRegister(const char[], const char[], const char[], PetscErrorCode (*)(TS));
133 extern PetscErrorCode   TSRegisterAll(const char[]);
134 extern PetscErrorCode   TSRegisterDestroy(void);
135 
136 /*MC
137   TSRegisterDynamic - Adds a creation method to the TS package.
138 
139   Synopsis:
140   PetscErrorCode TSRegisterDynamic(const char *name, const char *path, const char *func_name, PetscErrorCode (*create_func)(TS))
141 
142   Not Collective
143 
144   Input Parameters:
145 + name        - The name of a new user-defined creation routine
146 . path        - The path (either absolute or relative) of the library containing this routine
147 . func_name   - The name of the creation routine
148 - create_func - The creation routine itself
149 
150   Notes:
151   TSRegisterDynamic() may be called multiple times to add several user-defined tses.
152 
153   If dynamic libraries are used, then the fourth input argument (create_func) is ignored.
154 
155   Sample usage:
156 .vb
157   TSRegisterDynamic("my_ts", "/home/username/my_lib/lib/libO/solaris/libmy.a", "MyTSCreate", MyTSCreate);
158 .ve
159 
160   Then, your ts type can be chosen with the procedural interface via
161 .vb
162     TS ts;
163     TSCreate(MPI_Comm, &ts);
164     TSSetType(ts, "my_ts")
165 .ve
166   or at runtime via the option
167 .vb
168     -ts_type my_ts
169 .ve
170 
171   Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
172         If your function is not being put into a shared library then use TSRegister() instead
173 
174   Level: advanced
175 
176 .keywords: TS, register
177 .seealso: TSRegisterAll(), TSRegisterDestroy()
178 M*/
179 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
180 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,0)
181 #else
182 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,d)
183 #endif
184 
185 extern PetscErrorCode   TSGetSNES(TS,SNES*);
186 extern PetscErrorCode   TSGetKSP(TS,KSP*);
187 
188 extern PetscErrorCode   TSView(TS,PetscViewer);
189 extern PetscErrorCode   TSViewFromOptions(TS,const char[]);
190 
191 extern PetscErrorCode   TSSetApplicationContext(TS,void *);
192 extern PetscErrorCode   TSGetApplicationContext(TS,void **);
193 
194 extern PetscErrorCode   TSMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG *);
195 extern PetscErrorCode   TSMonitorLG(TS,PetscInt,PetscReal,Vec,void *);
196 extern PetscErrorCode   TSMonitorLGDestroy(PetscDrawLG);
197 
198 /*S
199    TSGLAdapt - Abstract object that manages time-step adaptivity
200 
201    Level: beginner
202 
203 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType
204 S*/
205 typedef struct _p_TSGLAdapt *TSGLAdapt;
206 
207 /*E
208     TSGLAdaptType - String with the name of TSGLAdapt scheme or the creation function
209        with an optional dynamic library name, for example
210        http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate()
211 
212    Level: beginner
213 
214 .seealso: TSGLAdaptSetType(), TS
215 E*/
216 #define TSGLAdaptType  char*
217 #define TSGLADAPT_NONE "none"
218 #define TSGLADAPT_SIZE "size"
219 #define TSGLADAPT_BOTH "both"
220 
221 /*MC
222    TSGLAdaptRegisterDynamic - adds a TSGLAdapt implementation
223 
224    Synopsis:
225    PetscErrorCode TSGLAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
226 
227    Not Collective
228 
229    Input Parameters:
230 +  name_scheme - name of user-defined adaptivity scheme
231 .  path - path (either absolute or relative) the library containing this scheme
232 .  name_create - name of routine to create method context
233 -  routine_create - routine to create method context
234 
235    Notes:
236    TSGLAdaptRegisterDynamic() may be called multiple times to add several user-defined families.
237 
238    If dynamic libraries are used, then the fourth input argument (routine_create)
239    is ignored.
240 
241    Sample usage:
242 .vb
243    TSGLAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
244                             "MySchemeCreate",MySchemeCreate);
245 .ve
246 
247    Then, your scheme can be chosen with the procedural interface via
248 $     TSGLAdaptSetType(ts,"my_scheme")
249    or at runtime via the option
250 $     -ts_adapt_type my_scheme
251 
252    Level: advanced
253 
254    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
255           and others of the form ${any_environmental_variable} occuring in pathname will be
256           replaced with appropriate values.
257 
258 .keywords: TSGLAdapt, register
259 
260 .seealso: TSGLAdaptRegisterAll()
261 M*/
262 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
263 #  define TSGLAdaptRegisterDynamic(a,b,c,d)  TSGLAdaptRegister(a,b,c,0)
264 #else
265 #  define TSGLAdaptRegisterDynamic(a,b,c,d)  TSGLAdaptRegister(a,b,c,d)
266 #endif
267 
268 extern PetscErrorCode  TSGLAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSGLAdapt));
269 extern PetscErrorCode  TSGLAdaptRegisterAll(const char[]);
270 extern PetscErrorCode  TSGLAdaptRegisterDestroy(void);
271 extern PetscErrorCode  TSGLAdaptInitializePackage(const char[]);
272 extern PetscErrorCode  TSGLAdaptFinalizePackage(void);
273 extern PetscErrorCode  TSGLAdaptCreate(MPI_Comm,TSGLAdapt*);
274 extern PetscErrorCode  TSGLAdaptSetType(TSGLAdapt,const TSGLAdaptType);
275 extern PetscErrorCode  TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]);
276 extern PetscErrorCode  TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
277 extern PetscErrorCode  TSGLAdaptView(TSGLAdapt,PetscViewer);
278 extern PetscErrorCode  TSGLAdaptSetFromOptions(TSGLAdapt);
279 extern PetscErrorCode  TSGLAdaptDestroy(TSGLAdapt);
280 
281 /*E
282     TSGLAcceptType - String with the name of TSGLAccept scheme or the function
283        with an optional dynamic library name, for example
284        http://www.mcs.anl.gov/petsc/lib.a:mytsglaccept()
285 
286    Level: beginner
287 
288 .seealso: TSGLSetAcceptType(), TS
289 E*/
290 #define TSGLAcceptType  char*
291 #define TSGLACCEPT_ALWAYS "always"
292 
293 typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *);
294 extern PetscErrorCode  TSGLAcceptRegister(const char[],const char[],const char[],TSGLAcceptFunction);
295 
296 /*MC
297    TSGLAcceptRegisterDynamic - adds a TSGL acceptance scheme
298 
299    Synopsis:
300    PetscErrorCode TSGLAcceptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
301 
302    Not Collective
303 
304    Input Parameters:
305 +  name_scheme - name of user-defined acceptance scheme
306 .  path - path (either absolute or relative) the library containing this scheme
307 .  name_create - name of routine to create method context
308 -  routine_create - routine to create method context
309 
310    Notes:
311    TSGLAcceptRegisterDynamic() may be called multiple times to add several user-defined families.
312 
313    If dynamic libraries are used, then the fourth input argument (routine_create)
314    is ignored.
315 
316    Sample usage:
317 .vb
318    TSGLAcceptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
319                              "MySchemeCreate",MySchemeCreate);
320 .ve
321 
322    Then, your scheme can be chosen with the procedural interface via
323 $     TSGLSetAcceptType(ts,"my_scheme")
324    or at runtime via the option
325 $     -ts_gl_accept_type my_scheme
326 
327    Level: advanced
328 
329    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
330           and others of the form ${any_environmental_variable} occuring in pathname will be
331           replaced with appropriate values.
332 
333 .keywords: TSGL, TSGLAcceptType, register
334 
335 .seealso: TSGLRegisterAll()
336 M*/
337 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
338 #  define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,0)
339 #else
340 #  define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,d)
341 #endif
342 
343 /*E
344   TSGLType - family of time integration method within the General Linear class
345 
346   Level: beginner
347 
348 .seealso: TSGLSetType(), TSGLRegister()
349 E*/
350 #define TSGLType char*
351 #define TSGL_IRKS   "irks"
352 
353 /*MC
354    TSGLRegisterDynamic - adds a TSGL implementation
355 
356    Synopsis:
357    PetscErrorCode TSGLRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
358 
359    Not Collective
360 
361    Input Parameters:
362 +  name_scheme - name of user-defined general linear scheme
363 .  path - path (either absolute or relative) the library containing this scheme
364 .  name_create - name of routine to create method context
365 -  routine_create - routine to create method context
366 
367    Notes:
368    TSGLRegisterDynamic() may be called multiple times to add several user-defined families.
369 
370    If dynamic libraries are used, then the fourth input argument (routine_create)
371    is ignored.
372 
373    Sample usage:
374 .vb
375    TSGLRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
376                        "MySchemeCreate",MySchemeCreate);
377 .ve
378 
379    Then, your scheme can be chosen with the procedural interface via
380 $     TSGLSetType(ts,"my_scheme")
381    or at runtime via the option
382 $     -ts_gl_type my_scheme
383 
384    Level: advanced
385 
386    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
387           and others of the form ${any_environmental_variable} occuring in pathname will be
388           replaced with appropriate values.
389 
390 .keywords: TSGL, register
391 
392 .seealso: TSGLRegisterAll()
393 M*/
394 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
395 #  define TSGLRegisterDynamic(a,b,c,d)       TSGLRegister(a,b,c,0)
396 #else
397 #  define TSGLRegisterDynamic(a,b,c,d)       TSGLRegister(a,b,c,d)
398 #endif
399 
400 extern PetscErrorCode   TSGLRegister(const char[],const char[],const char[],PetscErrorCode(*)(TS));
401 extern PetscErrorCode  TSGLRegisterAll(const char[]);
402 extern PetscErrorCode  TSGLRegisterDestroy(void);
403 extern PetscErrorCode  TSGLInitializePackage(const char[]);
404 extern PetscErrorCode  TSGLFinalizePackage(void);
405 extern PetscErrorCode  TSGLSetType(TS,const TSGLType);
406 extern PetscErrorCode  TSGLGetAdapt(TS,TSGLAdapt*);
407 extern PetscErrorCode  TSGLSetAcceptType(TS,const TSGLAcceptType);
408 
409 /*
410        PETSc interface to Sundials
411 */
412 #ifdef PETSC_HAVE_SUNDIALS
413 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
414 extern const char *TSSundialsLmmTypes[];
415 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
416 extern const char *TSSundialsGramSchmidtTypes[];
417 extern PetscErrorCode   TSSundialsSetType(TS,TSSundialsLmmType);
418 extern PetscErrorCode   TSSundialsGetPC(TS,PC*);
419 extern PetscErrorCode   TSSundialsSetTolerance(TS,PetscReal,PetscReal);
420 extern PetscErrorCode   TSSundialsSetMinTimeStep(TS,PetscReal);
421 extern PetscErrorCode   TSSundialsSetMaxTimeStep(TS,PetscReal);
422 extern PetscErrorCode   TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
423 extern PetscErrorCode   TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
424 extern PetscErrorCode   TSSundialsSetGMRESRestart(TS,PetscInt);
425 extern PetscErrorCode   TSSundialsSetLinearTolerance(TS,PetscReal);
426 extern PetscErrorCode   TSSundialsSetExactFinalTime(TS,PetscBool );
427 extern PetscErrorCode   TSSundialsMonitorInternalSteps(TS,PetscBool );
428 extern PetscErrorCode   TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
429 #endif
430 
431 extern PetscErrorCode   TSRKSetTolerance(TS,PetscReal);
432 
433 extern PetscErrorCode  TSThetaSetTheta(TS,PetscReal);
434 extern PetscErrorCode  TSThetaGetTheta(TS,PetscReal*);
435 
436 extern PetscErrorCode  TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*);
437 extern PetscErrorCode  TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*);
438 extern PetscErrorCode  TSAlphaSetRadius(TS,PetscReal);
439 extern PetscErrorCode  TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal);
440 extern PetscErrorCode  TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*);
441 
442 extern PetscErrorCode  TSSetDM(TS,DM);
443 extern PetscErrorCode  TSGetDM(TS,DM*);
444 
445 extern PetscErrorCode  SNESTSFormFunction(SNES,Vec,Vec,void*);
446 extern PetscErrorCode  SNESTSFormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
447 
448 PETSC_EXTERN_CXX_END
449 #endif
450