Crocoddyl
action.hpp
1 // BSD 3-Clause License
3 //
4 // Copyright (C) 2019-2024, LAAS-CNRS, New York University,
5 // Max Planck Gesellschaft, University of Edinburgh,
6 // Heriot-Watt University
7 // Copyright note valid unless otherwise stated in individual files.
8 // All rights reserved.
10 
11 #ifndef CROCODDYL_CORE_NUMDIFF_ACTION_HPP_
12 #define CROCODDYL_CORE_NUMDIFF_ACTION_HPP_
13 
14 #include <vector>
15 
16 #include "crocoddyl/core/action-base.hpp"
17 #include "crocoddyl/core/fwd.hpp"
18 
19 namespace crocoddyl {
20 
45 template <typename _Scalar>
47  public:
48  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
49 
50  typedef _Scalar Scalar;
55  typedef typename MathBaseTpl<Scalar>::VectorXs VectorXs;
56  typedef typename MathBaseTpl<Scalar>::MatrixXs MatrixXs;
57 
66  explicit ActionModelNumDiffTpl(boost::shared_ptr<Base> model,
67  bool with_gauss_approx = false);
68  virtual ~ActionModelNumDiffTpl();
69 
73  virtual void calc(const boost::shared_ptr<ActionDataAbstract>& data,
74  const Eigen::Ref<const VectorXs>& x,
75  const Eigen::Ref<const VectorXs>& u);
76 
81  virtual void calc(const boost::shared_ptr<ActionDataAbstract>& data,
82  const Eigen::Ref<const VectorXs>& x);
83 
87  virtual void calcDiff(const boost::shared_ptr<ActionDataAbstract>& data,
88  const Eigen::Ref<const VectorXs>& x,
89  const Eigen::Ref<const VectorXs>& u);
90 
95  virtual void calcDiff(const boost::shared_ptr<ActionDataAbstract>& data,
96  const Eigen::Ref<const VectorXs>& x);
97 
101  virtual boost::shared_ptr<ActionDataAbstract> createData();
102 
106  virtual void quasiStatic(const boost::shared_ptr<ActionDataAbstract>& data,
107  Eigen::Ref<VectorXs> u,
108  const Eigen::Ref<const VectorXs>& x,
109  const std::size_t maxiter = 100,
110  const Scalar tol = Scalar(1e-9));
111 
115  const boost::shared_ptr<Base>& get_model() const;
116 
121  const Scalar get_disturbance() const;
122 
127  void set_disturbance(const Scalar disturbance);
128 
133 
139  virtual void print(std::ostream& os) const;
140 
141  protected:
144  using Base::nr_;
145  using Base::nu_;
146  using Base::state_;
147  using Base::u_lb_;
148  using Base::u_ub_;
149 
150  private:
162  void assertStableStateFD(const Eigen::Ref<const VectorXs>& x);
163 
164  boost::shared_ptr<Base> model_;
166  Scalar e_jac_;
168  Scalar e_hess_;
170  bool with_gauss_approx_;
172 };
173 
174 template <typename _Scalar>
175 struct ActionDataNumDiffTpl : public ActionDataAbstractTpl<_Scalar> {
176  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
177 
178  typedef _Scalar Scalar;
181  typedef typename MathBaseTpl<Scalar>::VectorXs VectorXs;
182  typedef typename MathBaseTpl<Scalar>::MatrixXs MatrixXs;
183 
190  template <template <typename Scalar> class Model>
191  explicit ActionDataNumDiffTpl(Model<Scalar>* const model)
192  : Base(model),
193  Rx(model->get_model()->get_nr(),
194  model->get_model()->get_state()->get_ndx()),
195  Ru(model->get_model()->get_nr(), model->get_model()->get_nu()),
196  dx(model->get_model()->get_state()->get_ndx()),
197  du(model->get_model()->get_nu()),
198  xp(model->get_model()->get_state()->get_nx()) {
199  Rx.setZero();
200  Ru.setZero();
201  dx.setZero();
202  du.setZero();
203  xp.setZero();
204 
205  const std::size_t ndx = model->get_model()->get_state()->get_ndx();
206  const std::size_t nu = model->get_model()->get_nu();
207  data_0 = model->get_model()->createData();
208  for (std::size_t i = 0; i < ndx; ++i) {
209  data_x.push_back(model->get_model()->createData());
210  }
211  for (std::size_t i = 0; i < nu; ++i) {
212  data_u.push_back(model->get_model()->createData());
213  }
214  }
215 
216  using Base::cost;
217  using Base::Fu;
218  using Base::Fx;
219  using Base::Lu;
220  using Base::Luu;
221  using Base::Lx;
222  using Base::Lxu;
223  using Base::Lxx;
224  using Base::r;
225  using Base::xnext;
226 
227  Scalar x_norm;
228  Scalar
230  Scalar
232  Scalar xh_hess;
234  Scalar uh_hess;
236  Scalar xh_hess_pow2;
237  Scalar uh_hess_pow2;
238  Scalar xuh_hess_pow2;
239  MatrixXs Rx;
240  MatrixXs Ru;
241  VectorXs dx;
242  VectorXs du;
243  VectorXs xp;
245  boost::shared_ptr<Base> data_0;
246  std::vector<boost::shared_ptr<Base> >
248  std::vector<boost::shared_ptr<Base> >
250 };
251 
252 } // namespace crocoddyl
253 
254 /* --- Details -------------------------------------------------------------- */
255 /* --- Details -------------------------------------------------------------- */
256 /* --- Details -------------------------------------------------------------- */
257 #include "crocoddyl/core/numdiff/action.hxx"
258 
259 #endif // CROCODDYL_CORE_NUMDIFF_ACTION_HPP_
Abstract class for action model.
Definition: action-base.hpp:95
VectorXs u_lb_
Lower control limits.
VectorXs u_ub_
Upper control limits.
boost::shared_ptr< StateAbstract > state_
Model of the state.
std::size_t nu_
Control dimension.
std::size_t nr_
Dimension of the cost residual.
This class computes the numerical differentiation of an action model.
Definition: action.hpp:46
const Scalar get_disturbance() const
Return the disturbance constant used in the numerical differentiation routine.
bool get_with_gauss_approx()
Identify if the Gauss approximation is going to be used or not.
void set_disturbance(const Scalar disturbance)
Modify the disturbance constant used in the numerical differentiation routine.
virtual void print(std::ostream &os) const
Print relevant information of the diff-action numdiff model.
ActionModelNumDiffTpl(boost::shared_ptr< Base > model, bool with_gauss_approx=false)
Initialize the numdiff action model.
virtual void calc(const boost::shared_ptr< ActionDataAbstract > &data, const Eigen::Ref< const VectorXs > &x, const Eigen::Ref< const VectorXs > &u)
Compute the next state and cost value.
virtual void calcDiff(const boost::shared_ptr< ActionDataAbstract > &data, const Eigen::Ref< const VectorXs > &x, const Eigen::Ref< const VectorXs > &u)
Compute the derivatives of the dynamics and cost functions.
const boost::shared_ptr< Base > & get_model() const
Return the acton model that we use to numerical differentiate.
virtual void calc(const boost::shared_ptr< ActionDataAbstract > &data, const Eigen::Ref< const VectorXs > &x)
virtual void calcDiff(const boost::shared_ptr< ActionDataAbstract > &data, const Eigen::Ref< const VectorXs > &x)
virtual boost::shared_ptr< ActionDataAbstract > createData()
Create the action data.
virtual void quasiStatic(const boost::shared_ptr< ActionDataAbstract > &data, Eigen::Ref< VectorXs > u, const Eigen::Ref< const VectorXs > &x, const std::size_t maxiter=100, const Scalar tol=Scalar(1e-9))
Computes the quasic static commands.
VectorXs xnext
evolution state
MatrixXs Fx
Jacobian of the dynamics w.r.t. the state .
MatrixXs Fu
Jacobian of the dynamics w.r.t. the control .
MatrixXs Luu
Hessian of the cost w.r.t. the control .
VectorXs Lx
Jacobian of the cost w.r.t. the state .
MatrixXs Lxx
Hessian of the cost w.r.t. the state .
VectorXs Lu
Jacobian of the cost w.r.t. the control .
VectorXs r
Cost residual.
std::vector< boost::shared_ptr< Base > > data_u
The temporary data associated with the control variation.
Definition: action.hpp:249
Scalar x_norm
Norm of the state vector.
Definition: action.hpp:227
MatrixXs Ru
Cost residual jacobian: .
Definition: action.hpp:240
MatrixXs Rx
Cost residual jacobian: .
Definition: action.hpp:239
boost::shared_ptr< Base > data_0
The data that contains the final results.
Definition: action.hpp:245
Scalar uh_jac
Disturbance value used for computing .
Definition: action.hpp:231
Scalar xh_jac
Disturbance value used for computing .
Definition: action.hpp:229
VectorXs du
Control disturbance.
Definition: action.hpp:242
std::vector< boost::shared_ptr< Base > > data_x
The temporary data associated with the state variation.
Definition: action.hpp:247
ActionDataNumDiffTpl(Model< Scalar > *const model)
Initialize the numdiff action data.
Definition: action.hpp:191
VectorXs dx
State disturbance.
Definition: action.hpp:241