timings-cholesky.cpp
Go to the documentation of this file.
1 //
2 // Copyright (c) 2018 CNRS
3 //
4 
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"
16 
17 #include <iostream>
18 
19 #include "pinocchio/utils/timer.hpp"
20 
21 
22 int main(int argc, const char ** argv)
23 {
24  using namespace Eigen;
25  using namespace pinocchio;
26 
28  #ifdef NDEBUG
29  const int NBT = 1000*100;
30  #else
31  const int NBT = 1;
32  std::cout << "(the time score in debug mode is not relevant) " << std::endl;
33  #endif
34 
36 
37  std::string filename = PINOCCHIO_MODEL_DIR + std::string("/simple_humanoid.urdf");
38  if(argc>1) filename = argv[1];
39  if( filename == "HS")
41  else
43  std::cout << "nq = " << model.nq << std::endl;
44 
45  pinocchio::Data data(model);
46  VectorXd qmax = Eigen::VectorXd::Ones(model.nq);
47 
48  MatrixXd A(model.nv,model.nv), B(model.nv,model.nv);
49  A.setZero(); B.setRandom();
50 
54  for(size_t i=0;i<NBT;++i)
55  {
56  qs[i] = randomConfiguration(model,-qmax,qmax);
57  lhs[i] = Eigen::VectorXd::Zero(model.nv);
58  rhs[i] = Eigen::VectorXd::Random(model.nv);
59  }
60 
61  double total = 0;
62  SMOOTH(NBT)
63  {
64  crba(model,data,qs[_smooth]);
65  timer.tic();
66  cholesky::decompose(model,data);
67  total += timer.toc(timer.DEFAULT_UNIT);
68  }
69  std::cout << "Cholesky = \t" << (total/NBT)
70  << " " << timer.unitName(timer.DEFAULT_UNIT) <<std::endl;
71 
72  total = 0;
73  Eigen::LDLT<Eigen::MatrixXd> Mldlt(data.M);
74  SMOOTH(NBT)
75  {
76  crba(model,data,qs[_smooth]);
77  data.M.triangularView<Eigen::StrictlyLower>()
78  = data.M.transpose().triangularView<Eigen::StrictlyLower>();
79  timer.tic();
80  Mldlt.compute(data.M);
81  total += timer.toc(timer.DEFAULT_UNIT);
82  }
83  std::cout << "Dense Eigen Cholesky = \t" << (total/NBT)
84  << " " << timer.unitName(timer.DEFAULT_UNIT) <<std::endl;
85 
86  timer.tic();
87  SMOOTH(NBT)
88  {
89  cholesky::solve(model,data,rhs[_smooth]);
90  }
91  std::cout << "Cholesky solve vector = \t\t"; timer.toc(std::cout,NBT);
92 
93  timer.tic();
94  SMOOTH(NBT)
95  {
96  cholesky::UDUtv(model,data,rhs[_smooth]);
97  }
98  std::cout << "UDUtv = \t\t"; timer.toc(std::cout,NBT);
99 
100  MatrixXd Minv(model.nv,model.nv); Minv.setZero();
101  timer.tic();
102  SMOOTH(NBT)
103  {
104  cholesky::computeMinv(model,data,Minv);
105  }
106  std::cout << "Minv from cholesky = \t\t"; timer.toc(std::cout,NBT);
107 
108  timer.tic();
109  SMOOTH(NBT)
110  {
111  cholesky::solve(model,data,Minv.col(10));
112  }
113  std::cout << "Cholesky solve column = \t\t"; timer.toc(std::cout,NBT);
114 
115  timer.tic();
116  SMOOTH(NBT)
117  {
118  lhs[_smooth].noalias() = Minv*rhs[_smooth];
119  }
120  std::cout << "Minv*v = \t\t"; timer.toc(std::cout,NBT);
121 
122  timer.tic();
123  SMOOTH(NBT)
124  {
125  A.noalias() = Minv*B;
126  }
127  std::cout << "A = Minv*B = \t\t"; timer.toc(std::cout,NBT);
128 
129  data.M.triangularView<Eigen::StrictlyLower>()
130  = data.M.transpose().triangularView<Eigen::StrictlyLower>();
131  timer.tic();
132  SMOOTH(NBT)
133  {
134  A.noalias() = data.M.inverse();
135  }
136  std::cout << "M.inverse() = \t\t"; timer.toc(std::cout,NBT);
137 
138  timer.tic();
139  SMOOTH(NBT)
140  {
141  computeMinverse(model,data,qs[_smooth]);
142  }
143  std::cout << "computeMinverse = \t\t"; timer.toc(std::cout,NBT);
144 
145  return 0;
146 }
double toc()
Definition: timer.hpp:68
Unit DEFAULT_UNIT
Definition: timer.hpp:50
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)
Definition: timer.hpp:52
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
void tic()
Definition: timer.hpp:63
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)
#define SMOOTH(s)
Definition: timer.hpp:38
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.
Definition: timings.cpp:28
B
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 ...
filename
JointModelFreeFlyerTpl< double > JointModelFreeFlyer
A
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
Definition: conversions.cpp:14
int nq
Dimension of the configuration vector representation.


pinocchio
Author(s):
autogenerated on Fri Jun 23 2023 02:38:33