xref: /petsc/include/petscts.h (revision b93fea0e7692c4b741e57ced863a893e4886d68d)
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 
2176bdecfbSBarry Smith /*J
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
2976bdecfbSBarry Smith J*/
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"
43e27a552bSJed Brown #define TSROSW            "rosw"
4482bf6240SBarry Smith 
45435da068SBarry Smith /*E
46435da068SBarry Smith     TSProblemType - Determines the type of problem this TS object is to be used to solve
47435da068SBarry Smith 
48435da068SBarry Smith    Level: beginner
49435da068SBarry Smith 
50435da068SBarry Smith .seealso: TSCreate()
51435da068SBarry Smith E*/
5219bcc07fSBarry Smith typedef enum {TS_LINEAR,TS_NONLINEAR} TSProblemType;
53818ad0c1SBarry Smith 
54*b93fea0eSJed Brown /*E
55*b93fea0eSJed Brown    TSConvergedReason - reason a TS method has converged or not
56*b93fea0eSJed Brown 
57*b93fea0eSJed Brown    Level: beginner
58*b93fea0eSJed Brown 
59*b93fea0eSJed Brown    Developer Notes: this must match finclude/petscts.h
60*b93fea0eSJed Brown 
61*b93fea0eSJed Brown    Each reason has its own manual page.
62*b93fea0eSJed Brown 
63*b93fea0eSJed Brown .seealso: TSGetConvergedReason()
64*b93fea0eSJed Brown E*/
65193ac0bcSJed Brown typedef enum {
66193ac0bcSJed Brown   TS_CONVERGED_ITERATING      = 0,
67193ac0bcSJed Brown   TS_CONVERGED_TIME           = 1,
68193ac0bcSJed Brown   TS_CONVERGED_ITS            = 2,
69193ac0bcSJed Brown   TS_DIVERGED_NONLINEAR_SOLVE = -1,
70193ac0bcSJed Brown   TS_DIVERGED_STEP_REJECTED   = -2
71193ac0bcSJed Brown } TSConvergedReason;
72193ac0bcSJed Brown extern const char *const*TSConvergedReasons;
73193ac0bcSJed Brown 
74*b93fea0eSJed Brown /*MC
75*b93fea0eSJed Brown    TS_CONVERGED_ITERATING - this only occurs if TSGetConvergedReason() is called during the TSSolve()
76*b93fea0eSJed Brown 
77*b93fea0eSJed Brown    Level: beginner
78*b93fea0eSJed Brown 
79*b93fea0eSJed Brown .seealso: TSSolve(), TSConvergedReason(), TSGetAdapt()
80*b93fea0eSJed Brown M*/
81*b93fea0eSJed Brown 
82*b93fea0eSJed Brown /*MC
83*b93fea0eSJed Brown    TS_CONVERGED_TIME - the final time was reached
84*b93fea0eSJed Brown 
85*b93fea0eSJed Brown    Level: beginner
86*b93fea0eSJed Brown 
87*b93fea0eSJed Brown .seealso: TSSolve(), TSConvergedReason(), TSGetAdapt(), TSSetDuration()
88*b93fea0eSJed Brown M*/
89*b93fea0eSJed Brown 
90*b93fea0eSJed Brown /*MC
91*b93fea0eSJed Brown    TS_CONVERGED_ITS - the maximum number of iterations was reached prior to the final time
92*b93fea0eSJed Brown 
93*b93fea0eSJed Brown    Level: beginner
94*b93fea0eSJed Brown 
95*b93fea0eSJed Brown .seealso: TSSolve(), TSConvergedReason(), TSGetAdapt(), TSSetDuration()
96*b93fea0eSJed Brown M*/
97*b93fea0eSJed Brown 
98*b93fea0eSJed Brown /*MC
99*b93fea0eSJed Brown    TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed
100*b93fea0eSJed Brown 
101*b93fea0eSJed Brown    Level: beginner
102*b93fea0eSJed Brown 
103*b93fea0eSJed Brown .seealso: TSSolve(), TSConvergedReason(), TSGetAdapt(), TSGetSNES(), SNESGetConvergedReason()
104*b93fea0eSJed Brown M*/
105*b93fea0eSJed Brown 
106*b93fea0eSJed Brown /*MC
107*b93fea0eSJed Brown    TS_DIVERGED_STEP_REJECTED - too many steps were rejected
108*b93fea0eSJed Brown 
109*b93fea0eSJed Brown    Level: beginner
110*b93fea0eSJed Brown 
111*b93fea0eSJed Brown .seealso: TSSolve(), TSConvergedReason(), TSGetAdapt()
112*b93fea0eSJed Brown M*/
113*b93fea0eSJed Brown 
114000e7ae3SMatthew Knepley /* Logging support */
1157087cfbeSBarry Smith extern PetscClassId  TS_CLASSID;
116000e7ae3SMatthew Knepley 
1177087cfbeSBarry Smith extern PetscErrorCode   TSInitializePackage(const char[]);
1188ba1e511SMatthew Knepley 
1197087cfbeSBarry Smith extern PetscErrorCode   TSCreate(MPI_Comm,TS*);
120fcfd50ebSBarry Smith extern PetscErrorCode   TSDestroy(TS*);
121818ad0c1SBarry Smith 
1227087cfbeSBarry Smith extern PetscErrorCode   TSSetProblemType(TS,TSProblemType);
1237087cfbeSBarry Smith extern PetscErrorCode   TSGetProblemType(TS,TSProblemType*);
124228d79bcSJed Brown extern PetscErrorCode   TSMonitor(TS,PetscInt,PetscReal,Vec);
125c2efdce3SBarry Smith extern PetscErrorCode   TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**));
1267087cfbeSBarry Smith extern PetscErrorCode   TSMonitorCancel(TS);
127818ad0c1SBarry Smith 
1287087cfbeSBarry Smith extern PetscErrorCode   TSSetOptionsPrefix(TS,const char[]);
1297087cfbeSBarry Smith extern PetscErrorCode   TSAppendOptionsPrefix(TS,const char[]);
1307087cfbeSBarry Smith extern PetscErrorCode   TSGetOptionsPrefix(TS,const char *[]);
1317087cfbeSBarry Smith extern PetscErrorCode   TSSetFromOptions(TS);
1327087cfbeSBarry Smith extern PetscErrorCode   TSSetUp(TS);
133277b19d0SLisandro Dalcin extern PetscErrorCode   TSReset(TS);
134818ad0c1SBarry Smith 
1357087cfbeSBarry Smith extern PetscErrorCode   TSSetSolution(TS,Vec);
1367087cfbeSBarry Smith extern PetscErrorCode   TSGetSolution(TS,Vec*);
137818ad0c1SBarry Smith 
1387087cfbeSBarry Smith extern PetscErrorCode   TSSetDuration(TS,PetscInt,PetscReal);
1397087cfbeSBarry Smith extern PetscErrorCode   TSGetDuration(TS,PetscInt*,PetscReal*);
140a43b19c4SJed Brown extern PetscErrorCode   TSSetExactFinalTime(TS,PetscBool);
141818ad0c1SBarry Smith 
1427087cfbeSBarry Smith extern PetscErrorCode   TSMonitorDefault(TS,PetscInt,PetscReal,Vec,void*);
1437087cfbeSBarry Smith extern PetscErrorCode   TSMonitorSolution(TS,PetscInt,PetscReal,Vec,void*);
1446083293cSBarry Smith extern PetscErrorCode   TSMonitorSolutionCreate(TS,PetscViewer,void**);
1456083293cSBarry Smith extern PetscErrorCode   TSMonitorSolutionDestroy(void**);
146fb1732b5SBarry Smith extern PetscErrorCode   TSMonitorSolutionBinary(TS,PetscInt,PetscReal,Vec,void*);
147ed81e22dSJed Brown extern PetscErrorCode   TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*);
148ed81e22dSJed Brown extern PetscErrorCode   TSMonitorSolutionVTKDestroy(void*);
149fb1732b5SBarry Smith 
150193ac0bcSJed Brown extern PetscErrorCode   TSStep(TS);
1511c3436cfSJed Brown extern PetscErrorCode   TSEvaluateStep(TS,PetscInt,Vec,PetscBool*);
1525a3a76d0SJed Brown extern PetscErrorCode   TSSolve(TS,Vec,PetscReal*);
153193ac0bcSJed Brown extern PetscErrorCode   TSGetConvergedReason(TS,TSConvergedReason*);
1549f67acb7SJed Brown extern PetscErrorCode   TSGetNonlinearSolveIterations(TS,PetscInt*);
1559f67acb7SJed Brown extern PetscErrorCode   TSGetLinearSolveIterations(TS,PetscInt*);
156818ad0c1SBarry Smith 
1577087cfbeSBarry Smith extern PetscErrorCode   TSSetInitialTimeStep(TS,PetscReal,PetscReal);
1587087cfbeSBarry Smith extern PetscErrorCode   TSGetTimeStep(TS,PetscReal*);
1597087cfbeSBarry Smith extern PetscErrorCode   TSGetTime(TS,PetscReal*);
1607087cfbeSBarry Smith extern PetscErrorCode   TSSetTime(TS,PetscReal);
1617087cfbeSBarry Smith extern PetscErrorCode   TSGetTimeStepNumber(TS,PetscInt*);
1627087cfbeSBarry Smith extern PetscErrorCode   TSSetTimeStep(TS,PetscReal);
163818ad0c1SBarry Smith 
1642d4b5aceSLisandro Dalcin typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*);
1652d4b5aceSLisandro Dalcin typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
166089b2837SJed Brown extern PetscErrorCode   TSSetRHSFunction(TS,Vec,TSRHSFunction,void*);
167089b2837SJed Brown extern PetscErrorCode   TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**);
1682d4b5aceSLisandro Dalcin extern PetscErrorCode   TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*);
169089b2837SJed Brown extern PetscErrorCode   TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**);
170818ad0c1SBarry Smith 
171316643e7SJed Brown typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*);
172316643e7SJed Brown typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*);
173089b2837SJed Brown extern PetscErrorCode   TSSetIFunction(TS,Vec,TSIFunction,void*);
174089b2837SJed Brown extern PetscErrorCode   TSGetIFunction(TS,Vec*,TSIFunction*,void**);
1757087cfbeSBarry Smith extern PetscErrorCode   TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*);
1767087cfbeSBarry Smith extern PetscErrorCode   TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**);
177316643e7SJed Brown 
1780e4ef248SJed Brown extern PetscErrorCode   TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*);
1790e4ef248SJed Brown extern PetscErrorCode   TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
1800026cea9SSean Farley extern PetscErrorCode   TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*);
1810026cea9SSean Farley extern PetscErrorCode   TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*);
1827087cfbeSBarry Smith extern PetscErrorCode   TSDefaultComputeJacobianColor(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
1837087cfbeSBarry Smith extern PetscErrorCode   TSDefaultComputeJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
184e34be4c2SBarry Smith 
1857087cfbeSBarry Smith extern PetscErrorCode   TSSetPreStep(TS, PetscErrorCode (*)(TS));
1867087cfbeSBarry Smith extern PetscErrorCode   TSSetPostStep(TS, PetscErrorCode (*)(TS));
1877087cfbeSBarry Smith extern PetscErrorCode   TSPreStep(TS);
1887087cfbeSBarry Smith extern PetscErrorCode   TSPostStep(TS);
189cd652676SJed Brown extern PetscErrorCode   TSSetRetainStages(TS,PetscBool);
190cd652676SJed Brown extern PetscErrorCode   TSInterpolate(TS,PetscReal,Vec);
1911c3436cfSJed Brown extern PetscErrorCode   TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec);
1921c3436cfSJed Brown extern PetscErrorCode   TSErrorNormWRMS(TS,Vec,PetscReal*);
1938d59e960SJed Brown extern PetscErrorCode   TSSetCFLTimeLocal(TS,PetscReal);
1948d59e960SJed Brown extern PetscErrorCode   TSGetCFLTime(TS,PetscReal*);
195000e7ae3SMatthew Knepley 
1967087cfbeSBarry Smith extern PetscErrorCode   TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*);
1977087cfbeSBarry Smith extern PetscErrorCode   TSPseudoDefaultTimeStep(TS,PetscReal*,void*);
1987087cfbeSBarry Smith extern PetscErrorCode   TSPseudoComputeTimeStep(TS,PetscReal *);
199ad24ef9bSSatish Balay extern PetscErrorCode   TSPseudoSetMaxTimeStep(TS,PetscReal);
200d8345c25SBarry Smith 
2017087cfbeSBarry Smith extern PetscErrorCode   TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*);
2027087cfbeSBarry Smith extern PetscErrorCode   TSPseudoDefaultVerifyTimeStep(TS,Vec,void*,PetscReal*,PetscBool *);
2037087cfbeSBarry Smith extern PetscErrorCode   TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *);
2047087cfbeSBarry Smith extern PetscErrorCode   TSPseudoSetTimeStepIncrement(TS,PetscReal);
2057087cfbeSBarry Smith extern PetscErrorCode   TSPseudoIncrementDtFromInitialDt(TS);
20621c89e3eSBarry Smith 
2077087cfbeSBarry Smith extern PetscErrorCode   TSPythonSetType(TS,const char[]);
2081d6018f0SLisandro Dalcin 
2097087cfbeSBarry Smith extern PetscErrorCode   TSComputeRHSFunction(TS,PetscReal,Vec,Vec);
2107087cfbeSBarry Smith extern PetscErrorCode   TSComputeRHSJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*);
211214bc6a2SJed Brown extern PetscErrorCode   TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool);
212214bc6a2SJed Brown extern PetscErrorCode   TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,PetscBool);
213818ad0c1SBarry Smith 
214d6ebe24aSShri Abhyankar extern PetscErrorCode   TSVISetVariableBounds(TS,Vec,Vec);
215d6ebe24aSShri Abhyankar 
216bdad233fSMatthew Knepley /* Dynamic creation and loading functions */
217b0a32e0cSBarry Smith extern PetscFList TSList;
218ace3abfcSBarry Smith extern PetscBool  TSRegisterAllCalled;
2197087cfbeSBarry Smith extern PetscErrorCode   TSGetType(TS,const TSType*);
2207087cfbeSBarry Smith extern PetscErrorCode   TSSetType(TS,const TSType);
2217087cfbeSBarry Smith extern PetscErrorCode   TSRegister(const char[], const char[], const char[], PetscErrorCode (*)(TS));
2227087cfbeSBarry Smith extern PetscErrorCode   TSRegisterAll(const char[]);
2237087cfbeSBarry Smith extern PetscErrorCode   TSRegisterDestroy(void);
22430de9b25SBarry Smith 
22530de9b25SBarry Smith /*MC
22630de9b25SBarry Smith   TSRegisterDynamic - Adds a creation method to the TS package.
22730de9b25SBarry Smith 
22830de9b25SBarry Smith   Synopsis:
2291890ba74SBarry Smith   PetscErrorCode TSRegisterDynamic(const char *name, const char *path, const char *func_name, PetscErrorCode (*create_func)(TS))
23030de9b25SBarry Smith 
23130de9b25SBarry Smith   Not Collective
23230de9b25SBarry Smith 
23330de9b25SBarry Smith   Input Parameters:
23430de9b25SBarry Smith + name        - The name of a new user-defined creation routine
23530de9b25SBarry Smith . path        - The path (either absolute or relative) of the library containing this routine
23630de9b25SBarry Smith . func_name   - The name of the creation routine
23730de9b25SBarry Smith - create_func - The creation routine itself
23830de9b25SBarry Smith 
23930de9b25SBarry Smith   Notes:
24030de9b25SBarry Smith   TSRegisterDynamic() may be called multiple times to add several user-defined tses.
24130de9b25SBarry Smith 
24230de9b25SBarry Smith   If dynamic libraries are used, then the fourth input argument (create_func) is ignored.
24330de9b25SBarry Smith 
24430de9b25SBarry Smith   Sample usage:
24530de9b25SBarry Smith .vb
24630de9b25SBarry Smith   TSRegisterDynamic("my_ts", "/home/username/my_lib/lib/libO/solaris/libmy.a", "MyTSCreate", MyTSCreate);
24730de9b25SBarry Smith .ve
24830de9b25SBarry Smith 
24930de9b25SBarry Smith   Then, your ts type can be chosen with the procedural interface via
25030de9b25SBarry Smith .vb
25165f2c13aSKai Germaschewski     TS ts;
25265f2c13aSKai Germaschewski     TSCreate(MPI_Comm, &ts);
25365f2c13aSKai Germaschewski     TSSetType(ts, "my_ts")
25430de9b25SBarry Smith .ve
25530de9b25SBarry Smith   or at runtime via the option
25630de9b25SBarry Smith .vb
25730de9b25SBarry Smith     -ts_type my_ts
25830de9b25SBarry Smith .ve
25930de9b25SBarry Smith 
260ab901514SBarry Smith   Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
26130de9b25SBarry Smith         If your function is not being put into a shared library then use TSRegister() instead
26230de9b25SBarry Smith 
26330de9b25SBarry Smith   Level: advanced
26430de9b25SBarry Smith 
26530de9b25SBarry Smith .keywords: TS, register
26630de9b25SBarry Smith .seealso: TSRegisterAll(), TSRegisterDestroy()
26730de9b25SBarry Smith M*/
268aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
269f1af5d2fSBarry Smith #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,0)
2706df38c32SLois Curfman McInnes #else
271f1af5d2fSBarry Smith #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,d)
2726df38c32SLois Curfman McInnes #endif
2736df38c32SLois Curfman McInnes 
2747087cfbeSBarry Smith extern PetscErrorCode   TSGetSNES(TS,SNES*);
2757087cfbeSBarry Smith extern PetscErrorCode   TSGetKSP(TS,KSP*);
276818ad0c1SBarry Smith 
2777087cfbeSBarry Smith extern PetscErrorCode   TSView(TS,PetscViewer);
27821c89e3eSBarry Smith 
2797087cfbeSBarry Smith extern PetscErrorCode   TSSetApplicationContext(TS,void *);
280e71120c6SJed Brown extern PetscErrorCode   TSGetApplicationContext(TS,void *);
28121c89e3eSBarry Smith 
2827087cfbeSBarry Smith extern PetscErrorCode   TSMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG *);
2837087cfbeSBarry Smith extern PetscErrorCode   TSMonitorLG(TS,PetscInt,PetscReal,Vec,void *);
284fcfd50ebSBarry Smith extern PetscErrorCode   TSMonitorLGDestroy(PetscDrawLG*);
2853914022bSBarry Smith 
28676bdecfbSBarry Smith /*J
287815f1ad5SJed Brown    TSSSPType - string with the name of TSSSP scheme.
288815f1ad5SJed Brown 
289815f1ad5SJed Brown    Level: beginner
290815f1ad5SJed Brown 
291815f1ad5SJed Brown .seealso: TSSSPSetType(), TS
29276bdecfbSBarry Smith J*/
293815f1ad5SJed Brown #define TSSSPType char*
294815f1ad5SJed Brown #define TSSSPRKS2  "rks2"
295815f1ad5SJed Brown #define TSSSPRKS3  "rks3"
296815f1ad5SJed Brown #define TSSSPRK104 "rk104"
297815f1ad5SJed Brown 
298815f1ad5SJed Brown extern PetscErrorCode TSSSPSetType(TS,const TSSSPType);
299815f1ad5SJed Brown extern PetscErrorCode TSSSPGetType(TS,const TSSSPType*);
300815f1ad5SJed Brown extern PetscErrorCode TSSSPSetNumStages(TS,PetscInt);
301815f1ad5SJed Brown extern PetscErrorCode TSSSPGetNumStages(TS,PetscInt*);
302815f1ad5SJed Brown 
303ed657a08SJed Brown /*S
30484df9cb4SJed Brown    TSAdapt - Abstract object that manages time-step adaptivity
30584df9cb4SJed Brown 
30684df9cb4SJed Brown    Level: beginner
30784df9cb4SJed Brown 
30884df9cb4SJed Brown .seealso: TS, TSAdaptCreate(), TSAdaptType
30984df9cb4SJed Brown S*/
31084df9cb4SJed Brown typedef struct _p_TSAdapt *TSAdapt;
31184df9cb4SJed Brown 
31284df9cb4SJed Brown /*E
31384df9cb4SJed Brown     TSAdaptType - String with the name of TSAdapt scheme or the creation function
31484df9cb4SJed Brown        with an optional dynamic library name, for example
31584df9cb4SJed Brown        http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate()
31684df9cb4SJed Brown 
31784df9cb4SJed Brown    Level: beginner
31884df9cb4SJed Brown 
31984df9cb4SJed Brown .seealso: TSAdaptSetType(), TS
32084df9cb4SJed Brown E*/
32184df9cb4SJed Brown #define TSAdaptType  char*
32284df9cb4SJed Brown #define TSADAPTBASIC "basic"
323b0f836d7SJed Brown #define TSADAPTNONE  "none"
3248d59e960SJed Brown #define TSADAPTCFL   "cfl"
32584df9cb4SJed Brown 
32684df9cb4SJed Brown /*MC
32784df9cb4SJed Brown    TSAdaptRegisterDynamic - adds a TSAdapt implementation
32884df9cb4SJed Brown 
32984df9cb4SJed Brown    Synopsis:
33084df9cb4SJed Brown    PetscErrorCode TSAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
33184df9cb4SJed Brown 
33284df9cb4SJed Brown    Not Collective
33384df9cb4SJed Brown 
33484df9cb4SJed Brown    Input Parameters:
33584df9cb4SJed Brown +  name_scheme - name of user-defined adaptivity scheme
33684df9cb4SJed Brown .  path - path (either absolute or relative) the library containing this scheme
33784df9cb4SJed Brown .  name_create - name of routine to create method context
33884df9cb4SJed Brown -  routine_create - routine to create method context
33984df9cb4SJed Brown 
34084df9cb4SJed Brown    Notes:
34184df9cb4SJed Brown    TSAdaptRegisterDynamic() may be called multiple times to add several user-defined families.
34284df9cb4SJed Brown 
34384df9cb4SJed Brown    If dynamic libraries are used, then the fourth input argument (routine_create)
34484df9cb4SJed Brown    is ignored.
34584df9cb4SJed Brown 
34684df9cb4SJed Brown    Sample usage:
34784df9cb4SJed Brown .vb
34884df9cb4SJed Brown    TSAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
34984df9cb4SJed Brown                             "MySchemeCreate",MySchemeCreate);
35084df9cb4SJed Brown .ve
35184df9cb4SJed Brown 
35284df9cb4SJed Brown    Then, your scheme can be chosen with the procedural interface via
35384df9cb4SJed Brown $     TSAdaptSetType(ts,"my_scheme")
35484df9cb4SJed Brown    or at runtime via the option
35584df9cb4SJed Brown $     -ts_adapt_type my_scheme
35684df9cb4SJed Brown 
35784df9cb4SJed Brown    Level: advanced
35884df9cb4SJed Brown 
35984df9cb4SJed Brown    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
36084df9cb4SJed Brown           and others of the form ${any_environmental_variable} occuring in pathname will be
36184df9cb4SJed Brown           replaced with appropriate values.
36284df9cb4SJed Brown 
36384df9cb4SJed Brown .keywords: TSAdapt, register
36484df9cb4SJed Brown 
36584df9cb4SJed Brown .seealso: TSAdaptRegisterAll()
36684df9cb4SJed Brown M*/
36784df9cb4SJed Brown #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
36884df9cb4SJed Brown #  define TSAdaptRegisterDynamic(a,b,c,d)  TSAdaptRegister(a,b,c,0)
36984df9cb4SJed Brown #else
37084df9cb4SJed Brown #  define TSAdaptRegisterDynamic(a,b,c,d)  TSAdaptRegister(a,b,c,d)
37184df9cb4SJed Brown #endif
37284df9cb4SJed Brown 
37384df9cb4SJed Brown extern PetscErrorCode TSGetAdapt(TS,TSAdapt*);
37484df9cb4SJed Brown extern PetscErrorCode TSAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSAdapt));
37584df9cb4SJed Brown extern PetscErrorCode TSAdaptRegisterAll(const char[]);
37684df9cb4SJed Brown extern PetscErrorCode TSAdaptRegisterDestroy(void);
37784df9cb4SJed Brown extern PetscErrorCode TSAdaptInitializePackage(const char[]);
37884df9cb4SJed Brown extern PetscErrorCode TSAdaptFinalizePackage(void);
37984df9cb4SJed Brown extern PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*);
38084df9cb4SJed Brown extern PetscErrorCode TSAdaptSetType(TSAdapt,const TSAdaptType);
38184df9cb4SJed Brown extern PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]);
38284df9cb4SJed Brown extern PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
38384df9cb4SJed Brown extern PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool);
38484df9cb4SJed Brown extern PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**);
38584df9cb4SJed Brown extern PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*);
38684df9cb4SJed Brown extern PetscErrorCode TSAdaptView(TSAdapt,PetscViewer);
38784df9cb4SJed Brown extern PetscErrorCode TSAdaptSetFromOptions(TSAdapt);
38884df9cb4SJed Brown extern PetscErrorCode TSAdaptDestroy(TSAdapt*);
3891c3436cfSJed Brown extern PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool);
3901c3436cfSJed Brown extern PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal);
39184df9cb4SJed Brown 
39284df9cb4SJed Brown /*S
393ed657a08SJed Brown    TSGLAdapt - Abstract object that manages time-step adaptivity
394ed657a08SJed Brown 
395ed657a08SJed Brown    Level: beginner
396ed657a08SJed Brown 
39784df9cb4SJed Brown    Developer Notes:
39884df9cb4SJed Brown    This functionality should be replaced by the TSAdapt.
39984df9cb4SJed Brown 
400ed657a08SJed Brown .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType
401ed657a08SJed Brown S*/
402ed657a08SJed Brown typedef struct _p_TSGLAdapt *TSGLAdapt;
403ed657a08SJed Brown 
40476bdecfbSBarry Smith /*J
405ed657a08SJed Brown     TSGLAdaptType - String with the name of TSGLAdapt scheme or the creation function
406ed657a08SJed Brown        with an optional dynamic library name, for example
407ed657a08SJed Brown        http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate()
408ed657a08SJed Brown 
409ed657a08SJed Brown    Level: beginner
410ed657a08SJed Brown 
411ed657a08SJed Brown .seealso: TSGLAdaptSetType(), TS
41276bdecfbSBarry Smith J*/
413ed657a08SJed Brown #define TSGLAdaptType  char*
414ed657a08SJed Brown #define TSGLADAPT_NONE "none"
415ed657a08SJed Brown #define TSGLADAPT_SIZE "size"
416ed657a08SJed Brown #define TSGLADAPT_BOTH "both"
417ed657a08SJed Brown 
418ed657a08SJed Brown /*MC
419ed657a08SJed Brown    TSGLAdaptRegisterDynamic - adds a TSGLAdapt implementation
420ed657a08SJed Brown 
421ed657a08SJed Brown    Synopsis:
4221890ba74SBarry Smith    PetscErrorCode TSGLAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
423ed657a08SJed Brown 
424ed657a08SJed Brown    Not Collective
425ed657a08SJed Brown 
426ed657a08SJed Brown    Input Parameters:
427ed657a08SJed Brown +  name_scheme - name of user-defined adaptivity scheme
428ed657a08SJed Brown .  path - path (either absolute or relative) the library containing this scheme
429ed657a08SJed Brown .  name_create - name of routine to create method context
430ed657a08SJed Brown -  routine_create - routine to create method context
431ed657a08SJed Brown 
432ed657a08SJed Brown    Notes:
433ed657a08SJed Brown    TSGLAdaptRegisterDynamic() may be called multiple times to add several user-defined families.
434ed657a08SJed Brown 
435ed657a08SJed Brown    If dynamic libraries are used, then the fourth input argument (routine_create)
436ed657a08SJed Brown    is ignored.
437ed657a08SJed Brown 
438ed657a08SJed Brown    Sample usage:
439ed657a08SJed Brown .vb
440ed657a08SJed Brown    TSGLAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
441ed657a08SJed Brown                             "MySchemeCreate",MySchemeCreate);
442ed657a08SJed Brown .ve
443ed657a08SJed Brown 
444ed657a08SJed Brown    Then, your scheme can be chosen with the procedural interface via
445ed657a08SJed Brown $     TSGLAdaptSetType(ts,"my_scheme")
446ed657a08SJed Brown    or at runtime via the option
447ed657a08SJed Brown $     -ts_adapt_type my_scheme
448ed657a08SJed Brown 
449ed657a08SJed Brown    Level: advanced
450ed657a08SJed Brown 
451ed657a08SJed Brown    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
452ed657a08SJed Brown           and others of the form ${any_environmental_variable} occuring in pathname will be
453ed657a08SJed Brown           replaced with appropriate values.
454ed657a08SJed Brown 
455ed657a08SJed Brown .keywords: TSGLAdapt, register
456ed657a08SJed Brown 
457ed657a08SJed Brown .seealso: TSGLAdaptRegisterAll()
458ed657a08SJed Brown M*/
459ed657a08SJed Brown #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
460ed657a08SJed Brown #  define TSGLAdaptRegisterDynamic(a,b,c,d)  TSGLAdaptRegister(a,b,c,0)
461ed657a08SJed Brown #else
462ed657a08SJed Brown #  define TSGLAdaptRegisterDynamic(a,b,c,d)  TSGLAdaptRegister(a,b,c,d)
463ed657a08SJed Brown #endif
464ed657a08SJed Brown 
4657087cfbeSBarry Smith extern PetscErrorCode  TSGLAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSGLAdapt));
4667087cfbeSBarry Smith extern PetscErrorCode  TSGLAdaptRegisterAll(const char[]);
4677087cfbeSBarry Smith extern PetscErrorCode  TSGLAdaptRegisterDestroy(void);
4687087cfbeSBarry Smith extern PetscErrorCode  TSGLAdaptInitializePackage(const char[]);
4697087cfbeSBarry Smith extern PetscErrorCode  TSGLAdaptFinalizePackage(void);
4707087cfbeSBarry Smith extern PetscErrorCode  TSGLAdaptCreate(MPI_Comm,TSGLAdapt*);
4717087cfbeSBarry Smith extern PetscErrorCode  TSGLAdaptSetType(TSGLAdapt,const TSGLAdaptType);
4727087cfbeSBarry Smith extern PetscErrorCode  TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]);
4737087cfbeSBarry Smith extern PetscErrorCode  TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
4747087cfbeSBarry Smith extern PetscErrorCode  TSGLAdaptView(TSGLAdapt,PetscViewer);
4757087cfbeSBarry Smith extern PetscErrorCode  TSGLAdaptSetFromOptions(TSGLAdapt);
476fcfd50ebSBarry Smith extern PetscErrorCode  TSGLAdaptDestroy(TSGLAdapt*);
477ed657a08SJed Brown 
47876bdecfbSBarry Smith /*J
479ed657a08SJed Brown     TSGLAcceptType - String with the name of TSGLAccept scheme or the function
480ed657a08SJed Brown        with an optional dynamic library name, for example
481ed657a08SJed Brown        http://www.mcs.anl.gov/petsc/lib.a:mytsglaccept()
482ed657a08SJed Brown 
483ed657a08SJed Brown    Level: beginner
484ed657a08SJed Brown 
485ed657a08SJed Brown .seealso: TSGLSetAcceptType(), TS
48676bdecfbSBarry Smith J*/
487ed657a08SJed Brown #define TSGLAcceptType  char*
488ed657a08SJed Brown #define TSGLACCEPT_ALWAYS "always"
489ed657a08SJed Brown 
490ace3abfcSBarry Smith typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *);
4917087cfbeSBarry Smith extern PetscErrorCode  TSGLAcceptRegister(const char[],const char[],const char[],TSGLAcceptFunction);
492ed657a08SJed Brown 
493ed657a08SJed Brown /*MC
494ed657a08SJed Brown    TSGLAcceptRegisterDynamic - adds a TSGL acceptance scheme
495ed657a08SJed Brown 
496ed657a08SJed Brown    Synopsis:
4971890ba74SBarry Smith    PetscErrorCode TSGLAcceptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
498ed657a08SJed Brown 
499ed657a08SJed Brown    Not Collective
500ed657a08SJed Brown 
501ed657a08SJed Brown    Input Parameters:
502ed657a08SJed Brown +  name_scheme - name of user-defined acceptance scheme
503ed657a08SJed Brown .  path - path (either absolute or relative) the library containing this scheme
504ed657a08SJed Brown .  name_create - name of routine to create method context
505ed657a08SJed Brown -  routine_create - routine to create method context
506ed657a08SJed Brown 
507ed657a08SJed Brown    Notes:
508ed657a08SJed Brown    TSGLAcceptRegisterDynamic() may be called multiple times to add several user-defined families.
509ed657a08SJed Brown 
510ed657a08SJed Brown    If dynamic libraries are used, then the fourth input argument (routine_create)
511ed657a08SJed Brown    is ignored.
512ed657a08SJed Brown 
513ed657a08SJed Brown    Sample usage:
514ed657a08SJed Brown .vb
515ed657a08SJed Brown    TSGLAcceptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
516ed657a08SJed Brown                              "MySchemeCreate",MySchemeCreate);
517ed657a08SJed Brown .ve
518ed657a08SJed Brown 
519ed657a08SJed Brown    Then, your scheme can be chosen with the procedural interface via
520ed657a08SJed Brown $     TSGLSetAcceptType(ts,"my_scheme")
521ed657a08SJed Brown    or at runtime via the option
522ed657a08SJed Brown $     -ts_gl_accept_type my_scheme
523ed657a08SJed Brown 
524ed657a08SJed Brown    Level: advanced
525ed657a08SJed Brown 
526ed657a08SJed Brown    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
527ed657a08SJed Brown           and others of the form ${any_environmental_variable} occuring in pathname will be
528ed657a08SJed Brown           replaced with appropriate values.
529ed657a08SJed Brown 
530ed657a08SJed Brown .keywords: TSGL, TSGLAcceptType, register
531ed657a08SJed Brown 
532ed657a08SJed Brown .seealso: TSGLRegisterAll()
533ed657a08SJed Brown M*/
534ed657a08SJed Brown #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
535ed657a08SJed Brown #  define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,0)
536ed657a08SJed Brown #else
537ed657a08SJed Brown #  define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,d)
538ed657a08SJed Brown #endif
539ed657a08SJed Brown 
54076bdecfbSBarry Smith /*J
54118b56cb9SJed Brown   TSGLType - family of time integration method within the General Linear class
54218b56cb9SJed Brown 
54318b56cb9SJed Brown   Level: beginner
54418b56cb9SJed Brown 
54518b56cb9SJed Brown .seealso: TSGLSetType(), TSGLRegister()
54676bdecfbSBarry Smith J*/
54718b56cb9SJed Brown #define TSGLType char*
548ed657a08SJed Brown #define TSGL_IRKS   "irks"
54918b56cb9SJed Brown 
550ed657a08SJed Brown /*MC
551ed657a08SJed Brown    TSGLRegisterDynamic - adds a TSGL implementation
55218b56cb9SJed Brown 
553ed657a08SJed Brown    Synopsis:
5541890ba74SBarry Smith    PetscErrorCode TSGLRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
555ed657a08SJed Brown 
556ed657a08SJed Brown    Not Collective
557ed657a08SJed Brown 
558ed657a08SJed Brown    Input Parameters:
559ed657a08SJed Brown +  name_scheme - name of user-defined general linear scheme
560ed657a08SJed Brown .  path - path (either absolute or relative) the library containing this scheme
561ed657a08SJed Brown .  name_create - name of routine to create method context
562ed657a08SJed Brown -  routine_create - routine to create method context
563ed657a08SJed Brown 
564ed657a08SJed Brown    Notes:
565ed657a08SJed Brown    TSGLRegisterDynamic() may be called multiple times to add several user-defined families.
566ed657a08SJed Brown 
567ed657a08SJed Brown    If dynamic libraries are used, then the fourth input argument (routine_create)
568ed657a08SJed Brown    is ignored.
569ed657a08SJed Brown 
570ed657a08SJed Brown    Sample usage:
571ed657a08SJed Brown .vb
572ed657a08SJed Brown    TSGLRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
573ed657a08SJed Brown                        "MySchemeCreate",MySchemeCreate);
574ed657a08SJed Brown .ve
575ed657a08SJed Brown 
576ed657a08SJed Brown    Then, your scheme can be chosen with the procedural interface via
577ed657a08SJed Brown $     TSGLSetType(ts,"my_scheme")
578ed657a08SJed Brown    or at runtime via the option
579ed657a08SJed Brown $     -ts_gl_type my_scheme
580ed657a08SJed Brown 
581ed657a08SJed Brown    Level: advanced
582ed657a08SJed Brown 
583ed657a08SJed Brown    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
584ed657a08SJed Brown           and others of the form ${any_environmental_variable} occuring in pathname will be
585ed657a08SJed Brown           replaced with appropriate values.
586ed657a08SJed Brown 
587ed657a08SJed Brown .keywords: TSGL, register
588ed657a08SJed Brown 
589ed657a08SJed Brown .seealso: TSGLRegisterAll()
590ed657a08SJed Brown M*/
591ed657a08SJed Brown #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
592ed657a08SJed Brown #  define TSGLRegisterDynamic(a,b,c,d)       TSGLRegister(a,b,c,0)
593ed657a08SJed Brown #else
594ed657a08SJed Brown #  define TSGLRegisterDynamic(a,b,c,d)       TSGLRegister(a,b,c,d)
595ed657a08SJed Brown #endif
596ed657a08SJed Brown 
5977087cfbeSBarry Smith extern PetscErrorCode  TSGLRegister(const char[],const char[],const char[],PetscErrorCode(*)(TS));
5987087cfbeSBarry Smith extern PetscErrorCode  TSGLRegisterAll(const char[]);
5997087cfbeSBarry Smith extern PetscErrorCode  TSGLRegisterDestroy(void);
6007087cfbeSBarry Smith extern PetscErrorCode  TSGLInitializePackage(const char[]);
6017087cfbeSBarry Smith extern PetscErrorCode  TSGLFinalizePackage(void);
6027087cfbeSBarry Smith extern PetscErrorCode  TSGLSetType(TS,const TSGLType);
6037087cfbeSBarry Smith extern PetscErrorCode  TSGLGetAdapt(TS,TSGLAdapt*);
6047087cfbeSBarry Smith extern PetscErrorCode  TSGLSetAcceptType(TS,const TSGLAcceptType);
60518b56cb9SJed Brown 
6068a381b04SJed Brown #define TSARKIMEXType char*
6078a381b04SJed Brown #define TSARKIMEX2D     "2d"
608a3a57f36SJed Brown #define TSARKIMEX2E     "2e"
6096cf0794eSJed Brown #define TSARKIMEXPRSSP2 "prssp2"
6108a381b04SJed Brown #define TSARKIMEX3      "3"
6116cf0794eSJed Brown #define TSARKIMEXBPR3   "bpr3"
6126cf0794eSJed Brown #define TSARKIMEXARS443 "ars443"
6138a381b04SJed Brown #define TSARKIMEX4      "4"
6148a381b04SJed Brown #define TSARKIMEX5      "5"
6158a381b04SJed Brown extern PetscErrorCode TSARKIMEXGetType(TS ts,const TSARKIMEXType*);
6168a381b04SJed Brown extern PetscErrorCode TSARKIMEXSetType(TS ts,const TSARKIMEXType);
6174cc180ffSJed Brown extern PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool);
618108c343cSJed Brown extern PetscErrorCode TSARKIMEXRegister(const TSARKIMEXType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[],const PetscReal[]);
6198a381b04SJed Brown extern PetscErrorCode TSARKIMEXFinalizePackage(void);
6208a381b04SJed Brown extern PetscErrorCode TSARKIMEXInitializePackage(const char path[]);
6218a381b04SJed Brown extern PetscErrorCode TSARKIMEXRegisterDestroy(void);
6228a381b04SJed Brown extern PetscErrorCode TSARKIMEXRegisterAll(void);
6238a381b04SJed Brown 
624e27a552bSJed Brown #define TSRosWType char*
62561692a83SJed Brown #define TSROSW2M          "2m"
62661692a83SJed Brown #define TSROSW2P          "2p"
627fe7e6d57SJed Brown #define TSROSWRA3PW       "ra3pw"
628fe7e6d57SJed Brown #define TSROSWRA34PW2     "ra34pw2"
629ef3c5b88SJed Brown #define TSROSWRODAS3      "rodas3"
630ef3c5b88SJed Brown #define TSROSWSANDU3      "sandu3"
631961f28d0SJed Brown #define TSROSWASSP3P3S1C  "assp3p3s1c"
632961f28d0SJed Brown #define TSROSWLASSP3P4S2C "lassp3p4s2c"
633961f28d0SJed Brown #define TSROSWLLSSP3P3S2C "llssp3p3s2c"
634753f8adbSEmil Constantinescu #define TSROSWARK3        "ark3"
635b1c69cc3SEmil Constantinescu 
636e27a552bSJed Brown extern PetscErrorCode TSRosWGetType(TS ts,const TSRosWType*);
637e27a552bSJed Brown extern PetscErrorCode TSRosWSetType(TS ts,const TSRosWType);
63861692a83SJed Brown extern PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool);
639fe7e6d57SJed Brown extern PetscErrorCode TSRosWRegister(const TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[]);
640e27a552bSJed Brown extern PetscErrorCode TSRosWFinalizePackage(void);
641e27a552bSJed Brown extern PetscErrorCode TSRosWInitializePackage(const char path[]);
642e27a552bSJed Brown extern PetscErrorCode TSRosWRegisterDestroy(void);
643e27a552bSJed Brown extern PetscErrorCode TSRosWRegisterAll(void);
644e27a552bSJed Brown 
64583e2fdc7SBarry Smith /*
6460f3b3ca1SHong Zhang        PETSc interface to Sundials
64783e2fdc7SBarry Smith */
6486dc17ff5SMatthew Knepley #ifdef PETSC_HAVE_SUNDIALS
6496fadb2cdSHong Zhang typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
6506fadb2cdSHong Zhang extern const char *TSSundialsLmmTypes[];
6516fadb2cdSHong Zhang typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
652a04cf4d8SBarry Smith extern const char *TSSundialsGramSchmidtTypes[];
6537087cfbeSBarry Smith extern PetscErrorCode   TSSundialsSetType(TS,TSSundialsLmmType);
6547087cfbeSBarry Smith extern PetscErrorCode   TSSundialsGetPC(TS,PC*);
6557087cfbeSBarry Smith extern PetscErrorCode   TSSundialsSetTolerance(TS,PetscReal,PetscReal);
6567087cfbeSBarry Smith extern PetscErrorCode   TSSundialsSetMinTimeStep(TS,PetscReal);
6577087cfbeSBarry Smith extern PetscErrorCode   TSSundialsSetMaxTimeStep(TS,PetscReal);
6587087cfbeSBarry Smith extern PetscErrorCode   TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
6597087cfbeSBarry Smith extern PetscErrorCode   TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
6607087cfbeSBarry Smith extern PetscErrorCode   TSSundialsSetGMRESRestart(TS,PetscInt);
6617087cfbeSBarry Smith extern PetscErrorCode   TSSundialsSetLinearTolerance(TS,PetscReal);
6627087cfbeSBarry Smith extern PetscErrorCode   TSSundialsMonitorInternalSteps(TS,PetscBool );
6637087cfbeSBarry Smith extern PetscErrorCode   TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
66454c84890SJed Brown extern PetscErrorCode   TSSundialsSetMaxl(TS,PetscInt);
6656dc17ff5SMatthew Knepley #endif
66683e2fdc7SBarry Smith 
6677087cfbeSBarry Smith extern PetscErrorCode   TSRKSetTolerance(TS,PetscReal);
668262ff7bbSBarry Smith 
6697087cfbeSBarry Smith extern PetscErrorCode  TSThetaSetTheta(TS,PetscReal);
6707087cfbeSBarry Smith extern PetscErrorCode  TSThetaGetTheta(TS,PetscReal*);
67126f2ff8fSLisandro Dalcin extern PetscErrorCode  TSThetaGetEndpoint(TS,PetscBool*);
672eb284becSJed Brown extern PetscErrorCode  TSThetaSetEndpoint(TS,PetscBool);
6730de4c49aSJed Brown 
6747087cfbeSBarry Smith extern PetscErrorCode  TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*);
6757087cfbeSBarry Smith extern PetscErrorCode  TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*);
6767087cfbeSBarry Smith extern PetscErrorCode  TSAlphaSetRadius(TS,PetscReal);
6777087cfbeSBarry Smith extern PetscErrorCode  TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal);
6787087cfbeSBarry Smith extern PetscErrorCode  TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*);
67988df8a41SLisandro Dalcin 
6807087cfbeSBarry Smith extern PetscErrorCode  TSSetDM(TS,DM);
6817087cfbeSBarry Smith extern PetscErrorCode  TSGetDM(TS,DM*);
6826c699258SBarry Smith 
6837087cfbeSBarry Smith extern PetscErrorCode  SNESTSFormFunction(SNES,Vec,Vec,void*);
6847087cfbeSBarry Smith extern PetscErrorCode  SNESTSFormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
6850f5c6efeSJed Brown 
686e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
687818ad0c1SBarry Smith #endif
688