cppadcg-algo.cpp
Go to the documentation of this file.
1 //
2 // Copyright (c) 2018-2019 CNRS INRIA
3 //
4 
5 #include "pinocchio/codegen/cppadcg.hpp"
6 
7 #include "pinocchio/multibody/model.hpp"
8 #include "pinocchio/multibody/data.hpp"
9 
10 #include "pinocchio/algorithm/kinematics.hpp"
11 #include "pinocchio/algorithm/jacobian.hpp"
12 #include "pinocchio/algorithm/crba.hpp"
13 #include "pinocchio/algorithm/rnea.hpp"
14 #include "pinocchio/algorithm/aba.hpp"
15 #include "pinocchio/algorithm/joint-configuration.hpp"
16 
17 #include "pinocchio/parsers/sample-models.hpp"
18 
19 #include <iostream>
20 
21 #include <boost/test/unit_test.hpp>
22 #include <boost/utility/binary.hpp>
23 
24 BOOST_AUTO_TEST_SUITE(BOOST_TEST_MODULE)
25 
26  BOOST_AUTO_TEST_CASE(test_crba_code_generation)
27  {
28  typedef double Scalar;
29  typedef CppAD::cg::CG<Scalar> CGScalar;
30  typedef CppAD::AD<CGScalar> ADScalar;
31 
32  typedef Eigen::Matrix<ADScalar,Eigen::Dynamic,1> ADVector;
33 
35  typedef Model::Data Data;
36 
37  typedef pinocchio::ModelTpl<ADScalar> ADModel;
38  typedef ADModel::Data ADData;
39 
40  Model model;
42  model.lowerPositionLimit.head<3>().fill(-1.);
43  model.upperPositionLimit.head<3>().fill(1.);
44  Data data(model);
45 
46  ADModel ad_model = model.cast<ADScalar>();
47  ADData ad_data(ad_model);
48 
49  // Sample random configuration
50  typedef Model::ConfigVectorType ConfigVectorType;
51  typedef Model::TangentVectorType TangentVectorType;
52  ConfigVectorType q(model.nq);
54 
55  TangentVectorType v(TangentVectorType::Random(model.nv));
56  TangentVectorType a(TangentVectorType::Random(model.nv));
57 
58  typedef ADModel::ConfigVectorType ADConfigVectorType;
59  typedef ADModel::TangentVectorType ADTangentVectorType;
60 
61  ADConfigVectorType ad_q = q.cast<ADScalar>();
62  ADTangentVectorType ad_v = v.cast<ADScalar>();
63  ADTangentVectorType ad_a = a.cast<ADScalar>();
64 
65  ADTangentVectorType & X = ad_a;
66  CppAD::Independent(X);
67 
68  pinocchio::rnea(ad_model,ad_data,ad_q,ad_v,ad_a);
69  ADVector Y(model.nv); Y = ad_data.tau;
70 
71  CppAD::ADFun<CGScalar> fun(X,Y);
72 
73  // generates source code
74  CppAD::cg::ModelCSourceGen<Scalar> cgen(fun, "rnea");
75  cgen.setCreateJacobian(true);
76  cgen.setCreateForwardZero(true);
77  cgen.setCreateForwardOne(true);
78  cgen.setCreateReverseOne(true);
79  cgen.setCreateReverseTwo(true);
80  CppAD::cg::ModelLibraryCSourceGen<Scalar> libcgen(cgen);
81 
82  // compile source code
83  CppAD::cg::DynamicModelLibraryProcessor<Scalar> p(libcgen);
84 
85  CppAD::cg::GccCompiler<Scalar> compiler;
86  std::unique_ptr<CppAD::cg::DynamicLib<Scalar>> dynamicLib = p.createDynamicLibrary(compiler);
87 
88  // save to files (not really required)
89  CppAD::cg::SaveFilesModelLibraryProcessor<Scalar> p2(libcgen);
90  p2.saveSources();
91 
92  // use the generated code
93  std::unique_ptr<CppAD::cg::GenericModel<Scalar> > rnea_generated = dynamicLib->model("rnea");
94 
95  CPPAD_TESTVECTOR(Scalar) x((size_t)model.nv);
96  Eigen::Map<TangentVectorType>(x.data(),model.nv,1) = a;
97 
98  CPPAD_TESTVECTOR(Scalar) tau = rnea_generated->ForwardZero(x);
99 
100  Eigen::Map<TangentVectorType> tau_map(tau.data(),model.nv,1);
101  Data::TangentVectorType tau_ref = pinocchio::rnea(model,data,q,v,a);
102  BOOST_CHECK(tau_map.isApprox(tau_ref));
103 
104  pinocchio::crba(model,data,q);
105  data.M.triangularView<Eigen::StrictlyLower>()
106  = data.M.transpose().triangularView<Eigen::StrictlyLower>();
107 
108  CPPAD_TESTVECTOR(Scalar) dtau_da = rnea_generated->Jacobian(x);
109  Eigen::Map<PINOCCHIO_EIGEN_PLAIN_ROW_MAJOR_TYPE(Data::MatrixXs)> M_map(dtau_da.data(),model.nv,model.nv);
110  BOOST_CHECK(M_map.isApprox(data.M));
111  }
112 
113 BOOST_AUTO_TEST_SUITE_END()
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...
ModelTpl< double > Model
q
void rnea(const int num_threads, ModelPoolTpl< Scalar, Options, JointCollectionTpl > &pool, const Eigen::MatrixBase< ConfigVectorPool > &q, const Eigen::MatrixBase< TangentVectorPool1 > &v, const Eigen::MatrixBase< TangentVectorPool2 > &a, const Eigen::MatrixBase< TangentVectorPool3 > &tau)
The Recursive Newton-Euler algorithm. It computes the inverse dynamics, aka the joint torques accordi...
v
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.
data
SE3::Scalar Scalar
Definition: conversions.cpp:13
DataTpl< double > Data
list a
x
— Training
Definition: continuous.py:157
BOOST_AUTO_TEST_CASE(test_crba_code_generation)
void humanoidRandom(ModelTpl< Scalar, Options, JointCollectionTpl > &model, bool usingFF=true)
Create a humanoid kinematic tree with 6-DOF limbs and random joint placements.
void fill(Eigen::Ref< MatType > mat, const typename MatType::Scalar &value)


pinocchio
Author(s):
autogenerated on Tue Jun 1 2021 02:45:02