5 #include "pinocchio/spatial/fwd.hpp" 6 #include "pinocchio/spatial/se3.hpp" 7 #include "pinocchio/multibody/visitor.hpp" 8 #include "pinocchio/multibody/model.hpp" 9 #include "pinocchio/multibody/data.hpp" 10 #include "pinocchio/algorithm/joint-configuration.hpp" 11 #include "pinocchio/algorithm/crba.hpp" 12 #include "pinocchio/algorithm/aba.hpp" 13 #include "pinocchio/algorithm/cholesky.hpp" 14 #include "pinocchio/parsers/urdf.hpp" 15 #include "pinocchio/parsers/sample-models.hpp" 19 #include "pinocchio/utils/timer.hpp" 22 int main(
int argc,
const char ** argv)
24 using namespace Eigen;
29 const int NBT = 1000*100;
32 std::cout <<
"(the time score in debug mode is not relevant) " << std::endl;
38 if(argc>1) filename = argv[1];
43 std::cout <<
"nq = " << model.
nq << std::endl;
46 VectorXd qmax = Eigen::VectorXd::Ones(model.
nq);
48 MatrixXd
A(model.
nv,model.
nv),
B(model.
nv,model.
nv);
49 A.setZero(); B.setRandom();
54 for(
size_t i=0;
i<NBT;++
i)
57 lhs[
i] = Eigen::VectorXd::Zero(model.
nv);
58 rhs[
i] = Eigen::VectorXd::Random(model.
nv);
64 crba(model,data,
qs[_smooth]);
69 std::cout <<
"Cholesky = \t" << (total/NBT)
73 Eigen::LDLT<Eigen::MatrixXd> Mldlt(data.
M);
76 crba(model,data,
qs[_smooth]);
77 data.
M.triangularView<Eigen::StrictlyLower>()
78 = data.
M.transpose().triangularView<Eigen::StrictlyLower>();
80 Mldlt.compute(data.
M);
83 std::cout <<
"Dense Eigen Cholesky = \t" << (total/NBT)
91 std::cout <<
"Cholesky solve vector = \t\t"; timer.
toc(std::cout,NBT);
98 std::cout <<
"UDUtv = \t\t"; timer.
toc(std::cout,NBT);
100 MatrixXd Minv(model.
nv,model.
nv); Minv.setZero();
106 std::cout <<
"Minv from cholesky = \t\t"; timer.
toc(std::cout,NBT);
113 std::cout <<
"Cholesky solve column = \t\t"; timer.
toc(std::cout,NBT);
118 lhs[_smooth].noalias() = Minv*rhs[_smooth];
120 std::cout <<
"Minv*v = \t\t"; timer.
toc(std::cout,NBT);
125 A.noalias() = Minv*B;
127 std::cout <<
"A = Minv*B = \t\t"; timer.
toc(std::cout,NBT);
129 data.
M.triangularView<Eigen::StrictlyLower>()
130 = data.
M.transpose().triangularView<Eigen::StrictlyLower>();
134 A.noalias() = data.
M.inverse();
136 std::cout <<
"M.inverse() = \t\t"; timer.
toc(std::cout,NBT);
143 std::cout <<
"computeMinverse = \t\t"; timer.
toc(std::cout,NBT);
const DataTpl< Scalar, Options, JointCollectionTpl >::MatrixXs & crba(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q)
Computes the upper triangular part of the joint space inertia matrix M by using the Composite Rigid B...
Mat & computeMinv(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< Mat > &Minv)
Computes the inverse of the joint space inertia matrix M from its Cholesky factorization.
const DataTpl< Scalar, Options, JointCollectionTpl >::RowMatrixXs & computeMinverse(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q)
Computes the inverse of the joint space inertia matrix using Articulated Body formulation.
const DataTpl< Scalar, Options, JointCollectionTpl >::MatrixXs & decompose(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data)
Compute the Cholesky decomposition of the joint space inertia matrix M contained in data...
static std::string unitName(Unit u)
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.
#define PINOCCHIO_MODEL_DIR
int main(int argc, const char **argv)
MatrixXs M
The joint space inertia matrix (a square matrix of dim model.nv).
#define PINOCCHIO_ALIGNED_STD_VECTOR(Type)
ModelTpl< Scalar, Options, JointCollectionTpl > & buildModel(const std::string &filename, const typename ModelTpl< Scalar, Options, JointCollectionTpl >::JointModel &rootJoint, ModelTpl< Scalar, Options, JointCollectionTpl > &model, const bool verbose=false)
Build the model from a URDF file with a particular joint as root of the model tree inside the model g...
Main pinocchio namespace.
Mat & UDUtv(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< Mat > &m)
Performs the multiplication by using the Cholesky decomposition of M stored in data.
int nv
Dimension of the velocity vector space.
Mat & solve(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< Mat > &y)
Return the solution of using the Cholesky decomposition stored in data given the entry ...
JointModelFreeFlyerTpl< double > JointModelFreeFlyer
void humanoidRandom(ModelTpl< Scalar, Options, JointCollectionTpl > &model, bool usingFF=true)
Create a humanoid kinematic tree with 6-DOF limbs and random joint placements.
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > VectorXd
int nq
Dimension of the configuration vector representation.