crocoddyl 1.9.0
Contact RObot COntrol by Differential DYnamic programming Library (Crocoddyl)
 
Loading...
Searching...
No Matches
control-base.hpp
1
2// 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
18namespace crocoddyl {
19
39template <typename _Scalar>
41 public:
42 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
43
44 typedef _Scalar Scalar;
47 typedef typename MathBase::VectorXs VectorXs;
48 typedef typename MathBase::MatrixXs MatrixXs;
49
56 ControlParametrizationModelAbstractTpl(const std::size_t nw, const std::size_t nu);
58
66 virtual void calc(const boost::shared_ptr<ControlParametrizationDataAbstract>& data, const Scalar t,
67 const Eigen::Ref<const VectorXs>& u) const = 0;
68
78 virtual void calcDiff(const boost::shared_ptr<ControlParametrizationDataAbstract>& data, const Scalar t,
79 const Eigen::Ref<const VectorXs>& u) const = 0;
80
86 virtual boost::shared_ptr<ControlParametrizationDataAbstract> createData();
87
95 virtual void params(const boost::shared_ptr<ControlParametrizationDataAbstract>& data, const Scalar t,
96 const Eigen::Ref<const VectorXs>& w) const = 0;
97
106 virtual void convertBounds(const Eigen::Ref<const VectorXs>& w_lb, const Eigen::Ref<const VectorXs>& w_ub,
107 Eigen::Ref<VectorXs> u_lb, Eigen::Ref<VectorXs> u_ub) const = 0;
108
120 virtual void multiplyByJacobian(const boost::shared_ptr<ControlParametrizationDataAbstract>& data,
121 const Eigen::Ref<const MatrixXs>& A, Eigen::Ref<MatrixXs> out,
122 const AssignmentOp = setto) const = 0;
123
124 virtual MatrixXs multiplyByJacobian_J(const boost::shared_ptr<ControlParametrizationDataAbstract>& data,
125 const Eigen::Ref<const MatrixXs>& A, const AssignmentOp = setto) const;
126
139 virtual void multiplyJacobianTransposeBy(const boost::shared_ptr<ControlParametrizationDataAbstract>& data,
140 const Eigen::Ref<const MatrixXs>& A, Eigen::Ref<MatrixXs> out,
141 const AssignmentOp = setto) const = 0;
142
143 virtual MatrixXs multiplyJacobianTransposeBy_J(const boost::shared_ptr<ControlParametrizationDataAbstract>& data,
144 const Eigen::Ref<const MatrixXs>& A,
145 const AssignmentOp = setto) const;
146
150 virtual bool checkData(const boost::shared_ptr<ControlParametrizationDataAbstract>& data);
151
155 std::size_t get_nw() const;
156
160 std::size_t get_nu() const;
161
162 protected:
163 std::size_t nw_;
164 std::size_t nu_;
165};
166
167template <typename _Scalar>
169 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
170
171 typedef _Scalar Scalar;
173 typedef typename MathBase::VectorXs VectorXs;
174 typedef typename MathBase::MatrixXs MatrixXs;
175
176 template <template <typename Scalar> class Model>
177 explicit ControlParametrizationDataAbstractTpl(Model<Scalar>* const model)
178 : w(model->get_nw()), u(model->get_nu()), dw_du(model->get_nw(), model->get_nu()) {
179 w.setZero();
180 u.setZero();
181 dw_du.setZero();
182 }
184
185 VectorXs w;
186 VectorXs u;
187 MatrixXs dw_du;
188};
189
190} // namespace crocoddyl
191
192/* --- Details -------------------------------------------------------------- */
193/* --- Details -------------------------------------------------------------- */
194/* --- Details -------------------------------------------------------------- */
195#include "crocoddyl/core/control-base.hxx"
196
197#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 boost::shared_ptr< ControlParametrizationDataAbstract > createData()
Create the control-parametrization data.
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.
std::size_t get_nu() const
Return the dimension of control parameters.
VectorXs w
value of the differential control
MatrixXs dw_du
Jacobian of the differential control with respect to the parameters.
VectorXs u
value of the control parameters