|
DOLFIN
DOLFIN C++ interface
|
Wrapper class to SUNDIALS CVODE. More...
#include <CVode.h>
Public Types | |
| enum | LMM { cv_bdf = CV_BDF, cv_adams = CV_ADAMS } |
| enum | ITER { cv_functional = CV_FUNCTIONAL, cv_newton = CV_NEWTON } |
Public Member Functions | |
| CVode (LMM cv_lmm, ITER cv_iter) | |
| virtual | ~CVode () |
| Destructor. | |
| void | init (std::shared_ptr< GenericVector > u0, double atol, double rtol, long int mxsteps=0) |
| double | step (double dt) |
| double | get_time () const |
| void | set_time (double t0) |
| virtual void | derivs (double t, std::shared_ptr< GenericVector > u, std::shared_ptr< GenericVector > udot) |
| virtual int | jacobian (std::shared_ptr< const GenericVector > v, std::shared_ptr< GenericVector > Jv, double t, std::shared_ptr< const GenericVector > y, std::shared_ptr< const GenericVector > fy) |
| virtual int | jacobian_setup (double t, std::shared_ptr< GenericVector > Jv, std::shared_ptr< GenericVector > y) |
| virtual int | psolve (double tn, std::shared_ptr< const GenericVector >y, std::shared_ptr< const GenericVector > fy, std::shared_ptr< const GenericVector > r, std::shared_ptr< GenericVector > z, double gamma, double delta, int lr) |
| std::map< std::string, double > | statistics () |
Wrapper class to SUNDIALS CVODE.
| CVode::CVode | ( | LMM | cv_lmm, |
| ITER | cv_iter | ||
| ) |
Constructor
| cv_lmm | linear multistep method |
| cv_iter | iteration type |
|
virtual |
Overloaded function for time derivatives of u at time t. Given the vector u, at time t, provide the time derivative udot.
| t | time |
| u | input vector of values u |
| udot | output vector containing computed derivative of u at time t |
| double CVode::get_time | ( | ) | const |
Get current time
| void CVode::init | ( | std::shared_ptr< GenericVector > | u0, |
| double | atol, | ||
| double | rtol, | ||
| long int | mxsteps = 0 |
||
| ) |
Initialise CVode
| u0 | Input vector |
| atol | absolute tolerance |
| rtol | relative tolerance |
| mxsteps | maximum number of steps |
|
virtual |
Given the values (t, y, fy, v), compute Jv = (df/dy)v
| v | vector to be multiplied by the Jacobian df/dy |
| Jv | output vector of (df/dy)*v |
| t | current value of the independent variable. |
| y | current value of the ependent variable. |
| fy | vector f(t,y) |
|
virtual |
User-defined setup function called once per Newton iteration. Data structures for usage by the Jacobian function can be setup here
| t | current value of the independent variable |
| Jv | current value of the dependent variable vector, namely the predicted value of y(t). |
| y | vector f(t,y). |
|
virtual |
Overloaded preconditioner solver function
| tn | current value of the independent variable. |
| y | current value of the dependent variable vector. |
| fy | vector f(t,y) |
| r | right-hand side vector of the linear system. |
| z | output vector computed by PrecSolve. |
| gamma | scalar appearing in the Newton matrix. |
| delta | input tolerance if an iterative method is used. |
| lr | input flag indicating whether to use left or right preconditioner. |
Overloaded preconditioner solver function
| void CVode::set_time | ( | double | t0 | ) |
Set the current time
| t0 | current time |
| std::map< std::string, double > CVode::statistics | ( | ) |
Return statistics
| double CVode::step | ( | double | dt | ) |
Advance time by timestep dt
| dt | timestep |
1.8.11