expose-contact-dynamics.cpp
Go to the documentation of this file.
1 //
2 // Copyright (c) 2015-2020 CNRS, INRIA
3 //
4 
6 #include "pinocchio/algorithm/contact-dynamics.hpp"
7 
8 namespace pinocchio
9 {
10  namespace python
11  {
12 
14  Data & data,
15  const Eigen::VectorXd & q,
16  const Eigen::VectorXd & v,
17  const Eigen::VectorXd & tau,
18  const Eigen::MatrixXd & J,
19  const Eigen::VectorXd & gamma,
20  const double inv_damping = 0.0)
21  {
22  return forwardDynamics(model, data, q, v, tau, J, gamma, inv_damping);
23  }
24 
25  BOOST_PYTHON_FUNCTION_OVERLOADS(forwardDynamics_overloads, forwardDynamics_proxy, 7, 8)
26 
28  Data & data,
29  const Eigen::VectorXd & tau,
30  const Eigen::MatrixXd & J,
31  const Eigen::VectorXd & gamma,
32  const double inv_damping = 0.0)
33  {
34  return forwardDynamics(model, data, tau, J, gamma, inv_damping);
35  }
36 
37  BOOST_PYTHON_FUNCTION_OVERLOADS(forwardDynamics_overloads_no_q, forwardDynamics_proxy_no_q, 5, 6)
38 
40  Data & data,
41  const Eigen::VectorXd & q,
42  const Eigen::VectorXd & v_before,
43  const Eigen::MatrixXd & J,
44  const double r_coeff = 0.,
45  const double inv_damping = 0.)
46  {
47  return impulseDynamics(model, data, q, v_before, J, r_coeff, inv_damping);
48  }
49 
50  BOOST_PYTHON_FUNCTION_OVERLOADS(impulseDynamics_overloads, impulseDynamics_proxy, 5, 7)
51 
53  Data & data,
54  const Eigen::VectorXd & v_before,
55  const Eigen::MatrixXd & J,
56  const double r_coeff = 0.,
57  const double inv_damping = 0.)
58  {
59  return impulseDynamics(model, data, v_before, J, r_coeff, inv_damping);
60  }
61 
62  BOOST_PYTHON_FUNCTION_OVERLOADS(impulseDynamics_overloads_no_q, impulseDynamics_proxy_no_q, 4, 6)
63 
65  Data & data,
66  const Eigen::VectorXd & q,
67  const Eigen::MatrixXd & J,
68  const double mu = 0)
69  {
70  Eigen::MatrixXd KKTMatrix_inv(model.nv+J.rows(), model.nv+J.rows());
71  computeKKTContactDynamicMatrixInverse(model, data, q, J, KKTMatrix_inv, mu);
72  return KKTMatrix_inv;
73  }
74 
75  BOOST_PYTHON_FUNCTION_OVERLOADS(computeKKTContactDynamicMatrixInverse_overload,
77 
78  static const Eigen::MatrixXd getKKTContactDynamicMatrixInverse_proxy(const Model & model,
79  Data & data,
80  const Eigen::MatrixXd & J)
81  {
82  Eigen::MatrixXd KKTMatrix_inv(model.nv+J.rows(), model.nv+J.rows());
83  getKKTContactDynamicMatrixInverse(model, data, J, KKTMatrix_inv);
84  return KKTMatrix_inv;
85  }
86 
88  {
89  using namespace Eigen;
90 
91  bp::def("forwardDynamics",
93  forwardDynamics_overloads(
94  bp::args("Model","Data",
95  "Joint configuration q (size Model::nq)",
96  "Joint velocity v (size Model::nv)",
97  "Joint torque tau (size Model::nv)",
98  "Contact Jacobian J (size nb_constraint * Model::nv)",
99  "Contact drift gamma (size nb_constraint)",
100  "(double) Damping factor for cholesky decomposition of JMinvJt. Set to zero if constraints are full rank."),
101  "Solves the forward dynamics problem with contacts, puts the result in Data::ddq and return it. The contact forces are stored in data.lambda_c."
102  " Internally, pinocchio.computeAllTerms is called."
103  ));
104 
105  bp::def("forwardDynamics",
107  forwardDynamics_overloads_no_q(
108  bp::args("Model","Data",
109  "Joint torque tau (size Model::nv)",
110  "Contact Jacobian J (size nb_constraint * Model::nv)",
111  "Contact drift gamma (size nb_constraint)",
112  "(double) Damping factor for cholesky decomposition of JMinvJt. Set to zero if constraints are full rank."),
113  "Solves the forward dynamics problem with contacts, puts the result in Data::ddq and return it. The contact forces are stored in data.lambda_c."
114  " Assumes pinocchio.computeAllTerms has been called."
115  ));
116 
117  bp::def("impulseDynamics",
119  impulseDynamics_overloads(
120  bp::args("Model","Data",
121  "Joint configuration q (size Model::nq)",
122  "Joint velocity before impact v_before (size Model::nv)",
123  "Contact Jacobian J (size nb_constraint * Model::nv)",
124  "Coefficient of restitution r_coeff (0 = rigid impact; 1 = fully elastic impact)",
125  "Damping factor when J is rank deficient."
126  ),
127  "Solves the impact dynamics problem with contacts, store the result in Data::dq_after and return it. The contact impulses are stored in data.impulse_c."
128  " Internally, pinocchio.crba is called."
129  ));
130 
131  bp::def("impulseDynamics",
133  impulseDynamics_overloads_no_q(
134  bp::args("Model","Data",
135  "Joint velocity before impact v_before (size Model::nv)",
136  "Contact Jacobian J (size nb_constraint * Model::nv)",
137  "Coefficient of restitution r_coeff (0 = rigid impact; 1 = fully elastic impact)",
138  "Damping factor when J is rank deficient."),
139  "Solves the impact dynamics problem with contacts, store the result in Data::dq_after and return it. The contact impulses are stored in data.impulse_c."
140  " Assumes pinocchio.crba has been called."
141  ));
142 
143  bp::def("computeKKTContactDynamicMatrixInverse",
145  computeKKTContactDynamicMatrixInverse_overload(bp::args("model","data","q","J","damping"),
146  "Computes the inverse of the constraint matrix [[M J^T], [J 0]]."));
147 
148  bp::def("getKKTContactDynamicMatrixInverse",
149  getKKTContactDynamicMatrixInverse_proxy,
150  bp::args("Model","Data",
151  "Contact Jacobian J(size nb_constraint * Model::nv)"),
152  "Computes the inverse of the constraint matrix [[M JT], [J 0]]. forward/impulseDynamics must be called first. The jacobian should be the same that was provided to forward/impulseDynamics.");
153  }
154 
155  } // namespace python
156 } // namespace pinocchio
static const Eigen::VectorXd impulseDynamics_proxy(const Model &model, Data &data, const Eigen::VectorXd &q, const Eigen::VectorXd &v_before, const Eigen::MatrixXd &J, const double r_coeff=0., const double inv_damping=0.)
static const Eigen::VectorXd forwardDynamics_proxy(const Model &model, Data &data, const Eigen::VectorXd &q, const Eigen::VectorXd &v, const Eigen::VectorXd &tau, const Eigen::MatrixXd &J, const Eigen::VectorXd &gamma, const double inv_damping=0.0)
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...
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...
BOOST_PYTHON_FUNCTION_OVERLOADS(computeKKTContactDynamicMatrixInverse_overload, computeKKTContactDynamicMatrixInverse_proxy, 4, 5) static const Eigen
data
static Eigen::MatrixXd computeKKTContactDynamicMatrixInverse_proxy(const Model &model, Data &data, const Eigen::VectorXd &q, const Eigen::MatrixXd &J, const double mu=0)
JointCollectionTpl const Eigen::MatrixBase< ConfigVectorType > & q
static const Eigen::VectorXd impulseDynamics_proxy_no_q(const Model &model, Data &data, const Eigen::VectorXd &v_before, const Eigen::MatrixXd &J, const double r_coeff=0., const double inv_damping=0.)
static const Eigen::VectorXd forwardDynamics_proxy_no_q(const Model &model, Data &data, const Eigen::VectorXd &tau, const Eigen::MatrixXd &J, const Eigen::VectorXd &gamma, const double inv_damping=0.0)
Main pinocchio namespace.
Definition: timings.cpp:30
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...
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > VectorXd
Definition: conversions.cpp:14
JointCollectionTpl & model


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