Go to the documentation of this file.
18 #include <boost/test/unit_test.hpp>
19 #include <boost/utility/binary.hpp>
21 template<
typename Jo
intModel>
25 const std::string & parent_name,
26 const std::string &
name,
28 bool setRandomLimits =
true)
31 typedef typename JointModel::ConfigVector_t CV;
32 typedef typename JointModel::TangentVector_t TV;
39 TV::Random() + TV::Constant(1), TV::Random() + TV::Constant(1),
40 CV::Random() - CV::Constant(1), CV::Random() + CV::Constant(1));
44 model.addJointFrame(idx);
50 BOOST_AUTO_TEST_SUITE(BOOST_TEST_MODULE)
58 model.lowerPositionLimit.head<3>().
fill(-1.);
59 model.upperPositionLimit.head<3>().
fill(1.);
61 Eigen::VectorXd
v = Eigen::VectorXd::Ones(
model.nv);
64 data_ref.
M.triangularView<Eigen::StrictlyLower>() =
65 data_ref.
M.transpose().triangularView<Eigen::StrictlyLower>();
67 data_ref.Ycrb[0] = data_ref.liMi[1].act(data_ref.Ycrb[1]);
71 data_ref_other.
M.triangularView<Eigen::StrictlyLower>() =
72 data_ref_other.
M.transpose().triangularView<Eigen::StrictlyLower>();
73 BOOST_CHECK(data_ref_other.
M.isApprox(data_ref.
M));
75 pinocchio::SE3 cMo(pinocchio::SE3::Matrix3::Identity(), -data_ref.Ycrb[0].lever());
76 BOOST_CHECK(data_ref.Ycrb[0].isApprox(data_ref_other.oYcrb[0]));
80 BOOST_CHECK(
data.oYcrb[0].matrix().isApprox(data_ref.Ycrb[0].matrix(), 1e-12));
83 BOOST_CHECK(
data.Ig.matrix().isApprox(Ig_ref.matrix(), 1e-12));
89 cM1.
inverse().toActionMatrix().transpose() * data_ref.
M.topRows<6>());
90 BOOST_CHECK(
data.Ag.isApprox(Ag_ref, 1e-12));
99 model.lowerPositionLimit.head<3>().
fill(-1.);
100 model.upperPositionLimit.head<3>().
fill(1.);
102 Eigen::VectorXd
v = Eigen::VectorXd::Zero(
model.nv);
107 BOOST_CHECK(
data.J.isApprox(data_ref.
J));
108 BOOST_CHECK(
data.Ag.isApprox(data_ref.
Ag));
112 BOOST_CHECK(
data.J.isApprox(data_ref.
J));
123 model.lowerPositionLimit.head<7>().
fill(-1.);
124 model.upperPositionLimit.head<7>().
fill(1.);
128 Eigen::VectorXd
v = Eigen::VectorXd::Random(
model.nv);
129 Eigen::VectorXd
a = Eigen::VectorXd::Random(
model.nv);
131 const Eigen::VectorXd g =
rnea(
model, data_ref,
q, 0 *
v, 0 *
a);
135 data_ref.
M.triangularView<Eigen::StrictlyLower>() =
136 data_ref.
M.transpose().triangularView<Eigen::StrictlyLower>();
139 data_ref.Ycrb[0] = data_ref.liMi[1].act(data_ref.Ycrb[1]);
140 cMo.translation() = -data_ref.Ycrb[0].lever();
142 SE3 oM1(data_ref.liMi[1]);
144 Data::Matrix6x Ag_ref(cM1.toDualActionMatrix() * data_ref.
M.topRows<6>());
146 Force hdot_ref(cM1.act(
Force(data_ref.
tau.head<6>() - g.head<6>())));
150 BOOST_CHECK(
data.Ag.isApprox(Ag_ref));
151 BOOST_CHECK(
data.Ag.isApprox(data_ref.
Ag));
152 BOOST_CHECK(
data.hg.isApprox(data_ref.
hg));
155 BOOST_CHECK(data_ref.com[0].isApprox(data_ref.Ycrb[0].lever()));
156 BOOST_CHECK(data_ref.vcom[0].isApprox(
data.hg.linear() / data_ref.
M(0, 0)));
157 BOOST_CHECK(data_ref.vcom[0].isApprox(
data.vcom[0]));
158 BOOST_CHECK(data_ref.acom[0].isApprox(hdot_ref.linear() / data_ref.
M(0, 0)));
161 BOOST_CHECK(hdot.isApprox(hdot_ref));
164 BOOST_CHECK(
data.dAg.isZero());
171 const double alpha = 1e-8;
172 Eigen::VectorXd q_plus(
model.nq);
181 for (
size_t i = 1;
i < (size_t)
model.njoints; ++
i)
183 Inertia::Matrix6 dYcrb = (data_ref_plus.oMi[
i].act(data_ref_plus.Ycrb[
i]).matrix()
184 - data_ref.oMi[
i].act(data_ref.Ycrb[
i]).matrix())
186 BOOST_CHECK(
data.doYcrb[
i].isApprox(dYcrb, sqrt(
alpha)));
198 data_ref.oMi[1].toDualActionMatrix() * data_ref.
M.topRows<6>();
200 const double alpha = 1e-8;
201 Eigen::VectorXd q_plus(
model.nq);
206 const Data::Matrix6x Ag_plus_ref = oMc_ref_plus.toDualActionMatrix() * data_ref.
Ag;
209 data_ref.oMi[1].toDualActionMatrix() * data_ref.
M.topRows<6>();
217 BOOST_CHECK(oMc.isApprox(oMc_ref));
218 BOOST_CHECK(dAg_ref.isApprox(dAg_ref_from_M, sqrt(
alpha)));
227 const double alpha = 1e-8;
228 Eigen::VectorXd q_plus(
model.nq);
235 BOOST_CHECK(
data.dAg.isApprox(dAg_ref, sqrt(
alpha)));
240 std::vector<Data::Matrix6x> dAgdq((
size_t)
model.nv, Data::Matrix6x::Zero(6,
model.nv));
242 Eigen::VectorXd v_fd(Eigen::VectorXd::Zero(
model.nv));
248 const Force hg0 = oMc_ref.act(data_fd.
hg);
252 const double alpha = 1e-8;
253 Eigen::VectorXd q_plus(
model.nq);
255 for (
int k = 0; k <
model.nv; ++k)
262 Ag_fd = oMc_fd.toDualActionMatrix() * data_fd.
Ag;
263 hg_fd = oMc_fd.act(data_fd.
hg);
264 dAgdq[(size_t)k] = (Ag_fd - Ag0) /
alpha;
265 dhdq.col(k) = (hg_fd - hg0).toVector() /
alpha;
271 for (
int k = 0; k <
model.nv; ++k)
273 dAg_ref.col(k) = dAgdq[(size_t)k] *
v;
276 BOOST_CHECK(dhdq.isApprox(dAg_ref, sqrt(
alpha)));
287 model.lowerPositionLimit.head<3>().
fill(-1.);
288 model.upperPositionLimit.head<3>().
fill(1.);
290 Eigen::VectorXd
v = Eigen::VectorXd::Zero(
model.nv);
295 BOOST_CHECK(
data.J.isApprox(data_ref.
J));
296 BOOST_CHECK(
data.dJ.isApprox(data_ref.
dJ));
297 BOOST_CHECK(
data.Ag.isApprox(data_ref.
Ag));
298 BOOST_CHECK(
data.dAg.isApprox(data_ref.
dAg));
302 BOOST_CHECK(
data.J.isApprox(data_ref.
J));
303 BOOST_CHECK(
data.dJ.isApprox(data_ref.
dJ));
314 model.lowerPositionLimit.head<7>().
fill(-1.);
315 model.upperPositionLimit.head<7>().
fill(1.);
319 Eigen::VectorXd
v = Eigen::VectorXd::Random(
model.nv);
320 Eigen::VectorXd
a = Eigen::VectorXd::Random(
model.nv);
327 BOOST_CHECK(
data.mass[0] == data_ref.
mass[0]);
328 BOOST_CHECK(
data.com[0].isApprox(data_ref.com[0]));
329 BOOST_CHECK(
data.hg.isApprox(data_ref.
hg));
330 for (
size_t k = 1; k < (size_t)
model.njoints; ++k)
332 BOOST_CHECK(
data.mass[k] == data_ref.
mass[k]);
333 BOOST_CHECK(
data.com[k].isApprox(data_ref.com[k]));
334 BOOST_CHECK(
data.v[k].isApprox(data_ref.v[k]));
341 BOOST_CHECK(data_fk1.hg.isApprox(
data.hg));
344 model.gravity.setZero();
349 BOOST_CHECK(
data.mass[0] == data_ref.
mass[0]);
350 BOOST_CHECK(
data.com[0].isApprox(data_ref.com[0]));
351 BOOST_CHECK(
data.hg.isApprox(data_ref.
hg));
352 BOOST_CHECK(
data.dhg.isApprox(hgdot));
353 for (
size_t k = 1; k < (size_t)
model.njoints; ++k)
355 BOOST_CHECK(
data.mass[k] == data_ref.
mass[k]);
356 BOOST_CHECK(
data.com[k].isApprox(data_ref.com[k]));
357 BOOST_CHECK(
data.v[k].isApprox(data_ref.v[k]));
358 BOOST_CHECK(
data.a[k].isApprox(data_ref.a_gf[k]));
359 BOOST_CHECK(
data.f[k].isApprox(data_ref.f[k]));
366 BOOST_CHECK(data_fk2.hg.isApprox(
data.hg));
367 BOOST_CHECK(data_fk2.dhg.isApprox(
data.dhg));
371 const double eps = 1e-8;
372 Eigen::VectorXd v_plus =
v +
eps *
a;
384 BOOST_CHECK(
data.dhg.isApprox(dhg_ref, sqrt(
eps)));
387 BOOST_AUTO_TEST_SUITE_END()
MatrixXs M
The joint space inertia matrix (a square matrix of dim model.nv).
void addJointAndBody(Model &model, const JointModelBase< D > &jmodel, const Model::JointIndex parent_id, const SE3 &joint_placement, const std::string &name, const Inertia &Y)
traits< SE3Tpl >::Vector3 Vector3
HPP_FCL_DLLAPI Halfspace transform(const Halfspace &a, const Transform3f &tf)
SE3Tpl inverse() const
aXb = bXa.inverse()
Matrix6x dAg
Centroidal Momentum Matrix Time Variation.
std::vector< Scalar > mass
Vector of subtree mass. In other words, mass[j] is the mass of the subtree supported by joint ....
void forwardKinematics(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q)
Update the joint placements according to the current joint configuration.
Matrix6x dJ
Derivative of the Jacobian with respect to the time.
std::string name(const LieGroupGenericTpl< LieGroupCollection > &lg)
Visit a LieGroupVariant to get the name of it.
static void addJointAndBody(pinocchio::Model &model, const pinocchio::JointModelBase< JointModel > &joint, const std::string &parent_name, const std::string &name, const pinocchio::SE3 placement=pinocchio::SE3::Random(), bool setRandomLimits=true)
JointModelSphericalTpl< context::Scalar > JointModelSpherical
void humanoidRandom(ModelTpl< Scalar, Options, JointCollectionTpl > &model, bool usingFF=true)
Create a humanoid kinematic tree with 6-DOF limbs and random joint placements.
const DataTpl< Scalar, Options, JointCollectionTpl >::MatrixXs & crba(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Convention convention=Convention::LOCAL)
Computes the upper triangular part of the joint space inertia matrix M by using the Composite Rigid B...
const DataTpl< Scalar, Options, JointCollectionTpl >::Matrix6x & computeJointJacobiansTimeVariation(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType > &v)
Computes the full model Jacobian variations with respect to time. It corresponds to dJ/dt which depen...
const DataTpl< Scalar, Options, JointCollectionTpl >::Matrix6x & computeJointJacobians(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q)
Computes the full model Jacobian, i.e. the stack of all motion subspace expressed in the world frame....
const DataTpl< Scalar, Options, JointCollectionTpl >::Matrix6x & dccrba(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType > &v)
Computes the time derivative of the Centroidal Momentum Matrix according to the current configuration...
void randomConfiguration(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const Eigen::MatrixBase< ConfigVectorIn1 > &lowerLimits, const Eigen::MatrixBase< ConfigVectorIn2 > &upperLimits, const Eigen::MatrixBase< ReturnType > &qout)
Generate a configuration vector uniformly sampled among provided limits.
const DataTpl< Scalar, Options, JointCollectionTpl >::Matrix6x & computeCentroidalMapTimeVariation(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType > &v)
Computes the Centroidal Momentum Matrix time derivative.
ConstLinearRef translation() const
pinocchio::JointIndex JointIndex
const MotionDense< Motion2 > const SE3Tpl< SE3Scalar, SE3Options > & placement
const DataTpl< Scalar, Options, JointCollectionTpl >::Matrix6x & computeCentroidalMap(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q)
Computes the Centroidal Momentum Matrix.
void integrate(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType > &v, const Eigen::MatrixBase< ReturnType > &qout)
Integrate a configuration vector for the specified model for a tangent vector during one unit time.
TangentVectorType tau
Vector of joint torques (dim model.nv).
const DataTpl< Scalar, Options, JointCollectionTpl >::Matrix6x & ccrba(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType > &v)
Computes the Centroidal Momentum Matrix, the Composite Ridig Body Inertia as well as the centroidal m...
static InertiaTpl Random()
const DataTpl< Scalar, Options, JointCollectionTpl >::Force & computeCentroidalMomentumTimeVariation(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data)
Computes the Centroidal momemtum and its time derivatives, a.k.a. the total momenta of the system and...
Matrix6x Ag
Centroidal Momentum Matrix.
JointCollectionTpl const Eigen::MatrixBase< ConfigVectorType > & q
JointCollectionTpl const Eigen::MatrixBase< ConfigVectorType > const Eigen::MatrixBase< TangentVectorType > & v
Matrix6x J
Jacobian of joint placements.
const DataTpl< Scalar, Options, JointCollectionTpl >::TangentVectorType & rnea(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType1 > &v, const Eigen::MatrixBase< TangentVectorType2 > &a)
The Recursive Newton-Euler algorithm. It computes the inverse dynamics, aka the joint torques accordi...
Force hg
Centroidal momentum quantity.
Eigen::Matrix< Scalar, 6, Eigen::Dynamic, Options > Matrix6x
The 6d jacobian type (temporary)
const DataTpl< Scalar, Options, JointCollectionTpl >::Force & computeCentroidalMomentum(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data)
Computes the Centroidal momentum, a.k.a. the total momenta of the system expressed around the center ...
BOOST_AUTO_TEST_CASE(test_ccrba)
JointCollectionTpl & model
Main pinocchio namespace.
bool isApprox(const Box &s1, const Box &s2, const FCL_REAL tol)
const DataTpl< Scalar, Options, JointCollectionTpl >::Vector3 & centerOfMass(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const bool computeSubtreeComs=true)
Computes the center of mass position of a given model according to a particular joint configuration....
pinocchio
Author(s):
autogenerated on Thu Dec 19 2024 03:41:26