5 #include "pinocchio/spatial/se3.hpp" 6 #include "pinocchio/multibody/model.hpp" 7 #include "pinocchio/multibody/data.hpp" 8 #include "pinocchio/algorithm/jacobian.hpp" 9 #include "pinocchio/algorithm/contact-dynamics.hpp" 10 #include "pinocchio/algorithm/joint-configuration.hpp" 11 #include "pinocchio/parsers/sample-models.hpp" 12 #include "pinocchio/utils/timer.hpp" 16 #include <boost/test/unit_test.hpp> 17 #include <boost/utility/binary.hpp> 19 BOOST_AUTO_TEST_SUITE ( BOOST_TEST_MODULE )
23 using namespace Eigen;
38 const std::string RF =
"rleg6_joint";
39 const std::string LF =
"lleg6_joint";
48 Eigen::MatrixXd
J (12, model.
nv);
50 J.topRows<6> () = J_RF;
51 J.bottomRows<6> () = J_LF;
55 Eigen::MatrixXd
H(J.transpose());
58 data.
M.triangularView<Eigen::StrictlyLower>() = data.
M.transpose().triangularView<Eigen::StrictlyLower>();
60 MatrixXd Minv (data.
M.inverse());
61 MatrixXd JMinvJt (J * Minv * J.transpose());
63 Eigen::MatrixXd G_ref(J.transpose());
65 for(
int k=0;k<model.
nv;++k) G_ref.row(k) /= sqrt(data.
D[k]);
66 Eigen::MatrixXd H_ref(G_ref.transpose() * G_ref);
67 BOOST_CHECK(H_ref.isApprox(JMinvJt,1e-12));
69 VectorXd lambda_ref = -JMinvJt.inverse() * (J*Minv*(tau - data.
nle) + gamma);
70 BOOST_CHECK(data.
lambda_c.isApprox(lambda_ref, 1e-12));
72 VectorXd a_ref = Minv*(tau - data.
nle + J.transpose()*lambda_ref);
74 Eigen::VectorXd dynamics_residual_ref (data.
M * a_ref + data.
nle - tau - J.transpose()*lambda_ref);
75 BOOST_CHECK(dynamics_residual_ref.norm() <= 1e-11);
78 BOOST_CHECK(constraint_residual.norm() <= 1e-12);
81 BOOST_CHECK(dynamics_residual.norm() <= 1e-12);
87 using namespace Eigen;
98 const std::string RF =
"rleg6_joint";
99 const std::string LF =
"lleg6_joint";
108 Eigen::MatrixXd
J(12, model.
nv);
110 J.topRows<6>() = J_RF;
111 J.bottomRows<6>() = J_LF;
115 data_ref.
M.triangularView<Eigen::StrictlyLower>() = data_ref.
M.transpose().triangularView<Eigen::StrictlyLower>();
117 Eigen::MatrixXd MJtJ(model.
nv+12, model.
nv+12);
118 MJtJ << data_ref.
M, J.transpose(),
119 J, Eigen::MatrixXd::Zero(12, 12);
121 Eigen::MatrixXd KKTMatrix_inv(model.
nv+12, model.
nv+12);
124 BOOST_CHECK(KKTMatrix_inv.isApprox(MJtJ.inverse()));
129 using namespace Eigen;
143 const std::string RF =
"rleg6_joint";
144 const std::string LF =
"lleg6_joint";
153 Eigen::MatrixXd
J (12, model.
nv);
155 J.topRows<6> () = J_RF;
156 J.bottomRows<6> () = J_LF;
160 Eigen::MatrixXd
H(J.transpose());
164 data.
M.triangularView<Eigen::StrictlyLower>() = data.
M.transpose().triangularView<Eigen::StrictlyLower>();
166 Eigen::MatrixXd MJtJ(model.
nv+12, model.
nv+12);
167 MJtJ << data.
M, J.transpose(),
168 J, Eigen::MatrixXd::Zero(12, 12);
170 Eigen::MatrixXd KKTMatrix_inv(model.
nv+12, model.
nv+12);
173 BOOST_CHECK(KKTMatrix_inv.isApprox(MJtJ.inverse()));
179 data.
M.triangularView<Eigen::StrictlyLower>() = data.
M.transpose().triangularView<Eigen::StrictlyLower>();
180 MJtJ << data.
M, J.transpose(),
181 J, Eigen::MatrixXd::Zero(12, 12);
185 BOOST_CHECK(KKTMatrix_inv.isApprox(MJtJ.inverse()));
190 using namespace Eigen;
205 const std::string RF =
"rleg6_joint";
211 Eigen::MatrixXd
J (12, model.
nv);
213 J.topRows<6> () = J_RF;
214 J.bottomRows<6> () = J_RF;
222 Eigen::MatrixXd
H(J.transpose());
223 data.
M.triangularView<Eigen::StrictlyLower>() =
224 data.
M.transpose().triangularView<Eigen::StrictlyLower>();
226 MatrixXd Minv (data.
M.inverse());
227 MatrixXd JMinvJt (J * Minv * J.transpose());
230 Eigen::MatrixXd G_ref(J.transpose());
232 for(
int k=0;k<model.
nv;++k) G_ref.row(k) /= sqrt(data.
D[k]);
233 Eigen::MatrixXd H_ref(G_ref.transpose() * G_ref);
234 BOOST_CHECK(H_ref.isApprox(JMinvJt,1e-12));
239 BOOST_CHECK(constraint_residual.norm() <= 1e-9);
240 BOOST_CHECK(dynamics_residual.norm() <= 1e-12);
245 using namespace Eigen;
259 const std::string RF =
"rleg6_joint";
260 const std::string LF =
"lleg6_joint";
269 Eigen::MatrixXd
J (12, model.
nv);
271 J.topRows<6> () = J_RF;
272 J.bottomRows<6> () = J_LF;
276 Eigen::MatrixXd
H(J.transpose());
279 data.
M.triangularView<Eigen::StrictlyLower>() = data.
M.transpose().triangularView<Eigen::StrictlyLower>();
281 MatrixXd Minv (data.
M.inverse());
282 MatrixXd JMinvJt (J * Minv * J.transpose());
284 Eigen::MatrixXd G_ref(J.transpose());
286 for(
int k=0;k<model.
nv;++k) G_ref.row(k) /= sqrt(data.
D[k]);
287 Eigen::MatrixXd H_ref(G_ref.transpose() * G_ref);
288 BOOST_CHECK(H_ref.isApprox(JMinvJt,1e-12));
290 VectorXd lambda_ref = JMinvJt.inverse() * (-r_coeff * J * v_before - J * v_before);
291 BOOST_CHECK(data.
impulse_c.isApprox(lambda_ref, 1e-12));
293 VectorXd v_after_ref = Minv*(data.
M * v_before + J.transpose()*lambda_ref);
296 BOOST_CHECK(constraint_residual.norm() <= 1e-12);
299 BOOST_CHECK(dynamics_residual.norm() <= 1e-12);
305 using namespace Eigen;
313 #ifdef _INTENSE_TESTING_ 314 const size_t NBT = 1000*1000;
316 const size_t NBT = 100;
320 const size_t NBT = 1;
321 std::cout <<
"(the time score in debug mode is not relevant) " ;
322 #endif // ifndef NDEBUG 332 const std::string RF =
"rleg6_joint";
333 const std::string LF =
"lleg6_joint";
342 Eigen::MatrixXd
J (12, model.
nv);
343 J.topRows<6> () = J_RF;
344 J.bottomRows<6> () = J_LF;
358 timer.
toc(std::cout,NBT);
362 BOOST_AUTO_TEST_SUITE_END ()
void getKKTContactDynamicMatrixInverse(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConstraintMatrixType > &J, const Eigen::MatrixBase< KKTMatrixType > &KKTMatrix_inv)
Computes the inverse of the KKT matrix for dynamics with contact constraints. It computes the followi...
const DataTpl< Scalar, Options, JointCollectionTpl >::TangentVectorType & impulseDynamics(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType > &v_before, const Eigen::MatrixBase< ConstraintMatrixType > &J, const Scalar r_coeff=0., const Scalar inv_damping=0.)
Compute the impulse dynamics with contact constraints. Internally, pinocchio::crba is called...
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...
void computeKKTContactDynamicMatrixInverse(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< ConstraintMatrixType > &J, const Eigen::MatrixBase< KKTMatrixType > &KKTMatrix_inv, const Scalar &inv_damping=0.)
Computes the inverse of the KKT matrix for dynamics with contact constraints. It computes the followi...
ConfigVectorType lowerPositionLimit
Lower joint configuration limit.
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.
VectorXs D
Diagonal of the joint space intertia matrix obtained by a Cholesky Decomposition. ...
void getJointJacobian(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const DataTpl< Scalar, Options, JointCollectionTpl > &data, const typename ModelTpl< Scalar, Options, JointCollectionTpl >::JointIndex jointId, const ReferenceFrame rf, const Eigen::MatrixBase< Matrix6Like > &J)
Computes the Jacobian of a specific joint frame expressed either in the world (rf = WORLD) frame or i...
VectorXs nle
Vector of Non Linear Effects (dim model.nv). It corresponds to concatenation of the Coriolis...
The LOCAL frame convention corresponds to the frame directly attached to the moving part (Joint...
MatrixXs M
The joint space inertia matrix (a square matrix of dim model.nv).
TangentVectorType ddq
The joint accelerations computed from ABA.
VectorXs lambda_c
Lagrange Multipliers corresponding to the contact forces in pinocchio::forwardDynamics.
Main pinocchio namespace.
JointIndex getJointId(const std::string &name) const
Return the index of a joint given by its name.
const DataTpl< Scalar, Options, JointCollectionTpl >::Matrix6x & computeJointJacobians(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q)
Computes the full model Jacobian, i.e. the stack of all motion subspace expressed in the world frame...
TangentVectorType dq_after
Generalized velocity after impact.
int nv
Dimension of the velocity vector space.
Mat & Uiv(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< Mat > &v)
Perform the pivot inversion using the Cholesky decomposition stored in data and acting in place...
const DataTpl< Scalar, Options, JointCollectionTpl >::TangentVectorType & forwardDynamics(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< ConstraintMatrixType > &J, const Eigen::MatrixBase< DriftVectorType > &gamma, const Scalar inv_damping=0.)
Compute the forward dynamics with contact constraints. Internally, pinocchio::computeAllTerms is call...
void normalize(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const Eigen::MatrixBase< ConfigVectorType > &qout)
Normalize a configuration vector.
ConfigVectorType upperPositionLimit
Upper joint configuration limit.
VectorXs impulse_c
Lagrange Multipliers corresponding to the contact impulses in pinocchio::impulseDynamics.
void humanoidRandom(ModelTpl< Scalar, Options, JointCollectionTpl > &model, bool usingFF=true)
Create a humanoid kinematic tree with 6-DOF limbs and random joint placements.
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > VectorXd
Eigen::Matrix< Scalar, 6, Eigen::Dynamic, Options > Matrix6x
The 6d jacobian type (temporary)
int nq
Dimension of the configuration vector representation.