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