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