24 #ifndef _DYNAMICALSYSTEM_H_ 25 #define _DYNAMICALSYSTEM_H_ 31 #include <eigen3/Eigen/Core> 106 const Eigen::Ref<const Eigen::VectorXd>& x,
107 Eigen::Ref<Eigen::VectorXd> xd
120 virtual void analyticalSolution(
const Eigen::VectorXd& ts, Eigen::MatrixXd& xs, Eigen::MatrixXd& xds)
const = 0;
130 virtual void integrateStart(Eigen::Ref<Eigen::VectorXd> x, Eigen::Ref<Eigen::VectorXd> xd)
const;
139 void integrateStart(
const Eigen::VectorXd& x_init, Eigen::Ref<Eigen::VectorXd> x, Eigen::Ref<Eigen::VectorXd> xd);
152 virtual void integrateStep(
double dt,
const Eigen::Ref<const Eigen::VectorXd> x, Eigen::Ref<Eigen::VectorXd> x_updated, Eigen::Ref<Eigen::VectorXd> xd_updated)
const;
177 virtual std::string
toString(
void)
const = 0;
196 integration_method_ = integration_method;
203 inline int dim(
void)
const {
229 inline double tau(
void)
const {
return tau_; }
252 initial_state=initial_state_;
259 assert(initial_state.size()==dim_orig_);
275 attractor_state=attractor_state_;
282 assert(attractor_state.size()==dim_orig_);
290 inline std::string
name(
void)
const {
return name_; }
320 void integrateStepEuler(
double dt,
const Eigen::Ref<const Eigen::VectorXd> x, Eigen::Ref<Eigen::VectorXd> x_updated, Eigen::Ref<Eigen::VectorXd> xd_updated)
const;
330 void integrateStepRungeKutta(
double dt,
const Eigen::Ref<const Eigen::VectorXd> x, Eigen::Ref<Eigen::VectorXd> x_updated, Eigen::Ref<Eigen::VectorXd> xd_updated)
const;
351 Eigen::VectorXd initial_state_;
359 Eigen::VectorXd attractor_state_;
386 template<
class Archive>
387 void serialize(Archive & ar,
const unsigned int version)
389 ar & BOOST_SERIALIZATION_NVP(dim_);
390 ar & BOOST_SERIALIZATION_NVP(dim_orig_);
393 ar & BOOST_SERIALIZATION_NVP(tau_);
394 ar & BOOST_SERIALIZATION_NVP(initial_state_);
395 ar & BOOST_SERIALIZATION_NVP(attractor_state_);
396 ar & BOOST_SERIALIZATION_NVP(name_);
397 ar & BOOST_SERIALIZATION_NVP(integration_method_);
404 #include <boost/serialization/assume_abstract.hpp> 408 #include <boost/serialization/export.hpp> 412 #endif // _DYNAMICALSYSTEM_H_ double tau(void) const
Accessor function for the time constant.
void set_dim(int dim)
Set the dimensionality of the dynamical system, i.e.
virtual void analyticalSolution(const Eigen::VectorXd &ts, Eigen::MatrixXd &xs, Eigen::MatrixXd &xds) const =0
Return analytical solution of the system at certain times.
IntegrationMethod
The possible integration methods that can be used.
virtual void set_attractor_state(const Eigen::Ref< const Eigen::VectorXd > &attractor_state)
Mutator function for the attractor state of the dynamical system.
int dim_orig(void) const
Get the dimensionality of the dynamical system, i.e.
BOOST_CLASS_IMPLEMENTATION(DmpBbo::DynamicalSystem, boost::serialization::object_serializable)
Don't add version information to archives.
void set_integration_method(IntegrationMethod integration_method)
Choose the integration method.
void initial_state(Eigen::VectorXd &initial_state) const
Accessor function for the initial state of the dynamical system.
Eigen::VectorXd initial_state(void) const
Accessor function for the initial state of the dynamical system.
virtual void integrateStep(double dt, const Eigen::Ref< const Eigen::VectorXd > x, Eigen::Ref< Eigen::VectorXd > x_updated, Eigen::Ref< Eigen::VectorXd > xd_updated) const
Integrate the system one time step.
DynamicalSystem(void)
Default constructor.
virtual void set_initial_state(const Eigen::VectorXd &initial_state)
Mutator function for the initial state of the dynamical system.
int dim(void) const
Get the dimensionality of the dynamical system, i.e.
virtual DynamicalSystem * clone(void) const =0
Return a pointer to a deep copy of the DynamicalSystem object.
virtual void differentialEquation(const Eigen::Ref< const Eigen::VectorXd > &x, Eigen::Ref< Eigen::VectorXd > xd) const =0
The differential equation which defines the system.
BOOST_SERIALIZATION_ASSUME_ABSTRACT(DmpBbo::DynamicalSystem)
Don't add version information to archives.
std::string name(void) const
Accessor function for the name of the dynamical system.
virtual ~DynamicalSystem(void)
Destructor.
Interface for implementing dynamical systems.
Header file for adding real-time debugging using several macros.
virtual std::string toString(void) const =0
Returns a string representation of the object.
virtual void integrateStart(Eigen::Ref< Eigen::VectorXd > x, Eigen::Ref< Eigen::VectorXd > xd) const
Start integrating the system.
friend class boost::serialization::access
Give boost serialization access to private members.
virtual void set_name(std::string name)
Mutator function for the name of the dynamical system.
Eigen::VectorXd attractor_state(void) const
Accessor function for the attractor state of the dynamical system.
void attractor_state(Eigen::VectorXd &attractor_state) const
Accessor function for the attractor state of the dynamical system.
virtual void set_tau(double tau)
Mutator function for the time constant.
Header file for serialization of Eigen matrices.
friend std::ostream & operator<<(std::ostream &output, const DynamicalSystem &dyn_sys)
Write a DynamicalSystem to an output stream.