A class for integrating ODEs with CVODES. More...
#include <ODE.h>
A class for integrating ODEs with CVODES.
Configuration Parameters:
Parameter Key | Type | Default Value | Description |
---|---|---|---|
"NonlinearSolver" | ptree | - | A child tree of options for the nonlinear solver. See ODEPiece::NonlinearSolverOptions. |
"LinearSolver" | ptree | - | A child tree of options for the linear solver. See ODEPiece::LinearSolverOptions. |
"Integrator" | ptree | - | A child tree of options for the nonlinear solver See ODEPiece::IntegratorOptions. |
Classes | |
struct | IntegratorOptions |
struct | Interface |
struct | LinearSolverOptions |
struct | NonlinearSolverOptions |
Public Member Functions | |
ODE (std::shared_ptr< ModPiece > const &rhsIn, Eigen::VectorXd const &evalTimesIn, boost::property_tree::ptree opts=boost::property_tree::ptree()) | |
ODE (std::shared_ptr< ModPiece > const &rhsIn, Eigen::VectorXd const &evalTimesIn, bool isAuto, boost::property_tree::ptree opts=boost::property_tree::ptree()) | |
virtual | ~ODE ()=default |
Public Member Functions inherited from muq::Modeling::ModPiece | |
ModPiece (std::vector< int > const &inputSizes, std::vector< int > const &outputSizes) | |
ModPiece (Eigen::VectorXi const &inputSizes, Eigen::VectorXi const &outputSizes) | |
virtual | ~ModPiece ()=default |
virtual double | GetRunTime (const std::string &method="Evaluate") const override |
Get the average run time for one of the implemented methods. More... | |
virtual void | ResetCallTime () override |
Resets the number of call and times. More... | |
virtual unsigned long int | GetNumCalls (const std::string &method="Evaluate") const override |
get the number of times one of the implemented methods has been called. More... | |
virtual std::vector< Eigen::VectorXd > const & | Evaluate (std::vector< Eigen::VectorXd > const &input) |
Evaluate the ModPiece. More... | |
virtual std::vector< Eigen::VectorXd > const & | Evaluate (ref_vector< Eigen::VectorXd > const &input) |
VARIADIC_TO_REFVECTOR (Evaluate, Eigen::VectorXd, std::vector< Eigen::VectorXd > const &) | |
virtual Eigen::VectorXd const & | Gradient (unsigned int const outputDimWrt, unsigned int const inputDimWrt, std::vector< Eigen::VectorXd > const &input, Eigen::VectorXd const &sensitivity) |
Compute the Gradient \(J^Tv\). More... | |
virtual Eigen::VectorXd const & | Gradient (unsigned int const outputDimWrt, unsigned int const inputDimWrt, ref_vector< Eigen::VectorXd > const &input, Eigen::VectorXd const &sensitivity) |
Eigen::VectorXd const & | Gradient (unsigned int outWrt, unsigned int inWrt, Eigen::VectorXd const &last, Eigen::VectorXd const &sens) |
template<typename... Args> | |
Eigen::VectorXd const & | Gradient (unsigned int wrtOut, unsigned int wrtIn, Args const &... args) |
Eigen::VectorXd const & | ApplyHessian (unsigned int outWrt, unsigned int inWrt1, unsigned int inWrt2, Eigen::VectorXd const &last, Eigen::VectorXd const &sens, Eigen::VectorXd const &vec) |
template<typename... Args> | |
Eigen::VectorXd const & | ApplyHessian (unsigned int wrtOut, unsigned int wrtIn1, unsigned int wrtIn2, Args const &... args) |
virtual Eigen::MatrixXd const & | Jacobian (unsigned int const outputDimWrt, unsigned int const inputDimWrt, std::vector< Eigen::VectorXd > const &input) |
Compute the Jacobian of this ModPiece. More... | |
virtual Eigen::MatrixXd const & | Jacobian (unsigned int const outputDimWrt, unsigned int const inputDimWrt, ref_vector< Eigen::VectorXd > const &input) |
template<typename... Args> | |
Eigen::MatrixXd const & | Jacobian (unsigned int outWrt, unsigned int inWrt, Args const &... args) |
template<typename... Args> | |
Eigen::MatrixXd | JacobianByFD (unsigned int outWrt, unsigned int inWrt, Args const &... args) |
template<typename... Args> | |
Eigen::MatrixXd | ApplyJacobianByFD (unsigned int outWrt, unsigned int inWrt, Args const &... args) |
virtual Eigen::VectorXd const & | ApplyJacobian (unsigned int const outputDimWrt, unsigned int const inputDimWrt, std::vector< Eigen::VectorXd > const &input, Eigen::VectorXd const &vec) |
Apply the Jacobian of this ModPiece to a vector. More... | |
virtual Eigen::VectorXd const & | ApplyJacobian (unsigned int const outputDimWrt, unsigned int const inputDimWrt, ref_vector< Eigen::VectorXd > const &input, Eigen::VectorXd const &vec) |
Eigen::VectorXd const & | ApplyJacobian (unsigned int outWrt, unsigned int inWrt, Eigen::VectorXd const &last, Eigen::VectorXd const &sens) |
template<typename... Args> | |
Eigen::VectorXd const & | ApplyJacobian (unsigned int wrtOut, unsigned int wrtIn, Args const &... args) |
virtual Eigen::VectorXd | GradientByFD (unsigned int const outputDimWrt, unsigned int const inputDimWrt, std::vector< Eigen::VectorXd > const &input, Eigen::VectorXd const &sensitivity) |
virtual Eigen::VectorXd | GradientByFD (unsigned int const outputDimWrt, unsigned int const inputDimWrt, ref_vector< Eigen::VectorXd > const &input, Eigen::VectorXd const &sensitivity) |
virtual Eigen::MatrixXd | JacobianByFD (unsigned int const outputDimWrt, unsigned int const inputDimWrt, std::vector< Eigen::VectorXd > const &input) |
virtual Eigen::MatrixXd | JacobianByFD (unsigned int const outputDimWrt, unsigned int const inputDimWrt, ref_vector< Eigen::VectorXd > const &input) |
virtual Eigen::VectorXd | ApplyJacobianByFD (unsigned int const outputDimWrt, unsigned int const inputDimWrt, std::vector< Eigen::VectorXd > const &input, Eigen::VectorXd const &vec) |
virtual Eigen::VectorXd | ApplyJacobianByFD (unsigned int const outputDimWrt, unsigned int const inputDimWrt, ref_vector< Eigen::VectorXd > const &input, Eigen::VectorXd const &vec) |
virtual Eigen::VectorXd const & | ApplyHessian (unsigned int const outWrt, unsigned int const inWrt1, unsigned int const inWrt2, std::vector< Eigen::VectorXd > const &input, Eigen::VectorXd const &sens, Eigen::VectorXd const &vec) |
virtual Eigen::VectorXd const & | ApplyHessian (unsigned int const outWrt, unsigned int const inWrt1, unsigned int const inWrt2, ref_vector< Eigen::VectorXd > const &input, Eigen::VectorXd const &sens, Eigen::VectorXd const &vec) |
virtual Eigen::VectorXd | ApplyHessianByFD (unsigned int const outWrt, unsigned int const inWrt1, unsigned int const inWrt2, std::vector< Eigen::VectorXd > const &input, Eigen::VectorXd const &sens, Eigen::VectorXd const &vec) |
virtual Eigen::VectorXd | ApplyHessianByFD (unsigned int const outWrt, unsigned int const inWrt1, unsigned int const inWrt2, ref_vector< Eigen::VectorXd > const &input, Eigen::VectorXd const &sens, Eigen::VectorXd const &vec) |
void | EnableCache () |
void | DisableCache () |
bool | CacheStatus () const |
virtual void | SetWarnLevel (unsigned int newLevel) |
Public Member Functions inherited from muq::Modeling::WorkPiece | |
WorkPiece () | |
Create a muq::Modeling::WorkPiece with no fixed number of inputs and outputs and variable input/output types. More... | |
WorkPiece (int const num, WorkPiece::Fix const fix=WorkPiece::Fix::Inputs) | |
Create a muq::Modeling::WorkPiece with either a fixed number of inputs or outputs and variable input/output types. More... | |
WorkPiece (int const numIns, int const numOuts) | |
Create a muq::Modeling::WorkPiece with a fixed number of inputs and outputs but variable input/output types. More... | |
WorkPiece (std::vector< std::string > const &types, WorkPiece::Fix const fix=WorkPiece::Fix::Inputs) | |
Create a muq::Modeling::WorkPiece with either a fixed number of inputs with specified types or a fixed number of outputs with specified types. More... | |
WorkPiece (std::map< unsigned int, std::string > const &types, WorkPiece::Fix const fix=WorkPiece::Fix::Inputs) | |
Create a muq::Modeling::WorkPiece where either some of the inputs have specified types or some of the outputs have specified types. More... | |
WorkPiece (std::map< unsigned int, std::string > const &types, int const num, WorkPiece::Fix const fixTypes=WorkPiece::Fix::Inputs, WorkPiece::Fix const fixNum=WorkPiece::Fix::Inputs) | |
Create a muq::Modeling::WorkPiece where either some of the inputs have specified types or some of the outputs have specified types and either the number of inputs or the number of outputs is fixed. More... | |
WorkPiece (std::vector< std::string > const &types, int const num) | |
Create a muq::Modeling::WorkPiece with a fixed number of inputs with specified types and a fixed number of outputs (of uknown type) More... | |
WorkPiece (int const num, std::vector< std::string > const &types) | |
Create a muq::Modeling::WorkPiece with a fixed number of outputs with specified types and a fixed number of inputs (of uknown type) More... | |
WorkPiece (std::map< unsigned int, std::string > const &inTypes, int const numIns, int const numOuts) | |
Create a muq::Modeling::WorkPiece where some of the inputs are known and we know the input and output numbers. More... | |
WorkPiece (int const numIns, std::map< unsigned int, std::string > const &outTypes, int const numOuts) | |
Create a muq::Modeling::WorkPiece where some of the outputs are known and we know the input and output numbers. More... | |
WorkPiece (std::vector< std::string > const &inTypes, std::vector< std::string > const &outTypes) | |
Create a muq::Modeling::WorkPiece with a fixed number of inputs and outputs with specified types. More... | |
WorkPiece (std::map< unsigned int, std::string > const &inTypes, std::vector< std::string > const &outTypes) | |
Create a muq::Modeling::WorkPiece where some of the inputs are known and all of the outputs have specified types. More... | |
WorkPiece (std::map< unsigned int, std::string > const &inTypes, int const num, std::vector< std::string > const &outTypes) | |
Create a muq::Modeling::WorkPiece where some of the inputs are known with a known number of inputs and all of the outputs have specified types. More... | |
WorkPiece (std::vector< std::string > const &inTypes, std::map< unsigned int, std::string > const &outTypes) | |
Create a muq::Modeling::WorkPiece where some of the outputs and all of the inputs have specified types. More... | |
WorkPiece (std::vector< std::string > const &inTypes, std::map< unsigned int, std::string > const &outTypes, int const num) | |
Create a muq::Modeling::WorkPiece where some of the outputs with a known number of outputs and all of the inputs have specified types. More... | |
WorkPiece (std::map< unsigned int, std::string > const &inTypes, std::map< unsigned int, std::string > const &outTypes) | |
Create a muq::Mdoeling::WorkPiece where some of the inputs and some of the outputs have specified types. More... | |
WorkPiece (std::map< unsigned int, std::string > const &inTypes, int const numIn, std::map< unsigned int, std::string > const &outTypes) | |
Create a muq::Mdoeling::WorkPiece where some of the inputs and some of the outputs have specified types with a fixed number of inputs. More... | |
WorkPiece (std::map< unsigned int, std::string > const &inTypes, std::map< unsigned int, std::string > const &outTypes, int const numOut) | |
Create a muq::Mdoeling::WorkPiece where some of the inputs and some of the outputs have specified types with a fixed number of outputs. More... | |
WorkPiece (std::map< unsigned int, std::string > const &inTypes, int const numIn, std::map< unsigned int, std::string > const &outTypes, int const numOut) | |
Create a muq::Mdoeling::WorkPiece where some of the inputs and some of the outputs have specified types with a fixed number of inputs and outputs. More... | |
virtual | ~WorkPiece () |
Default destructor. More... | |
std::vector< boost::any > const & | Evaluate (std::vector< boost::any > const &ins) |
Evaluate this muq::Modeling::WorkPiece. More... | |
std::vector< boost::any > const & | Evaluate (ref_vector< boost::any > const &ins) |
Evaluate this muq::Modeling::WorkPiece using references to the inputs. More... | |
std::vector< boost::any > const & | Evaluate () |
Evaluate this muq::Modeling::WorkPiece in the case that there are no inputs. More... | |
template<typename... Args> | |
std::vector< boost::any > const & | Evaluate (Args... args) |
Evalaute this muq::Modeling::WorkPiece using multiple arguments. More... | |
std::string const & | Name () |
Get the (unique) name of this work piece. More... | |
void | SetName (std::string const &newName) |
Set the name of this work piece. More... | |
std::string | InputType (unsigned int inputNum, bool const demangle=true) const |
Get the input type (if we know it) for a specific input. More... | |
int | InputSize (unsigned int inputNum) const |
Get the length of a vector valued input with fixed size. More... | |
void | SetInputSize (unsigned int inputNum, int newSize) |
std::string | OutputType (unsigned int outputNum, bool const demangle=true) const |
Get the output type (if we know it) for a specific output. More... | |
std::map< unsigned int, std::string > | OutputTypes () const |
Get the output types. More... | |
std::map< unsigned int, std::string > | InputTypes () const |
Get the input types. More... | |
unsigned int | ID () const |
Get the unique ID number. More... | |
Public Attributes | |
const bool | isAutonomous |
Does the RHS explicitly depend on time? More... | |
const double | startTime |
The time when this simulation starts – either set in the opts ptree or set to the minimum of evalTimes. More... | |
Public Attributes inherited from muq::Modeling::ModPiece | |
const Eigen::VectorXi | inputSizes |
const Eigen::VectorXi | outputSizes |
Public Attributes inherited from muq::Modeling::WorkPiece | |
int | numInputs |
The number of inputs. More... | |
int | numOutputs |
The number of outputs. More... | |
Additional Inherited Members | |
Static Public Member Functions inherited from muq::Modeling::WorkPiece | |
static ref_vector< const boost::any > | ToRefVector (std::vector< boost::any > const &anyVec) |
Create vector of references from a vector of boost::any's. More... | |
static ref_vector< const Eigen::VectorXd > | ToRefVector (std::vector< Eigen::VectorXd > const &anyVec) |
|
inline |
ODE::ODE | ( | std::shared_ptr< ModPiece > const & | rhsIn, |
Eigen::VectorXd const & | evalTimesIn, | ||
bool | isAuto, | ||
boost::property_tree::ptree | opts = boost::property_tree::ptree() |
||
) |
Definition at line 24 of file ODE.cpp.
References intOpts, linOpts, nonlinOpts, muq::Modeling::ODE::NonlinearSolverOptions::SetOptions(), muq::Modeling::ODE::LinearSolverOptions::SetOptions(), and muq::Modeling::ODE::IntegratorOptions::SetOptions().
|
virtualdefault |
|
overrideprotectedvirtual |
Reimplemented from muq::Modeling::ModPiece.
Definition at line 822 of file ODE.cpp.
References muq::Modeling::ODE::IntegratorOptions::abstol, muq::Modeling::ModPiece::cacheEnabled, muq::Modeling::ModPiece::cacheInput, CheckFlag(), CreateLinearSolver(), CreateNonlinearSolver(), nlohmann::detail::dtoa_impl::e, evalTimes, intOpts, isAutonomous, muq::Modeling::ModPiece::jacobianAction, muq::Modeling::ODE::IntegratorOptions::maxNumSteps, muq::Modeling::ODE::NonlinearSolverOptions::method, muq::Modeling::ODE::IntegratorOptions::method, nonlinOpts, muq::Modeling::ModPiece::outputs, muq::Modeling::ODE::IntegratorOptions::reltol, muq::Modeling::ODE::Interface::RHS(), rhs, and startTime.
|
staticprivate |
Checks the return value produced by CVODES. Throws and exception if error occured.
Definition at line 1212 of file ODE.cpp.
Referenced by ApplyJacobianImpl(), CreateAdjointLinearSolver(), CreateAdjointNonlinearSolver(), CreateLinearSolver(), CreateNonlinearSolver(), EvaluateImpl(), GradientImpl(), and JacobianImpl().
|
private |
Definition at line 218 of file ODE.cpp.
References CheckFlag(), linOpts, muq::Modeling::ODE::LinearSolverOptions::maxIters, muq::Modeling::ODE::NonlinearSolverOptions::method, muq::Modeling::ODE::LinearSolverOptions::method, and nonlinOpts.
Referenced by GradientImpl().
|
private |
Definition at line 273 of file ODE.cpp.
References CheckFlag(), muq::Modeling::ODE::NonlinearSolverOptions::maxIters, muq::Modeling::ODE::NonlinearSolverOptions::method, and nonlinOpts.
Referenced by GradientImpl().
|
private |
Creates a linear solver and, if necessary, a dense matrix for storing the Jacobian. If the nonlinear solver is fixed point iteration, then no linear solver is necessary and this function returns an empty linear solver. If no Jacobian is needed, either because a fixed point nonlinear solver is specified or because an iterative linear solver is specified, then NULL is returned. The options read by the constructor and stored in linOpts
are used to construct the solver.
Definition at line 133 of file ODE.cpp.
References CheckFlag(), linOpts, muq::Modeling::ODE::LinearSolverOptions::maxIters, muq::Modeling::ODE::NonlinearSolverOptions::method, muq::Modeling::ODE::LinearSolverOptions::method, nonlinOpts, muq::Modeling::ODE::Interface::RHSJacobian(), and muq::Modeling::ODE::Interface::RHSJacobianAction().
Referenced by ApplyJacobianImpl(), EvaluateImpl(), GradientImpl(), and JacobianImpl().
|
private |
Definition at line 189 of file ODE.cpp.
References CheckFlag(), muq::Modeling::ODE::NonlinearSolverOptions::maxIters, muq::Modeling::ODE::NonlinearSolverOptions::method, and nonlinOpts.
Referenced by ApplyJacobianImpl(), EvaluateImpl(), GradientImpl(), and JacobianImpl().
|
overrideprotectedvirtual |
Integrates the ODE forward in time.
Implements muq::Modeling::ModPiece.
Definition at line 302 of file ODE.cpp.
References muq::Modeling::ODE::IntegratorOptions::abstol, CheckFlag(), CreateLinearSolver(), CreateNonlinearSolver(), nlohmann::detail::dtoa_impl::e, evalTimes, intOpts, isAutonomous, muq::Modeling::ODE::IntegratorOptions::maxNumSteps, muq::Modeling::ODE::IntegratorOptions::method, muq::Modeling::ModPiece::outputs, muq::Modeling::ModPiece::outputSizes, muq::Modeling::ODE::IntegratorOptions::reltol, muq::Modeling::ODE::Interface::RHS(), rhs, and startTime.
|
staticprivate |
|
staticprivate |
|
staticprivate |
|
overrideprotectedvirtual |
Reimplemented from muq::Modeling::ModPiece.
Definition at line 399 of file ODE.cpp.
References muq::Modeling::ODE::IntegratorOptions::abstol, muq::Modeling::ModPiece::cacheEnabled, muq::Modeling::ModPiece::cacheInput, CheckFlag(), muq::Modeling::ODE::IntegratorOptions::checkPtGap, CreateAdjointLinearSolver(), CreateAdjointNonlinearSolver(), CreateLinearSolver(), CreateNonlinearSolver(), nlohmann::detail::dtoa_impl::e, evalTimes, muq::Modeling::ModPiece::gradient, muq::Modeling::ModPiece::inputSizes, intOpts, isAutonomous, muq::Modeling::ODE::IntegratorOptions::maxNumSteps, muq::Modeling::ODE::IntegratorOptions::method, muq::Modeling::ModPiece::outputs, muq::Modeling::ModPiece::outputSizes, muq::Modeling::ODE::IntegratorOptions::reltol, muq::Modeling::ODE::Interface::RHS(), rhs, and startTime.
|
overrideprotectedvirtual |
Reimplemented from muq::Modeling::ModPiece.
Definition at line 638 of file ODE.cpp.
References muq::Modeling::ODE::IntegratorOptions::abstol, muq::Modeling::ModPiece::cacheEnabled, muq::Modeling::ModPiece::cacheInput, CheckFlag(), CreateLinearSolver(), CreateNonlinearSolver(), nlohmann::detail::dtoa_impl::e, evalTimes, muq::Modeling::ModPiece::inputSizes, intOpts, isAutonomous, muq::Modeling::ModPiece::jacobian, muq::Modeling::ODE::IntegratorOptions::maxNumSteps, muq::Modeling::ODE::NonlinearSolverOptions::method, muq::Modeling::ODE::IntegratorOptions::method, nonlinOpts, muq::Modeling::ModPiece::outputs, muq::Modeling::ModPiece::outputSizes, muq::Modeling::ODE::IntegratorOptions::reltol, muq::Modeling::ODE::Interface::RHS(), rhs, muq::Modeling::ODE::IntegratorOptions::sensAbsTol, muq::Modeling::ODE::IntegratorOptions::sensRelTol, and startTime.
|
protected |
Definition at line 189 of file ODE.h.
Referenced by ApplyJacobianImpl(), EvaluateImpl(), GetOutputSizes(), GradientImpl(), and JacobianImpl().
|
protected |
Definition at line 186 of file ODE.h.
Referenced by ApplyJacobianImpl(), EvaluateImpl(), GradientImpl(), JacobianImpl(), and ODE().
const bool muq::Modeling::ODE::isAutonomous |
Does the RHS explicitly depend on time?
Definition at line 48 of file ODE.h.
Referenced by ApplyJacobianImpl(), EvaluateImpl(), GradientImpl(), and JacobianImpl().
|
protected |
Definition at line 185 of file ODE.h.
Referenced by CreateAdjointLinearSolver(), CreateLinearSolver(), and ODE().
|
protected |
Definition at line 184 of file ODE.h.
Referenced by ApplyJacobianImpl(), CreateAdjointLinearSolver(), CreateAdjointNonlinearSolver(), CreateLinearSolver(), CreateNonlinearSolver(), JacobianImpl(), and ODE().
|
protected |
Definition at line 192 of file ODE.h.
Referenced by ApplyJacobianImpl(), EvaluateImpl(), GradientImpl(), and JacobianImpl().
const double muq::Modeling::ODE::startTime |
The time when this simulation starts – either set in the opts ptree or set to the minimum of evalTimes.
Definition at line 51 of file ODE.h.
Referenced by ApplyJacobianImpl(), EvaluateImpl(), GetStartTime(), GradientImpl(), and JacobianImpl().