pinocchio  UNKNOWN
centroidal.hxx
1 //
2 // Copyright (c) 2015-2018 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_centroidal_hxx__
19 #define __se3_centroidal_hxx__
20 
21 #include "pinocchio/multibody/visitor.hpp"
22 #include "pinocchio/spatial/act-on-set.hpp"
23 #include "pinocchio/algorithm/kinematics.hpp"
24 #include "pinocchio/algorithm/check.hpp"
25 
27 
28 namespace se3
29 {
30 
31  struct CcrbaForwardStep : public fusion::JointVisitor<CcrbaForwardStep>
32  {
33  typedef boost::fusion::vector< const se3::Model &,
34  se3::Data &,
35  const Eigen::VectorXd &
36  > ArgsType;
37 
38  JOINT_VISITOR_INIT(CcrbaForwardStep);
39 
40  template<typename JointModel>
41  static void algo(const se3::JointModelBase<JointModel> & jmodel,
43  const se3::Model & model,
44  se3::Data & data,
45  const Eigen::VectorXd & q)
46  {
47  const Model::JointIndex & i = (Model::JointIndex) jmodel.id();
48  const Model::Index & parent = model.parents[i];
49 
50  jmodel.calc(jdata.derived(),q);
51 
52  data.liMi[i] = model.jointPlacements[i]*jdata.M();
53  data.Ycrb[i] = model.inertias[i];
54 
55  if (parent>0) data.oMi[i] = data.oMi[parent]*data.liMi[i];
56  else data.oMi[i] = data.liMi[i];
57  }
58 
59  }; // struct CcrbaForwardStep
60 
61  struct CcrbaBackwardStep : public fusion::JointVisitor<CcrbaBackwardStep>
62  {
63  typedef boost::fusion::vector< const se3::Model &,
64  se3::Data &
65  > ArgsType;
66 
67  JOINT_VISITOR_INIT(CcrbaBackwardStep);
68 
69  template<typename JointModel>
70  static void algo(const se3::JointModelBase<JointModel> & jmodel,
72  const se3::Model & model,
73  se3::Data & data)
74  {
75  typedef typename SizeDepType<JointModel::NV>::template ColsReturn<Data::Matrix6x>::Type ColsBlock;
76 
77  const Model::JointIndex & i = (Model::JointIndex) jmodel.id();
78  const Model::Index & parent = model.parents[i];
79 
80  data.Ycrb[parent] += data.liMi[i].act(data.Ycrb[i]);
81 
82  jdata.U() = data.Ycrb[i] * jdata.S();
83 
84  ColsBlock jF
85  = data.Ag.middleCols <JointModel::NV> (jmodel.idx_v());
86  // = data.Ag.middleCols(jmodel.idx_v(), jmodel.nv());
87 
88  forceSet::se3Action(data.oMi[i],jdata.U(),jF);
89  }
90 
91  }; // struct CcrbaBackwardStep
92 
93  inline const Data::Matrix6x &
94  ccrba(const Model & model, Data & data,
95  const Eigen::VectorXd & q,
96  const Eigen::VectorXd & v)
97  {
98  assert(model.check(data) && "data is not consistent with model.");
99  typedef Eigen::Block<Data::Matrix6x,3,-1> Block3x;
100 
101  forwardKinematics(model, data, q);
102  data.Ycrb[0].setZero();
103  for(Model::Index i=1;i<(Model::Index)(model.njoints);++i )
104  data.Ycrb[i] = model.inertias[i];
105 
106 
107  for(Model::Index i=(Model::Index)(model.njoints-1);i>0;--i)
108  {
109  CcrbaBackwardStep::run(model.joints[i],data.joints[i],
110  CcrbaBackwardStep::ArgsType(model,data));
111  }
112 
113  // Express the centroidal map around the center of mass
114  data.com[0] = data.Ycrb[0].lever();
115 
116  const Block3x Ag_lin = data.Ag.middleRows<3>(Force::LINEAR);
117  Block3x Ag_ang = data.Ag.middleRows<3>(Force::ANGULAR);
118  for (long i = 0; i<model.nv; ++i)
119  Ag_ang.col(i) += Ag_lin.col(i).cross(data.com[0]);
120 
121  data.hg.toVector().noalias() = data.Ag*v;
122 
123  data.Ig.mass() = data.Ycrb[0].mass();
124  data.Ig.lever().setZero();
125  data.Ig.inertia() = data.Ycrb[0].inertia();
126 
127  return data.Ag;
128  }
129 
130  struct DCcrbaForwardStep : public fusion::JointVisitor<DCcrbaForwardStep>
131  {
132  typedef boost::fusion::vector< const se3::Model &,
133  se3::Data &,
134  const Eigen::VectorXd &,
135  const Eigen::VectorXd &
136  > ArgsType;
137 
138  JOINT_VISITOR_INIT(DCcrbaForwardStep);
139 
140  template<typename JointModel>
141  static void algo(const se3::JointModelBase<JointModel> & jmodel,
143  const se3::Model & model,
144  se3::Data & data,
145  const Eigen::VectorXd & q,
146  const Eigen::VectorXd & v)
147  {
148  const Model::JointIndex & i = (Model::JointIndex) jmodel.id();
149  const Model::Index & parent = model.parents[i];
150 
151  jmodel.calc(jdata.derived(),q,v);
152 
153  data.liMi[i] = model.jointPlacements[i]*jdata.M();
154  data.Ycrb[i] = model.inertias[i];
155 
156  data.v[i] = jdata.v();
157 
158  if (parent>0)
159  {
160  data.oMi[i] = data.oMi[parent]*data.liMi[i];
161  data.v[i] += data.liMi[i].actInv(data.v[parent]);
162  }
163  else data.oMi[i] = data.liMi[i];
164  }
165 
166  }; // struct DCcrbaForwardStep
167 
168  struct DCcrbaBackwardStep : public fusion::JointVisitor<DCcrbaBackwardStep>
169  {
170  typedef boost::fusion::vector< const se3::Model &,
171  se3::Data &
172  > ArgsType;
173 
174  JOINT_VISITOR_INIT(DCcrbaBackwardStep);
175 
176  template<typename JointModel>
177  static void algo(const se3::JointModelBase<JointModel> & jmodel,
179  const se3::Model & model,
180  se3::Data & data)
181  {
182  typedef typename SizeDepType<JointModel::NV>::template ColsReturn<Data::Matrix6x>::Type ColsBlock;
183 
184  const Model::JointIndex & i = (Model::JointIndex) jmodel.id();
185  const Model::Index & parent = model.parents[i];
186  const Inertia & Y = data.oYcrb[i];
187  const Inertia::Matrix6 & doYcrb = data.doYcrb[i];
188 
189  ColsBlock J_cols = jmodel.jointCols(data.J);
190  J_cols = data.oMi[i].act(jdata.S());
191 
192  ColsBlock dJ_cols = jmodel.jointCols(data.dJ);
193  motionSet::motionAction(data.ov[i],J_cols,dJ_cols);
194 
195  data.oYcrb[parent] += Y;
196  if(parent > 0)
197  data.doYcrb[parent] += doYcrb;
198 
199  // Calc Ag
200  ColsBlock Ag_cols = jmodel.jointCols(data.Ag);
201  motionSet::inertiaAction(Y,J_cols,Ag_cols);
202 
203  // Calc dAg = Ivx + vxI
204  ColsBlock dAg_cols = jmodel.jointCols(data.dAg);
205  motionSet::inertiaAction(Y,dJ_cols,dAg_cols);
206  dAg_cols += doYcrb * J_cols;
207  }
208 
209  }; // struct DCcrbaBackwardStep
210 
211  inline const Data::Matrix6x &
212  dccrba(const Model & model, Data & data,
213  const Eigen::VectorXd & q,
214  const Eigen::VectorXd & v)
215  {
216  assert(model.check(data) && "data is not consistent with model.");
217  typedef Eigen::Block <Data::Matrix6x,3,-1> Block3x;
218 
219  forwardKinematics(model,data,q,v);
220  data.oYcrb[0].setZero();
221  for(Model::Index i=1;i<(Model::Index)(model.njoints);++i)
222  {
223  data.oYcrb[i] = data.oMi[i].act(model.inertias[i]);
224  data.ov[i] = data.oMi[i].act(data.v[i]); // v_i expressed in the world frame
225  data.doYcrb[i] = data.oYcrb[i].variation(data.ov[i]);
226  }
227 
228  for(Model::Index i=(Model::Index)(model.njoints-1);i>0;--i)
229  {
230  DCcrbaBackwardStep::run(model.joints[i],data.joints[i],
231  DCcrbaBackwardStep::ArgsType(model,data));
232  }
233 
234  // Express the centroidal map around the center of mass
235  data.com[0] = data.oYcrb[0].lever();
236 
237  const Block3x Ag_lin = data.Ag.middleRows<3> (Force::LINEAR);
238  Block3x Ag_ang = data.Ag.middleRows<3> (Force::ANGULAR);
239  for (long i = 0; i<model.nv; ++i)
240  Ag_ang.col(i) += Ag_lin.col(i).cross(data.com[0]);
241 
242  data.hg.toVector().noalias() = data.Ag*v;
243  data.vcom[0].noalias() = data.hg.linear()/data.oYcrb[0].mass();
244 
245  const Block3x dAg_lin = data.dAg.middleRows<3>(Force::LINEAR);
246  Block3x dAg_ang = data.dAg.middleRows<3>(Force::ANGULAR);
247  for (long i = 0; i<model.nv; ++i)
248  dAg_ang.col(i) += dAg_lin.col(i).cross(data.com[0]);
249 
250  data.Ig.mass() = data.oYcrb[0].mass();
251  data.Ig.lever().setZero();
252  data.Ig.inertia() = data.oYcrb[0].inertia();
253 
254  return data.dAg;
255  }
256 
257  // --- CHECKER ---------------------------------------------------------------
258  // --- CHECKER ---------------------------------------------------------------
259  // --- CHECKER ---------------------------------------------------------------
260 
261 
262 } // namespace se3
263 
265 
266 #endif // ifndef __se3_centroidal_hxx__
267 
const Data::Matrix6x & dccrba(const Model &model, Data &data, const Eigen::VectorXd &q, const Eigen::VectorXd &v)
Computes the time derivative of the Centroidal Momentum Matrix according to the current configuration...
Definition: centroidal.hxx:212
Matrix6x Ag
Centroidal Momentum Matrix.
Definition: data.hpp:173
Inertia Ig
Centroidal Composite Rigid Body Inertia.
Definition: data.hpp:187
container::aligned_vector< Eigen::Vector3d > vcom
Vector of subtree center of mass linear velocities expressed in the root joint of the subtree...
Definition: data.hpp:250
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
container::aligned_vector< Inertia > oYcrb
Inertia quantities expressed in the world frame.
Definition: data.hpp:149
EIGEN_MAKE_ALIGNED_OPERATOR_NEW typedef Eigen::Matrix< double, 6, Eigen::Dynamic > Matrix6x
The 6d jacobian type (temporary)
Definition: data.hpp:45
JointModelVector joints
Model of joint i, encapsulated in a JointModelAccessor.
Definition: model.hpp:70
container::aligned_vector< Inertia > Ycrb
Vector of sub-tree composite rigid body inertias, i.e. the apparent inertia of the subtree supported ...
Definition: data.hpp:110
Matrix6x J
Jacobian of joint placements.
Definition: data.hpp:217
Matrix6x dAg
Centroidal Momentum Matrix Time Variation.
Definition: data.hpp:178
container::aligned_vector< Inertia::Matrix6 > doYcrb
Time variation of the inertia quantities expressed in the world frame.
Definition: data.hpp:152
Force hg
Centroidal momentum quantity.
Definition: data.hpp:183
ConstLinearType linear() const
Return the linear part of the force vector.
Definition: force-base.hpp:60
std::vector< JointIndex > parents
Joint parent of joint i, denoted li (li==parents[i]).
Definition: model.hpp:73
container::aligned_vector< Motion > ov
Vector of joint velocities expressed at the origin.
Definition: data.hpp:70
Matrix6x dJ
Derivative of the Jacobian with respect to the time.
Definition: data.hpp:220
container::aligned_vector< Eigen::Vector3d > com
Vector of subtree center of mass positions expressed in the root joint of the subtree. In other words, com[j] is the CoM position of the subtree supported by joint and expressed in the joint frame . The element com[0] corresponds to the center of mass position of the whole model and expressed in the global frame.
Definition: data.hpp:247
JointDataVector joints
Vector of se3::JointData associated to the se3::JointModel stored in model, encapsulated in JointData...
Definition: data.hpp:55
void forwardKinematics(const Model &model, Data &data, const Eigen::VectorXd &q)
Update the joint placements according to the current joint configuration.
Definition: kinematics.hxx:106
container::aligned_vector< SE3 > liMi
Vector of relative joint placements (wrt the body parent).
Definition: data.hpp:90
const Data::Matrix6x & ccrba(const Model &model, Data &data, const Eigen::VectorXd &q, const Eigen::VectorXd &v)
Computes the Centroidal Momentum Matrix, the Composite Ridig Body Inertia as well as the centroidal m...
Definition: centroidal.hxx:94
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
ToVectorConstReturnType toVector() const
Return the force as an Eigen vector.
Definition: force-base.hpp:99
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