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