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