5 #ifndef __pinocchio_codegen_code_generator_algo_hpp__ 6 #define __pinocchio_codegen_code_generator_algo_hpp__ 8 #include "pinocchio/codegen/code-generator-base.hpp" 10 #include "pinocchio/algorithm/joint-configuration.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::ADConfigVectorType ADConfigVectorType;
27 typedef typename Base::ADTangentVectorType ADTangentVectorType;
28 typedef typename Base::MatrixXs MatrixXs;
29 typedef typename Base::VectorXs VectorXs;
36 ad_q = ADConfigVectorType(model.
nq); ad_q =
neutral(ad_model);
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;
99 MatrixXs dtau_dq, dtau_dv, dtau_da;
103 using Base::ad_model;
114 ADConfigVectorType ad_q, ad_q_plus;
115 ADTangentVectorType ad_dq, ad_v, ad_a;
118 template<
typename _Scalar>
122 typedef typename Base::Scalar Scalar;
125 typedef typename Base::ADConfigVectorType ADConfigVectorType;
126 typedef typename Base::ADTangentVectorType ADTangentVectorType;
127 typedef typename Base::MatrixXs MatrixXs;
128 typedef typename Base::VectorXs VectorXs;
135 ad_q = ADConfigVectorType(model.
nq); ad_q =
neutral(ad_model);
136 ad_v = ADTangentVectorType(model.
nv); ad_v.setZero();
137 ad_tau = ADTangentVectorType(model.
nv); ad_tau.setZero();
141 da_dq = MatrixXs::Zero(model.
nv,model.
nq);
142 da_dv = MatrixXs::Zero(model.
nv,model.
nv);
143 da_dtau = MatrixXs::Zero(model.
nv,model.
nv);
148 CppAD::Independent(ad_X);
150 Eigen::DenseIndex it = 0;
151 ad_q = ad_X.segment(it,ad_model.nq); it += ad_model.nq;
152 ad_v = ad_X.segment(it,ad_model.nv); it += ad_model.nv;
153 ad_tau = ad_X.segment(it,ad_model.nv); it += ad_model.nv;
158 ad_fun.Dependent(ad_X,ad_Y);
159 ad_fun.optimize(
"no_compare_op");
162 using Base::evalFunction;
163 template<
typename ConfigVectorType,
typename TangentVector1,
typename TangentVector2>
164 void evalFunction(
const Eigen::MatrixBase<ConfigVectorType> & q,
165 const Eigen::MatrixBase<TangentVector1> & v,
166 const Eigen::MatrixBase<TangentVector2> & tau)
169 Eigen::DenseIndex it = 0;
170 x.segment(it,ad_model.nq) = q; it += ad_model.nq;
171 x.segment(it,ad_model.nv) = v; it += ad_model.nv;
172 x.segment(it,ad_model.nv) = tau; it += ad_model.nv;
178 using Base::evalJacobian;
179 template<
typename ConfigVectorType,
typename TangentVector1,
typename TangentVector2>
180 void evalJacobian(
const Eigen::MatrixBase<ConfigVectorType> & q,
181 const Eigen::MatrixBase<TangentVector1> & v,
182 const Eigen::MatrixBase<TangentVector2> & tau)
185 Eigen::DenseIndex it = 0;
186 x.segment(it,ad_model.nq) = q; it += ad_model.nq;
187 x.segment(it,ad_model.nv) = v; it += ad_model.nv;
188 x.segment(it,ad_model.nv) = tau; it += ad_model.nv;
193 da_dq = Base::jac.middleCols(it,ad_model.nq); it += ad_model.nq;
194 da_dv = Base::jac.middleCols(it,ad_model.nv); it += ad_model.nv;
195 da_dtau = Base::jac.middleCols(it,ad_model.nv); it += ad_model.nv;
198 MatrixXs da_dq,da_dv,da_dtau;
202 using Base::ad_model;
213 ADConfigVectorType ad_q, ad_q_plus;
214 ADTangentVectorType ad_dq, ad_v, ad_tau;
217 template<
typename _Scalar>
221 typedef typename Base::Scalar Scalar;
224 typedef typename Base::ADConfigVectorType ADConfigVectorType;
225 typedef typename Base::ADTangentVectorType ADTangentVectorType;
226 typedef typename Base::MatrixXs MatrixXs;
227 typedef typename Base::VectorXs VectorXs;
234 ad_q = ADConfigVectorType(model.
nq); ad_q =
neutral(ad_model);
238 M = MatrixXs::Zero(model.
nv,model.
nq);
245 CppAD::Independent(ad_X);
247 Eigen::DenseIndex it = 0;
248 ad_q = ad_X.segment(it,ad_model.nq); it += ad_model.nq;
251 Eigen::DenseIndex it_Y = 0;
253 for(Eigen::DenseIndex i = 0; i < ad_model.nv; ++i)
255 for(Eigen::DenseIndex j = i; j < ad_model.nv; ++j)
257 ad_Y[it_Y++] = ad_data.M(i,j);
263 ad_fun.Dependent(ad_X,ad_Y);
264 ad_fun.optimize(
"no_compare_op");
267 template<
typename ConfigVectorType>
268 void evalFunction(
const Eigen::MatrixBase<ConfigVectorType> & q)
271 Eigen::DenseIndex it = 0;
272 x.segment(it,ad_model.nq) = q; it += ad_model.nq;
274 Base::evalFunction(x);
277 Eigen::DenseIndex it_Y = 0;
278 for(Eigen::DenseIndex i = 0; i < ad_model.nv; ++i)
280 for(Eigen::DenseIndex j = i; j < ad_model.nv; ++j)
282 M(i,j) = Base::y[it_Y++];
293 using Base::ad_model;
303 ADConfigVectorType ad_q;
306 template<
typename _Scalar>
310 typedef typename Base::Scalar Scalar;
313 typedef typename Base::ADConfigVectorType ADConfigVectorType;
314 typedef typename Base::ADTangentVectorType ADTangentVectorType;
315 typedef typename Base::MatrixXs MatrixXs;
316 typedef typename Base::VectorXs VectorXs;
323 ad_q = ADConfigVectorType(model.
nq); ad_q =
neutral(ad_model);
327 Minv = MatrixXs::Zero(model.
nv,model.
nq);
334 CppAD::Independent(ad_X);
336 Eigen::DenseIndex it = 0;
337 ad_q = ad_X.segment(it,ad_model.nq); it += ad_model.nq;
340 Eigen::DenseIndex it_Y = 0;
341 for(Eigen::DenseIndex i = 0; i < ad_model.nv; ++i)
343 for(Eigen::DenseIndex j = i; j < ad_model.nv; ++j)
345 ad_Y[it_Y++] = ad_data.Minv(i,j);
351 ad_fun.Dependent(ad_X,ad_Y);
352 ad_fun.optimize(
"no_compare_op");
355 template<
typename ConfigVectorType>
356 void evalFunction(
const Eigen::MatrixBase<ConfigVectorType> & q)
359 Eigen::DenseIndex it = 0;
360 x.segment(it,ad_model.nq) = q; it += ad_model.nq;
362 Base::evalFunction(x);
365 Eigen::DenseIndex it_Y = 0;
366 for(Eigen::DenseIndex i = 0; i < ad_model.nv; ++i)
368 for(Eigen::DenseIndex j = i; j < ad_model.nv; ++j)
370 Minv(i,j) = Base::y[it_Y++];
379 using Base::ad_model;
389 ADConfigVectorType ad_q;
392 template<
typename _Scalar>
396 typedef typename Base::Scalar Scalar;
399 typedef typename Base::ADConfigVectorType ADConfigVectorType;
400 typedef typename Base::ADTangentVectorType ADTangentVectorType;
401 typedef typename Base::MatrixXs MatrixXs;
402 typedef typename PINOCCHIO_EIGEN_PLAIN_ROW_MAJOR_TYPE(MatrixXs) RowMatrixXs;
403 typedef typename Base::VectorXs VectorXs;
406 typedef typename ADData::MatrixXs ADMatrixXs;
407 typedef typename PINOCCHIO_EIGEN_PLAIN_ROW_MAJOR_TYPE(ADMatrixXs) RowADMatrixXs;
411 const std::string &
library_name =
"cg_partial_rnea_eval")
414 ad_q = ADConfigVectorType(model.
nq); ad_q =
neutral(ad_model);
415 ad_v = ADTangentVectorType(model.
nv); ad_v.setZero();
416 ad_a = ADTangentVectorType(model.
nv); ad_a.setZero();
421 ad_dtau_dq = ADMatrixXs::Zero(model.
nv,model.
nv);
422 ad_dtau_dv = ADMatrixXs::Zero(model.
nv,model.
nv);
423 ad_dtau_da = ADMatrixXs::Zero(model.
nv,model.
nv);
425 dtau_dq = MatrixXs::Zero(model.
nv,model.
nv);
426 dtau_dv = MatrixXs::Zero(model.
nv,model.
nv);
427 dtau_da = MatrixXs::Zero(model.
nv,model.
nv);
434 CppAD::Independent(ad_X);
436 Eigen::DenseIndex it = 0;
437 ad_q = ad_X.segment(it,ad_model.nq); it += ad_model.nq;
438 ad_v = ad_X.segment(it,ad_model.nv); it += ad_model.nv;
439 ad_a = ad_X.segment(it,ad_model.nv); it += ad_model.nv;
443 ad_dtau_dq,ad_dtau_dv,ad_dtau_da);
447 Eigen::DenseIndex it_Y = 0;
448 Eigen::Map<RowADMatrixXs>(ad_Y.data()+it_Y,ad_model.nv,ad_model.nv) = ad_dtau_dq;
449 it_Y += ad_model.nv*ad_model.nv;
450 Eigen::Map<RowADMatrixXs>(ad_Y.data()+it_Y,ad_model.nv,ad_model.nv) = ad_dtau_dv;
451 it_Y += ad_model.nv*ad_model.nv;
452 Eigen::Map<RowADMatrixXs>(ad_Y.data()+it_Y,ad_model.nv,ad_model.nv) = ad_dtau_da;
453 it_Y += ad_model.nv*ad_model.nv;
455 ad_fun.Dependent(ad_X,ad_Y);
456 ad_fun.optimize(
"no_compare_op");
459 template<
typename ConfigVectorType,
typename TangentVector1,
typename TangentVector2>
460 void evalFunction(
const Eigen::MatrixBase<ConfigVectorType> & q,
461 const Eigen::MatrixBase<TangentVector1> & v,
462 const Eigen::MatrixBase<TangentVector2> & a)
465 Eigen::DenseIndex it_x = 0;
466 x.segment(it_x,ad_model.nq) = q; it_x += ad_model.nq;
467 x.segment(it_x,ad_model.nq) = v; it_x += ad_model.nv;
468 x.segment(it_x,ad_model.nq) = a; it_x += ad_model.nv;
470 Base::evalFunction(x);
473 Eigen::DenseIndex it_y = 0;
474 dtau_dq = Eigen::Map<RowMatrixXs>(Base::y.data()+it_y,ad_model.nv,ad_model.nv);
475 it_y += ad_model.nv*ad_model.nv;
476 dtau_dv = Eigen::Map<RowMatrixXs>(Base::y.data()+it_y,ad_model.nv,ad_model.nv);
477 it_y += ad_model.nv*ad_model.nv;
478 dtau_da = Eigen::Map<RowMatrixXs>(Base::y.data()+it_y,ad_model.nv,ad_model.nv);
479 it_y += ad_model.nv*ad_model.nv;
485 using Base::ad_model;
493 VectorXs partial_derivatives;
494 ADMatrixXs ad_dtau_dq, ad_dtau_dv, ad_dtau_da;
495 MatrixXs dtau_dq, dtau_dv, dtau_da;
497 ADConfigVectorType ad_q;
498 ADTangentVectorType ad_v, ad_a;
501 template<
typename _Scalar>
505 typedef typename Base::Scalar Scalar;
508 typedef typename Base::ADConfigVectorType ADConfigVectorType;
509 typedef typename Base::ADTangentVectorType ADTangentVectorType;
510 typedef typename Base::MatrixXs MatrixXs;
511 typedef typename PINOCCHIO_EIGEN_PLAIN_ROW_MAJOR_TYPE(MatrixXs) RowMatrixXs;
512 typedef typename Base::VectorXs VectorXs;
515 typedef typename ADData::MatrixXs ADMatrixXs;
516 typedef typename PINOCCHIO_EIGEN_PLAIN_ROW_MAJOR_TYPE(ADMatrixXs) RowADMatrixXs;
520 const std::string &
library_name =
"cg_partial_aba_eval")
523 ad_q = ADConfigVectorType(model.
nq); ad_q =
neutral(ad_model);
524 ad_v = ADTangentVectorType(model.
nv); ad_v.setZero();
525 ad_tau = ADTangentVectorType(model.
nv); ad_tau.setZero();
530 ad_dddq_dq = ADMatrixXs::Zero(model.
nv,model.
nv);
531 ad_dddq_dv = ADMatrixXs::Zero(model.
nv,model.
nv);
532 ad_dddq_dtau = ADMatrixXs::Zero(model.
nv,model.
nv);
534 dddq_dq = MatrixXs::Zero(model.
nv,model.
nv);
535 dddq_dv = MatrixXs::Zero(model.
nv,model.
nv);
536 dddq_dtau = MatrixXs::Zero(model.
nv,model.
nv);
543 CppAD::Independent(ad_X);
545 Eigen::DenseIndex it = 0;
546 ad_q = ad_X.segment(it,ad_model.nq); it += ad_model.nq;
547 ad_v = ad_X.segment(it,ad_model.nv); it += ad_model.nv;
548 ad_tau = ad_X.segment(it,ad_model.nv); it += ad_model.nv;
552 ad_dddq_dq,ad_dddq_dv,ad_dddq_dtau);
556 Eigen::DenseIndex it_Y = 0;
557 Eigen::Map<RowADMatrixXs>(ad_Y.data()+it_Y,ad_model.nv,ad_model.nv) = ad_dddq_dq;
558 it_Y += ad_model.nv*ad_model.nv;
559 Eigen::Map<RowADMatrixXs>(ad_Y.data()+it_Y,ad_model.nv,ad_model.nv) = ad_dddq_dv;
560 it_Y += ad_model.nv*ad_model.nv;
561 Eigen::Map<RowADMatrixXs>(ad_Y.data()+it_Y,ad_model.nv,ad_model.nv) = ad_dddq_dtau;
562 it_Y += ad_model.nv*ad_model.nv;
564 ad_fun.Dependent(ad_X,ad_Y);
565 ad_fun.optimize(
"no_compare_op");
568 template<
typename ConfigVectorType,
typename TangentVector1,
typename TangentVector2>
569 void evalFunction(
const Eigen::MatrixBase<ConfigVectorType> & q,
570 const Eigen::MatrixBase<TangentVector1> & v,
571 const Eigen::MatrixBase<TangentVector2> & tau)
574 Eigen::DenseIndex it_x = 0;
575 x.segment(it_x,ad_model.nq) = q; it_x += ad_model.nq;
576 x.segment(it_x,ad_model.nq) = v; it_x += ad_model.nv;
577 x.segment(it_x,ad_model.nq) = tau; it_x += ad_model.nv;
579 Base::evalFunction(x);
582 Eigen::DenseIndex it_y = 0;
583 dddq_dq = Eigen::Map<RowMatrixXs>(Base::y.data()+it_y,ad_model.nv,ad_model.nv);
584 it_y += ad_model.nv*ad_model.nv;
585 dddq_dv = Eigen::Map<RowMatrixXs>(Base::y.data()+it_y,ad_model.nv,ad_model.nv);
586 it_y += ad_model.nv*ad_model.nv;
587 dddq_dtau = Eigen::Map<RowMatrixXs>(Base::y.data()+it_y,ad_model.nv,ad_model.nv);
588 it_y += ad_model.nv*ad_model.nv;
594 using Base::ad_model;
602 VectorXs partial_derivatives;
603 ADMatrixXs ad_dddq_dq, ad_dddq_dv, ad_dddq_dtau;
604 MatrixXs dddq_dq, dddq_dv, dddq_dtau;
606 ADConfigVectorType ad_q;
607 ADTangentVectorType ad_v, ad_tau;
610 template<
typename _Scalar>
614 typedef typename Base::Scalar Scalar;
617 typedef typename Base::ADConfigVectorType ADConfigVectorType;
618 typedef typename Base::ADTangentVectorType ADTangentVectorType;
619 typedef typename Base::MatrixXs MatrixXs;
620 typedef typename Base::VectorXs VectorXs;
627 ad_q = ADConfigVectorType(model.
nq); ad_q =
neutral(ad_model);
628 ad_v = ADTangentVectorType(model.
nv); ad_v.setZero();
635 CppAD::Independent(ad_X);
637 Eigen::DenseIndex it = 0;
638 ad_q = ad_X.segment(it,ad_model.nq); it += ad_model.nq;
639 ad_v = ad_X.segment(it,ad_model.nv);
642 ad_fun.Dependent(ad_X,ad_Y);
643 ad_fun.optimize(
"no_compare_op");
646 using Base::evalFunction;
647 template<
typename ConfigVectorType1,
typename TangentVector,
typename ConfigVectorType2>
648 void evalFunction(
const Eigen::MatrixBase<ConfigVectorType1> & q,
649 const Eigen::MatrixBase<TangentVector> & v,
650 const Eigen::MatrixBase<ConfigVectorType2> & qout)
653 Eigen::DenseIndex it = 0;
654 x.segment(it,ad_model.nq) = q; it += ad_model.nq;
655 x.segment(it,ad_model.nv) = v;
658 PINOCCHIO_EIGEN_CONST_CAST(ConfigVectorType2,qout) = Base::y;
663 using Base::ad_model;
672 ADConfigVectorType ad_q;
673 ADTangentVectorType ad_v;
677 template<
typename _Scalar>
681 typedef typename Base::Scalar Scalar;
684 typedef typename Base::ADConfigVectorType ADConfigVectorType;
685 typedef typename Base::ADTangentVectorType ADTangentVectorType;
686 typedef typename Base::MatrixXs MatrixXs;
687 typedef typename Base::VectorXs VectorXs;
691 const std::string &
library_name =
"cg_difference_eval")
694 ad_q0 = ADConfigVectorType(model.
nq); ad_q0 =
neutral(ad_model);
695 ad_q1 = ADConfigVectorType(model.
nq); ad_q1 =
neutral(ad_model);
702 CppAD::Independent(ad_X);
704 Eigen::DenseIndex it = 0;
705 ad_q0 = ad_X.segment(it,ad_model.nq); it += ad_model.nq;
706 ad_q1 = ad_X.segment(it,ad_model.nq);
709 ad_fun.Dependent(ad_X,ad_Y);
710 ad_fun.optimize(
"no_compare_op");
713 using Base::evalFunction;
714 template<
typename ConfigVectorType1,
typename ConfigVectorType2,
typename TangentVector>
715 void evalFunction(
const Eigen::MatrixBase<ConfigVectorType1> & q0,
716 const Eigen::MatrixBase<ConfigVectorType2> & q1,
717 const Eigen::MatrixBase<TangentVector> & v)
720 Eigen::DenseIndex it = 0;
721 x.segment(it,ad_model.nq) = q0; it += ad_model.nq;
722 x.segment(it,ad_model.nq) = q1;
725 PINOCCHIO_EIGEN_CONST_CAST(TangentVector,v) = Base::y;
730 using Base::ad_model;
739 ADConfigVectorType ad_q0;
740 ADConfigVectorType ad_q1;
744 template<
typename _Scalar>
748 typedef typename Base::Scalar Scalar;
751 typedef typename Base::ADConfigVectorType ADConfigVectorType;
752 typedef typename Base::ADTangentVectorType ADTangentVectorType;
753 typedef typename Base::ADVectorXs ADVectorXs;
754 typedef typename Base::ADMatrixXs ADMatrixXs;
755 typedef typename Base::MatrixXs MatrixXs;
756 typedef typename Base::VectorXs VectorXs;
760 const std::string &
library_name =
"cg_dDifference_eval")
765 ad_J0 = ADMatrixXs::Zero(ad_model.nv, ad_model.nv);
766 ad_J1 = ADMatrixXs::Zero(ad_model.nv, ad_model.nv);
772 CppAD::Independent(ad_X);
774 Eigen::DenseIndex it = 0;
775 ad_q0 = ad_X.segment(it,ad_model.nq); it += ad_model.nq;
776 ad_q1 = ad_X.segment(it,ad_model.nq);
780 Eigen::Map<ADMatrixXs>(ad_Y.data(), 2*ad_model.nv, ad_model.nv).topRows(ad_model.nv) = ad_J0;
781 Eigen::Map<ADMatrixXs>(ad_Y.data(), 2*ad_model.nv, ad_model.nv).bottomRows(ad_model.nv) = ad_J1;
782 ad_fun.Dependent(ad_X,ad_Y);
783 ad_fun.optimize(
"no_compare_op");
786 using Base::evalFunction;
787 template<
typename ConfigVectorType1,
typename ConfigVectorType2,
typename JacobianMatrix>
788 void evalFunction(
const Eigen::MatrixBase<ConfigVectorType1> & q0,
789 const Eigen::MatrixBase<ConfigVectorType2> & q1,
790 const Eigen::MatrixBase<JacobianMatrix> & J,
794 Eigen::DenseIndex it = 0;
795 x.segment(it,ad_model.nq) = q0; it += ad_model.nq;
796 x.segment(it,ad_model.nq) = q1;
801 case pinocchio::ARG0:
802 PINOCCHIO_EIGEN_CONST_CAST(JacobianMatrix,J)
803 = Eigen::Map<MatrixXs>(Base::y.data(), 2*ad_model.nv, ad_model.nv).topRows(ad_model.nv);
805 case pinocchio::ARG1:
806 PINOCCHIO_EIGEN_CONST_CAST(JacobianMatrix,J)
807 = Eigen::Map<MatrixXs>(Base::y.data(), 2*ad_model.nv, ad_model.nv).bottomRows(ad_model.nv);
810 assert(
false &&
"Wrong argument");
817 using Base::ad_model;
824 ADConfigVectorType ad_q0;
825 ADConfigVectorType ad_q1;
832 #endif // ifndef __pinocchio_codegen_code_generator_algo_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...
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...
Eigen::DenseIndex getOutputDimension() const
Dimension of the output vector.
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.
Eigen::DenseIndex getInputDimension() const
Dimension of the input vector.
void integrate(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType > &v, const Eigen::MatrixBase< ReturnType > &qout)
Integrate a configuration vector for the specified model for a tangent vector during one unit time...
void dDifference(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const Eigen::MatrixBase< ConfigVector1 > &q0, const Eigen::MatrixBase< ConfigVector2 > &q1, const Eigen::MatrixBase< JacobianMatrix > &J, const ArgumentPosition arg)
Computes the Jacobian of a small variation of the configuration vector into the tangent space at iden...
void difference(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const Eigen::MatrixBase< ConfigVectorIn1 > &q0, const Eigen::MatrixBase< ConfigVectorIn2 > &q1, const Eigen::MatrixBase< ReturnType > &dvout)
Compute the tangent vector that must be integrated during one unit time to go from q0 to q1...
void buildMap()
build the mapping Y = f(X)
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 partial derivatives of the Recursive Newton Euler Algorithms with respect to the joint c...
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)
ArgumentPosition
Argument position. Used as template parameter to refer to an argument.
void buildMap()
build the mapping Y = f(X)
void buildMap()
build the mapping Y = f(X)
Main pinocchio namespace.
TangentVectorType tau
Vector of joint torques (dim model.nv).
void neutral(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const Eigen::MatrixBase< ReturnType > &qout)
Return the neutral configuration element related to the model configuration space.
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.