expose-contact-dynamics.cpp
Go to the documentation of this file.
1 //
2 // Copyright (c) 2020-2021 INRIA
3 //
4 
7 
9 
10 namespace pinocchio
11 {
12  namespace python
13  {
14 
16  const context::Model & model,
18  const context::VectorXs & q,
19  const context::VectorXs & v,
20  const context::VectorXs & tau,
21  const context::MatrixXs & J,
22  const context::VectorXs & gamma,
24  {
25 
28  return forwardDynamics(model, data, q, v, tau, J, gamma, inv_damping);
30  }
31 
33  const context::Model & model,
35  const context::VectorXs & tau,
36  const context::MatrixXs & J,
37  const context::VectorXs & gamma,
39  {
40 
43  return forwardDynamics(model, data, tau, J, gamma, inv_damping);
45  }
46 
48  const context::Model & model,
50  const context::VectorXs & q,
51  const context::VectorXs & v_before,
52  const context::MatrixXs & J,
55  {
56 
59  return impulseDynamics(model, data, q, v_before, J, r_coeff, inv_damping);
61  }
62 
64  const context::Model & model,
66  const context::VectorXs & v_before,
67  const context::MatrixXs & J,
70  {
71 
74  return impulseDynamics(model, data, v_before, J, r_coeff, inv_damping);
76  }
77 
79  const context::Model & model,
81  const context::VectorXs & q,
82  const context::MatrixXs & J,
84  {
85  context::MatrixXs KKTMatrix_inv(model.nv + J.rows(), model.nv + J.rows());
87  return KKTMatrix_inv;
88  }
89 
90  static const context::MatrixXs getKKTContactDynamicMatrixInverse_proxy(
91  const context::Model & model, context::Data & data, const context::MatrixXs & J)
92  {
93  context::MatrixXs MJtJ_inv(model.nv + J.rows(), model.nv + J.rows());
94 
99 
100  return MJtJ_inv;
101  }
102 
104  {
105  using namespace Eigen;
106 
107  bp::def(
108  "forwardDynamics", &forwardDynamics_proxy,
109  (bp::arg("model"), bp::arg("data"), bp::arg("q"), bp::arg("v"), bp::arg("tau"),
110  bp::arg("constraint_jacobian"), bp::arg("constraint_drift"), bp::arg("damping") = 0),
111  "Solves the constrained dynamics problem with contacts, puts the result in "
112  "context::Data::ddq and return it. The contact forces are stored in data.lambda_c.\n"
113  "Note: internally, pinocchio.computeAllTerms is called.",
115 
116  bp::def(
117  "forwardDynamics", &forwardDynamics_proxy_no_q,
118  (bp::arg("model"), bp::arg("data"), bp::arg("tau"), bp::arg("constraint_jacobian"),
119  bp::arg("constraint_drift"), bp::arg("damping") = 0),
120  "Solves the forward dynamics problem with contacts, puts the result in "
121  "context::Data::ddq and return it. The contact forces are stored in data.lambda_c.\n"
122  "Note: this function assumes that pinocchio.computeAllTerms has been called first.",
124 
125  bp::def(
126  "impulseDynamics", &impulseDynamics_proxy,
127  (bp::arg("model"), bp::arg("data"), bp::arg("q"), bp::arg("v_before"),
128  bp::arg("constraint_jacobian"), bp::arg("restitution_coefficient") = 0,
129  bp::arg("damping") = 0),
130  "Solves the impact dynamics problem with contacts, store the result in "
131  "context::Data::dq_after and return it. The contact impulses are stored in "
132  "data.impulse_c.\n"
133  "Note: internally, pinocchio.crba is called.",
135 
136  bp::def(
137  "impulseDynamics", &impulseDynamics_proxy_no_q,
138  (bp::arg("model"), bp::arg("data"), bp::arg("v_before"), bp::arg("constraint_jacobian"),
139  bp::arg("restitution_coefficient") = 0, bp::arg("damping") = 0),
140  "Solves the impact dynamics problem with contacts, store the result in "
141  "context::Data::dq_after and return it. The contact impulses are stored in "
142  "data.impulse_c.\n"
143  "Note: this function assumes that pinocchio.crba has been called first.",
145 
146  bp::def(
147  "computeKKTContactDynamicMatrixInverse", computeKKTContactDynamicMatrixInverse_proxy,
148  (bp::arg("model"), bp::arg("data"), bp::arg("q"), bp::arg("constraint_jacobian"),
149  bp::arg("damping") = 0),
150  "Computes the inverse of the constraint matrix [[M J^T], [J 0]].",
152 
153  bp::def(
154  "getKKTContactDynamicMatrixInverse", getKKTContactDynamicMatrixInverse_proxy,
155  bp::args("model", "data", "constraint_jacobian"),
156  "Computes the inverse of the constraint matrix [[M Jt], [J 0]].\n forwardDynamics or "
157  "impulseDynamics must have been called first.\n"
158  "Note: the constraint Jacobian should be the same that was provided to "
159  "forwardDynamics or impulseDynamics.",
161  }
162 
163  } // namespace python
164 } // namespace pinocchio
pinocchio::getKKTContactDynamicMatrixInverse
PINOCCHIO_DEPRECATED 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.
Eigen
pinocchio::python::getKKTContactDynamicMatrixInverse_proxy
static const context::MatrixXs getKKTContactDynamicMatrixInverse_proxy(const context::Model &model, context::Data &data, const context::MatrixXs &J)
Definition: expose-contact-dynamics.cpp:90
pinocchio::python::forwardDynamics_proxy_no_q
static const context::VectorXs forwardDynamics_proxy_no_q(const context::Model &model, context::Data &data, const context::VectorXs &tau, const context::MatrixXs &J, const context::VectorXs &gamma, const context::Scalar inv_damping=context::Scalar(0.0))
Definition: expose-contact-dynamics.cpp:32
pinocchio::DataTpl
Definition: context/generic.hpp:25
pinocchio::python::exposeContactDynamics
void exposeContactDynamics()
Definition: expose-contact-dynamics.cpp:103
pinocchio::python::impulseDynamics_proxy_no_q
static const context::VectorXs impulseDynamics_proxy_no_q(const context::Model &model, context::Data &data, const context::VectorXs &v_before, const context::MatrixXs &J, const context::Scalar r_coeff=context::Scalar(0.), const context::Scalar inv_damping=context::Scalar(0.))
Definition: expose-contact-dynamics.cpp:63
setup.data
data
Definition: cmake/cython/setup.in.py:48
pinocchio::python::impulseDynamics_proxy
static const context::VectorXs impulseDynamics_proxy(const context::Model &model, context::Data &data, const context::VectorXs &q, const context::VectorXs &v_before, const context::MatrixXs &J, const context::Scalar r_coeff=context::Scalar(0.), const context::Scalar inv_damping=context::Scalar(0.))
Definition: expose-contact-dynamics.cpp:47
PINOCCHIO_COMPILER_DIAGNOSTIC_POP
#define PINOCCHIO_COMPILER_DIAGNOSTIC_POP
Definition: include/pinocchio/macros.hpp:134
inverse-kinematics-3d.J
J
Definition: inverse-kinematics-3d.py:28
pinocchio::computeKKTContactDynamicMatrixInverse
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...
pinocchio::python::VectorXs
context::VectorXs VectorXs
Definition: admm-solver.cpp:30
bindings_dynamics.r_coeff
float r_coeff
Definition: bindings_dynamics.py:10
forward-dynamics-derivatives.tau
tau
Definition: forward-dynamics-derivatives.py:25
algorithms.hpp
python
pinocchio::python::v
const Vector3Like & v
Definition: bindings/python/spatial/explog.hpp:66
pinocchio::python::computeKKTContactDynamicMatrixInverse_proxy
static context::MatrixXs computeKKTContactDynamicMatrixInverse_proxy(const context::Model &model, context::Data &data, const context::VectorXs &q, const context::MatrixXs &J, const context::Scalar mu=context::Scalar(0))
Definition: expose-contact-dynamics.cpp:78
pinocchio::forwardDynamics
const PINOCCHIO_DEPRECATED 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...
bindings_dynamics.inv_damping
float inv_damping
Definition: bindings_dynamics.py:11
PINOCCHIO_COMPILER_DIAGNOSTIC_IGNORED_DEPRECECATED_DECLARATIONS
#define PINOCCHIO_COMPILER_DIAGNOSTIC_IGNORED_DEPRECECATED_DECLARATIONS
Definition: include/pinocchio/macros.hpp:135
PINOCCHIO_COMPILER_DIAGNOSTIC_PUSH
#define PINOCCHIO_COMPILER_DIAGNOSTIC_PUSH
macros for pragma push/pop/ignore deprecated warnings
Definition: include/pinocchio/macros.hpp:133
pinocchio::q
JointCollectionTpl const Eigen::MatrixBase< ConfigVectorType > & q
Definition: joint-configuration.hpp:1083
pinocchio::python::mimic_not_supported_function
Definition: model-checker.hpp:22
mu
double mu
Definition: delassus.cpp:58
pinocchio::python::forwardDynamics_proxy
static const context::VectorXs forwardDynamics_proxy(const context::Model &model, context::Data &data, const context::VectorXs &q, const context::VectorXs &v, const context::VectorXs &tau, const context::MatrixXs &J, const context::VectorXs &gamma, const context::Scalar inv_damping=context::Scalar(0.0))
Definition: expose-contact-dynamics.cpp:15
pinocchio::impulseDynamics
const PINOCCHIO_DEPRECATED 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.
contact-dynamics.hpp
model-checker.hpp
pinocchio::ModelTpl
Definition: context/generic.hpp:20
pinocchio::python::context::Scalar
PINOCCHIO_PYTHON_SCALAR_TYPE Scalar
Definition: bindings/python/context/generic.hpp:37
pinocchio::model
JointCollectionTpl & model
Definition: joint-configuration.hpp:1082
pinocchio
Main pinocchio namespace.
Definition: timings.cpp:33


pinocchio
Author(s):
autogenerated on Thu Apr 10 2025 02:42:18