18 #ifndef __se3_kinematics_derivatives_hxx__ 19 #define __se3_kinematics_derivatives_hxx__ 21 #include "pinocchio/multibody/visitor.hpp" 22 #include "pinocchio/algorithm/check.hpp" 29 typedef boost::fusion::vector<
const se3::Model &,
31 const Eigen::VectorXd &,
32 const Eigen::VectorXd &,
33 const Eigen::VectorXd &
38 template<
typename Jo
intModel>
41 const se3::Model & model,
43 const Eigen::VectorXd & q,
44 const Eigen::VectorXd & v,
45 const Eigen::VectorXd & a)
47 const Model::JointIndex & i = jmodel.id();
48 const Model::JointIndex & parent = model.
parents[i];
55 jmodel.calc(jdata.derived(),q,v);
62 oMi = data.
oMi[parent]*data.
liMi[i];
63 vi += data.
liMi[i].actInv(data.
v[parent]);
68 ai = jdata.S() * jmodel.jointVelocitySelector(a) + jdata.c() + (vi ^ jdata.v());
70 ai += data.
liMi[i].actInv(data.
a[parent]);
73 ColsBlock dJcols = jmodel.jointCols(data.
dJ);
74 ColsBlock Jcols = jmodel.jointCols(data.
J);
76 Jcols = oMi.
act(jdata.S());
78 motionSet::motionAction(ov,Jcols,dJcols);
86 const Eigen::VectorXd & q,
87 const Eigen::VectorXd & v,
88 const Eigen::VectorXd & a)
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.");
98 for(Model::JointIndex i=1; i<(Model::JointIndex) model.
njoints; ++i)
100 ForwardKinematicsDerivativesForwardStep::run(model.
joints[i],data.
joints[i],
101 ForwardKinematicsDerivativesForwardStep::ArgsType(model,data,q,v,a));
105 template<ReferenceFrame rf>
108 typedef boost::fusion::vector<
const se3::Model &,
118 template<
typename Jo
intModel>
120 const se3::Model & model,
123 const Motion & vlast,
124 Data::Matrix6x & v_partial_dq,
125 Data::Matrix6x & v_partial_dv)
127 const Model::JointIndex & i = jmodel.id();
128 const Model::JointIndex & parent = model.
parents[i];
129 Motion & vtmp = data.
ov[0];
132 ColsBlock Jcols = jmodel.jointCols(data.
J);
135 ColsBlock v_partial_dv_cols = jmodel.jointCols(v_partial_dv);
137 v_partial_dv_cols = Jcols;
139 motionSet::se3ActionInverse(oMlast,Jcols,v_partial_dv_cols);
142 ColsBlock v_partial_dq_cols = jmodel.jointCols(v_partial_dq);
146 vtmp = data.
ov[parent] - vlast;
149 motionSet::motionAction(vtmp,Jcols,v_partial_dq_cols);
155 vtmp = oMlast.
actInv(data.
ov[parent]);
156 motionSet::motionAction(vtmp,v_partial_dv_cols,v_partial_dq_cols);
165 template<ReferenceFrame rf>
168 const Model::JointIndex jointId,
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.");
176 const SE3 & oMlast = data.
oMi[jointId];
177 const Motion & vlast = data.
ov[jointId];
179 for(Model::JointIndex i = jointId; i > 0; i = model.
parents[i])
182 typename JointVelocityDerivativesBackwardStep<rf>::ArgsType(model,data,
189 data.
ov[0].setZero();
192 template<ReferenceFrame rf>
196 typedef boost::fusion::vector<
const se3::Model &,
198 const Model::JointIndex,
207 template<
typename Jo
intModel>
209 const se3::Model & model,
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)
217 const Model::JointIndex & i = jmodel.id();
218 const Model::JointIndex & parent = model.
parents[i];
222 const SE3 & oMlast = data.
oMi[jointId];
223 const Motion & vlast = data.
ov[jointId];
224 const Motion & alast = data.
oa[jointId];
227 ColsBlock dJcols = jmodel.jointCols(data.
dJ);
228 ColsBlock Jcols = jmodel.jointCols(data.
J);
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);
237 a_partial_da_cols = Jcols;
239 motionSet::se3ActionInverse(oMlast,Jcols,a_partial_da_cols);
245 vtmp = data.
ov[parent] - vlast;
250 motionSet::motionAction(vtmp,Jcols,v_partial_dq_cols);
252 a_partial_dv_cols = v_partial_dq_cols + dJcols;
259 vtmp = oMlast.
actInv(data.
ov[parent]);
260 motionSet::motionAction(vtmp,a_partial_da_cols,v_partial_dq_cols);
264 vtmp -= data.
v[jointId];
266 vtmp = -data.
v[jointId];
268 motionSet::motionAction(vtmp,a_partial_da_cols,a_partial_dv_cols);
269 motionSet::se3ActionInverse<ADDTO>(oMlast,dJcols,a_partial_dv_cols);
276 atmp = data.
oa[parent] - alast;
279 motionSet::motionAction(atmp,Jcols,a_partial_dq_cols);
282 motionSet::motionAction<ADDTO>(vtmp,dJcols,a_partial_dq_cols);
288 atmp = oMlast.
actInv(data.
oa[parent]);
289 motionSet::motionAction(atmp,a_partial_da_cols,a_partial_dq_cols);
292 motionSet::motionAction<ADDTO>(vtmp,v_partial_dq_cols,a_partial_dq_cols);
300 template<ReferenceFrame rf>
303 const Model::JointIndex jointId,
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.");
315 for(Model::JointIndex i = jointId; i > 0; i = model.
parents[i])
318 typename JointAccelerationDerivativesBackwardStep<rf>::ArgsType(model,data,
328 data.
ov[0].setZero();
329 data.
oa[0].setZero();
334 #endif // ifndef __se3_kinematics_derivatives_hxx__ container::aligned_vector< Motion > a
Vector of joint accelerations expressed at the centers of the joints.
int nq
Dimension of the configuration vector representation.
int nv
Dimension of the velocity vector space.
container::aligned_vector< SE3 > jointPlacements
Placement (SE3) of the input of joint i regarding to the parent joint output li.
EIGEN_MAKE_ALIGNED_OPERATOR_NEW typedef Eigen::Matrix< double, 6, Eigen::Dynamic > Matrix6x
The 6d jacobian type (temporary)
JointModelVector joints
Model of joint i, encapsulated in a JointModelAccessor.
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.
container::aligned_vector< Motion > oa
Vector of joint accelerations expressed at the origin.
std::vector< JointIndex > parents
Joint parent of joint i, denoted li (li==parents[i]).
container::aligned_vector< Motion > ov
Vector of joint velocities expressed at the origin.
Matrix6x dJ
Derivative of the Jacobian with respect to the time.
JointDataVector joints
Vector of se3::JointData associated to the se3::JointModel stored in model, encapsulated in JointData...
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).
container::aligned_vector< SE3 > oMi
Vector of absolute joint placements (wrt the world).
bool check(const AlgorithmCheckerBase< D > &checker=AlgorithmCheckerBase< D >()) const
internal::SE3GroupAction< D >::ReturnType act(const D &d) const
ay = aXb.act(by)
container::aligned_vector< Motion > v
Vector of joint velocities expressed at the centers of the joints.
int njoints
Number of joints.
internal::SE3GroupAction< D >::ReturnType actInv(const D &d) const
by = aXb.actInv(ay)
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...