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