5 #ifndef __pinocchio_utils_code_generator_algo_hpp__ 6 #define __pinocchio_utils_code_generator_algo_hpp__ 8 #ifdef PINOCCHIO_WITH_CPPADCG_SUPPORT 10 #include "pinocchio/codegen/code-generator-base.hpp" 11 #include "pinocchio/algorithm/crba.hpp" 12 #include "pinocchio/algorithm/rnea.hpp" 13 #include "pinocchio/algorithm/aba.hpp" 14 #include "pinocchio/algorithm/rnea-derivatives.hpp" 15 #include "pinocchio/algorithm/aba-derivatives.hpp" 19 template<
typename _Scalar>
23 typedef typename Base::Scalar Scalar;
26 typedef typename Base::ADCongigVectorType ADCongigVectorType;
27 typedef typename Base::ADTangentVectorType ADTangentVectorType;
28 typedef typename Base::MatrixXs MatrixXs;
29 typedef typename Base::VectorXs VectorXs;
37 ad_v = ADTangentVectorType(model.
nv); ad_v.setZero();
38 ad_a = ADTangentVectorType(model.
nv); ad_a.setZero();
42 dtau_dq = MatrixXs::Zero(model.
nv,model.
nq);
43 dtau_dv = MatrixXs::Zero(model.
nv,model.
nv);
44 dtau_da = MatrixXs::Zero(model.
nv,model.
nv);
49 CppAD::Independent(ad_X);
51 Eigen::DenseIndex it = 0;
52 ad_q = ad_X.segment(it,ad_model.
nq); it += ad_model.
nq;
53 ad_v = ad_X.segment(it,ad_model.
nv); it += ad_model.
nv;
54 ad_a = ad_X.segment(it,ad_model.
nv); it += ad_model.
nv;
60 ad_fun.Dependent(ad_X,ad_Y);
61 ad_fun.optimize(
"no_compare_op");
64 using Base::evalFunction;
65 template<
typename ConfigVectorType,
typename TangentVector1,
typename TangentVector2>
66 void evalFunction(
const Eigen::MatrixBase<ConfigVectorType> & q,
67 const Eigen::MatrixBase<TangentVector1> & v,
68 const Eigen::MatrixBase<TangentVector2> & a)
71 Eigen::DenseIndex it = 0;
72 x.segment(it,ad_model.
nq) = q; it += ad_model.
nq;
73 x.segment(it,ad_model.
nv) = v; it += ad_model.
nv;
74 x.segment(it,ad_model.
nv) = a; it += ad_model.
nv;
80 using Base::evalJacobian;
81 template<
typename ConfigVectorType,
typename TangentVector1,
typename TangentVector2>
82 void evalJacobian(
const Eigen::MatrixBase<ConfigVectorType> & q,
83 const Eigen::MatrixBase<TangentVector1> & v,
84 const Eigen::MatrixBase<TangentVector2> & a)
87 Eigen::DenseIndex it = 0;
88 x.segment(it,ad_model.
nq) = q; it += ad_model.
nq;
89 x.segment(it,ad_model.
nv) = v; it += ad_model.
nv;
90 x.segment(it,ad_model.
nv) = a; it += ad_model.
nv;
94 dtau_dq = Base::jac.middleCols(it,ad_model.
nq); it += ad_model.
nq;
95 dtau_dv = Base::jac.middleCols(it,ad_model.
nv); it += ad_model.
nv;
96 dtau_da = Base::jac.middleCols(it,ad_model.
nv); it += ad_model.
nv;
101 using Base::ad_model;
111 MatrixXs dtau_dq, dtau_dv, dtau_da;
113 ADCongigVectorType ad_q, ad_q_plus;
114 ADTangentVectorType ad_dq, ad_v, ad_a;
117 template<
typename _Scalar>
121 typedef typename Base::Scalar Scalar;
124 typedef typename Base::ADCongigVectorType ADCongigVectorType;
125 typedef typename Base::ADTangentVectorType ADTangentVectorType;
126 typedef typename Base::MatrixXs MatrixXs;
127 typedef typename Base::VectorXs VectorXs;
134 ad_q = ADCongigVectorType(model.
nq); ad_q = ad_model.neutralConfiguration;
135 ad_v = ADTangentVectorType(model.
nv); ad_v.setZero();
136 ad_tau = ADTangentVectorType(model.
nv); ad_tau.setZero();
140 da_dq = MatrixXs::Zero(model.
nv,model.
nq);
141 da_dv = MatrixXs::Zero(model.
nv,model.
nv);
142 da_dtau = MatrixXs::Zero(model.
nv,model.
nv);
147 CppAD::Independent(ad_X);
149 Eigen::DenseIndex it = 0;
150 ad_q = ad_X.segment(it,ad_model.nq); it += ad_model.nq;
151 ad_v = ad_X.segment(it,ad_model.nv); it += ad_model.nv;
152 ad_tau = ad_X.segment(it,ad_model.nv); it += ad_model.nv;
157 ad_fun.Dependent(ad_X,ad_Y);
158 ad_fun.optimize(
"no_compare_op");
161 using Base::evalFunction;
162 template<
typename ConfigVectorType,
typename TangentVector1,
typename TangentVector2>
163 void evalFunction(
const Eigen::MatrixBase<ConfigVectorType> & q,
164 const Eigen::MatrixBase<TangentVector1> & v,
165 const Eigen::MatrixBase<TangentVector2> & tau)
168 Eigen::DenseIndex it = 0;
169 x.segment(it,ad_model.nq) = q; it += ad_model.nq;
170 x.segment(it,ad_model.nv) = v; it += ad_model.nv;
171 x.segment(it,ad_model.nv) = tau; it += ad_model.nv;
177 using Base::evalJacobian;
178 template<
typename ConfigVectorType,
typename TangentVector1,
typename TangentVector2>
179 void evalJacobian(
const Eigen::MatrixBase<ConfigVectorType> & q,
180 const Eigen::MatrixBase<TangentVector1> & v,
181 const Eigen::MatrixBase<TangentVector2> & tau)
184 Eigen::DenseIndex it = 0;
185 x.segment(it,ad_model.nq) = q; it += ad_model.nq;
186 x.segment(it,ad_model.nv) = v; it += ad_model.nv;
187 x.segment(it,ad_model.nv) = tau; it += ad_model.nv;
192 da_dq = Base::jac.middleCols(it,ad_model.nq); it += ad_model.nq;
193 da_dv = Base::jac.middleCols(it,ad_model.nv); it += ad_model.nv;
194 da_dtau = Base::jac.middleCols(it,ad_model.nv); it += ad_model.nv;
199 using Base::ad_model;
209 MatrixXs da_dq,da_dv,da_dtau;
211 ADCongigVectorType ad_q, ad_q_plus;
212 ADTangentVectorType ad_dq, ad_v, ad_tau;
215 template<
typename _Scalar>
219 typedef typename Base::Scalar Scalar;
222 typedef typename Base::ADCongigVectorType ADCongigVectorType;
223 typedef typename Base::ADTangentVectorType ADTangentVectorType;
224 typedef typename Base::MatrixXs MatrixXs;
225 typedef typename Base::VectorXs VectorXs;
232 ad_q = ADCongigVectorType(model.
nq); ad_q = ad_model.neutralConfiguration;
236 M = MatrixXs::Zero(model.
nv,model.
nq);
243 CppAD::Independent(ad_X);
245 Eigen::DenseIndex it = 0;
246 ad_q = ad_X.segment(it,ad_model.nq); it += ad_model.nq;
249 Eigen::DenseIndex it_Y = 0;
251 for(Eigen::DenseIndex i = 0; i < ad_model.nv; ++i)
253 for(Eigen::DenseIndex j = i; j < ad_model.nv; ++j)
255 ad_Y[it_Y++] = ad_data.M(i,j);
261 ad_fun.Dependent(ad_X,ad_Y);
262 ad_fun.optimize(
"no_compare_op");
265 template<
typename ConfigVectorType>
266 void evalFunction(
const Eigen::MatrixBase<ConfigVectorType> & q)
269 Eigen::DenseIndex it = 0;
270 x.segment(it,ad_model.nq) = q; it += ad_model.nq;
272 Base::evalFunction(x);
275 Eigen::DenseIndex it_Y = 0;
276 for(Eigen::DenseIndex i = 0; i < ad_model.nv; ++i)
278 for(Eigen::DenseIndex j = i; j < ad_model.nv; ++j)
280 M(i,j) = Base::y[it_Y++];
289 using Base::ad_model;
300 ADCongigVectorType ad_q;
303 template<
typename _Scalar>
307 typedef typename Base::Scalar Scalar;
310 typedef typename Base::ADCongigVectorType ADCongigVectorType;
311 typedef typename Base::ADTangentVectorType ADTangentVectorType;
312 typedef typename Base::MatrixXs MatrixXs;
313 typedef typename Base::VectorXs VectorXs;
320 ad_q = ADCongigVectorType(model.
nq); ad_q = ad_model.neutralConfiguration;
324 Minv = MatrixXs::Zero(model.
nv,model.
nq);
331 CppAD::Independent(ad_X);
333 Eigen::DenseIndex it = 0;
334 ad_q = ad_X.segment(it,ad_model.nq); it += ad_model.nq;
337 Eigen::DenseIndex it_Y = 0;
338 for(Eigen::DenseIndex i = 0; i < ad_model.nv; ++i)
340 for(Eigen::DenseIndex j = i; j < ad_model.nv; ++j)
342 ad_Y[it_Y++] = ad_data.Minv(i,j);
348 ad_fun.Dependent(ad_X,ad_Y);
349 ad_fun.optimize(
"no_compare_op");
352 template<
typename ConfigVectorType>
353 void evalFunction(
const Eigen::MatrixBase<ConfigVectorType> & q)
356 Eigen::DenseIndex it = 0;
357 x.segment(it,ad_model.nq) = q; it += ad_model.nq;
359 Base::evalFunction(x);
362 Eigen::DenseIndex it_Y = 0;
363 for(Eigen::DenseIndex i = 0; i < ad_model.nv; ++i)
365 for(Eigen::DenseIndex j = i; j < ad_model.nv; ++j)
367 Minv(i,j) = Base::y[it_Y++];
374 using Base::ad_model;
385 ADCongigVectorType ad_q;
388 template<
typename _Scalar>
392 typedef typename Base::Scalar Scalar;
395 typedef typename Base::ADCongigVectorType ADCongigVectorType;
396 typedef typename Base::ADTangentVectorType ADTangentVectorType;
397 typedef typename Base::MatrixXs MatrixXs;
398 typedef typename PINOCCHIO_EIGEN_PLAIN_ROW_MAJOR_TYPE(MatrixXs) RowMatrixXs;
399 typedef typename Base::VectorXs VectorXs;
402 typedef typename ADData::MatrixXs ADMatrixXs;
403 typedef typename PINOCCHIO_EIGEN_PLAIN_ROW_MAJOR_TYPE(ADMatrixXs) RowADMatrixXs;
407 const std::string &
library_name =
"cg_partial_rnea_eval")
410 ad_q = ADCongigVectorType(model.
nq); ad_q = ad_model.neutralConfiguration;
411 ad_v = ADTangentVectorType(model.
nv); ad_v.setZero();
412 ad_a = ADTangentVectorType(model.
nv); ad_a.setZero();
417 ad_dtau_dq = ADMatrixXs::Zero(model.
nv,model.
nv);
418 ad_dtau_dv = ADMatrixXs::Zero(model.
nv,model.
nv);
419 ad_dtau_da = ADMatrixXs::Zero(model.
nv,model.
nv);
421 dtau_dq = MatrixXs::Zero(model.
nv,model.
nv);
422 dtau_dv = MatrixXs::Zero(model.
nv,model.
nv);
423 dtau_da = MatrixXs::Zero(model.
nv,model.
nv);
430 CppAD::Independent(ad_X);
432 Eigen::DenseIndex it = 0;
433 ad_q = ad_X.segment(it,ad_model.nq); it += ad_model.nq;
434 ad_v = ad_X.segment(it,ad_model.nv); it += ad_model.nv;
435 ad_a = ad_X.segment(it,ad_model.nv); it += ad_model.nv;
439 ad_dtau_dq,ad_dtau_dv,ad_dtau_da);
443 Eigen::DenseIndex it_Y = 0;
444 Eigen::Map<RowADMatrixXs>(ad_Y.data()+it_Y,ad_model.nv,ad_model.nv) = ad_dtau_dq;
445 it_Y += ad_model.nv*ad_model.nv;
446 Eigen::Map<RowADMatrixXs>(ad_Y.data()+it_Y,ad_model.nv,ad_model.nv) = ad_dtau_dv;
447 it_Y += ad_model.nv*ad_model.nv;
448 Eigen::Map<RowADMatrixXs>(ad_Y.data()+it_Y,ad_model.nv,ad_model.nv) = ad_dtau_da;
449 it_Y += ad_model.nv*ad_model.nv;
451 ad_fun.Dependent(ad_X,ad_Y);
452 ad_fun.optimize(
"no_compare_op");
455 template<
typename ConfigVectorType,
typename TangentVector1,
typename TangentVector2>
456 void evalFunction(
const Eigen::MatrixBase<ConfigVectorType> & q,
457 const Eigen::MatrixBase<TangentVector1> & v,
458 const Eigen::MatrixBase<TangentVector2> & a)
461 Eigen::DenseIndex it_x = 0;
462 x.segment(it_x,ad_model.nq) = q; it_x += ad_model.nq;
463 x.segment(it_x,ad_model.nq) = v; it_x += ad_model.nv;
464 x.segment(it_x,ad_model.nq) = a; it_x += ad_model.nv;
466 Base::evalFunction(x);
469 Eigen::DenseIndex it_y = 0;
470 dtau_dq = Eigen::Map<RowMatrixXs>(Base::y.data()+it_y,ad_model.nv,ad_model.nv);
471 it_y += ad_model.nv*ad_model.nv;
472 dtau_dv = Eigen::Map<RowMatrixXs>(Base::y.data()+it_y,ad_model.nv,ad_model.nv);
473 it_y += ad_model.nv*ad_model.nv;
474 dtau_da = Eigen::Map<RowMatrixXs>(Base::y.data()+it_y,ad_model.nv,ad_model.nv);
475 it_y += ad_model.nv*ad_model.nv;
481 using Base::ad_model;
489 VectorXs partial_derivatives;
490 ADMatrixXs ad_dtau_dq, ad_dtau_dv, ad_dtau_da;
491 MatrixXs dtau_dq, dtau_dv, dtau_da;
493 ADCongigVectorType ad_q;
494 ADTangentVectorType ad_v, ad_a;
497 template<
typename _Scalar>
501 typedef typename Base::Scalar Scalar;
504 typedef typename Base::ADCongigVectorType ADCongigVectorType;
505 typedef typename Base::ADTangentVectorType ADTangentVectorType;
506 typedef typename Base::MatrixXs MatrixXs;
507 typedef typename PINOCCHIO_EIGEN_PLAIN_ROW_MAJOR_TYPE(MatrixXs) RowMatrixXs;
508 typedef typename Base::VectorXs VectorXs;
511 typedef typename ADData::MatrixXs ADMatrixXs;
512 typedef typename PINOCCHIO_EIGEN_PLAIN_ROW_MAJOR_TYPE(ADMatrixXs) RowADMatrixXs;
516 const std::string &
library_name =
"cg_partial_aba_eval")
519 ad_q = ADCongigVectorType(model.
nq); ad_q = ad_model.neutralConfiguration;
520 ad_v = ADTangentVectorType(model.
nv); ad_v.setZero();
521 ad_tau = ADTangentVectorType(model.
nv); ad_tau.setZero();
526 ad_dddq_dq = ADMatrixXs::Zero(model.
nv,model.
nv);
527 ad_dddq_dv = ADMatrixXs::Zero(model.
nv,model.
nv);
528 ad_dddq_dtau = ADMatrixXs::Zero(model.
nv,model.
nv);
530 dddq_dq = MatrixXs::Zero(model.
nv,model.
nv);
531 dddq_dv = MatrixXs::Zero(model.
nv,model.
nv);
532 dddq_dtau = MatrixXs::Zero(model.
nv,model.
nv);
539 CppAD::Independent(ad_X);
541 Eigen::DenseIndex it = 0;
542 ad_q = ad_X.segment(it,ad_model.nq); it += ad_model.nq;
543 ad_v = ad_X.segment(it,ad_model.nv); it += ad_model.nv;
544 ad_tau = ad_X.segment(it,ad_model.nv); it += ad_model.nv;
548 ad_dddq_dq,ad_dddq_dv,ad_dddq_dtau);
552 Eigen::DenseIndex it_Y = 0;
553 Eigen::Map<RowADMatrixXs>(ad_Y.data()+it_Y,ad_model.nv,ad_model.nv) = ad_dddq_dq;
554 it_Y += ad_model.nv*ad_model.nv;
555 Eigen::Map<RowADMatrixXs>(ad_Y.data()+it_Y,ad_model.nv,ad_model.nv) = ad_dddq_dv;
556 it_Y += ad_model.nv*ad_model.nv;
557 Eigen::Map<RowADMatrixXs>(ad_Y.data()+it_Y,ad_model.nv,ad_model.nv) = ad_dddq_dtau;
558 it_Y += ad_model.nv*ad_model.nv;
560 ad_fun.Dependent(ad_X,ad_Y);
561 ad_fun.optimize(
"no_compare_op");
564 template<
typename ConfigVectorType,
typename TangentVector1,
typename TangentVector2>
565 void evalFunction(
const Eigen::MatrixBase<ConfigVectorType> & q,
566 const Eigen::MatrixBase<TangentVector1> & v,
567 const Eigen::MatrixBase<TangentVector2> & tau)
570 Eigen::DenseIndex it_x = 0;
571 x.segment(it_x,ad_model.nq) = q; it_x += ad_model.nq;
572 x.segment(it_x,ad_model.nq) = v; it_x += ad_model.nv;
573 x.segment(it_x,ad_model.nq) = tau; it_x += ad_model.nv;
575 Base::evalFunction(x);
578 Eigen::DenseIndex it_y = 0;
579 dddq_dq = Eigen::Map<RowMatrixXs>(Base::y.data()+it_y,ad_model.nv,ad_model.nv);
580 it_y += ad_model.nv*ad_model.nv;
581 dddq_dv = Eigen::Map<RowMatrixXs>(Base::y.data()+it_y,ad_model.nv,ad_model.nv);
582 it_y += ad_model.nv*ad_model.nv;
583 dddq_dtau = Eigen::Map<RowMatrixXs>(Base::y.data()+it_y,ad_model.nv,ad_model.nv);
584 it_y += ad_model.nv*ad_model.nv;
590 using Base::ad_model;
598 VectorXs partial_derivatives;
599 ADMatrixXs ad_dddq_dq, ad_dddq_dv, ad_dddq_dtau;
600 MatrixXs dddq_dq, dddq_dv, dddq_dtau;
602 ADCongigVectorType ad_q;
603 ADTangentVectorType ad_v, ad_tau;
608 #endif // ifdef PINOCCHIO_WITH_CPPADCG_SUPPORT 610 #endif // ifndef __pinocchio_utils_code_generator_base_hpp__ const DataTpl< Scalar, Options, JointCollectionTpl >::TangentVectorType & aba(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType1 > &v, const Eigen::MatrixBase< TangentVectorType2 > &tau)
The Articulated-Body algorithm. It computes the forward dynamics, aka the joint accelerations given t...
PINOCCHIO_DEPRECATED ConfigVectorType neutralConfiguration
Vector of joint's neutral configurations.
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...
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 std::string function_name
Name of the function.
void buildMap()
build the mapping Y = f(X)
void buildMap()
build the mapping Y = f(X)
void buildMap()
build the mapping Y = f(X)
int nv
Dimension of the velocity vector space.
void computeRNEADerivatives(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType1 > &v, const Eigen::MatrixBase< TangentVectorType2 > &a, const Eigen::MatrixBase< MatrixType1 > &rnea_partial_dq, const Eigen::MatrixBase< MatrixType2 > &rnea_partial_dv, const Eigen::MatrixBase< MatrixType3 > &rnea_partial_da)
Computes the derivatives of the Recursive Newton Euler Algorithms with respect to the joint configura...
const std::string library_name
Name of the library.
bool build_jacobian
Options to build or not the Jacobian of he function.
void buildMap()
build the mapping Y = f(X)
Eigen::DenseIndex getOutputDimension() const
Dimension of the output vector.
Eigen::DenseIndex getInputDimension() const
Dimension of the input vector.
Main pinocchio namespace.
TangentVectorType tau
Vector of joint torques (dim model.nv).
const DataTpl< Scalar, Options, JointCollectionTpl >::TangentVectorType & rnea(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType1 > &v, const Eigen::MatrixBase< TangentVectorType2 > &a)
The Recursive Newton-Euler algorithm. It computes the inverse dynamics, aka the joint torques accordi...
void buildMap()
build the mapping Y = f(X)
void computeABADerivatives(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType1 > &v, const Eigen::MatrixBase< TangentVectorType2 > &tau, const Eigen::MatrixBase< MatrixType1 > &aba_partial_dq, const Eigen::MatrixBase< MatrixType2 > &aba_partial_dv, const Eigen::MatrixBase< MatrixType3 > &aba_partial_dtau)
The derivatives of the Articulated-Body algorithm.
void buildMap()
build the mapping Y = f(X)
int nq
Dimension of the configuration vector representation.