pinocchio  UNKNOWN
kinematics-derivatives.hxx
1 //
2 // Copyright (c) 2017-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_kinematics_derivatives_hxx__
19 #define __se3_kinematics_derivatives_hxx__
20 
21 #include "pinocchio/multibody/visitor.hpp"
22 #include "pinocchio/algorithm/check.hpp"
23 
24 namespace se3
25 {
26 
27  struct ForwardKinematicsDerivativesForwardStep : public fusion::JointVisitor<ForwardKinematicsDerivativesForwardStep>
28  {
29  typedef boost::fusion::vector<const se3::Model &,
30  se3::Data &,
31  const Eigen::VectorXd &,
32  const Eigen::VectorXd &,
33  const Eigen::VectorXd &
34  > ArgsType;
35 
36  JOINT_VISITOR_INIT(ForwardKinematicsDerivativesForwardStep);
37 
38  template<typename JointModel>
39  static void algo(const se3::JointModelBase<JointModel> & jmodel,
41  const se3::Model & model,
42  se3::Data & data,
43  const Eigen::VectorXd & q,
44  const Eigen::VectorXd & v,
45  const Eigen::VectorXd & a)
46  {
47  const Model::JointIndex & i = jmodel.id();
48  const Model::JointIndex & parent = model.parents[i];
49  SE3 & oMi = data.oMi[i];
50  Motion & vi = data.v[i];
51  Motion & ai = data.a[i];
52  Motion & ov = data.ov[i];
53  Motion & oa = data.oa[i];
54 
55  jmodel.calc(jdata.derived(),q,v);
56 
57  data.liMi[i] = model.jointPlacements[i]*jdata.M();
58  vi = jdata.v();
59 
60  if(parent>0)
61  {
62  oMi = data.oMi[parent]*data.liMi[i];
63  vi += data.liMi[i].actInv(data.v[parent]);
64  }
65  else
66  oMi = data.liMi[i];
67 
68  ai = jdata.S() * jmodel.jointVelocitySelector(a) + jdata.c() + (vi ^ jdata.v());
69  if(parent>0)
70  ai += data.liMi[i].actInv(data.a[parent]);
71 
72  typedef typename SizeDepType<JointModel::NV>::template ColsReturn<Data::Matrix6x>::Type ColsBlock;
73  ColsBlock dJcols = jmodel.jointCols(data.dJ);
74  ColsBlock Jcols = jmodel.jointCols(data.J);
75 
76  Jcols = oMi.act(jdata.S());
77  ov = oMi.act(vi); // Spatial velocity of joint i expressed in the global frame o
78  motionSet::motionAction(ov,Jcols,dJcols);
79  oa = oMi.act(ai); // Spatial acceleration of joint i expressed in the global frame o
80  }
81 
82  };
83 
84  inline void
86  const Eigen::VectorXd & q,
87  const Eigen::VectorXd & v,
88  const Eigen::VectorXd & a)
89  {
90  assert(q.size() == model.nq && "The configuration vector is not of right size");
91  assert(v.size() == model.nv && "The velocity vector is not of right size");
92  assert(a.size() == model.nv && "The acceleration vector is not of right size");
93  assert(model.check(data) && "data is not consistent with model.");
94 
95  data.v[0].setZero();
96  data.a[0].setZero();
97 
98  for(Model::JointIndex i=1; i<(Model::JointIndex) model.njoints; ++i)
99  {
100  ForwardKinematicsDerivativesForwardStep::run(model.joints[i],data.joints[i],
101  ForwardKinematicsDerivativesForwardStep::ArgsType(model,data,q,v,a));
102  }
103  }
104 
105  template<ReferenceFrame rf>
106  struct JointVelocityDerivativesBackwardStep : public fusion::JointModelVisitor< JointVelocityDerivativesBackwardStep<rf> >
107  {
108  typedef boost::fusion::vector<const se3::Model &,
109  se3::Data &,
110  const SE3 &,
111  const Motion &,
112  Data::Matrix6x &,
113  Data::Matrix6x &
114  > ArgsType;
115 
116  JOINT_MODEL_VISITOR_INIT(JointVelocityDerivativesBackwardStep);
117 
118  template<typename JointModel>
119  static void algo(const se3::JointModelBase<JointModel> & jmodel,
120  const se3::Model & model,
121  se3::Data & data,
122  const SE3 & oMlast,
123  const Motion & vlast,
124  Data::Matrix6x & v_partial_dq,
125  Data::Matrix6x & v_partial_dv)
126  {
127  const Model::JointIndex & i = jmodel.id();
128  const Model::JointIndex & parent = model.parents[i];
129  Motion & vtmp = data.ov[0]; // Temporary variable
130 
131  typedef typename SizeDepType<JointModel::NV>::template ColsReturn<Data::Matrix6x>::Type ColsBlock;
132  ColsBlock Jcols = jmodel.jointCols(data.J);
133 
134  // dvec/dv
135  ColsBlock v_partial_dv_cols = jmodel.jointCols(v_partial_dv);
136  if(rf == WORLD)
137  v_partial_dv_cols = Jcols;
138  else
139  motionSet::se3ActionInverse(oMlast,Jcols,v_partial_dv_cols);
140 
141  // dvec/dq
142  ColsBlock v_partial_dq_cols = jmodel.jointCols(v_partial_dq);
143  if(rf == WORLD)
144  {
145  if(parent > 0)
146  vtmp = data.ov[parent] - vlast;
147  else
148  vtmp = -vlast;
149  motionSet::motionAction(vtmp,Jcols,v_partial_dq_cols);
150  }
151  else
152  {
153  if(parent > 0)
154  {
155  vtmp = oMlast.actInv(data.ov[parent]);
156  motionSet::motionAction(vtmp,v_partial_dv_cols,v_partial_dq_cols);
157  }
158  }
159 
160 
161  }
162 
163  };
164 
165  template<ReferenceFrame rf>
166  inline void getJointVelocityDerivatives(const Model & model,
167  Data & data,
168  const Model::JointIndex jointId,
169  Data::Matrix6x & v_partial_dq,
170  Data::Matrix6x & v_partial_dv)
171  {
172  assert( v_partial_dq.cols() == model.nv );
173  assert( v_partial_dv.cols() == model.nv );
174  assert(model.check(data) && "data is not consistent with model.");
175 
176  const SE3 & oMlast = data.oMi[jointId];
177  const Motion & vlast = data.ov[jointId];
178 
179  for(Model::JointIndex i = jointId; i > 0; i = model.parents[i])
180  {
182  typename JointVelocityDerivativesBackwardStep<rf>::ArgsType(model,data,
183  oMlast,vlast,
184  v_partial_dq,
185  v_partial_dv));
186  }
187 
188  // Set back ov[0] to a zero value
189  data.ov[0].setZero();
190  }
191 
192  template<ReferenceFrame rf>
194  : public fusion::JointModelVisitor< JointAccelerationDerivativesBackwardStep<rf> >
195  {
196  typedef boost::fusion::vector<const se3::Model &,
197  se3::Data &,
198  const Model::JointIndex,
199  Data::Matrix6x &,
200  Data::Matrix6x &,
201  Data::Matrix6x &,
202  Data::Matrix6x &
203  > ArgsType;
204 
205  JOINT_MODEL_VISITOR_INIT(JointAccelerationDerivativesBackwardStep);
206 
207  template<typename JointModel>
208  static void algo(const se3::JointModelBase<JointModel> & jmodel,
209  const se3::Model & model,
210  se3::Data & data,
211  const Model::JointIndex jointId,
212  Data::Matrix6x & v_partial_dq,
213  Data::Matrix6x & a_partial_dq,
214  Data::Matrix6x & a_partial_dv,
215  Data::Matrix6x & a_partial_da)
216  {
217  const Model::JointIndex & i = jmodel.id();
218  const Model::JointIndex & parent = model.parents[i];
219  Motion & vtmp = data.ov[0]; // Temporary variable
220  Motion & atmp = data.oa[0]; // Temporary variable
221 
222  const SE3 & oMlast = data.oMi[jointId];
223  const Motion & vlast = data.ov[jointId];
224  const Motion & alast = data.oa[jointId];
225 
226  typedef typename SizeDepType<JointModel::NV>::template ColsReturn<Data::Matrix6x>::Type ColsBlock;
227  ColsBlock dJcols = jmodel.jointCols(data.dJ);
228  ColsBlock Jcols = jmodel.jointCols(data.J);
229 
230  ColsBlock v_partial_dq_cols = jmodel.jointCols(v_partial_dq);
231  ColsBlock a_partial_dq_cols = jmodel.jointCols(a_partial_dq);
232  ColsBlock a_partial_dv_cols = jmodel.jointCols(a_partial_dv);
233  ColsBlock a_partial_da_cols = jmodel.jointCols(a_partial_da);
234 
235  // dacc/da
236  if(rf == WORLD)
237  a_partial_da_cols = Jcols;
238  else
239  motionSet::se3ActionInverse(oMlast,Jcols,a_partial_da_cols);
240 
241  // dacc/dv
242  if(rf == WORLD)
243  {
244  if(parent > 0)
245  vtmp = data.ov[parent] - vlast;
246  else
247  vtmp = -vlast;
248 
250  motionSet::motionAction(vtmp,Jcols,v_partial_dq_cols);
251 
252  a_partial_dv_cols = v_partial_dq_cols + dJcols;
253  }
254  else
255  {
257  if(parent > 0)
258  {
259  vtmp = oMlast.actInv(data.ov[parent]);
260  motionSet::motionAction(vtmp,a_partial_da_cols,v_partial_dq_cols);
261  }
262 
263  if(parent > 0)
264  vtmp -= data.v[jointId];
265  else
266  vtmp = -data.v[jointId];
267 
268  motionSet::motionAction(vtmp,a_partial_da_cols,a_partial_dv_cols);
269  motionSet::se3ActionInverse<ADDTO>(oMlast,dJcols,a_partial_dv_cols);
270  }
271 
272  // dacc/dq
273  if(rf == WORLD)
274  {
275  if(parent > 0)
276  atmp = data.oa[parent] - alast;
277  else
278  atmp = -alast;
279  motionSet::motionAction(atmp,Jcols,a_partial_dq_cols);
280 
281  if(parent >0)
282  motionSet::motionAction<ADDTO>(vtmp,dJcols,a_partial_dq_cols);
283  }
284  else
285  {
286  if(parent > 0)
287  {
288  atmp = oMlast.actInv(data.oa[parent]);
289  motionSet::motionAction(atmp,a_partial_da_cols,a_partial_dq_cols);
290  }
291 
292  motionSet::motionAction<ADDTO>(vtmp,v_partial_dq_cols,a_partial_dq_cols);
293  }
294 
295 
296  }
297 
298  };
299 
300  template<ReferenceFrame rf>
301  inline void getJointAccelerationDerivatives(const Model & model,
302  Data & data,
303  const Model::JointIndex jointId,
304  Data::Matrix6x & v_partial_dq,
305  Data::Matrix6x & a_partial_dq,
306  Data::Matrix6x & a_partial_dv,
307  Data::Matrix6x & a_partial_da)
308  {
309  assert( v_partial_dq.cols() == model.nv );
310  assert( a_partial_dq.cols() == model.nv );
311  assert( a_partial_dv.cols() == model.nv );
312  assert( a_partial_da.cols() == model.nv );
313  assert(model.check(data) && "data is not consistent with model.");
314 
315  for(Model::JointIndex i = jointId; i > 0; i = model.parents[i])
316  {
318  typename JointAccelerationDerivativesBackwardStep<rf>::ArgsType(model,data,
319  jointId,
320  v_partial_dq,
321  a_partial_dq,
322  a_partial_dv,
323  a_partial_da));
324 
325  }
326 
327  // Set Zero to temporary spatial variables
328  data.ov[0].setZero();
329  data.oa[0].setZero();
330  }
331 
332 } // namespace se3
333 
334 #endif // ifndef __se3_kinematics_derivatives_hxx__
335 
container::aligned_vector< Motion > a
Vector of joint accelerations expressed at the centers of the joints.
Definition: data.hpp:58
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
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
void getJointVelocityDerivatives(const Model &model, Data &data, const Model::JointIndex jointId, Data::Matrix6x &v_partial_dq, Data::Matrix6x &v_partial_dv)
Computes the partial derivaties of the spatial velocity of a given with respect to the joint configur...
static void algo(const se3::JointModelBase< JointModel > &jmodel, const se3::Model &model, se3::Data &data, const Model::JointIndex jointId, Data::Matrix6x &v_partial_dq, Data::Matrix6x &a_partial_dq, Data::Matrix6x &a_partial_dv, Data::Matrix6x &a_partial_da)
Matrix6x J
Jacobian of joint placements.
Definition: data.hpp:217
container::aligned_vector< Motion > oa
Vector of joint accelerations expressed at the origin.
Definition: data.hpp:61
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
JointDataVector joints
Vector of se3::JointData associated to the se3::JointModel stored in model, encapsulated in JointData...
Definition: data.hpp:55
void getJointAccelerationDerivatives(const Model &model, Data &data, const Model::JointIndex jointId, Data::Matrix6x &v_partial_dq, Data::Matrix6x &a_partial_dq, Data::Matrix6x &a_partial_dv, Data::Matrix6x &a_partial_da)
Computes the partial derivaties of the spatial acceleration of a given with respect to the joint conf...
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
internal::SE3GroupAction< D >::ReturnType act(const D &d) const
ay = aXb.act(by)
Definition: se3.hpp:94
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
internal::SE3GroupAction< D >::ReturnType actInv(const D &d) const
by = aXb.actInv(ay)
Definition: se3.hpp:100
void computeForwardKinematicsDerivatives(const Model &model, Data &data, const Eigen::VectorXd &q, const Eigen::VectorXd &v, const Eigen::VectorXd &a)
Computes all the terms required to compute the derivatives of the placement, spatial velocity and acc...