Go to the documentation of this file.
13 #ifndef EIGEN_RUNTIME_NO_MALLOC
14 #define EIGEN_RUNTIME_NO_MALLOC
30 #include <boost/test/unit_test.hpp>
31 #include <boost/utility/binary.hpp>
33 template<
typename Jo
intModel>
37 const std::string & parent_name,
38 const std::string &
name,
40 bool setRandomLimits =
true)
43 typedef typename JointModel::ConfigVector_t CV;
44 typedef typename JointModel::TangentVector_t TV;
51 TV::Random() + TV::Constant(1), TV::Random() + TV::Constant(1),
52 CV::Random() - CV::Constant(1), CV::Random() + CV::Constant(1));
56 model.addJointFrame(idx);
62 BOOST_AUTO_TEST_SUITE(BOOST_TEST_MODULE)
71 #ifdef _INTENSE_TESTING_
72 const size_t NBT = 1000 * 1000;
74 const size_t NBT = 10;
77 Eigen::VectorXd
q = Eigen::VectorXd::Zero(
model.nq);
85 timer.
toc(std::cout, NBT);
89 using namespace Eigen;
93 Eigen::VectorXd
q = Eigen::VectorXd::Ones(
model.nq);
95 Eigen::VectorXd
v = Eigen::VectorXd::Zero(
model.nv);
96 Eigen::VectorXd
a = Eigen::VectorXd::Zero(
model.nv);
99 data.M.triangularView<Eigen::StrictlyLower>() =
100 data.M.transpose().triangularView<Eigen::StrictlyLower>();
111 BOOST_CHECK(
M.isApprox(
data.M, 1e-12));
113 std::cout <<
"(the time score in debug mode is not relevant) ";
115 q = Eigen::VectorXd::Zero(
model.nq);
123 timer.
toc(std::cout, NBT);
125 #endif // ifndef NDEBUG
134 model.lowerPositionLimit.head<7>().
fill(-1.);
135 model.upperPositionLimit.head<7>().
fill(1.);
139 Eigen::VectorXd
v(Eigen::VectorXd::Random(
model.nv));
142 data_ref.
M.triangularView<Eigen::StrictlyLower>() =
143 data_ref.
M.transpose().triangularView<Eigen::StrictlyLower>();
146 data.M.triangularView<Eigen::StrictlyLower>() =
147 data.M.transpose().triangularView<Eigen::StrictlyLower>();
149 BOOST_CHECK(
data.M.isApprox(data_ref.
M));
153 BOOST_CHECK(
data.Ag.isApprox(data_ref.
Ag));
154 BOOST_CHECK(
data.J.isApprox(data_ref.
J));
162 BOOST_CHECK(data2.Ycrb[0] ==
data.Ycrb[0]);
172 BOOST_CHECK(
model == model_ref);
176 model.armature = Eigen::VectorXd::Random(
model.nv) + Eigen::VectorXd::Constant(
model.nv, 1.);
180 data_ref.
M.triangularView<Eigen::StrictlyLower>() =
181 data_ref.
M.transpose().triangularView<Eigen::StrictlyLower>();
182 data_ref.
M.diagonal() +=
model.armature;
185 data.M.triangularView<Eigen::StrictlyLower>() =
186 data.M.transpose().triangularView<Eigen::StrictlyLower>();
188 BOOST_CHECK(
data.M.isApprox(data_ref.
M));
201 SE3::Random(),
"revolute_unaligned");
205 Eigen::internal::set_is_malloc_allowed(
false);
207 Eigen::internal::set_is_malloc_allowed(
true);
212 BOOST_AUTO_TEST_SUITE_END()
MatrixXs M
The joint space inertia matrix (a square matrix of dim model.nv).
std::string name(const LieGroupGenericTpl< LieGroupCollection > &lg)
Visit a LieGroupVariant to get the name of it.
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 & 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....
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.
pinocchio::JointIndex JointIndex
const MotionDense< Motion2 > const SE3Tpl< SE3Scalar, SE3Options > & placement
BOOST_AUTO_TEST_CASE(test_crba)
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()
Matrix6x Ag
Centroidal Momentum Matrix.
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)
MotionTpl< Scalar, Options > bias(const JointDataTpl< Scalar, Options, JointCollectionTpl > &jdata)
Visit a JointDataTpl through JointBiasVisitor to get the joint bias as a dense motion.
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...
void neutral(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const Eigen::MatrixBase< ReturnType > &qout)
Return the neutral configuration element related to the model configuration space.
void normalize(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const Eigen::MatrixBase< ConfigVectorType > &qout)
Normalize a configuration vector.
JointCollectionTpl & model
Main pinocchio namespace.
pinocchio
Author(s):
autogenerated on Thu Dec 19 2024 03:41:28