18 #ifndef __se3_centroidal_hxx__ 19 #define __se3_centroidal_hxx__ 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" 33 typedef boost::fusion::vector<
const se3::Model &,
35 const Eigen::VectorXd &
40 template<
typename Jo
intModel>
43 const se3::Model & model,
45 const Eigen::VectorXd & q)
47 const Model::JointIndex & i = (Model::JointIndex) jmodel.id();
48 const Model::Index & parent = model.
parents[i];
50 jmodel.calc(jdata.derived(),q);
55 if (parent>0) data.
oMi[i] = data.
oMi[parent]*data.
liMi[i];
56 else data.
oMi[i] = data.
liMi[i];
63 typedef boost::fusion::vector<
const se3::Model &,
69 template<
typename Jo
intModel>
72 const se3::Model & model,
77 const Model::JointIndex & i = (Model::JointIndex) jmodel.id();
78 const Model::Index & parent = model.
parents[i];
82 jdata.U() = data.
Ycrb[i] * jdata.S();
85 = data.
Ag.middleCols <JointModel::NV> (jmodel.idx_v());
88 forceSet::se3Action(data.
oMi[i],jdata.U(),jF);
95 const Eigen::VectorXd & q,
96 const Eigen::VectorXd & v)
98 assert(model.
check(data) &&
"data is not consistent with model.");
102 data.
Ycrb[0].setZero();
103 for(Model::Index i=1;i<(Model::Index)(model.
njoints);++i )
107 for(Model::Index i=(Model::Index)(model.
njoints-1);i>0;--i)
110 CcrbaBackwardStep::ArgsType(model,data));
114 data.
com[0] = data.
Ycrb[0].lever();
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]);
123 data.
Ig.mass() = data.
Ycrb[0].mass();
124 data.
Ig.lever().setZero();
125 data.
Ig.inertia() = data.
Ycrb[0].inertia();
132 typedef boost::fusion::vector<
const se3::Model &,
134 const Eigen::VectorXd &,
135 const Eigen::VectorXd &
140 template<
typename Jo
intModel>
143 const se3::Model & model,
145 const Eigen::VectorXd & q,
146 const Eigen::VectorXd & v)
148 const Model::JointIndex & i = (Model::JointIndex) jmodel.id();
149 const Model::Index & parent = model.
parents[i];
151 jmodel.calc(jdata.derived(),q,v);
156 data.
v[i] = jdata.v();
161 data.
v[i] += data.
liMi[i].actInv(data.
v[parent]);
163 else data.
oMi[i] = data.
liMi[i];
170 typedef boost::fusion::vector<
const se3::Model &,
176 template<
typename Jo
intModel>
179 const se3::Model & model,
184 const Model::JointIndex & i = (Model::JointIndex) jmodel.id();
185 const Model::Index & parent = model.
parents[i];
187 const Inertia::Matrix6 & doYcrb = data.
doYcrb[i];
189 ColsBlock J_cols = jmodel.jointCols(data.
J);
190 J_cols = data.
oMi[i].act(jdata.S());
192 ColsBlock dJ_cols = jmodel.jointCols(data.
dJ);
193 motionSet::motionAction(data.
ov[i],J_cols,dJ_cols);
195 data.
oYcrb[parent] += Y;
197 data.
doYcrb[parent] += doYcrb;
200 ColsBlock Ag_cols = jmodel.jointCols(data.
Ag);
201 motionSet::inertiaAction(Y,J_cols,Ag_cols);
204 ColsBlock dAg_cols = jmodel.jointCols(data.
dAg);
205 motionSet::inertiaAction(Y,dJ_cols,dAg_cols);
206 dAg_cols += doYcrb * J_cols;
213 const Eigen::VectorXd & q,
214 const Eigen::VectorXd & v)
216 assert(model.
check(data) &&
"data is not consistent with model.");
220 data.
oYcrb[0].setZero();
221 for(Model::Index i=1;i<(Model::Index)(model.
njoints);++i)
224 data.
ov[i] = data.
oMi[i].act(data.
v[i]);
228 for(Model::Index i=(Model::Index)(model.
njoints-1);i>0;--i)
231 DCcrbaBackwardStep::ArgsType(model,data));
235 data.
com[0] = data.
oYcrb[0].lever();
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]);
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]);
250 data.
Ig.mass() = data.
oYcrb[0].mass();
251 data.
Ig.lever().setZero();
252 data.
Ig.inertia() = data.
oYcrb[0].inertia();
266 #endif // ifndef __se3_centroidal_hxx__ 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...
Matrix6x Ag
Centroidal Momentum Matrix.
Inertia Ig
Centroidal Composite Rigid Body Inertia.
container::aligned_vector< Eigen::Vector3d > vcom
Vector of subtree center of mass linear velocities expressed in the root joint of the subtree...
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.
container::aligned_vector< Inertia > oYcrb
Inertia quantities expressed in the world frame.
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.
container::aligned_vector< Inertia > Ycrb
Vector of sub-tree composite rigid body inertias, i.e. the apparent inertia of the subtree supported ...
Matrix6x J
Jacobian of joint placements.
Matrix6x dAg
Centroidal Momentum Matrix Time Variation.
container::aligned_vector< Inertia::Matrix6 > doYcrb
Time variation of the inertia quantities expressed in the world frame.
Force hg
Centroidal momentum quantity.
ConstLinearType linear() const
Return the linear part of the force vector.
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.
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.
JointDataVector joints
Vector of se3::JointData associated to the se3::JointModel stored in model, encapsulated in JointData...
void forwardKinematics(const Model &model, Data &data, const Eigen::VectorXd &q)
Update the joint placements according to the current joint configuration.
container::aligned_vector< SE3 > liMi
Vector of relative joint placements (wrt the body parent).
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...
container::aligned_vector< SE3 > oMi
Vector of absolute joint placements (wrt the world).
bool check(const AlgorithmCheckerBase< D > &checker=AlgorithmCheckerBase< D >()) const
container::aligned_vector< Inertia > inertias
Spatial inertias of the body i expressed in the supporting joint frame i.
ToVectorConstReturnType toVector() const
Return the force as an Eigen vector.
container::aligned_vector< Motion > v
Vector of joint velocities expressed at the centers of the joints.
int njoints
Number of joints.