Crocoddyl
control-base.hpp
1 // BSD 3-Clause License
3 //
4 // Copyright (C) 2021, University of Edinburgh, University of Trento
5 // Copyright note valid unless otherwise stated in individual files.
6 // All rights reserved.
8 
9 #ifndef CROCODDYL_CORE_CONTROL_BASE_HPP_
10 #define CROCODDYL_CORE_CONTROL_BASE_HPP_
11 
12 #include <boost/shared_ptr.hpp>
13 
14 #include "crocoddyl/core/fwd.hpp"
15 #include "crocoddyl/core/mathbase.hpp"
16 #include "crocoddyl/core/utils/exception.hpp"
17 
18 namespace crocoddyl {
19 
42 template <typename _Scalar>
44  public:
45  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
46 
47  typedef _Scalar Scalar;
51  typedef typename MathBase::VectorXs VectorXs;
52  typedef typename MathBase::MatrixXs MatrixXs;
53 
61  const std::size_t nu);
63 
71  virtual void calc(
72  const boost::shared_ptr<ControlParametrizationDataAbstract>& data,
73  const Scalar t, const Eigen::Ref<const VectorXs>& u) const = 0;
74 
85  virtual void calcDiff(
86  const boost::shared_ptr<ControlParametrizationDataAbstract>& data,
87  const Scalar t, const Eigen::Ref<const VectorXs>& u) const = 0;
88 
94  virtual boost::shared_ptr<ControlParametrizationDataAbstract> createData();
95 
104  virtual void params(
105  const boost::shared_ptr<ControlParametrizationDataAbstract>& data,
106  const Scalar t, const Eigen::Ref<const VectorXs>& w) const = 0;
107 
117  virtual void convertBounds(const Eigen::Ref<const VectorXs>& w_lb,
118  const Eigen::Ref<const VectorXs>& w_ub,
119  Eigen::Ref<VectorXs> u_lb,
120  Eigen::Ref<VectorXs> u_ub) const = 0;
121 
135  virtual void multiplyByJacobian(
136  const boost::shared_ptr<ControlParametrizationDataAbstract>& data,
137  const Eigen::Ref<const MatrixXs>& A, Eigen::Ref<MatrixXs> out,
138  const AssignmentOp = setto) const = 0;
139 
140  virtual MatrixXs multiplyByJacobian_J(
141  const boost::shared_ptr<ControlParametrizationDataAbstract>& data,
142  const Eigen::Ref<const MatrixXs>& A, const AssignmentOp = setto) const;
143 
159  const boost::shared_ptr<ControlParametrizationDataAbstract>& data,
160  const Eigen::Ref<const MatrixXs>& A, Eigen::Ref<MatrixXs> out,
161  const AssignmentOp = setto) const = 0;
162 
163  virtual MatrixXs multiplyJacobianTransposeBy_J(
164  const boost::shared_ptr<ControlParametrizationDataAbstract>& data,
165  const Eigen::Ref<const MatrixXs>& A, const AssignmentOp = setto) const;
166 
170  virtual bool checkData(
171  const boost::shared_ptr<ControlParametrizationDataAbstract>& data);
172 
176  std::size_t get_nw() const;
177 
181  std::size_t get_nu() const;
182 
183  protected:
184  std::size_t nw_;
185  std::size_t nu_;
186 };
187 
188 template <typename _Scalar>
190  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
191 
192  typedef _Scalar Scalar;
194  typedef typename MathBase::VectorXs VectorXs;
195  typedef typename MathBase::MatrixXs MatrixXs;
196 
197  template <template <typename Scalar> class Model>
198  explicit ControlParametrizationDataAbstractTpl(Model<Scalar>* const model)
199  : w(model->get_nw()),
200  u(model->get_nu()),
201  dw_du(model->get_nw(), model->get_nu()) {
202  w.setZero();
203  u.setZero();
204  dw_du.setZero();
205  }
207 
208  VectorXs w;
209  VectorXs u;
210  MatrixXs dw_du;
212 };
213 
214 } // namespace crocoddyl
215 
216 /* --- Details -------------------------------------------------------------- */
217 /* --- Details -------------------------------------------------------------- */
218 /* --- Details -------------------------------------------------------------- */
219 #include "crocoddyl/core/control-base.hxx"
220 
221 #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.
virtual void multiplyJacobianTransposeBy(const boost::shared_ptr< ControlParametrizationDataAbstract > &data, const Eigen::Ref< const MatrixXs > &A, Eigen::Ref< MatrixXs > out, const AssignmentOp=setto) const =0
Compute the product between the transpose of the derivative of the control input with respect to the ...
ControlParametrizationModelAbstractTpl(const std::size_t nw, const std::size_t nu)
Initialize the control dimensions.
virtual bool checkData(const boost::shared_ptr< ControlParametrizationDataAbstract > &data)
Checks that a specific data belongs to this model.
virtual void multiplyByJacobian(const boost::shared_ptr< ControlParametrizationDataAbstract > &data, const Eigen::Ref< const MatrixXs > &A, Eigen::Ref< MatrixXs > out, const AssignmentOp=setto) const =0
Compute the product between the given matrix A and the derivative of the control input with respect t...
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.
VectorXs w
value of the differential control
VectorXs u
value of the control parameters