Crocoddyl
state-base.hpp
1 // BSD 3-Clause License
3 //
4 // Copyright (C) 2019-2025, LAAS-CNRS, University of Edinburgh, INRIA,
5 // Heriot-Watt University
6 // Copyright note valid unless otherwise stated in individual files.
7 // All rights reserved.
9 
10 #ifndef CROCODDYL_CORE_STATE_BASE_HPP_
11 #define CROCODDYL_CORE_STATE_BASE_HPP_
12 
13 #include "crocoddyl/core/fwd.hpp"
14 
15 namespace crocoddyl {
16 
17 enum Jcomponent { both = 0, first = 1, second = 2 };
18 
19 inline bool is_a_Jcomponent(Jcomponent firstsecond) {
20  return (firstsecond == first || firstsecond == second || firstsecond == both);
21 }
22 
23 class StateBase {
24  public:
25  virtual ~StateBase() = default;
26 
27  CROCODDYL_BASE_CAST(StateBase, StateAbstractTpl)
28 };
29 
47 template <typename _Scalar>
48 class StateAbstractTpl : public StateBase {
49  public:
50  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
51 
52  typedef _Scalar Scalar;
54  typedef typename MathBase::VectorXs VectorXs;
55  typedef typename MathBase::MatrixXs MatrixXs;
56 
63  StateAbstractTpl(const std::size_t nx, const std::size_t ndx);
65  virtual ~StateAbstractTpl();
66 
70  virtual VectorXs zero() const = 0;
71 
75  virtual VectorXs rand() const = 0;
76 
94  virtual void diff(const Eigen::Ref<const VectorXs>& x0,
95  const Eigen::Ref<const VectorXs>& x1,
96  Eigen::Ref<VectorXs> dxout) const = 0;
97 
114  virtual void integrate(const Eigen::Ref<const VectorXs>& x,
115  const Eigen::Ref<const VectorXs>& dx,
116  Eigen::Ref<VectorXs> xout) const = 0;
117 
159  virtual void Jdiff(const Eigen::Ref<const VectorXs>& x0,
160  const Eigen::Ref<const VectorXs>& x1,
161  Eigen::Ref<MatrixXs> Jfirst, Eigen::Ref<MatrixXs> Jsecond,
162  const Jcomponent firstsecond = both) const = 0;
163 
203  virtual void Jintegrate(const Eigen::Ref<const VectorXs>& x,
204  const Eigen::Ref<const VectorXs>& dx,
205  Eigen::Ref<MatrixXs> Jfirst,
206  Eigen::Ref<MatrixXs> Jsecond,
207  const Jcomponent firstsecond = both,
208  const AssignmentOp op = setto) const = 0;
209 
224  virtual void JintegrateTransport(const Eigen::Ref<const VectorXs>& x,
225  const Eigen::Ref<const VectorXs>& dx,
226  Eigen::Ref<MatrixXs> Jin,
227  const Jcomponent firstsecond) const = 0;
228 
237  VectorXs diff_dx(const Eigen::Ref<const VectorXs>& x0,
238  const Eigen::Ref<const VectorXs>& x1);
239 
247  VectorXs integrate_x(const Eigen::Ref<const VectorXs>& x,
248  const Eigen::Ref<const VectorXs>& dx);
249 
257  std::vector<MatrixXs> Jdiff_Js(const Eigen::Ref<const VectorXs>& x0,
258  const Eigen::Ref<const VectorXs>& x1,
259  const Jcomponent firstsecond = both);
260 
268  std::vector<MatrixXs> Jintegrate_Js(const Eigen::Ref<const VectorXs>& x,
269  const Eigen::Ref<const VectorXs>& dx,
270  const Jcomponent firstsecond = both);
271 
275  template <class Scalar>
276  friend std::ostream& operator<<(std::ostream& os,
277  const ActionModelAbstractTpl<Scalar>& model);
278 
284  virtual void print(std::ostream& os) const;
285 
289  std::size_t get_nx() const;
290 
294  std::size_t get_ndx() const;
295 
299  std::size_t get_nq() const;
300 
304  std::size_t get_nv() const;
305 
309  const VectorXs& get_lb() const;
310 
314  const VectorXs& get_ub() const;
315 
319  bool get_has_limits() const;
320 
324  void set_lb(const VectorXs& lb);
325 
329  void set_ub(const VectorXs& ub);
330 
331  protected:
332  void update_has_limits();
333 
334  std::size_t nx_;
335  std::size_t ndx_;
336  std::size_t nq_;
337  std::size_t nv_;
338  VectorXs lb_;
339  VectorXs ub_;
340  bool has_limits_;
341 };
342 
343 } // namespace crocoddyl
344 
345 /* --- Details -------------------------------------------------------------- */
346 /* --- Details -------------------------------------------------------------- */
347 /* --- Details -------------------------------------------------------------- */
348 #include "crocoddyl/core/state-base.hxx"
349 
350 CROCODDYL_DECLARE_EXTERN_TEMPLATE_CLASS(crocoddyl::StateAbstractTpl)
351 
352 #endif // CROCODDYL_CORE_STATE_BASE_HPP_
Abstract class for action model.
Definition: action-base.hpp:97
Abstract class for the state representation.
Definition: state-base.hpp:48
virtual void Jdiff(const Eigen::Ref< const VectorXs > &x0, const Eigen::Ref< const VectorXs > &x1, Eigen::Ref< MatrixXs > Jfirst, Eigen::Ref< MatrixXs > Jsecond, const Jcomponent firstsecond=both) const =0
Compute the Jacobian of the state manifold differentiation.
std::vector< MatrixXs > Jdiff_Js(const Eigen::Ref< const VectorXs > &x0, const Eigen::Ref< const VectorXs > &x1, const Jcomponent firstsecond=both)
std::vector< MatrixXs > Jintegrate_Js(const Eigen::Ref< const VectorXs > &x, const Eigen::Ref< const VectorXs > &dx, const Jcomponent firstsecond=both)
Compute the Jacobian of the state manifold integration.
VectorXs diff_dx(const Eigen::Ref< const VectorXs > &x0, const Eigen::Ref< const VectorXs > &x1)
Compute the state manifold differentiation.
virtual VectorXs rand() const =0
Generate a random state.
std::size_t nv_
Velocity dimension.
Definition: state-base.hpp:337
virtual void JintegrateTransport(const Eigen::Ref< const VectorXs > &x, const Eigen::Ref< const VectorXs > &dx, Eigen::Ref< MatrixXs > Jin, const Jcomponent firstsecond) const =0
Parallel transport from integrate(x, dx) to x.
std::size_t get_nq() const
Return the dimension of the configuration tuple.
virtual void print(std::ostream &os) const
Print relevant information of the state model.
std::size_t get_ndx() const
Return the dimension of the tangent space of the state manifold.
std::size_t get_nv() const
Return the dimension of tangent space of the configuration manifold.
const VectorXs & get_lb() const
Return the state lower bound.
std::size_t nx_
State dimension.
Definition: state-base.hpp:334
VectorXs integrate_x(const Eigen::Ref< const VectorXs > &x, const Eigen::Ref< const VectorXs > &dx)
Compute the state manifold integration.
const VectorXs & get_ub() const
Return the state upper bound.
bool has_limits_
Indicates whether any of the state limits is finite.
Definition: state-base.hpp:340
virtual void diff(const Eigen::Ref< const VectorXs > &x0, const Eigen::Ref< const VectorXs > &x1, Eigen::Ref< VectorXs > dxout) const =0
Compute the state manifold differentiation.
virtual VectorXs zero() const =0
Generate a zero state.
StateAbstractTpl(const std::size_t nx, const std::size_t ndx)
Initialize the state dimensions.
std::size_t nq_
Configuration dimension.
Definition: state-base.hpp:336
bool get_has_limits() const
Indicate if the state has defined limits.
VectorXs lb_
Lower state limits.
Definition: state-base.hpp:338
void set_ub(const VectorXs &ub)
Modify the state upper bound.
virtual void integrate(const Eigen::Ref< const VectorXs > &x, const Eigen::Ref< const VectorXs > &dx, Eigen::Ref< VectorXs > xout) const =0
Compute the state manifold integration.
friend std::ostream & operator<<(std::ostream &os, const ActionModelAbstractTpl< Scalar > &model)
Print information on the state model.
VectorXs ub_
Upper state limits.
Definition: state-base.hpp:339
virtual void Jintegrate(const Eigen::Ref< const VectorXs > &x, const Eigen::Ref< const VectorXs > &dx, Eigen::Ref< MatrixXs > Jfirst, Eigen::Ref< MatrixXs > Jsecond, const Jcomponent firstsecond=both, const AssignmentOp op=setto) const =0
Compute the Jacobian of the state manifold integration.
std::size_t get_nx() const
Return the dimension of the state tuple.
void set_lb(const VectorXs &lb)
Modify the state lower bound.
std::size_t ndx_
State rate dimension.
Definition: state-base.hpp:335