Chaste Release::3.1
AbstractIvpOdeSolver Class Reference

#include <AbstractIvpOdeSolver.hpp>

Inheritance diagram for AbstractIvpOdeSolver:
Collaboration diagram for AbstractIvpOdeSolver:

List of all members.

Public Member Functions

virtual OdeSolution Solve (AbstractOdeSystem *pAbstractOdeSystem, std::vector< double > &rYValues, double startTime, double endTime, double timeStep, double timeSampling)=0
virtual void Solve (AbstractOdeSystem *pAbstractOdeSystem, std::vector< double > &rYValues, double startTime, double endTime, double timeStep)=0
virtual void SolveAndUpdateStateVariable (AbstractOdeSystem *pAbstractOdeSystem, double startTime, double endTime, double timeStep)
bool StoppingEventOccurred ()
double GetStoppingTime ()
 AbstractIvpOdeSolver ()
virtual ~AbstractIvpOdeSolver ()

Protected Attributes

bool mStoppingEventOccurred
double mStoppingTime

Private Member Functions

template<class Archive >
void serialize (Archive &archive, const unsigned int version)

Friends

class boost::serialization::access

Detailed Description

Abstract initial value problem ODE solver class. Sets up variables and functions for a numerical solution technique for an initial value ODE problem.

Definition at line 52 of file AbstractIvpOdeSolver.hpp.


Constructor & Destructor Documentation

AbstractIvpOdeSolver::AbstractIvpOdeSolver ( )

Constructor.

Definition at line 42 of file AbstractIvpOdeSolver.cpp.

AbstractIvpOdeSolver::~AbstractIvpOdeSolver ( ) [virtual]

Virtual destructor since we have virtual methods.

Definition at line 48 of file AbstractIvpOdeSolver.cpp.


Member Function Documentation

double AbstractIvpOdeSolver::GetStoppingTime ( )

Get the stopping time for the solver.

Returns:
mStoppingTime.

Definition at line 70 of file AbstractIvpOdeSolver.cpp.

References mStoppingTime.

template<class Archive >
void AbstractIvpOdeSolver::serialize ( Archive &  archive,
const unsigned int  version 
) [inline, private]

Archive the member variables.

Parameters:
archivethe archive
versionthe current version of this class

Reimplemented in AbstractOneStepIvpOdeSolver, BackwardEulerIvpOdeSolver, CvodeAdaptor, EulerIvpOdeSolver, HeunIvpOdeSolver, MockEulerIvpOdeSolver, RungeKutta2IvpOdeSolver, RungeKutta4IvpOdeSolver, and RungeKuttaFehlbergIvpOdeSolver.

Definition at line 63 of file AbstractIvpOdeSolver.hpp.

References mStoppingEventOccurred, and mStoppingTime.

virtual OdeSolution AbstractIvpOdeSolver::Solve ( AbstractOdeSystem pAbstractOdeSystem,
std::vector< double > &  rYValues,
double  startTime,
double  endTime,
double  timeStep,
double  timeSampling 
) [pure virtual]

Solves a system of ODEs using a specified one-step ODE solver and returns the solution as an OdeSolution object.

Parameters:
pAbstractOdeSystempointer to the concrete ODE system to be solved
rYValuesa standard vector specifying the intial condition of each solution variable in the system (this can be the initial conditions vector stored in the ODE system)
startTimethe time at which the initial conditions are specified
endTimethe time to which the system should be solved and the solution returned
timeStepthe time interval to be used by the solver
timeSamplingthe interval at which to sample the solution to the ODE system
Returns:
OdeSolution is an object containing an integer of the number of equations, a stdAbstractOdeSystem::vector of times and a std::vector of std::vectors where each of those vectors contains the solution for one variable of the ODE system at those times.

Implemented in AbstractOneStepIvpOdeSolver, CvodeAdaptor, and RungeKuttaFehlbergIvpOdeSolver.

Referenced by SolveAndUpdateStateVariable().

virtual void AbstractIvpOdeSolver::Solve ( AbstractOdeSystem pAbstractOdeSystem,
std::vector< double > &  rYValues,
double  startTime,
double  endTime,
double  timeStep 
) [pure virtual]

Second version of Solve. Solves a system of ODEs using a specified one-step ODE solver. This method does not return the solution and therefore does not take in a sampling time. Instead, the mStateVariables component in the ODE system object is updated.

Parameters:
pAbstractOdeSystempointer to the concrete ODE system to be solved
rYValuesa standard vector specifying the intial condition of each solution variable in the system (this can be the initial conditions vector stored in the ODE system)
startTimethe time at which the initial conditions are specified
endTimethe time to which the system should be solved and the solution returned
timeStepthe time interval to be used by the solver

Implemented in AbstractOneStepIvpOdeSolver, CvodeAdaptor, and RungeKuttaFehlbergIvpOdeSolver.

void AbstractIvpOdeSolver::SolveAndUpdateStateVariable ( AbstractOdeSystem pAbstractOdeSystem,
double  startTime,
double  endTime,
double  timeStep 
) [virtual]

Solves a system of ODEs using a specified one-step ODE solver and update the state variables.

Parameters:
pAbstractOdeSystempointer to the concrete ODE system to be solved
startTimethe time at which the initial conditions are specified
endTimethe time to which the system should be solved and the solution returned
timeStepthe time interval to be used by the solver

Definition at line 52 of file AbstractIvpOdeSolver.cpp.

References EXCEPTION, AbstractUntemplatedParameterisedSystem::GetNumberOfStateVariables(), AbstractParameterisedSystem< VECTOR >::rGetStateVariables(), and Solve().

Referenced by AbstractOdeBasedContractionModel::RunAndUpdate().

bool AbstractIvpOdeSolver::StoppingEventOccurred ( )

Determine whether the solver quit due to the ODE's stopping event triggering

Definition at line 65 of file AbstractIvpOdeSolver.cpp.

References mStoppingEventOccurred.


Friends And Related Function Documentation

friend class boost::serialization::access [friend]

Member Data Documentation

If a stopping event occurred the time is stored here. (Only valid when mStoppingEventOccurred==true)

Definition at line 78 of file AbstractIvpOdeSolver.hpp.

Referenced by GetStoppingTime(), RungeKuttaFehlbergIvpOdeSolver::InternalSolve(), AbstractOneStepIvpOdeSolver::InternalSolve(), serialize(), CvodeAdaptor::Solve(), and AbstractOneStepIvpOdeSolver::Solve().


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