Actual source code: petscts.h
1: /*
2: User interface for the timestepping package. This is package
3: is for use in solving time-dependent PDEs.
4: */
7: #include petscsnes.h
10: /*S
11: TS - Abstract PETSc object that manages all time-steppers (ODE integrators)
13: Level: beginner
15: Concepts: ODE solvers
17: .seealso: TSCreate(), TSSetType(), TSType, SNES, KSP, PC
18: S*/
19: typedef struct _p_TS* TS;
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()
26: Level: beginner
28: .seealso: TSSetType(), TS
29: E*/
30: #define TSType char*
31: #define TS_EULER "euler"
32: #define TS_BEULER "beuler"
33: #define TS_PSEUDO "pseudo"
34: #define TS_CRANK_NICHOLSON "crank-nicholson"
35: #define TS_SUNDIALS "sundials"
36: #define TS_RUNGE_KUTTA "runge-kutta"
37: #define TS_PYTHON "python"
39: /*E
40: TSProblemType - Determines the type of problem this TS object is to be used to solve
42: Level: beginner
44: .seealso: TSCreate()
45: E*/
46: typedef enum {TS_LINEAR,TS_NONLINEAR} TSProblemType;
48: /* Logging support */
51: EXTERN PetscErrorCode TSInitializePackage(const char[]);
53: EXTERN PetscErrorCode TSCreate(MPI_Comm,TS*);
54: EXTERN PetscErrorCode TSDestroy(TS);
56: EXTERN PetscErrorCode TSSetProblemType(TS,TSProblemType);
57: EXTERN PetscErrorCode TSGetProblemType(TS,TSProblemType*);
58: EXTERN PetscErrorCode TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void*));
59: EXTERN PetscErrorCode TSMonitorCancel(TS);
61: EXTERN PetscErrorCode TSSetOptionsPrefix(TS,const char[]);
62: EXTERN PetscErrorCode TSAppendOptionsPrefix(TS,const char[]);
63: EXTERN PetscErrorCode TSGetOptionsPrefix(TS,const char *[]);
64: EXTERN PetscErrorCode TSSetFromOptions(TS);
65: EXTERN PetscErrorCode TSSetUp(TS);
67: EXTERN PetscErrorCode TSSetSolution(TS,Vec);
68: EXTERN PetscErrorCode TSGetSolution(TS,Vec*);
70: EXTERN PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal);
71: EXTERN PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*);
73: EXTERN PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,void*);
74: EXTERN PetscErrorCode TSMonitorSolution(TS,PetscInt,PetscReal,Vec,void*);
75: EXTERN PetscErrorCode TSStep(TS,PetscInt *,PetscReal*);
76: EXTERN PetscErrorCode TSSolve(TS,Vec);
79: EXTERN PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal);
80: EXTERN PetscErrorCode TSGetTimeStep(TS,PetscReal*);
81: EXTERN PetscErrorCode TSGetTime(TS,PetscReal*);
82: EXTERN PetscErrorCode TSSetTime(TS,PetscReal);
83: EXTERN PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*);
84: EXTERN PetscErrorCode TSSetTimeStep(TS,PetscReal);
86: EXTERN PetscErrorCode TSSetRHSFunction(TS,PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*),void*);
87: EXTERN PetscErrorCode TSSetMatrices(TS,Mat,PetscErrorCode (*)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),Mat,PetscErrorCode (*)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),MatStructure,void*);
88: EXTERN PetscErrorCode TSGetMatrices(TS,Mat*,Mat*,void**);
89: EXTERN PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,PetscErrorCode (*)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*),void*);
90: EXTERN PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,void**);
92: EXTERN PetscErrorCode TSDefaultComputeJacobianColor(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
93: EXTERN PetscErrorCode TSDefaultComputeJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
95: EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
96: EXTERN PetscErrorCode TSSetUpdate(TS, PetscErrorCode (*)(TS, PetscReal, PetscReal *));
97: EXTERN PetscErrorCode TSSetPostUpdate(TS, PetscErrorCode (*)(TS, PetscReal, PetscReal *));
98: EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
99: EXTERN PetscErrorCode TSDefaultPreStep(TS);
100: EXTERN PetscErrorCode TSDefaultUpdate(TS, PetscReal, PetscReal *);
101: EXTERN PetscErrorCode TSDefaultPostStep(TS);
103: EXTERN PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*);
104: EXTERN PetscErrorCode TSPseudoDefaultTimeStep(TS,PetscReal*,void*);
105: EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *);
107: EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscTruth*),void*);
108: EXTERN PetscErrorCode TSPseudoDefaultVerifyTimeStep(TS,Vec,void*,PetscReal*,PetscTruth*);
109: EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscTruth*);
110: EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal);
111: EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);
113: EXTERN PetscErrorCode TSPythonSetType(TS,const char[]);
115: EXTERN PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec);
116: EXTERN PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*);
118: /* Dynamic creation and loading functions */
121: EXTERN PetscErrorCode TSGetType(TS,const TSType*);
122: EXTERN PetscErrorCode TSSetType(TS,const TSType);
123: EXTERN PetscErrorCode TSRegister(const char[], const char[], const char[], PetscErrorCode (*)(TS));
124: EXTERN PetscErrorCode TSRegisterAll(const char[]);
125: EXTERN PetscErrorCode TSRegisterDestroy(void);
127: /*MC
128: TSRegisterDynamic - Adds a creation method to the TS package.
130: Synopsis:
131: PetscErrorCode TSRegisterDynamic(char *name, char *path, char *func_name, PetscErrorCode (*create_func)(TS))
133: Not Collective
135: Input Parameters:
136: + name - The name of a new user-defined creation routine
137: . path - The path (either absolute or relative) of the library containing this routine
138: . func_name - The name of the creation routine
139: - create_func - The creation routine itself
141: Notes:
142: TSRegisterDynamic() may be called multiple times to add several user-defined tses.
144: If dynamic libraries are used, then the fourth input argument (create_func) is ignored.
146: Sample usage:
147: .vb
148: TSRegisterDynamic("my_ts", "/home/username/my_lib/lib/libO/solaris/libmy.a", "MyTSCreate", MyTSCreate);
149: .ve
151: Then, your ts type can be chosen with the procedural interface via
152: .vb
153: TSCreate(MPI_Comm, TS *);
154: TSSetType(vec, "my_ts")
155: .ve
156: or at runtime via the option
157: .vb
158: -ts_type my_ts
159: .ve
161: Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
162: If your function is not being put into a shared library then use TSRegister() instead
164: Level: advanced
166: .keywords: TS, register
167: .seealso: TSRegisterAll(), TSRegisterDestroy()
168: M*/
169: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
170: #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,0)
171: #else
172: #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,d)
173: #endif
175: EXTERN PetscErrorCode TSGetSNES(TS,SNES*);
176: EXTERN PetscErrorCode TSGetKSP(TS,KSP*);
178: EXTERN PetscErrorCode TSView(TS,PetscViewer);
179: EXTERN PetscErrorCode TSViewFromOptions(TS,const char[]);
181: EXTERN PetscErrorCode TSSetApplicationContext(TS,void *);
182: EXTERN PetscErrorCode TSGetApplicationContext(TS,void **);
184: EXTERN PetscErrorCode TSMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG *);
185: EXTERN PetscErrorCode TSMonitorLG(TS,PetscInt,PetscReal,Vec,void *);
186: EXTERN PetscErrorCode TSMonitorLGDestroy(PetscDrawLG);
188: /*
189: PETSc interface to Sundials
190: */
191: #ifdef PETSC_HAVE_SUNDIALS
192: typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
194: typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
196: EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType);
197: EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*);
198: EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal);
199: EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
200: EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
201: EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt);
202: EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal);
203: EXTERN PetscErrorCode TSSundialsSetExactFinalTime(TS,PetscTruth);
204: EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
205: #endif
207: EXTERN PetscErrorCode TSRKSetTolerance(TS,PetscReal);
210: #endif