pinocchio  3.7.0
A fast and flexible implementation of Rigid Body Dynamics algorithms and their analytical derivatives
 
Loading...
Searching...
No Matches
code-generator-base.hpp
1//
2// Copyright (c) 2018-2020 CNRS INRIA
3//
4
5#ifndef __pinocchio_utils_code_generator_base_hpp__
6#define __pinocchio_utils_code_generator_base_hpp__
7
8#include "pinocchio/codegen/cppadcg.hpp"
9
10#include "pinocchio/multibody/model.hpp"
11#include "pinocchio/multibody/data.hpp"
12
13namespace pinocchio
14{
15
16 template<typename _Scalar>
18 {
19 typedef _Scalar Scalar;
20 typedef CppAD::cg::CG<Scalar> CGScalar;
21 typedef CppAD::AD<CGScalar> ADScalar;
22
23 enum
24 {
25 Options = 0
26 };
27
34
35 typedef Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic, Options> MatrixXs;
36 typedef Eigen::Matrix<Scalar, Eigen::Dynamic, 1, Options> VectorXs;
37 typedef Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic, Options | Eigen::RowMajor>
38 RowMatrixXs;
39 typedef Eigen::Matrix<ADScalar, Eigen::Dynamic, 1, Options> ADVectorXs;
40 typedef Eigen::Matrix<ADScalar, Eigen::Dynamic, Eigen::Dynamic, Options> ADMatrixXs;
41
42 typedef typename Model::ConfigVectorType ConfigVectorType;
43 typedef typename Model::TangentVectorType TangentVectorType;
44
45 typedef typename ADModel::ConfigVectorType ADConfigVectorType;
46 typedef typename ADModel::TangentVectorType ADTangentVectorType;
47
48 typedef CppAD::ADFun<CGScalar> ADFun;
49
51 const Model & model,
52 const Eigen::DenseIndex dim_input,
53 const Eigen::DenseIndex dim_output,
54 const std::string & function_name,
55 const std::string & library_name)
56 : ad_model(model.template cast<ADScalar>())
57 , ad_data(ad_model)
59 , library_name(library_name + "_" + model.name)
60 , build_forward(true)
61 , build_jacobian(true)
62 {
63 ad_X = ADVectorXs(dim_input);
64 ad_Y = ADVectorXs(dim_output);
65
66 y = VectorXs(ad_Y.size());
67
68 jac = RowMatrixXs(ad_Y.size(), ad_X.size());
69 }
70
71 virtual ~CodeGenBase()
72 {
73 }
74
76 virtual void buildMap() = 0;
77
78 void initLib()
79 {
80 buildMap();
81
82 // generates source code
83 cgen_ptr = std::unique_ptr<CppAD::cg::ModelCSourceGen<Scalar>>(
84 new CppAD::cg::ModelCSourceGen<Scalar>(ad_fun, function_name));
85 cgen_ptr->setCreateForwardZero(build_forward);
86 cgen_ptr->setCreateJacobian(build_jacobian);
87 libcgen_ptr = std::unique_ptr<CppAD::cg::ModelLibraryCSourceGen<Scalar>>(
88 new CppAD::cg::ModelLibraryCSourceGen<Scalar>(*cgen_ptr));
89
90 dynamicLibManager_ptr = std::unique_ptr<CppAD::cg::DynamicModelLibraryProcessor<Scalar>>(
91 new CppAD::cg::DynamicModelLibraryProcessor<Scalar>(*libcgen_ptr, library_name));
92 }
93
94 CppAD::cg::ModelCSourceGen<Scalar> & codeGenerator()
95 {
96 return *cgen_ptr;
97 }
98
99 void compileLib(
100 const std::string & gcc_path = "/usr/bin/gcc", const std::string & compile_options = "-Ofast")
101 {
102 CppAD::cg::GccCompiler<Scalar> compiler(gcc_path);
103 std::vector<std::string> compile_flags = compiler.getCompileFlags();
104 compile_flags[0] = compile_options;
105 compiler.setCompileFlags(compile_flags);
106 dynamicLibManager_ptr->createDynamicLibrary(compiler, false);
107 }
108
109 bool existLib() const
110 {
111 const std::string filename = dynamicLibManager_ptr->getLibraryName()
112 + CppAD::cg::system::SystemInfo<>::DYNAMIC_LIB_EXTENSION;
113 std::ifstream file(filename.c_str());
114 return file.good();
115 }
116
117 void compileAndLoadLib(const std::string & gcc_path)
118 {
119 compileLib(gcc_path);
120 loadLib(false);
121 }
122
123 void loadLib(
124 const bool generate_if_not_exist = true,
125 const std::string & gcc_path = "/usr/bin/gcc",
126 const std::string & compile_options = "-Ofast")
127 {
128 if (!existLib() && generate_if_not_exist)
129 compileLib(gcc_path, compile_options);
130
131 const auto it = dynamicLibManager_ptr->getOptions().find("dlOpenMode");
132 if (it == dynamicLibManager_ptr->getOptions().end())
133 {
134 dynamicLib_ptr.reset(new CppAD::cg::LinuxDynamicLib<Scalar>(
135 dynamicLibManager_ptr->getLibraryName()
136 + CppAD::cg::system::SystemInfo<>::DYNAMIC_LIB_EXTENSION));
137 }
138 else
139 {
140 int dlOpenMode = std::stoi(it->second);
141 dynamicLib_ptr.reset(new CppAD::cg::LinuxDynamicLib<Scalar>(
142 dynamicLibManager_ptr->getLibraryName()
143 + CppAD::cg::system::SystemInfo<>::DYNAMIC_LIB_EXTENSION,
144 dlOpenMode));
145 }
146
147 generatedFun_ptr = dynamicLib_ptr->model(function_name.c_str());
148 }
149
150 template<typename Vector>
151 void evalFunction(const Eigen::MatrixBase<Vector> & x)
152 {
153 assert(build_forward);
154
155 generatedFun_ptr->ForwardZero(PINOCCHIO_EIGEN_CONST_CAST(Vector, x), y);
156 }
157
158 template<typename Vector>
159 void evalJacobian(const Eigen::MatrixBase<Vector> & x)
160 {
161 assert(build_jacobian);
162
163 CppAD::cg::ArrayView<const Scalar> x_(
164 PINOCCHIO_EIGEN_CONST_CAST(Vector, x).data(), (size_t)x.size());
165 CppAD::cg::ArrayView<Scalar> jac_(jac.data(), (size_t)jac.size());
166 generatedFun_ptr->Jacobian(x_, jac_);
167 }
168
170 Eigen::DenseIndex getInputDimension() const
171 {
172 return ad_X.size();
173 }
175 Eigen::DenseIndex getOutputDimension() const
176 {
177 return ad_Y.size();
178 }
179
180 protected:
181 ADModel ad_model;
182 ADData ad_data;
183
185 const std::string function_name;
187 const std::string library_name;
188
191
194
195 ADVectorXs ad_X, ad_Y;
196 ADFun ad_fun;
197
198 ADConfigVectorType ad_q, ad_q_plus;
199 ADTangentVectorType ad_dq, ad_v, ad_a;
200
201 VectorXs y;
202 RowMatrixXs jac;
203
204 std::unique_ptr<CppAD::cg::ModelCSourceGen<Scalar>> cgen_ptr;
205 std::unique_ptr<CppAD::cg::ModelLibraryCSourceGen<Scalar>> libcgen_ptr;
206 std::unique_ptr<CppAD::cg::DynamicModelLibraryProcessor<Scalar>> dynamicLibManager_ptr;
207 std::unique_ptr<CppAD::cg::DynamicLib<Scalar>> dynamicLib_ptr;
208 std::unique_ptr<CppAD::cg::GenericModel<Scalar>> generatedFun_ptr;
209
210 }; // struct CodeGenBase
211
212} // namespace pinocchio
213
214#endif // ifndef __pinocchio_utils_code_generator_base_hpp__
Main pinocchio namespace.
Definition treeview.dox:11
Eigen::DenseIndex getOutputDimension() const
Dimension of the output vector.
Eigen::DenseIndex getInputDimension() const
Dimension of the input vector.
virtual void buildMap()=0
build the mapping Y = f(X)
bool build_jacobian
Options to build or not the Jacobian of he function.
const std::string function_name
Name of the function.
bool build_forward
Options to generate or not the source code for the evaluation function.
const std::string library_name
Name of the library.
VectorXs TangentVectorType
Dense vectorized version of a joint tangent vector (e.g. velocity, acceleration, etc)....
Definition model.hpp:95
VectorXs ConfigVectorType
Dense vectorized version of a joint configuration vector.
Definition model.hpp:87