Crocoddyl
control-base.hpp
1 // BSD 3-Clause License
3 //
4 // Copyright (C) 2021-2025, University of Edinburgh, University of Trento,
5 // Heriot-Watt University
6 // Copyright note valid unless otherwise stated in individual files.
7 // All rights reserved.
9 
10 #ifndef CROCODDYL_CORE_CONTROL_BASE_HPP_
11 #define CROCODDYL_CORE_CONTROL_BASE_HPP_
12 
13 #include <boost/shared_ptr.hpp>
14 
15 #include "crocoddyl/core/fwd.hpp"
16 #include "crocoddyl/core/mathbase.hpp"
17 #include "crocoddyl/core/utils/exception.hpp"
18 
19 namespace crocoddyl {
20 
43 template <typename _Scalar>
45  public:
46  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
47 
48  typedef _Scalar Scalar;
52  typedef typename MathBase::VectorXs VectorXs;
53  typedef typename MathBase::MatrixXs MatrixXs;
54 
62  const std::size_t nu);
64 
72  virtual void calc(
73  const boost::shared_ptr<ControlParametrizationDataAbstract>& data,
74  const Scalar t, const Eigen::Ref<const VectorXs>& u) const = 0;
75 
86  virtual void calcDiff(
87  const boost::shared_ptr<ControlParametrizationDataAbstract>& data,
88  const Scalar t, const Eigen::Ref<const VectorXs>& u) const = 0;
89 
95  virtual boost::shared_ptr<ControlParametrizationDataAbstract> createData();
96 
105  virtual void params(
106  const boost::shared_ptr<ControlParametrizationDataAbstract>& data,
107  const Scalar t, const Eigen::Ref<const VectorXs>& w) const = 0;
108 
118  virtual void convertBounds(const Eigen::Ref<const VectorXs>& w_lb,
119  const Eigen::Ref<const VectorXs>& w_ub,
120  Eigen::Ref<VectorXs> u_lb,
121  Eigen::Ref<VectorXs> u_ub) const = 0;
122 
136  virtual void multiplyByJacobian(
137  const boost::shared_ptr<ControlParametrizationDataAbstract>& data,
138  const Eigen::Ref<const MatrixXs>& A, Eigen::Ref<MatrixXs> out,
139  const AssignmentOp op = setto) const = 0;
140 
141  virtual MatrixXs multiplyByJacobian_J(
142  const boost::shared_ptr<ControlParametrizationDataAbstract>& data,
143  const Eigen::Ref<const MatrixXs>& A, const AssignmentOp op = setto) const;
144 
160  const boost::shared_ptr<ControlParametrizationDataAbstract>& data,
161  const Eigen::Ref<const MatrixXs>& A, Eigen::Ref<MatrixXs> out,
162  const AssignmentOp op = setto) const = 0;
163 
164  virtual MatrixXs multiplyJacobianTransposeBy_J(
165  const boost::shared_ptr<ControlParametrizationDataAbstract>& data,
166  const Eigen::Ref<const MatrixXs>& A, const AssignmentOp op = setto) const;
167 
171  virtual bool checkData(
172  const boost::shared_ptr<ControlParametrizationDataAbstract>& data);
173 
177  std::size_t get_nw() const;
178 
182  std::size_t get_nu() const;
183 
184  protected:
185  std::size_t nw_;
186  std::size_t nu_;
187 };
188 
189 template <typename _Scalar>
191  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
192 
193  typedef _Scalar Scalar;
195  typedef typename MathBase::VectorXs VectorXs;
196  typedef typename MathBase::MatrixXs MatrixXs;
197 
198  template <template <typename Scalar> class Model>
199  explicit ControlParametrizationDataAbstractTpl(Model<Scalar>* const model)
200  : w(model->get_nw()),
201  u(model->get_nu()),
202  dw_du(model->get_nw(), model->get_nu()) {
203  w.setZero();
204  u.setZero();
205  dw_du.setZero();
206  }
208 
209  VectorXs w;
210  VectorXs u;
211  MatrixXs dw_du;
213 };
214 
215 } // namespace crocoddyl
216 
217 /* --- Details -------------------------------------------------------------- */
218 /* --- Details -------------------------------------------------------------- */
219 /* --- Details -------------------------------------------------------------- */
220 #include "crocoddyl/core/control-base.hxx"
221 
222 #endif // CROCODDYL_CORE_CONTROL_BASE_HPP_
Abstract class for the control trajectory parametrization.
virtual void params(const boost::shared_ptr< ControlParametrizationDataAbstract > &data, const Scalar t, const Eigen::Ref< const VectorXs > &w) const =0
Update the control parameters u for a specified time t given the control input w.
ControlParametrizationModelAbstractTpl(const std::size_t nw, const std::size_t nu)
Initialize the control dimensions.
virtual void multiplyJacobianTransposeBy(const boost::shared_ptr< ControlParametrizationDataAbstract > &data, const Eigen::Ref< const MatrixXs > &A, Eigen::Ref< MatrixXs > out, const AssignmentOp op=setto) const =0
Compute the product between the transpose of the derivative of the control input with respect to the ...
virtual bool checkData(const boost::shared_ptr< ControlParametrizationDataAbstract > &data)
Checks that a specific data belongs to this model.
virtual void calc(const boost::shared_ptr< ControlParametrizationDataAbstract > &data, const Scalar t, const Eigen::Ref< const VectorXs > &u) const =0
Get the value of the control at the specified time.
virtual void calcDiff(const boost::shared_ptr< ControlParametrizationDataAbstract > &data, const Scalar t, const Eigen::Ref< const VectorXs > &u) const =0
Get the value of the Jacobian of the control with respect to the parameters.
std::size_t nu_
Control parameters dimension.
std::size_t get_nw() const
Return the dimension of the control inputs.
virtual void convertBounds(const Eigen::Ref< const VectorXs > &w_lb, const Eigen::Ref< const VectorXs > &w_ub, Eigen::Ref< VectorXs > u_lb, Eigen::Ref< VectorXs > u_ub) const =0
Convert the bounds on the control inputs w to bounds on the control parameters u.
virtual boost::shared_ptr< ControlParametrizationDataAbstract > createData()
Create the control-parametrization data.
std::size_t get_nu() const
Return the dimension of control parameters.
virtual void multiplyByJacobian(const boost::shared_ptr< ControlParametrizationDataAbstract > &data, const Eigen::Ref< const MatrixXs > &A, Eigen::Ref< MatrixXs > out, const AssignmentOp op=setto) const =0
Compute the product between the given matrix A and the derivative of the control input with respect t...
VectorXs w
value of the differential control
VectorXs u
value of the control parameters