25 #ifndef SRC_SYSTEMS_INCLUDE_CORBO_SYSTEMS_BENCHMARK_NONLINEAR_BENCHMARK_SYSTEMS_H_ 26 #define SRC_SYSTEMS_INCLUDE_CORBO_SYSTEMS_BENCHMARK_NONLINEAR_BENCHMARK_SYSTEMS_H_ 39 Ptr getInstance()
const override {
return std::make_shared<VanDerPolOscillator>(); }
44 bool isLinear()
const override {
return false; }
56 assert(x.size() == f.size() &&
"VanDerPolOscillator::dynamics(): x and f are not of the same size, do not forget to pre-allocate f.");
59 f[1] = -
_a * (x[0] * x[0] - 1) * x[1] - x[0] + u[0];
66 #ifdef MESSAGE_SUPPORT 68 void toMessage(messages::SystemDynamics& message)
const override 70 SystemDynamicsInterface::toMessage(message);
72 message.mutable_van_der_pol_oscillator()->set_a(
_a);
75 void fromMessage(
const messages::SystemDynamics& message, std::stringstream* issues)
override 77 SystemDynamicsInterface::fromMessage(message, issues);
79 _a = message.van_der_pol_oscillator().a();
95 Ptr getInstance()
const override {
return std::make_shared<DuffingOscillator>(); }
112 assert(x.size() == f.size() &&
"DuffingOscillator::dynamics(): x and f are not of the same size, do not forget to pre-allocate f.");
115 f[1] = -_damping * x[1] - _spring_alpha * x[0] - _spring_beta * x[0] * x[0] * x[0] + u[0];
122 _spring_alpha = spring_alpha;
123 _spring_beta = spring_beta;
126 #ifdef MESSAGE_SUPPORT 128 void toMessage(messages::SystemDynamics& message)
const override 130 SystemDynamicsInterface::toMessage(message);
132 message.mutable_duffing_oscillator()->set_damping(_damping);
133 message.mutable_duffing_oscillator()->set_spring_alpha(_spring_alpha);
134 message.mutable_duffing_oscillator()->set_spring_beta(_spring_beta);
137 void fromMessage(
const messages::SystemDynamics& message, std::stringstream* issues)
override 139 SystemDynamicsInterface::fromMessage(message, issues);
141 _damping = message.duffing_oscillator().damping();
142 _spring_alpha = message.duffing_oscillator().spring_alpha();
143 _spring_beta = message.duffing_oscillator().spring_beta();
149 double _spring_alpha = 1;
150 double _spring_beta = 1;
178 assert(x.size() == f.size() &&
"FreeSpaceRocket::dynamics(): x and f are not of the same size, do not forget to pre-allocate f.");
181 f[1] = (u[0] - 0.02 * x[1] * x[1]) / x[2];
182 f[2] = -0.01 * u[0] * u[0];
211 assert(x.size() == f.size() &&
"SimplePendulum::dynamics(): x and f are not of the same size, do not forget to pre-allocate f.");
214 f[1] = u[0] - _rho / (_m * _l * _l) * x[1] - _g / _l *
std::sin(x[0]);
218 void setParameters(
double mass,
double length,
double gravitation,
double friction)
226 #ifdef MESSAGE_SUPPORT 228 void toMessage(messages::SystemDynamics& message)
const override 230 SystemDynamicsInterface::toMessage(message);
232 messages::SimplePendulum* msg = message.mutable_simple_pendulum();
240 void fromMessage(
const messages::SystemDynamics& message, std::stringstream* issues)
override 242 SystemDynamicsInterface::fromMessage(message, issues);
244 const messages::SimplePendulum& msg = message.simple_pendulum();
268 Ptr getInstance()
const override {
return std::make_shared<MasslessPendulum>(); }
285 assert(x.size() == f.size() &&
"MasslessPendulum::dynamics(): x and f are not of the same size, do not forget to pre-allocate f.");
288 f[1] = u[0] - _omega0 *
std::sin(x[0]);
294 #ifdef MESSAGE_SUPPORT 296 void toMessage(messages::SystemDynamics& message)
const override 298 SystemDynamicsInterface::toMessage(message);
299 messages::MasslessPendulum* msg = message.mutable_massless_pendulum();
300 msg->set_omega0(_omega0);
303 void fromMessage(
const messages::SystemDynamics& message, std::stringstream* issues)
override 305 SystemDynamicsInterface::fromMessage(message, issues);
307 const messages::MasslessPendulum& msg = message.massless_pendulum();
308 _omega0 = msg.omega0();
313 double _omega0 = 1.0;
341 assert(x.size() == f.size() &&
"CartPoleSystem::dynamics(): x and f are not of the same size, do not forget to pre-allocate f.");
345 double sin_phi_phidot_sq =
std::sin(x[1]) * x[3] * x[3];
350 f[2] = (_l * _mp * sin_phi_phidot_sq + u[0] + _mp * _g *
std::cos(x[1]) *
std::sin(x[1])) / denum;
351 f[3] = -(_l * _mp *
std::cos(x[1]) * sin_phi_phidot_sq + u[0] *
std::cos(x[1]) + (_mp + _mc) * _g *
std::sin(x[1])) / (_l * denum);
363 #ifdef MESSAGE_SUPPORT 365 void toMessage(messages::SystemDynamics& message)
const override 367 SystemDynamicsInterface::toMessage(message);
369 messages::CartPole* msg = message.mutable_cart_pole();
377 void fromMessage(
const messages::SystemDynamics& message, std::stringstream* issues)
override 379 SystemDynamicsInterface::fromMessage(message, issues);
381 const messages::CartPole& msg = message.cart_pole();
430 assert(x.size() == f.size() &&
"ToyExample::dynamics(): x and f are not of the same size, do not forget to pre-allocate f.");
434 f[0] = x[1] + u[0] * (_mu + (1.0 - _mu) * x[0]);
435 f[1] = x[0] + u[0] * (_mu - 4.0 * (1.0 - _mu) * x[1]);
441 #ifdef MESSAGE_SUPPORT 443 void toMessage(messages::SystemDynamics& message)
const override 445 SystemDynamicsInterface::toMessage(message);
447 message.mutable_toy_example()->set_mu(_mu);
450 void fromMessage(
const messages::SystemDynamics& message, std::stringstream* issues)
override 452 SystemDynamicsInterface::fromMessage(message, issues);
454 _mu = message.toy_example().mu();
487 assert(x.size() == f.size() &&
"ArtsteinsCircle::dynamics(): x and f are not of the same size, do not forget to pre-allocate f.");
489 f[0] = (x[0] * x[0] - x[1] * x[1]) * u[0];
490 f[1] = 2 * x[0] * x[1] * u[0];
493 #ifdef MESSAGE_SUPPORT 495 void toMessage(messages::SystemDynamics& message)
const override 497 SystemDynamicsInterface::toMessage(message);
499 message.mutable_artsteins_circle();
502 void fromMessage(
const messages::SystemDynamics& message, std::stringstream* issues)
override 504 SystemDynamicsInterface::fromMessage(message, issues);
512 #endif // SRC_SYSTEMS_INCLUDE_CORBO_SYSTEMS_BENCHMARK_NONLINEAR_BENCHMARK_SYSTEMS_H_
int getStateDimension() const override
Return state dimension (x)
void setParameters(double damping, double spring_alpha, double spring_beta)
EIGEN_DEVICE_FUNC const CosReturnType cos() const
ToyExample()
Default constructor.
Ptr getInstance() const override
Return a newly created shared instance of the implemented class.
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC half pow(const half &a, const half &b)
bool isLinear() const override
Check if the system dynamics are linear.
Ptr getInstance() const override
Return a newly created shared instance of the implemented class.
int getStateDimension() const override
Return state dimension (x)
VanDerPolOscillator()
Default constructor.
bool isLinear() const override
Check if the system dynamics are linear.
void dynamics(const Eigen::Ref< const StateVector > &x, const Eigen::Ref< const ControlVector > &u, Eigen::Ref< StateVector > f) const override
Evaluate the system dynamics equation.
bool isLinear() const override
Check if the system dynamics are linear.
FreeSpaceRocket()
Default constructor.
int getInputDimension() const override
Return the plant input dimension (u)
void setParameters(double mu)
bool isContinuousTime() const override
Check if the system dynamics are defined in continuous-time.
Ptr getInstance() const override
Return a newly created shared instance of the implemented class.
bool isContinuousTime() const override
Check if the system dynamics are defined in continuous-time.
int getInputDimension() const override
Return the plant input dimension (u)
void setParameters(double mass, double length, double gravitation, double friction)
bool isContinuousTime() const override
Check if the system dynamics are defined in continuous-time.
void setParameter(double omega0)
int getStateDimension() const override
Return state dimension (x)
bool isContinuousTime() const override
Check if the system dynamics are defined in continuous-time.
int getStateDimension() const override
Return state dimension (x)
int getInputDimension() const override
Return the plant input dimension (u)
int getStateDimension() const override
Return state dimension (x)
bool isLinear() const override
Check if the system dynamics are linear.
void dynamics(const Eigen::Ref< const StateVector > &x, const Eigen::Ref< const ControlVector > &u, Eigen::Ref< StateVector > f) const override
Evaluate the system dynamics equation.
int getStateDimension() const override
Return state dimension (x)
void dynamics(const Eigen::Ref< const StateVector > &x, const Eigen::Ref< const ControlVector > &u, Eigen::Ref< StateVector > f) const override
Evaluate the system dynamics equation.
void dynamics(const Eigen::Ref< const StateVector > &x, const Eigen::Ref< const ControlVector > &u, Eigen::Ref< StateVector > f) const override
Evaluate the system dynamics equation.
int getStateDimension() const override
Return state dimension (x)
bool isContinuousTime() const override
Check if the system dynamics are defined in continuous-time.
bool isLinear() const override
Check if the system dynamics are linear.
int getInputDimension() const override
Return the plant input dimension (u)
Ptr getInstance() const override
Return a newly created shared instance of the implemented class.
CartPole()
Default constructor.
Ptr getInstance() const override
Return a newly created shared instance of the implemented class.
bool isLinear() const override
Check if the system dynamics are linear.
void dynamics(const Eigen::Ref< const StateVector > &x, const Eigen::Ref< const ControlVector > &u, Eigen::Ref< StateVector > f) const override
Evaluate the system dynamics equation.
void dynamics(const Eigen::Ref< const StateVector > &x, const Eigen::Ref< const ControlVector > &u, Eigen::Ref< StateVector > f) const override
Evaluate the system dynamics equation.
void setDampingCoefficient(double a)
int getInputDimension() const override
Return the plant input dimension (u)
A matrix or vector expression mapping an existing expression.
int getInputDimension() const override
Return the plant input dimension (u)
ArtsteinsCircle()
Default constructor.
Interface class for system dynamic equations.
int getInputDimension() const override
Return the plant input dimension (u)
bool isLinear() const override
Check if the system dynamics are linear.
void dynamics(const Eigen::Ref< const StateVector > &x, const Eigen::Ref< const ControlVector > &u, Eigen::Ref< StateVector > f) const override
Evaluate the system dynamics equation.
const double & getDampingCoefficient() const
void dynamics(const Eigen::Ref< const StateVector > &x, const Eigen::Ref< const ControlVector > &u, Eigen::Ref< StateVector > f) const override
Evaluate the system dynamics equation.
DuffingOscillator()
Default constructor.
bool isLinear() const override
Check if the system dynamics are linear.
Ptr getInstance() const override
Return a newly created shared instance of the implemented class.
bool isContinuousTime() const override
Check if the system dynamics are defined in continuous-time.
int getInputDimension() const override
Return the plant input dimension (u)
Ptr getInstance() const override
Return a newly created shared instance of the implemented class.
bool isContinuousTime() const override
Check if the system dynamics are defined in continuous-time.
bool isContinuousTime() const override
Check if the system dynamics are defined in continuous-time.
int getStateDimension() const override
Return state dimension (x)
Ptr getInstance() const override
Return a newly created shared instance of the implemented class.
std::shared_ptr< SystemDynamicsInterface > Ptr
#define FACTORY_REGISTER_SYSTEM_DYNAMICS(type)
MasslessPendulum()
Default constructor.
SimplePendulum()
Default constructor.
EIGEN_DEVICE_FUNC const SinReturnType sin() const