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