pinocchio  2.4.4
A fast and flexible implementation of Rigid Body Dynamics algorithms and their analytical derivatives
compute-all-terms.hpp
1 //
2 // Copyright (c) 2015-2020 CNRS INRIA
3 //
4 
5 #ifndef __pinocchio_compute_all_terms_hpp__
6 #define __pinocchio_compute_all_terms_hpp__
7 
8 #include "pinocchio/multibody/visitor.hpp"
9 #include "pinocchio/multibody/model.hpp"
10 #include "pinocchio/multibody/data.hpp"
11 #include "pinocchio/spatial/act-on-set.hpp"
12 #include "pinocchio/algorithm/center-of-mass.hpp"
13 #include "pinocchio/algorithm/energy.hpp"
14 #include "pinocchio/algorithm/check.hpp"
15 
16 namespace pinocchio
17 {
39  template<typename Scalar, int Options, template<typename,int> class JointCollectionTpl, typename ConfigVectorType, typename TangentVectorType>
40  inline void computeAllTerms(const ModelTpl<Scalar,Options,JointCollectionTpl> & model,
41  DataTpl<Scalar,Options,JointCollectionTpl> & data,
42  const Eigen::MatrixBase<ConfigVectorType> & q,
43  const Eigen::MatrixBase<TangentVectorType> & v);
44 
45 } // namespace pinocchio
46 
47 
50 namespace pinocchio
51 {
52  template<typename Scalar, int Options, template<typename,int> class JointCollectionTpl, typename ConfigVectorType, typename TangentVectorType>
54  : public fusion::JointUnaryVisitorBase< CATForwardStep<Scalar,Options,JointCollectionTpl,ConfigVectorType,TangentVectorType> >
55  {
58 
59  typedef boost::fusion::vector<const Model &,
60  Data &,
61  const ConfigVectorType &,
62  const TangentVectorType &
63  > ArgsType;
64 
65  template<typename JointModel>
66  static void algo(const JointModelBase<JointModel> & jmodel,
68  const Model & model,
69  Data & data,
70  const Eigen::MatrixBase<ConfigVectorType> & q,
71  const Eigen::MatrixBase<TangentVectorType> & v)
72  {
73  typedef typename Model::JointIndex JointIndex;
74  typedef typename Data::Inertia Inertia;
75 
76  const JointIndex & i = jmodel.id();
77  const JointIndex & parent = model.parents[i];
78 
79  jmodel.calc(jdata.derived(),q.derived(),v.derived());
80 
81  // CRBA
82  data.liMi[i] = model.jointPlacements[i]*jdata.M();
83  data.Ycrb[i] = model.inertias[i];
84 
85  // Jacobian + NLE
86  data.v[i] = jdata.v();
87 
88  if(parent>0)
89  {
90  data.oMi[i] = data.oMi[parent]*data.liMi[i];
91  data.v[i] += data.liMi[i].actInv(data.v[parent]);
92  }
93  else
94  data.oMi[i] = data.liMi[i];
95 
96  jmodel.jointCols(data.J) = data.oMi[i].act(jdata.S());
97 
98  data.a_gf[i] = data.a[i] = jdata.c() + (data.v[i] ^ jdata.v());
99  if (parent > 0)
100  data.a[i] += data.liMi[i].actInv(data.a[parent]);
101 
102  data.a_gf[i] += data.liMi[i].actInv(data.a_gf[parent]);
103 
104  data.f[i] = model.inertias[i]*data.a_gf[i] + model.inertias[i].vxiv(data.v[i]); // -f_ext
105 
106  // CoM
107  const Scalar & mass = model.inertias[i].mass();
108  const typename Inertia::Vector3 & lever = model.inertias[i].lever();
109 
110  data.com[i].noalias() = mass * lever;
111  data.mass[i] = mass;
112 
113  data.vcom[i].noalias() = mass * (data.v[i].angular().cross(lever) + data.v[i].linear());
114  }
115 
116  };
117 
118  template<typename Scalar, int Options, template<typename,int> class JointCollectionTpl>
120  : public fusion::JointUnaryVisitorBase <CATBackwardStep<Scalar,Options,JointCollectionTpl> >
121  {
124 
125  typedef boost::fusion::vector<const Model &,
126  Data &
127  > ArgsType;
128 
129  template<typename JointModel>
130  static void algo(const JointModelBase<JointModel> & jmodel,
132  const Model & model,
133  Data & data)
134  {
135  /*
136  * F[1:6,i] = Y*S
137  * M[i,SUBTREE] = S'*F[1:6,SUBTREE]
138  * if li>0
139  * Yli += liXi Yi
140  * F[1:6,SUBTREE] = liXi F[1:6,SUBTREE]
141  */
142  typedef typename Model::JointIndex JointIndex;
143  typedef typename Data::SE3 SE3;
144 
145  const JointIndex & i = jmodel.id();
146  const JointIndex & parent = model.parents[i];
147  const SE3 & oMi = data.oMi[i];
148 
149  /* F[1:6,i] = Y*S */
150  jmodel.jointCols(data.Fcrb[i]) = data.Ycrb[i] * jdata.S();
151 
152  /* M[i,SUBTREE] = S'*F[1:6,SUBTREE] */
153  data.M.block(jmodel.idx_v(),jmodel.idx_v(),jmodel.nv(),data.nvSubtree[i])
154  = jdata.S().transpose()*data.Fcrb[i].middleCols(jmodel.idx_v(),data.nvSubtree[i]);
155 
156 
157  jmodel.jointVelocitySelector(data.nle) = jdata.S().transpose()*data.f[i];
158  if(parent>0)
159  {
160  /* Yli += liXi Yi */
161  data.Ycrb[parent] += data.liMi[i].act(data.Ycrb[i]);
162 
163  /* F[1:6,SUBTREE] = liXi F[1:6,SUBTREE] */
164  Eigen::Block<typename Data::Matrix6x> jF
165  = data.Fcrb[parent].block(0,jmodel.idx_v(),6,data.nvSubtree[i]);
166  Eigen::Block<typename Data::Matrix6x> iF
167  = data.Fcrb[i].block(0,jmodel.idx_v(),6,data.nvSubtree[i]);
168  forceSet::se3Action(data.liMi[i], iF, jF);
169 
170  data.f[parent] += data.liMi[i].act(data.f[i]);
171  }
172 
173  // CoM
174  const SE3 & liMi = data.liMi[i];
175 
176  data.com[parent] += (liMi.rotation()*data.com[i]
177  + data.mass[i] * liMi.translation());
178 
179  typename SE3::Vector3 com_in_world(oMi.rotation() * data.com[i] + data.mass[i] * oMi.translation());
180 
181  data.vcom[parent] += liMi.rotation()*data.vcom[i];
182  data.mass[parent] += data.mass[i];
183 
184  typedef typename Data::Matrix6x Matrix6x;
185  typedef typename SizeDepType<JointModel::NV>::template ColsReturn<Matrix6x>::Type ColBlock;
186 
187  ColBlock Jcols = jmodel.jointCols(data.J);
188 
189  if( JointModel::NV==1 )
190  data.Jcom.col(jmodel.idx_v())
191  = data.mass[i] * Jcols.template topLeftCorner<3,1>()
192  - com_in_world.cross(Jcols.template bottomLeftCorner<3,1>()) ;
193  else
194  jmodel.jointCols(data.Jcom)
195  = data.mass[i] * Jcols.template topRows<3>()
196  - skew(com_in_world) * Jcols.template bottomRows<3>();
197 
198  data.com[i] /= data.mass[i];
199  data.vcom[i] /= data.mass[i];
200  }
201  };
202 
203  template<typename Scalar, int Options, template<typename,int> class JointCollectionTpl, typename ConfigVectorType, typename TangentVectorType>
206  const Eigen::MatrixBase<ConfigVectorType> & q,
207  const Eigen::MatrixBase<TangentVectorType> & v)
208  {
209  assert(model.check(data) && "data is not consistent with model.");
210  PINOCCHIO_CHECK_INPUT_ARGUMENT(q.size() == model.nq, "The configuration vector is not of right size");
211  PINOCCHIO_CHECK_INPUT_ARGUMENT(v.size() == model.nv, "The velocity vector is not of right size");
212 
214  typedef typename Model::JointIndex JointIndex;
215 
216  data.v[0].setZero();
217  data.a[0].setZero();
218  data.a_gf[0] = -model.gravity;
219 
220  data.mass[0] = 0;
221  data.com[0].setZero();
222  data.vcom[0].setZero();
223 
225  for(JointIndex i=1;i<(JointIndex) model.njoints;++i)
226  {
227  Pass1::run(model.joints[i],data.joints[i],
228  typename Pass1::ArgsType(model,data,q.derived(),v.derived()));
229  }
230 
232  for(JointIndex i=(JointIndex)(model.njoints-1);i>0;--i)
233  {
234  Pass2::run(model.joints[i],data.joints[i],
235  typename Pass2::ArgsType(model,data));
236  }
237 
238  // CoM
239  data.com[0] /= data.mass[0];
240  data.vcom[0] /= data.mass[0];
241 
242  // JCoM
243  data.Jcom /= data.mass[0];
244 
245  // Energy
246  computeKineticEnergy(model, data);
247  computePotentialEnergy(model, data);
248 
249  }
250 } // namespace pinocchio
252 
253 
254 #endif // ifndef __pinocchio_compute_all_terms_hpp__
255 
Scalar computePotentialEnergy(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data)
Computes the potential energy of the system, i.e. the potential energy linked to the gravity field...
JointDataVector joints
Vector of pinocchio::JointData associated to the pinocchio::JointModel stored in model, encapsulated in JointDataAccessor.
Definition: data.hpp:87
JointModelVector joints
Model of joint i, encapsulated in a JointModelAccessor.
Definition: model.hpp:96
Matrix3x Jcom
Jacobian of center of mass.
Definition: data.hpp:318
void computeAllTerms(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType > &v)
Computes efficiently all the terms needed for dynamic simulation. It is equivalent to the call at the...
int njoints
Number of joints.
Definition: model.hpp:81
int nv
Dimension of the velocity vector space.
Definition: model.hpp:78
Eigen::Matrix< Scalar, 6, Eigen::Dynamic, Options > Matrix6x
The 6d jacobian type (temporary)
Definition: data.hpp:71
MatrixXs M
The joint space inertia matrix (a square matrix of dim model.nv).
Definition: data.hpp:151
Motion gravity
Spatial gravity of the model.
Definition: model.hpp:151
VectorXs nle
Vector of Non Linear Effects (dim model.nv). It corresponds to concatenation of the Coriolis...
Definition: data.hpp:133
static void se3Action(const SE3Tpl< Scalar, Options > &m, const Eigen::MatrixBase< Mat > &iF, Eigen::MatrixBase< MatRet > const &jF)
SE3 action on a set of forces, represented by a 6xN matrix whose each column represent a spatial forc...
Base structure for Unary visitation of a JointModel. This structure provides runners to call the righ...
std::vector< Scalar > mass
Vector of subtree mass. In other words, mass[j] is the mass of the subtree supported by joint ...
Definition: data.hpp:314
Matrix6x J
Jacobian of joint placements.
Definition: data.hpp:275
Main pinocchio namespace.
Definition: treeview.dox:24
Scalar computeKineticEnergy(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data)
Computes the kinetic energy of the system. The result is accessible through data.kinetic_energy.
std::vector< JointIndex > parents
Joint parent of joint i, denoted li (li==parents[i]).
Definition: model.hpp:111
std::vector< int > nvSubtree
Dimension of the subtree motion space (for CRBA)
Definition: data.hpp:243
void skew(const Eigen::MatrixBase< Vector3 > &v, const Eigen::MatrixBase< Matrix3 > &M)
Computes the skew representation of a given 3d vector, i.e. the antisymmetric matrix representation o...
Definition: skew.hpp:21
bool check(const AlgorithmCheckerBase< D > &checker=AlgorithmCheckerBase< D >()) const
Check the validity of the attributes of Model with respect to the specification of some algorithms...
Definition: model.hpp:525
int nq
Dimension of the configuration vector representation.
Definition: model.hpp:75