xref: /petsc/include/petscts.h (revision 5a3a76d07e72bd733c23cdbf132d31c2a5ff3d88)
1818ad0c1SBarry Smith /*
2316643e7SJed Brown    User interface for the timestepping package. This package
3f64a0f93SLois Curfman McInnes    is for use in solving time-dependent PDEs.
4818ad0c1SBarry Smith */
50a835dfdSSatish Balay #if !defined(__PETSCTS_H)
60a835dfdSSatish Balay #define __PETSCTS_H
70a835dfdSSatish Balay #include "petscsnes.h"
8e9fa29b7SSatish Balay PETSC_EXTERN_CXX_BEGIN
9818ad0c1SBarry Smith 
10435da068SBarry Smith /*S
11435da068SBarry Smith      TS - Abstract PETSc object that manages all time-steppers (ODE integrators)
12435da068SBarry Smith 
13435da068SBarry Smith    Level: beginner
14435da068SBarry Smith 
15435da068SBarry Smith   Concepts: ODE solvers
16435da068SBarry Smith 
1794b7f48cSBarry Smith .seealso:  TSCreate(), TSSetType(), TSType, SNES, KSP, PC
18435da068SBarry Smith S*/
19f09e8eb9SSatish Balay typedef struct _p_TS* TS;
20435da068SBarry Smith 
21435da068SBarry Smith /*E
22435da068SBarry Smith     TSType - String with the name of a PETSc TS method or the creation function
23435da068SBarry Smith        with an optional dynamic library name, for example
24435da068SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:mytscreate()
25435da068SBarry Smith 
26435da068SBarry Smith    Level: beginner
27435da068SBarry Smith 
28435da068SBarry Smith .seealso: TSSetType(), TS
29435da068SBarry Smith E*/
30a313700dSBarry Smith #define TSType char*
319596e0b4SJed Brown #define TSEULER           "euler"
329596e0b4SJed Brown #define TSBEULER          "beuler"
339596e0b4SJed Brown #define TSPSEUDO          "pseudo"
344d91e141SJed Brown #define TSCN              "cn"
359596e0b4SJed Brown #define TSSUNDIALS        "sundials"
364d91e141SJed Brown #define TSRK              "rk"
379596e0b4SJed Brown #define TSPYTHON          "python"
389596e0b4SJed Brown #define TSTHETA           "theta"
3988df8a41SLisandro Dalcin #define TSALPHA           "alpha"
409596e0b4SJed Brown #define TSGL              "gl"
41ef7bb5aaSJed Brown #define TSSSP             "ssp"
428a381b04SJed Brown #define TSARKIMEX         "arkimex"
4382bf6240SBarry Smith 
44435da068SBarry Smith /*E
45435da068SBarry Smith     TSProblemType - Determines the type of problem this TS object is to be used to solve
46435da068SBarry Smith 
47435da068SBarry Smith    Level: beginner
48435da068SBarry Smith 
49435da068SBarry Smith .seealso: TSCreate()
50435da068SBarry Smith E*/
5119bcc07fSBarry Smith typedef enum {TS_LINEAR,TS_NONLINEAR} TSProblemType;
52818ad0c1SBarry Smith 
53193ac0bcSJed Brown typedef enum {
54193ac0bcSJed Brown   TS_CONVERGED_ITERATING      = 0,
55193ac0bcSJed Brown   TS_CONVERGED_TIME           = 1,
56193ac0bcSJed Brown   TS_CONVERGED_ITS            = 2,
57193ac0bcSJed Brown   TS_DIVERGED_NONLINEAR_SOLVE = -1,
58193ac0bcSJed Brown   TS_DIVERGED_STEP_REJECTED   = -2
59193ac0bcSJed Brown } TSConvergedReason;
60193ac0bcSJed Brown extern const char *const*TSConvergedReasons;
61193ac0bcSJed Brown 
62000e7ae3SMatthew Knepley /* Logging support */
637087cfbeSBarry Smith extern PetscClassId  TS_CLASSID;
64000e7ae3SMatthew Knepley 
657087cfbeSBarry Smith extern PetscErrorCode   TSInitializePackage(const char[]);
668ba1e511SMatthew Knepley 
677087cfbeSBarry Smith extern PetscErrorCode   TSCreate(MPI_Comm,TS*);
68fcfd50ebSBarry Smith extern PetscErrorCode   TSDestroy(TS*);
69818ad0c1SBarry Smith 
707087cfbeSBarry Smith extern PetscErrorCode   TSSetProblemType(TS,TSProblemType);
717087cfbeSBarry Smith extern PetscErrorCode   TSGetProblemType(TS,TSProblemType*);
72c2efdce3SBarry Smith extern PetscErrorCode   TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**));
737087cfbeSBarry Smith extern PetscErrorCode   TSMonitorCancel(TS);
74818ad0c1SBarry Smith 
757087cfbeSBarry Smith extern PetscErrorCode   TSSetOptionsPrefix(TS,const char[]);
767087cfbeSBarry Smith extern PetscErrorCode   TSAppendOptionsPrefix(TS,const char[]);
777087cfbeSBarry Smith extern PetscErrorCode   TSGetOptionsPrefix(TS,const char *[]);
787087cfbeSBarry Smith extern PetscErrorCode   TSSetFromOptions(TS);
797087cfbeSBarry Smith extern PetscErrorCode   TSSetUp(TS);
80277b19d0SLisandro Dalcin extern PetscErrorCode   TSReset(TS);
81818ad0c1SBarry Smith 
827087cfbeSBarry Smith extern PetscErrorCode   TSSetSolution(TS,Vec);
837087cfbeSBarry Smith extern PetscErrorCode   TSGetSolution(TS,Vec*);
84818ad0c1SBarry Smith 
857087cfbeSBarry Smith extern PetscErrorCode   TSSetDuration(TS,PetscInt,PetscReal);
867087cfbeSBarry Smith extern PetscErrorCode   TSGetDuration(TS,PetscInt*,PetscReal*);
87a43b19c4SJed Brown extern PetscErrorCode   TSSetExactFinalTime(TS,PetscBool);
88818ad0c1SBarry Smith 
897087cfbeSBarry Smith extern PetscErrorCode   TSMonitorDefault(TS,PetscInt,PetscReal,Vec,void*);
907087cfbeSBarry Smith extern PetscErrorCode   TSMonitorSolution(TS,PetscInt,PetscReal,Vec,void*);
91193ac0bcSJed Brown extern PetscErrorCode   TSStep(TS);
92*5a3a76d0SJed Brown extern PetscErrorCode   TSSolve(TS,Vec,PetscReal*);
93193ac0bcSJed Brown extern PetscErrorCode   TSGetConvergedReason(TS,TSConvergedReason*);
94818ad0c1SBarry Smith 
957087cfbeSBarry Smith extern PetscErrorCode   TSSetInitialTimeStep(TS,PetscReal,PetscReal);
967087cfbeSBarry Smith extern PetscErrorCode   TSGetTimeStep(TS,PetscReal*);
977087cfbeSBarry Smith extern PetscErrorCode   TSGetTime(TS,PetscReal*);
987087cfbeSBarry Smith extern PetscErrorCode   TSSetTime(TS,PetscReal);
997087cfbeSBarry Smith extern PetscErrorCode   TSGetTimeStepNumber(TS,PetscInt*);
1007087cfbeSBarry Smith extern PetscErrorCode   TSSetTimeStep(TS,PetscReal);
101818ad0c1SBarry Smith 
1022d4b5aceSLisandro Dalcin typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*);
1032d4b5aceSLisandro Dalcin typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
104089b2837SJed Brown extern PetscErrorCode   TSSetRHSFunction(TS,Vec,TSRHSFunction,void*);
105089b2837SJed Brown extern PetscErrorCode   TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**);
1062d4b5aceSLisandro Dalcin extern PetscErrorCode   TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*);
107089b2837SJed Brown extern PetscErrorCode   TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**);
108818ad0c1SBarry Smith 
109316643e7SJed Brown typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*);
110316643e7SJed Brown typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*);
111089b2837SJed Brown extern PetscErrorCode   TSSetIFunction(TS,Vec,TSIFunction,void*);
112089b2837SJed Brown extern PetscErrorCode   TSGetIFunction(TS,Vec*,TSIFunction*,void**);
1137087cfbeSBarry Smith extern PetscErrorCode   TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*);
1147087cfbeSBarry Smith extern PetscErrorCode   TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**);
115316643e7SJed Brown 
1160e4ef248SJed Brown extern PetscErrorCode   TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*);
1170e4ef248SJed Brown extern PetscErrorCode   TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
1180026cea9SSean Farley extern PetscErrorCode   TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*);
1190026cea9SSean Farley extern PetscErrorCode   TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*);
1207087cfbeSBarry Smith extern PetscErrorCode   TSDefaultComputeJacobianColor(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
1217087cfbeSBarry Smith extern PetscErrorCode   TSDefaultComputeJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
122e34be4c2SBarry Smith 
1237087cfbeSBarry Smith extern PetscErrorCode   TSSetPreStep(TS, PetscErrorCode (*)(TS));
1247087cfbeSBarry Smith extern PetscErrorCode   TSSetPostStep(TS, PetscErrorCode (*)(TS));
1257087cfbeSBarry Smith extern PetscErrorCode   TSPreStep(TS);
1267087cfbeSBarry Smith extern PetscErrorCode   TSPostStep(TS);
127cd652676SJed Brown extern PetscErrorCode   TSSetRetainStages(TS,PetscBool);
128cd652676SJed Brown extern PetscErrorCode   TSInterpolate(TS,PetscReal,Vec);
129000e7ae3SMatthew Knepley 
1307087cfbeSBarry Smith extern PetscErrorCode   TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*);
1317087cfbeSBarry Smith extern PetscErrorCode   TSPseudoDefaultTimeStep(TS,PetscReal*,void*);
1327087cfbeSBarry Smith extern PetscErrorCode   TSPseudoComputeTimeStep(TS,PetscReal *);
133d8345c25SBarry Smith 
1347087cfbeSBarry Smith extern PetscErrorCode   TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*);
1357087cfbeSBarry Smith extern PetscErrorCode   TSPseudoDefaultVerifyTimeStep(TS,Vec,void*,PetscReal*,PetscBool *);
1367087cfbeSBarry Smith extern PetscErrorCode   TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *);
1377087cfbeSBarry Smith extern PetscErrorCode   TSPseudoSetTimeStepIncrement(TS,PetscReal);
1387087cfbeSBarry Smith extern PetscErrorCode   TSPseudoIncrementDtFromInitialDt(TS);
13921c89e3eSBarry Smith 
1407087cfbeSBarry Smith extern PetscErrorCode   TSPythonSetType(TS,const char[]);
1411d6018f0SLisandro Dalcin 
1427087cfbeSBarry Smith extern PetscErrorCode   TSComputeRHSFunction(TS,PetscReal,Vec,Vec);
1437087cfbeSBarry Smith extern PetscErrorCode   TSComputeRHSJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*);
144214bc6a2SJed Brown extern PetscErrorCode   TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool);
145214bc6a2SJed Brown extern PetscErrorCode   TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,PetscBool);
146818ad0c1SBarry Smith 
147bdad233fSMatthew Knepley /* Dynamic creation and loading functions */
148b0a32e0cSBarry Smith extern PetscFList TSList;
149ace3abfcSBarry Smith extern PetscBool  TSRegisterAllCalled;
1507087cfbeSBarry Smith extern PetscErrorCode   TSGetType(TS,const TSType*);
1517087cfbeSBarry Smith extern PetscErrorCode   TSSetType(TS,const TSType);
1527087cfbeSBarry Smith extern PetscErrorCode   TSRegister(const char[], const char[], const char[], PetscErrorCode (*)(TS));
1537087cfbeSBarry Smith extern PetscErrorCode   TSRegisterAll(const char[]);
1547087cfbeSBarry Smith extern PetscErrorCode   TSRegisterDestroy(void);
15530de9b25SBarry Smith 
15630de9b25SBarry Smith /*MC
15730de9b25SBarry Smith   TSRegisterDynamic - Adds a creation method to the TS package.
15830de9b25SBarry Smith 
15930de9b25SBarry Smith   Synopsis:
1601890ba74SBarry Smith   PetscErrorCode TSRegisterDynamic(const char *name, const char *path, const char *func_name, PetscErrorCode (*create_func)(TS))
16130de9b25SBarry Smith 
16230de9b25SBarry Smith   Not Collective
16330de9b25SBarry Smith 
16430de9b25SBarry Smith   Input Parameters:
16530de9b25SBarry Smith + name        - The name of a new user-defined creation routine
16630de9b25SBarry Smith . path        - The path (either absolute or relative) of the library containing this routine
16730de9b25SBarry Smith . func_name   - The name of the creation routine
16830de9b25SBarry Smith - create_func - The creation routine itself
16930de9b25SBarry Smith 
17030de9b25SBarry Smith   Notes:
17130de9b25SBarry Smith   TSRegisterDynamic() may be called multiple times to add several user-defined tses.
17230de9b25SBarry Smith 
17330de9b25SBarry Smith   If dynamic libraries are used, then the fourth input argument (create_func) is ignored.
17430de9b25SBarry Smith 
17530de9b25SBarry Smith   Sample usage:
17630de9b25SBarry Smith .vb
17730de9b25SBarry Smith   TSRegisterDynamic("my_ts", "/home/username/my_lib/lib/libO/solaris/libmy.a", "MyTSCreate", MyTSCreate);
17830de9b25SBarry Smith .ve
17930de9b25SBarry Smith 
18030de9b25SBarry Smith   Then, your ts type can be chosen with the procedural interface via
18130de9b25SBarry Smith .vb
18265f2c13aSKai Germaschewski     TS ts;
18365f2c13aSKai Germaschewski     TSCreate(MPI_Comm, &ts);
18465f2c13aSKai Germaschewski     TSSetType(ts, "my_ts")
18530de9b25SBarry Smith .ve
18630de9b25SBarry Smith   or at runtime via the option
18730de9b25SBarry Smith .vb
18830de9b25SBarry Smith     -ts_type my_ts
18930de9b25SBarry Smith .ve
19030de9b25SBarry Smith 
191ab901514SBarry Smith   Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
19230de9b25SBarry Smith         If your function is not being put into a shared library then use TSRegister() instead
19330de9b25SBarry Smith 
19430de9b25SBarry Smith   Level: advanced
19530de9b25SBarry Smith 
19630de9b25SBarry Smith .keywords: TS, register
19730de9b25SBarry Smith .seealso: TSRegisterAll(), TSRegisterDestroy()
19830de9b25SBarry Smith M*/
199aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
200f1af5d2fSBarry Smith #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,0)
2016df38c32SLois Curfman McInnes #else
202f1af5d2fSBarry Smith #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,d)
2036df38c32SLois Curfman McInnes #endif
2046df38c32SLois Curfman McInnes 
2057087cfbeSBarry Smith extern PetscErrorCode   TSGetSNES(TS,SNES*);
2067087cfbeSBarry Smith extern PetscErrorCode   TSGetKSP(TS,KSP*);
207818ad0c1SBarry Smith 
2087087cfbeSBarry Smith extern PetscErrorCode   TSView(TS,PetscViewer);
20921c89e3eSBarry Smith 
2107087cfbeSBarry Smith extern PetscErrorCode   TSSetApplicationContext(TS,void *);
211e71120c6SJed Brown extern PetscErrorCode   TSGetApplicationContext(TS,void *);
21221c89e3eSBarry Smith 
2137087cfbeSBarry Smith extern PetscErrorCode   TSMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG *);
2147087cfbeSBarry Smith extern PetscErrorCode   TSMonitorLG(TS,PetscInt,PetscReal,Vec,void *);
215fcfd50ebSBarry Smith extern PetscErrorCode   TSMonitorLGDestroy(PetscDrawLG*);
2163914022bSBarry Smith 
217ed657a08SJed Brown /*S
218ed657a08SJed Brown    TSGLAdapt - Abstract object that manages time-step adaptivity
219ed657a08SJed Brown 
220ed657a08SJed Brown    Level: beginner
221ed657a08SJed Brown 
222ed657a08SJed Brown .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType
223ed657a08SJed Brown S*/
224ed657a08SJed Brown typedef struct _p_TSGLAdapt *TSGLAdapt;
225ed657a08SJed Brown 
226ed657a08SJed Brown /*E
227ed657a08SJed Brown     TSGLAdaptType - String with the name of TSGLAdapt scheme or the creation function
228ed657a08SJed Brown        with an optional dynamic library name, for example
229ed657a08SJed Brown        http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate()
230ed657a08SJed Brown 
231ed657a08SJed Brown    Level: beginner
232ed657a08SJed Brown 
233ed657a08SJed Brown .seealso: TSGLAdaptSetType(), TS
234ed657a08SJed Brown E*/
235ed657a08SJed Brown #define TSGLAdaptType  char*
236ed657a08SJed Brown #define TSGLADAPT_NONE "none"
237ed657a08SJed Brown #define TSGLADAPT_SIZE "size"
238ed657a08SJed Brown #define TSGLADAPT_BOTH "both"
239ed657a08SJed Brown 
240ed657a08SJed Brown /*MC
241ed657a08SJed Brown    TSGLAdaptRegisterDynamic - adds a TSGLAdapt implementation
242ed657a08SJed Brown 
243ed657a08SJed Brown    Synopsis:
2441890ba74SBarry Smith    PetscErrorCode TSGLAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
245ed657a08SJed Brown 
246ed657a08SJed Brown    Not Collective
247ed657a08SJed Brown 
248ed657a08SJed Brown    Input Parameters:
249ed657a08SJed Brown +  name_scheme - name of user-defined adaptivity scheme
250ed657a08SJed Brown .  path - path (either absolute or relative) the library containing this scheme
251ed657a08SJed Brown .  name_create - name of routine to create method context
252ed657a08SJed Brown -  routine_create - routine to create method context
253ed657a08SJed Brown 
254ed657a08SJed Brown    Notes:
255ed657a08SJed Brown    TSGLAdaptRegisterDynamic() may be called multiple times to add several user-defined families.
256ed657a08SJed Brown 
257ed657a08SJed Brown    If dynamic libraries are used, then the fourth input argument (routine_create)
258ed657a08SJed Brown    is ignored.
259ed657a08SJed Brown 
260ed657a08SJed Brown    Sample usage:
261ed657a08SJed Brown .vb
262ed657a08SJed Brown    TSGLAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
263ed657a08SJed Brown                             "MySchemeCreate",MySchemeCreate);
264ed657a08SJed Brown .ve
265ed657a08SJed Brown 
266ed657a08SJed Brown    Then, your scheme can be chosen with the procedural interface via
267ed657a08SJed Brown $     TSGLAdaptSetType(ts,"my_scheme")
268ed657a08SJed Brown    or at runtime via the option
269ed657a08SJed Brown $     -ts_adapt_type my_scheme
270ed657a08SJed Brown 
271ed657a08SJed Brown    Level: advanced
272ed657a08SJed Brown 
273ed657a08SJed Brown    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
274ed657a08SJed Brown           and others of the form ${any_environmental_variable} occuring in pathname will be
275ed657a08SJed Brown           replaced with appropriate values.
276ed657a08SJed Brown 
277ed657a08SJed Brown .keywords: TSGLAdapt, register
278ed657a08SJed Brown 
279ed657a08SJed Brown .seealso: TSGLAdaptRegisterAll()
280ed657a08SJed Brown M*/
281ed657a08SJed Brown #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
282ed657a08SJed Brown #  define TSGLAdaptRegisterDynamic(a,b,c,d)  TSGLAdaptRegister(a,b,c,0)
283ed657a08SJed Brown #else
284ed657a08SJed Brown #  define TSGLAdaptRegisterDynamic(a,b,c,d)  TSGLAdaptRegister(a,b,c,d)
285ed657a08SJed Brown #endif
286ed657a08SJed Brown 
2877087cfbeSBarry Smith extern PetscErrorCode  TSGLAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSGLAdapt));
2887087cfbeSBarry Smith extern PetscErrorCode  TSGLAdaptRegisterAll(const char[]);
2897087cfbeSBarry Smith extern PetscErrorCode  TSGLAdaptRegisterDestroy(void);
2907087cfbeSBarry Smith extern PetscErrorCode  TSGLAdaptInitializePackage(const char[]);
2917087cfbeSBarry Smith extern PetscErrorCode  TSGLAdaptFinalizePackage(void);
2927087cfbeSBarry Smith extern PetscErrorCode  TSGLAdaptCreate(MPI_Comm,TSGLAdapt*);
2937087cfbeSBarry Smith extern PetscErrorCode  TSGLAdaptSetType(TSGLAdapt,const TSGLAdaptType);
2947087cfbeSBarry Smith extern PetscErrorCode  TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]);
2957087cfbeSBarry Smith extern PetscErrorCode  TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
2967087cfbeSBarry Smith extern PetscErrorCode  TSGLAdaptView(TSGLAdapt,PetscViewer);
2977087cfbeSBarry Smith extern PetscErrorCode  TSGLAdaptSetFromOptions(TSGLAdapt);
298fcfd50ebSBarry Smith extern PetscErrorCode  TSGLAdaptDestroy(TSGLAdapt*);
299ed657a08SJed Brown 
300ed657a08SJed Brown /*E
301ed657a08SJed Brown     TSGLAcceptType - String with the name of TSGLAccept scheme or the function
302ed657a08SJed Brown        with an optional dynamic library name, for example
303ed657a08SJed Brown        http://www.mcs.anl.gov/petsc/lib.a:mytsglaccept()
304ed657a08SJed Brown 
305ed657a08SJed Brown    Level: beginner
306ed657a08SJed Brown 
307ed657a08SJed Brown .seealso: TSGLSetAcceptType(), TS
308ed657a08SJed Brown E*/
309ed657a08SJed Brown #define TSGLAcceptType  char*
310ed657a08SJed Brown #define TSGLACCEPT_ALWAYS "always"
311ed657a08SJed Brown 
312ace3abfcSBarry Smith typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *);
3137087cfbeSBarry Smith extern PetscErrorCode  TSGLAcceptRegister(const char[],const char[],const char[],TSGLAcceptFunction);
314ed657a08SJed Brown 
315ed657a08SJed Brown /*MC
316ed657a08SJed Brown    TSGLAcceptRegisterDynamic - adds a TSGL acceptance scheme
317ed657a08SJed Brown 
318ed657a08SJed Brown    Synopsis:
3191890ba74SBarry Smith    PetscErrorCode TSGLAcceptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
320ed657a08SJed Brown 
321ed657a08SJed Brown    Not Collective
322ed657a08SJed Brown 
323ed657a08SJed Brown    Input Parameters:
324ed657a08SJed Brown +  name_scheme - name of user-defined acceptance scheme
325ed657a08SJed Brown .  path - path (either absolute or relative) the library containing this scheme
326ed657a08SJed Brown .  name_create - name of routine to create method context
327ed657a08SJed Brown -  routine_create - routine to create method context
328ed657a08SJed Brown 
329ed657a08SJed Brown    Notes:
330ed657a08SJed Brown    TSGLAcceptRegisterDynamic() may be called multiple times to add several user-defined families.
331ed657a08SJed Brown 
332ed657a08SJed Brown    If dynamic libraries are used, then the fourth input argument (routine_create)
333ed657a08SJed Brown    is ignored.
334ed657a08SJed Brown 
335ed657a08SJed Brown    Sample usage:
336ed657a08SJed Brown .vb
337ed657a08SJed Brown    TSGLAcceptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
338ed657a08SJed Brown                              "MySchemeCreate",MySchemeCreate);
339ed657a08SJed Brown .ve
340ed657a08SJed Brown 
341ed657a08SJed Brown    Then, your scheme can be chosen with the procedural interface via
342ed657a08SJed Brown $     TSGLSetAcceptType(ts,"my_scheme")
343ed657a08SJed Brown    or at runtime via the option
344ed657a08SJed Brown $     -ts_gl_accept_type my_scheme
345ed657a08SJed Brown 
346ed657a08SJed Brown    Level: advanced
347ed657a08SJed Brown 
348ed657a08SJed Brown    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
349ed657a08SJed Brown           and others of the form ${any_environmental_variable} occuring in pathname will be
350ed657a08SJed Brown           replaced with appropriate values.
351ed657a08SJed Brown 
352ed657a08SJed Brown .keywords: TSGL, TSGLAcceptType, register
353ed657a08SJed Brown 
354ed657a08SJed Brown .seealso: TSGLRegisterAll()
355ed657a08SJed Brown M*/
356ed657a08SJed Brown #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
357ed657a08SJed Brown #  define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,0)
358ed657a08SJed Brown #else
359ed657a08SJed Brown #  define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,d)
360ed657a08SJed Brown #endif
361ed657a08SJed Brown 
36218b56cb9SJed Brown /*E
36318b56cb9SJed Brown   TSGLType - family of time integration method within the General Linear class
36418b56cb9SJed Brown 
36518b56cb9SJed Brown   Level: beginner
36618b56cb9SJed Brown 
36718b56cb9SJed Brown .seealso: TSGLSetType(), TSGLRegister()
36818b56cb9SJed Brown E*/
36918b56cb9SJed Brown #define TSGLType char*
370ed657a08SJed Brown #define TSGL_IRKS   "irks"
37118b56cb9SJed Brown 
372ed657a08SJed Brown /*MC
373ed657a08SJed Brown    TSGLRegisterDynamic - adds a TSGL implementation
37418b56cb9SJed Brown 
375ed657a08SJed Brown    Synopsis:
3761890ba74SBarry Smith    PetscErrorCode TSGLRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
377ed657a08SJed Brown 
378ed657a08SJed Brown    Not Collective
379ed657a08SJed Brown 
380ed657a08SJed Brown    Input Parameters:
381ed657a08SJed Brown +  name_scheme - name of user-defined general linear scheme
382ed657a08SJed Brown .  path - path (either absolute or relative) the library containing this scheme
383ed657a08SJed Brown .  name_create - name of routine to create method context
384ed657a08SJed Brown -  routine_create - routine to create method context
385ed657a08SJed Brown 
386ed657a08SJed Brown    Notes:
387ed657a08SJed Brown    TSGLRegisterDynamic() may be called multiple times to add several user-defined families.
388ed657a08SJed Brown 
389ed657a08SJed Brown    If dynamic libraries are used, then the fourth input argument (routine_create)
390ed657a08SJed Brown    is ignored.
391ed657a08SJed Brown 
392ed657a08SJed Brown    Sample usage:
393ed657a08SJed Brown .vb
394ed657a08SJed Brown    TSGLRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
395ed657a08SJed Brown                        "MySchemeCreate",MySchemeCreate);
396ed657a08SJed Brown .ve
397ed657a08SJed Brown 
398ed657a08SJed Brown    Then, your scheme can be chosen with the procedural interface via
399ed657a08SJed Brown $     TSGLSetType(ts,"my_scheme")
400ed657a08SJed Brown    or at runtime via the option
401ed657a08SJed Brown $     -ts_gl_type my_scheme
402ed657a08SJed Brown 
403ed657a08SJed Brown    Level: advanced
404ed657a08SJed Brown 
405ed657a08SJed Brown    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
406ed657a08SJed Brown           and others of the form ${any_environmental_variable} occuring in pathname will be
407ed657a08SJed Brown           replaced with appropriate values.
408ed657a08SJed Brown 
409ed657a08SJed Brown .keywords: TSGL, register
410ed657a08SJed Brown 
411ed657a08SJed Brown .seealso: TSGLRegisterAll()
412ed657a08SJed Brown M*/
413ed657a08SJed Brown #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
414ed657a08SJed Brown #  define TSGLRegisterDynamic(a,b,c,d)       TSGLRegister(a,b,c,0)
415ed657a08SJed Brown #else
416ed657a08SJed Brown #  define TSGLRegisterDynamic(a,b,c,d)       TSGLRegister(a,b,c,d)
417ed657a08SJed Brown #endif
418ed657a08SJed Brown 
4197087cfbeSBarry Smith extern PetscErrorCode  TSGLRegister(const char[],const char[],const char[],PetscErrorCode(*)(TS));
4207087cfbeSBarry Smith extern PetscErrorCode  TSGLRegisterAll(const char[]);
4217087cfbeSBarry Smith extern PetscErrorCode  TSGLRegisterDestroy(void);
4227087cfbeSBarry Smith extern PetscErrorCode  TSGLInitializePackage(const char[]);
4237087cfbeSBarry Smith extern PetscErrorCode  TSGLFinalizePackage(void);
4247087cfbeSBarry Smith extern PetscErrorCode  TSGLSetType(TS,const TSGLType);
4257087cfbeSBarry Smith extern PetscErrorCode  TSGLGetAdapt(TS,TSGLAdapt*);
4267087cfbeSBarry Smith extern PetscErrorCode  TSGLSetAcceptType(TS,const TSGLAcceptType);
42718b56cb9SJed Brown 
4288a381b04SJed Brown #define TSARKIMEXType char*
4298a381b04SJed Brown #define TSARKIMEX2D "2d"
430a3a57f36SJed Brown #define TSARKIMEX2E "2e"
4318a381b04SJed Brown #define TSARKIMEX3  "3"
4328a381b04SJed Brown #define TSARKIMEX4  "4"
4338a381b04SJed Brown #define TSARKIMEX5  "5"
4348a381b04SJed Brown extern PetscErrorCode TSARKIMEXGetType(TS ts,const TSARKIMEXType*);
4358a381b04SJed Brown extern PetscErrorCode TSARKIMEXSetType(TS ts,const TSARKIMEXType);
436cd652676SJed Brown extern PetscErrorCode TSARKIMEXRegister(const TSARKIMEXType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[],const PetscReal[]);
4378a381b04SJed Brown extern PetscErrorCode TSARKIMEXFinalizePackage(void);
4388a381b04SJed Brown extern PetscErrorCode TSARKIMEXInitializePackage(const char path[]);
4398a381b04SJed Brown extern PetscErrorCode TSARKIMEXRegisterDestroy(void);
4408a381b04SJed Brown extern PetscErrorCode TSARKIMEXRegisterAll(void);
4418a381b04SJed Brown 
44283e2fdc7SBarry Smith /*
4430f3b3ca1SHong Zhang        PETSc interface to Sundials
44483e2fdc7SBarry Smith */
4456dc17ff5SMatthew Knepley #ifdef PETSC_HAVE_SUNDIALS
4466fadb2cdSHong Zhang typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
4476fadb2cdSHong Zhang extern const char *TSSundialsLmmTypes[];
4486fadb2cdSHong Zhang typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
449a04cf4d8SBarry Smith extern const char *TSSundialsGramSchmidtTypes[];
4507087cfbeSBarry Smith extern PetscErrorCode   TSSundialsSetType(TS,TSSundialsLmmType);
4517087cfbeSBarry Smith extern PetscErrorCode   TSSundialsGetPC(TS,PC*);
4527087cfbeSBarry Smith extern PetscErrorCode   TSSundialsSetTolerance(TS,PetscReal,PetscReal);
4537087cfbeSBarry Smith extern PetscErrorCode   TSSundialsSetMinTimeStep(TS,PetscReal);
4547087cfbeSBarry Smith extern PetscErrorCode   TSSundialsSetMaxTimeStep(TS,PetscReal);
4557087cfbeSBarry Smith extern PetscErrorCode   TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
4567087cfbeSBarry Smith extern PetscErrorCode   TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
4577087cfbeSBarry Smith extern PetscErrorCode   TSSundialsSetGMRESRestart(TS,PetscInt);
4587087cfbeSBarry Smith extern PetscErrorCode   TSSundialsSetLinearTolerance(TS,PetscReal);
4597087cfbeSBarry Smith extern PetscErrorCode   TSSundialsMonitorInternalSteps(TS,PetscBool );
4607087cfbeSBarry Smith extern PetscErrorCode   TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
4616dc17ff5SMatthew Knepley #endif
46283e2fdc7SBarry Smith 
4637087cfbeSBarry Smith extern PetscErrorCode   TSRKSetTolerance(TS,PetscReal);
464262ff7bbSBarry Smith 
4657087cfbeSBarry Smith extern PetscErrorCode  TSThetaSetTheta(TS,PetscReal);
4667087cfbeSBarry Smith extern PetscErrorCode  TSThetaGetTheta(TS,PetscReal*);
467eb284becSJed Brown extern PetscErrorCode  TSThetaSetEndpoint(TS,PetscBool);
4680de4c49aSJed Brown 
4697087cfbeSBarry Smith extern PetscErrorCode  TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*);
4707087cfbeSBarry Smith extern PetscErrorCode  TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*);
4717087cfbeSBarry Smith extern PetscErrorCode  TSAlphaSetRadius(TS,PetscReal);
4727087cfbeSBarry Smith extern PetscErrorCode  TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal);
4737087cfbeSBarry Smith extern PetscErrorCode  TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*);
47488df8a41SLisandro Dalcin 
4757087cfbeSBarry Smith extern PetscErrorCode  TSSetDM(TS,DM);
4767087cfbeSBarry Smith extern PetscErrorCode  TSGetDM(TS,DM*);
4776c699258SBarry Smith 
4787087cfbeSBarry Smith extern PetscErrorCode  SNESTSFormFunction(SNES,Vec,Vec,void*);
4797087cfbeSBarry Smith extern PetscErrorCode  SNESTSFormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
4800f5c6efeSJed Brown 
481e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
482818ad0c1SBarry Smith #endif
483