casadi-crba.cpp
Go to the documentation of this file.
2 
4 
7 
8 int main(int /*argc*/, char ** /*argv*/)
9 {
10  using namespace pinocchio;
11 
12  typedef double Scalar;
13  typedef ::casadi::SX ADScalar;
14 
15  typedef ModelTpl<Scalar> Model;
16  typedef Model::Data Data;
17 
19  typedef ADModel::Data ADData;
20 
21  // Create a random humanoid model
22  Model model;
24  model.lowerPositionLimit.head<3>().fill(-1.);
25  model.upperPositionLimit.head<3>().fill(1.);
26  Data data(model);
27 
28  // Pick up random configuration, velocity and acceleration vectors.
29  Eigen::VectorXd q(model.nq);
31 
32  // Create CasADi model and data from model
33  typedef ADModel::ConfigVectorType ConfigVectorAD;
34  ADModel ad_model = model.cast<ADScalar>();
35  ADData ad_data(ad_model);
36 
37  // Create symbolic CasADi vectors
38  ::casadi::SX cs_q = ::casadi::SX::sym("q", model.nq);
39  ConfigVectorAD q_ad(model.nq);
40  q_ad = Eigen::Map<ConfigVectorAD>(static_cast<std::vector<ADScalar>>(cs_q).data(), model.nq, 1);
41 
42  // Build CasADi function
43  crba(ad_model, ad_data, q_ad, pinocchio::Convention::WORLD);
44  ad_data.M.triangularView<Eigen::StrictlyLower>() =
45  ad_data.M.transpose().triangularView<Eigen::StrictlyLower>();
46  ::casadi::SX M_ad(model.nv, model.nv);
47  for (Eigen::DenseIndex j = 0; j < model.nv; ++j)
48  {
49  for (Eigen::DenseIndex i = 0; i < model.nv; ++i)
50  {
51  M_ad(i, j) = ad_data.M(i, j);
52  }
53  }
54 
55  ::casadi::Function eval_crba("eval_crba", ::casadi::SXVector{cs_q}, ::casadi::SXVector{M_ad});
56 
57  // Evaluate CasADi expression with real value
58  std::vector<double> q_vec((size_t)model.nq);
59  Eigen::Map<Eigen::VectorXd>(q_vec.data(), model.nq, 1) = q;
60 
61  ::casadi::DM M_res = eval_crba(::casadi::DMVector{q_vec})[0];
62  Data::MatrixXs M_mat =
63  Eigen::Map<Data::MatrixXs>(static_cast<std::vector<double>>(M_res).data(), model.nv, model.nv);
64 
65  // Eval CRBA using classic Pinocchio model
67  data.M.triangularView<Eigen::StrictlyLower>() =
68  data.M.transpose().triangularView<Eigen::StrictlyLower>();
69 
70  // Print both results
71  std::cout << "pinocchio double:\n" << "\tM =\n" << data.M << std::endl;
72  std::cout << "pinocchio CasADi:\n" << "\tM =\n" << M_mat.transpose() << std::endl;
73 }
pinocchio::DataTpl::MatrixXs
Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::Dynamic, Options > MatrixXs
Definition: multibody/data.hpp:74
pinocchio::DataTpl
Definition: context/generic.hpp:25
pinocchio.casadi::sym
void sym(const Eigen::MatrixBase< MatrixDerived > &eig_mat, std::string const &name)
Definition: autodiff/casadi.hpp:230
pinocchio::Convention::WORLD
@ WORLD
inverse-kinematics.i
int i
Definition: inverse-kinematics.py:17
setup.data
data
Definition: cmake/cython/setup.in.py:48
pinocchio::crba
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...
pinocchio::buildModels::humanoidRandom
void humanoidRandom(ModelTpl< Scalar, Options, JointCollectionTpl > &model, bool usingFF=true, bool mimic=false)
Create a humanoid kinematic tree with 6-DOF limbs and random joint placements.
pinocchio::randomConfiguration
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.
Definition: joint-configuration.hpp:315
joint-configuration.hpp
pinocchio::Data
DataTpl< context::Scalar, context::Options > Data
Definition: multibody/fwd.hpp:34
pinocchio::q
JointCollectionTpl const Eigen::MatrixBase< ConfigVectorType > & q
Definition: joint-configuration.hpp:1083
ADScalar
CppAD::AD< Scalar > ADScalar
Definition: timings-cppad-jit.cpp:38
ADModel
pinocchio::ModelTpl< ADScalar > ADModel
Definition: timings-cppad-jit.cpp:51
fill
fill
casadi.hpp
pinocchio::ModelTpl< ADScalar >::ConfigVectorType
VectorXs ConfigVectorType
Dense vectorized version of a joint configuration vector.
Definition: multibody/model.hpp:87
sample-models.hpp
Scalar
double Scalar
Definition: timings-cppad-jit.cpp:37
pinocchio::ModelTpl< Scalar >
ADData
ADModel::Data ADData
Definition: timings-cppad-jit.cpp:52
crba.hpp
pinocchio::Model
ModelTpl< context::Scalar, context::Options > Model
Definition: multibody/fwd.hpp:33
pinocchio::model
JointCollectionTpl & model
Definition: joint-configuration.hpp:1082
pinocchio
Main pinocchio namespace.
Definition: timings.cpp:33
main
int main(int, char **)
Definition: casadi-crba.cpp:8


pinocchio
Author(s):
autogenerated on Wed May 28 2025 02:41:15