DOLFIN
DOLFIN C++ interface
Public Types | Public Member Functions | List of all members
dolfin::CVode Class Reference

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 ()
 

Detailed Description

Wrapper class to SUNDIALS CVODE.

Constructor & Destructor Documentation

◆ CVode()

CVode::CVode ( LMM  cv_lmm,
ITER  cv_iter 
)

Constructor

Parameters
cv_lmmlinear multistep method
cv_iteriteration type

Member Function Documentation

◆ derivs()

void CVode::derivs ( double  t,
std::shared_ptr< GenericVector u,
std::shared_ptr< GenericVector udot 
)
virtual

Overloaded function for time derivatives of u at time t. Given the vector u, at time t, provide the time derivative udot.

Parameters
ttime
uinput vector of values u
udotoutput vector containing computed derivative of u at time t

◆ get_time()

double CVode::get_time ( ) const

Get current time

Returns
double current time

◆ init()

void CVode::init ( std::shared_ptr< GenericVector u0,
double  atol,
double  rtol,
long int  mxsteps = 0 
)

Initialise CVode

Parameters
u0Input vector
atolabsolute tolerance
rtolrelative tolerance
mxstepsmaximum number of steps

◆ jacobian()

int CVode::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

Given the values (t, y, fy, v), compute Jv = (df/dy)v

Parameters
vvector to be multiplied by the Jacobian df/dy
Jvoutput vector of (df/dy)*v
tcurrent value of the independent variable.
ycurrent value of the ependent variable.
fyvector f(t,y)
Returns
success flag, 0 if successful

◆ jacobian_setup()

int CVode::jacobian_setup ( double  t,
std::shared_ptr< GenericVector Jv,
std::shared_ptr< GenericVector y 
)
virtual

User-defined setup function called once per Newton iteration. Data structures for usage by the Jacobian function can be setup here

Parameters
tcurrent value of the independent variable
Jvcurrent value of the dependent variable vector, namely the predicted value of y(t).
yvector f(t,y).
Returns
success flag, 0 if success

◆ psolve()

int CVode::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 
)
virtual

Overloaded preconditioner solver function

Parameters
tncurrent value of the independent variable.
ycurrent value of the dependent variable vector.
fyvector f(t,y)
rright-hand side vector of the linear system.
zoutput vector computed by PrecSolve.
gammascalar appearing in the Newton matrix.
deltainput tolerance if an iterative method is used.
lrinput flag indicating whether to use left or right preconditioner.
Returns
success flag, 0 if success

Overloaded preconditioner solver function

◆ set_time()

void CVode::set_time ( double  t0)

Set the current time

Parameters
t0current time

◆ statistics()

std::map< std::string, double > CVode::statistics ( )

Return statistics

Returns
map structure containing information stored in the CVode structure, ie. number of solver steps, RHS evaluations, current time

◆ step()

double CVode::step ( double  dt)

Advance time by timestep dt

Parameters
dttimestep
Returns
CVODE return flag

The documentation for this class was generated from the following files: