pinocchio  UNKNOWN
rnea.hxx
1 //
2 // Copyright (c) 2015-2017 CNRS
3 //
4 // This file is part of Pinocchio
5 // Pinocchio is free software: you can redistribute it
6 // and/or modify it under the terms of the GNU Lesser General Public
7 // License as published by the Free Software Foundation, either version
8 // 3 of the License, or (at your option) any later version.
9 //
10 // Pinocchio is distributed in the hope that it will be
11 // useful, but WITHOUT ANY WARRANTY; without even the implied warranty
12 // of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13 // General Lesser Public License for more details. You should have
14 // received a copy of the GNU Lesser General Public License along with
15 // Pinocchio If not, see
16 // <http://www.gnu.org/licenses/>.
17 
18 #ifndef __se3_rnea_hxx__
19 #define __se3_rnea_hxx__
20 
22 
23 #include "pinocchio/multibody/visitor.hpp"
24 #include "pinocchio/algorithm/check.hpp"
25 
26 namespace se3
27 {
28  struct RneaForwardStep : public fusion::JointVisitor<RneaForwardStep>
29  {
30  typedef boost::fusion::vector<const se3::Model &,
31  se3::Data &,
32  const Eigen::VectorXd &,
33  const Eigen::VectorXd &,
34  const Eigen::VectorXd &
35  > ArgsType;
36 
37  JOINT_VISITOR_INIT(RneaForwardStep);
38 
39  template<typename JointModel>
40  static void algo(const se3::JointModelBase<JointModel> & jmodel,
42  const se3::Model & model,
43  se3::Data & data,
44  const Eigen::VectorXd & q,
45  const Eigen::VectorXd & v,
46  const Eigen::VectorXd & a)
47  {
48  const Model::JointIndex & i = jmodel.id();
49  const Model::JointIndex & parent = model.parents[i];
50 
51  jmodel.calc(jdata.derived(),q,v);
52 
53  data.liMi[i] = model.jointPlacements[i]*jdata.M();
54 
55  data.v[i] = jdata.v();
56  if(parent>0) data.v[i] += data.liMi[i].actInv(data.v[parent]);
57 
58  data.a_gf[i] = jdata.S()*jmodel.jointVelocitySelector(a) + jdata.c() + (data.v[i] ^ jdata.v()) ;
59  data.a_gf[i] += data.liMi[i].actInv(data.a_gf[parent]);
60 
61  data.f[i] = model.inertias[i]*data.a_gf[i] + model.inertias[i].vxiv(data.v[i]); // -f_ext
62  }
63 
64  };
65 
66  struct RneaBackwardStep : public fusion::JointVisitor<RneaBackwardStep>
67  {
68  typedef boost::fusion::vector<const Model &,
69  Data &
70  > ArgsType;
71 
72  JOINT_VISITOR_INIT(RneaBackwardStep);
73 
74  template<typename JointModel>
75  static void algo(const JointModelBase<JointModel> & jmodel,
77  const Model & model,
78  Data & data)
79  {
80  const Model::JointIndex & i = jmodel.id();
81  const Model::JointIndex & parent = model.parents[i];
82 
83  jmodel.jointVelocitySelector(data.tau) = jdata.S().transpose()*data.f[i];
84  if(parent>0) data.f[parent] += data.liMi[i].act(data.f[i]);
85  }
86  };
87 
88  inline const Eigen::VectorXd&
89  rnea(const Model & model, Data& data,
90  const Eigen::VectorXd & q,
91  const Eigen::VectorXd & v,
92  const Eigen::VectorXd & a)
93  {
94  assert(model.check(data) && "data is not consistent with model.");
95 
96  data.v[0].setZero();
97  data.a_gf[0] = -model.gravity;
98 
99  for( Model::JointIndex i=1;i<(Model::JointIndex)model.njoints;++i )
100  {
101  RneaForwardStep::run(model.joints[i],data.joints[i],
102  RneaForwardStep::ArgsType(model,data,q,v,a));
103  }
104 
105  for( Model::JointIndex i=(Model::JointIndex)model.njoints-1;i>0;--i )
106  {
107  RneaBackwardStep::run(model.joints[i],data.joints[i],
108  RneaBackwardStep::ArgsType(model,data));
109  }
110 
111  return data.tau;
112  }
113 
114  inline const Eigen::VectorXd &
115  rnea(const Model & model, Data & data,
116  const Eigen::VectorXd & q,
117  const Eigen::VectorXd & v,
118  const Eigen::VectorXd & a,
120  {
121  assert(fext.size() == model.joints.size());
122  assert(model.check(data) && "data is not consistent with model.");
123 
124  data.v[0].setZero();
125  data.a_gf[0] = -model.gravity;
126 
127  for( Model::JointIndex i=1;i<(Model::JointIndex)model.njoints;++i )
128  {
129  RneaForwardStep::run(model.joints[i],data.joints[i],
130  RneaForwardStep::ArgsType(model,data,q,v,a));
131  data.f[i] -= fext[i];
132  }
133 
134  for( Model::JointIndex i=(Model::JointIndex)model.njoints-1;i>0;--i )
135  {
136  RneaBackwardStep::run(model.joints[i],data.joints[i],
137  RneaBackwardStep::ArgsType(model,data));
138  }
139 
140  return data.tau;
141  }
142 
143  struct NLEForwardStep : public fusion::JointVisitor<NLEForwardStep>
144  {
145  typedef boost::fusion::vector< const se3::Model &,
146  se3::Data &,
147  const Eigen::VectorXd &,
148  const Eigen::VectorXd &
149  > ArgsType;
150 
151  JOINT_VISITOR_INIT(NLEForwardStep);
152 
153  template<typename JointModel>
154  static void algo(const se3::JointModelBase<JointModel> & jmodel,
156  const se3::Model & model,
157  se3::Data & data,
158  const Eigen::VectorXd & q,
159  const Eigen::VectorXd & v)
160  {
161  const Model::JointIndex & i = jmodel.id();
162  const Model::JointIndex & parent = model.parents[i];
163 
164  jmodel.calc(jdata.derived(),q,v);
165 
166  data.liMi[i] = model.jointPlacements[i]*jdata.M();
167 
168  data.v[i] = jdata.v();
169  if(parent>0) data.v[i] += data.liMi[i].actInv(data.v[parent]);
170 
171  data.a_gf[i] = jdata.c() + (data.v[i] ^ jdata.v());
172  data.a_gf[i] += data.liMi[i].actInv(data.a_gf[parent]);
173 
174  data.f[i] = model.inertias[i]*data.a_gf[i] + model.inertias[i].vxiv(data.v[i]); // -f_ext
175  }
176 
177  };
178 
179  struct NLEBackwardStep : public fusion::JointVisitor<NLEBackwardStep>
180  {
181  typedef boost::fusion::vector<const Model &,
182  Data &
183  > ArgsType;
184 
185  JOINT_VISITOR_INIT(NLEBackwardStep);
186 
187  template<typename JointModel>
188  static void algo(const JointModelBase<JointModel> & jmodel,
190  const Model & model,
191  Data & data)
192  {
193  const Model::JointIndex & i = jmodel.id();
194  const Model::JointIndex & parent = model.parents[i];
195 
196  jmodel.jointVelocitySelector(data.nle) = jdata.S().transpose()*data.f[i];
197  if(parent>0) data.f[parent] += data.liMi[i].act(data.f[i]);
198  }
199  };
200 
201  inline const Eigen::VectorXd &
202  nonLinearEffects(const Model & model, Data & data,
203  const Eigen::VectorXd & q,
204  const Eigen::VectorXd & v)
205  {
206  assert(model.check(data) && "data is not consistent with model.");
207 
208  data.v[0].setZero ();
209  data.a_gf[0] = -model.gravity;
210 
211  for( size_t i=1;i<(size_t) model.njoints;++i )
212  {
213  NLEForwardStep::run(model.joints[i],data.joints[i],
214  NLEForwardStep::ArgsType(model,data,q,v));
215  }
216 
217  for( size_t i=(size_t) (model.njoints-1);i>0;--i )
218  {
219  NLEBackwardStep::run(model.joints[i],data.joints[i],
220  NLEBackwardStep::ArgsType(model,data));
221  }
222 
223  return data.nle;
224  }
225 
226  struct computeGeneralizedGravityForwardStep : public fusion::JointVisitor<computeGeneralizedGravityForwardStep>
227  {
228  typedef boost::fusion::vector< const se3::Model &,
229  se3::Data &,
230  const Eigen::VectorXd &
231  > ArgsType;
232 
233  JOINT_VISITOR_INIT(computeGeneralizedGravityForwardStep);
234 
235  template<typename JointModel>
236  static void algo(const se3::JointModelBase<JointModel> & jmodel,
238  const se3::Model & model,
239  se3::Data & data,
240  const Eigen::VectorXd & q)
241  {
242  const Model::JointIndex & i = jmodel.id();
243  const Model::JointIndex & parent = model.parents[i];
244 
245  jmodel.calc(jdata.derived(),q);
246 
247  data.liMi[i] = model.jointPlacements[i]*jdata.M();
248 
249  data.a_gf[i] = data.liMi[i].actInv(data.a_gf[(size_t) parent]);
250  data.f[i] = model.inertias[i]*data.a_gf[i];
251  }
252 
253  };
254 
255  struct computeGeneralizedGravityBackwardStep : public fusion::JointVisitor<computeGeneralizedGravityBackwardStep>
256  {
257  typedef boost::fusion::vector<const Model &,
258  Data &
259  > ArgsType;
260 
261  JOINT_VISITOR_INIT(computeGeneralizedGravityBackwardStep);
262 
263  template<typename JointModel>
264  static void algo(const JointModelBase<JointModel> & jmodel,
266  const Model & model,
267  Data & data)
268  {
269  const Model::JointIndex & i = jmodel.id();
270  const Model::JointIndex & parent = model.parents[i];
271 
272  jmodel.jointVelocitySelector(data.g) = jdata.S().transpose()*data.f[i];
273  if(parent>0) data.f[(size_t) parent] += data.liMi[i].act(data.f[i]);
274  }
275  };
276 
277  inline const Eigen::VectorXd &
278  computeGeneralizedGravity(const Model & model, Data & data,
279  const Eigen::VectorXd & q)
280  {
281  assert(model.check(data) && "data is not consistent with model.");
282 
283  data.a_gf[0] = -model.gravity;
284 
285  for(size_t i=1;i<(size_t) model.njoints;++i)
286  {
287  computeGeneralizedGravityForwardStep::run(model.joints[i],data.joints[i],
288  computeGeneralizedGravityForwardStep::ArgsType(model,data,q));
289  }
290 
291  for(size_t i=(size_t) (model.njoints-1);i>0;--i)
292  {
293  computeGeneralizedGravityBackwardStep::run(model.joints[i],data.joints[i],
294  computeGeneralizedGravityBackwardStep::ArgsType(model,data));
295  }
296 
297  return data.g;
298  }
299 
300  struct CoriolisMatrixForwardStep : public fusion::JointVisitor<CoriolisMatrixForwardStep>
301  {
302  typedef boost::fusion::vector<const se3::Model &,
303  se3::Data &,
304  const Eigen::VectorXd &,
305  const Eigen::VectorXd &
306  > ArgsType;
307 
308  JOINT_VISITOR_INIT(CoriolisMatrixForwardStep);
309 
310  template<typename JointModel>
311  static void algo(const se3::JointModelBase<JointModel> & jmodel,
313  const se3::Model & model,
314  se3::Data & data,
315  const Eigen::VectorXd & q,
316  const Eigen::VectorXd & v)
317  {
318  const Model::JointIndex & i = jmodel.id();
319  const Model::JointIndex & parent = model.parents[i];
320 
321  jmodel.calc(jdata.derived(),q,v);
322 
323  data.liMi[i] = model.jointPlacements[i]*jdata.M();
324  if(parent>0) data.oMi[i] = data.oMi[parent] * data.liMi[i];
325  else data.oMi[i] = data.liMi[i];
326 
327  // express quantities in the world frame
328  data.oYcrb[i] = data.oMi[i].act(model.inertias[i]);
329 
330  data.v[i] = jdata.v();
331  if(parent>0) data.v[i] += data.liMi[i].actInv(data.v[parent]);
332  data.ov[i] = data.oMi[i].act(data.v[i]);
333 
334  // computes S expressed at the world frame
335  typedef typename SizeDepType<JointModel::NV>::template ColsReturn<Data::Matrix6x>::Type ColsBlock;
336  ColsBlock Jcols = jmodel.jointCols(data.J);
337  Jcols = data.oMi[i].act(jdata.S()); // collection of S expressed at the world frame
338 
339  // computes vxS expressed at the world frame
340  ColsBlock dJcols = jmodel.jointCols(data.dJ);
341  motionSet::motionAction(data.ov[i],Jcols,dJcols);
342 
343  // computes vxI
344  Inertia::vxi(data.ov[i],data.oYcrb[i],data.vxI[i]);
345  }
346 
347  };
348 
349  struct CoriolisMatrixBackwardStep : public fusion::JointModelVisitor<CoriolisMatrixBackwardStep>
350  {
351  typedef boost::fusion::vector<const Model &,
352  Data &
353  > ArgsType;
354 
355  JOINT_MODEL_VISITOR_INIT(CoriolisMatrixBackwardStep);
356 
357  template<typename JointModel>
358  static void algo(const JointModelBase<JointModel> & jmodel,
359  const Model & model,
360  Data & data)
361  {
362  const Model::JointIndex & i = jmodel.id();
363  const Model::JointIndex & parent = model.parents[i];
364  Data::RowMatrix6 & M6tmpR = data.M6tmpR;
365 
366  typedef typename SizeDepType<JointModel::NV>::template ColsReturn<Data::Matrix6x>::Type ColsBlock;
367  ColsBlock dJcols = jmodel.jointCols(data.dJ);
368  ColsBlock Jcols = jmodel.jointCols(data.J);
369 
370  motionSet::inertiaAction(data.oYcrb[i],dJcols,jmodel.jointCols(data.dFdv));
371  jmodel.jointCols(data.dFdv) += data.vxI[i] * Jcols;
372 
373  data.C.block(jmodel.idx_v(),jmodel.idx_v(),jmodel.nv(),data.nvSubtree[i]).noalias()
374  = Jcols.transpose()*data.dFdv.middleCols(jmodel.idx_v(),data.nvSubtree[i]);
375 
376  lhsInertiaMult(data.oYcrb[i],Jcols.transpose(),M6tmpR.topRows(jmodel.nv()));
377  for(int j = data.parents_fromRow[(Model::Index)jmodel.idx_v()];j >= 0; j = data.parents_fromRow[(Model::Index)j])
378  data.C.middleRows(jmodel.idx_v(),jmodel.nv()).col(j).noalias() = M6tmpR.topRows(jmodel.nv()) * data.dJ.col(j);
379 
380  M6tmpR.topRows(jmodel.nv()).noalias() = Jcols.transpose() * data.vxI[i];
381  for(int j = data.parents_fromRow[(Model::Index)jmodel.idx_v()];j >= 0; j = data.parents_fromRow[(Model::Index)j])
382  data.C.middleRows(jmodel.idx_v(),jmodel.nv()).col(j) += M6tmpR.topRows(jmodel.nv()) * data.J.col(j);
383 
384  if(parent>0)
385  {
386  data.oYcrb[parent] += data.oYcrb[i];
387  data.vxI[parent] += data.vxI[i];
388  }
389 
390  }
391 
392  template<typename Min, typename Mout>
393  static void lhsInertiaMult(const Inertia & Y,
394  const Eigen::MatrixBase<Min> & J,
395  const Eigen::MatrixBase<Mout> & F)
396  {
397  Mout & F_ = const_cast<Mout &>(F.derived());
398  motionSet::inertiaAction(Y,J.derived().transpose(),F_.transpose());
399  }
400  };
401 
402  inline const Eigen::MatrixXd &
403  computeCoriolisMatrix(const Model & model, Data & data,
404  const Eigen::VectorXd & q,
405  const Eigen::VectorXd & v)
406  {
407  assert(model.check(data) && "data is not consistent with model.");
408  assert(q.size() == model.nq);
409  assert(v.size() == model.nv);
410 
411  for(size_t i=1;i<(size_t) model.njoints;++i)
412  {
413  CoriolisMatrixForwardStep::run(model.joints[i],data.joints[i],
414  CoriolisMatrixForwardStep::ArgsType(model,data,q,v));
415  }
416 
417  for(size_t i=(size_t) (model.njoints-1);i>0;--i)
418  {
419  CoriolisMatrixBackwardStep::run(model.joints[i],
420  CoriolisMatrixBackwardStep::ArgsType(model,data));
421  }
422 
423  return data.C;
424  }
425 
426 } // namespace se3
427 
429 
430 #endif // ifndef __se3_rnea_hxx__
Eigen::VectorXd tau
Vector of joint torques (dim model.nv).
Definition: data.hpp:93
int nq
Dimension of the configuration vector representation.
Definition: model.hpp:49
int nv
Dimension of the velocity vector space.
Definition: model.hpp:52
container::aligned_vector< SE3 > jointPlacements
Placement (SE3) of the input of joint i regarding to the parent joint output li.
Definition: model.hpp:67
std::vector< int > nvSubtree
Dimension of the subtree motion space (for CRBA)
Definition: data.hpp:195
static void vxi(const Motion &v, const InertiaTpl< double, _Options > &I, const Eigen::MatrixBase< M6 > &Iout)
Time variation operator. It computes the time derivative of an inertia I corresponding to the formula...
Definition: inertia.hpp:78
container::aligned_vector< Inertia > oYcrb
Inertia quantities expressed in the world frame.
Definition: data.hpp:149
const Eigen::VectorXd & rnea(const Model &model, Data &data, const Eigen::VectorXd &q, const Eigen::VectorXd &v, const Eigen::VectorXd &a)
The Recursive Newton-Euler algorithm. It computes the inverse dynamics, aka the joint torques accordi...
Definition: rnea.hxx:89
Eigen::VectorXd g
Vector of generalized gravity (dim model.nv).
Definition: data.hpp:103
JointModelVector joints
Model of joint i, encapsulated in a JointModelAccessor.
Definition: model.hpp:70
container::aligned_vector< Inertia::Matrix6 > vxI
Right variation of the inertia matrix.
Definition: data.hpp:143
Matrix6x J
Jacobian of joint placements.
Definition: data.hpp:217
Matrix6x dFdv
Variation of the forceset with respect to the joint velocity.
Definition: data.hpp:128
std::vector< JointIndex > parents
Joint parent of joint i, denoted li (li==parents[i]).
Definition: model.hpp:73
const Eigen::VectorXd & nonLinearEffects(const Model &model, Data &data, const Eigen::VectorXd &q, const Eigen::VectorXd &v)
Computes the non-linear effects (Corriolis, centrifual and gravitationnal effects), also called the bias terms of the Lagrangian dynamics:
Definition: rnea.hxx:202
container::aligned_vector< Motion > ov
Vector of joint velocities expressed at the origin.
Definition: data.hpp:70
Eigen::VectorXd nle
Vector of Non Linear Effects (dim model.nv). It corresponds to concatenation of the Coriolis...
Definition: data.hpp:98
Matrix6x dJ
Derivative of the Jacobian with respect to the time.
Definition: data.hpp:220
container::aligned_vector< Motion > a_gf
Vector of joint accelerations due to the gravity field.
Definition: data.hpp:64
JointDataVector joints
Vector of se3::JointData associated to the se3::JointModel stored in model, encapsulated in JointData...
Definition: data.hpp:55
Motion gravity
Spatial gravity of the model.
Definition: model.hpp:106
container::aligned_vector< SE3 > liMi
Vector of relative joint placements (wrt the body parent).
Definition: data.hpp:90
container::aligned_vector< SE3 > oMi
Vector of absolute joint placements (wrt the world).
Definition: data.hpp:87
bool check(const AlgorithmCheckerBase< D > &checker=AlgorithmCheckerBase< D >()) const
Definition: model.hpp:331
container::aligned_vector< Inertia > inertias
Spatial inertias of the body i expressed in the supporting joint frame i.
Definition: model.hpp:64
container::aligned_vector< Force > f
Vector of body forces expressed in the local frame of the joint. For each body, the force represents ...
Definition: data.hpp:74
const Eigen::MatrixXd & computeCoriolisMatrix(const Model &model, Data &data, const Eigen::VectorXd &q, const Eigen::VectorXd &v)
Computes the Coriolis Matrix of the Lagrangian dynamics:
Definition: rnea.hxx:403
RowMatrix6 M6tmpR
Temporary for derivative algorithms.
Definition: data.hpp:158
std::vector< int > parents_fromRow
First previous non-zero row in M (used in Cholesky Decomposition).
Definition: data.hpp:210
Specialization of an std::vector with an aligned allocator. This specialization might be used when th...
container::aligned_vector< Motion > v
Vector of joint velocities expressed at the centers of the joints.
Definition: data.hpp:67
int njoints
Number of joints.
Definition: model.hpp:55
const Eigen::VectorXd & computeGeneralizedGravity(const Model &model, Data &data, const Eigen::VectorXd &q)
Computes the generalized gravity contribution of the Lagrangian dynamics:
Definition: rnea.hxx:278
Eigen::MatrixXd C
The Coriolis matrix (a square matrix of dim model.nv).
Definition: data.hpp:122