robotoc
robotoc - efficient ROBOT Optimal Control solvers
Loading...
Searching...
No Matches
point_contact.hpp
Go to the documentation of this file.
1#ifndef ROBOTOC_POINT_CONTACT_HPP_
2#define ROBOTOC_POINT_CONTACT_HPP_
3
4#include "Eigen/Core"
5#include "pinocchio/multibody/model.hpp"
6#include "pinocchio/multibody/data.hpp"
7#include "pinocchio/container/aligned-vector.hpp"
8#include "pinocchio/spatial/force.hpp"
9
11#include "robotoc/robot/se3.hpp"
12
13#include <iostream>
14
15
16namespace robotoc {
17
23public:
24 using Matrix6xd = Eigen::Matrix<double, 6, Eigen::Dynamic>;
25
32 PointContact(const pinocchio::Model& model,
33 const ContactModelInfo& contact_model_info);
34
39
43 ~PointContact() = default;
44
48 PointContact(const PointContact&) = default;
49
54
58 PointContact(PointContact&&) noexcept = default;
59
63 PointContact& operator=(PointContact&&) noexcept = default;
64
72 const Eigen::Vector3d& contact_force,
73 pinocchio::container::aligned_vector<pinocchio::Force>& joint_forces) const;
74
86 template <typename VectorType1, typename VectorType2>
88 const pinocchio::Model& model, const pinocchio::Data& data,
89 const Eigen::MatrixBase<VectorType1>& desired_contact_position,
90 const Eigen::MatrixBase<VectorType2>& baumgarte_residual) const;
91
106 template <typename MatrixType1, typename MatrixType2, typename MatrixType3>
108 const pinocchio::Model& model, pinocchio::Data& data,
109 const Eigen::MatrixBase<MatrixType1>& baumgarte_partial_dq,
110 const Eigen::MatrixBase<MatrixType2>& baumgarte_partial_dv,
111 const Eigen::MatrixBase<MatrixType3>& baumgarte_partial_da);
112
121 template <typename VectorType>
123 const pinocchio::Model& model, const pinocchio::Data& data,
124 const Eigen::MatrixBase<VectorType>& velocity_residual) const;
125
137 template <typename MatrixType1, typename MatrixType2>
139 const pinocchio::Model& model, pinocchio::Data& data,
140 const Eigen::MatrixBase<MatrixType1>& velocity_partial_dq,
141 const Eigen::MatrixBase<MatrixType2>& velocity_partial_dv);
142
153 template <typename VectorType1, typename VectorType2>
155 const pinocchio::Model& model, const pinocchio::Data& data,
156 const Eigen::MatrixBase<VectorType1>& desired_contact_position,
157 const Eigen::MatrixBase<VectorType2>& position_residual) const;
158
169 template <typename MatrixType>
171 const pinocchio::Model& model, pinocchio::Data& data,
172 const Eigen::MatrixBase<MatrixType>& position_partial_dq);
173
181 const Eigen::Vector3d& contactPosition(const pinocchio::Data& data) const;
182
190 void setBaumgarteGains(const double baumgarte_position_gain,
191 const double baumgarte_velocity_gain);
192
197 int contactFrameId() const;
198
204 int parentJointId() const;
205
211
215 void disp(std::ostream& os) const;
216
217 friend std::ostream& operator<<(std::ostream& os,
218 const PointContact& point_contact);
219
220 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
221
222private:
223 ContactModelInfo info_;
224 int contact_frame_id_, parent_joint_id_, dimv_;
225 SE3 jXf_;
226 pinocchio::Motion v_frame_;
227 Eigen::Matrix3d v_linear_skew_, v_angular_skew_;
228 Matrix6xd J_frame_, frame_v_partial_dq_, frame_a_partial_dq_,
229 frame_a_partial_dv_, frame_a_partial_da_;
230};
231
232} // namespace robotoc
233
234#include "robotoc/robot/point_contact.hxx"
235
236#endif // ROBOTOC_POINT_CONTACT_HPP_
Kinematics and dynamic model of a point contact.
Definition: point_contact.hpp:22
void computeBaumgarteResidual(const pinocchio::Model &model, const pinocchio::Data &data, const Eigen::MatrixBase< VectorType1 > &desired_contact_position, const Eigen::MatrixBase< VectorType2 > &baumgarte_residual) const
Computes the residual of the contact constraints considered by the Baumgarte's stabilization method....
Definition: point_contact.hxx:17
int contactFrameId() const
Returns contact frame id, i.e., the index of the contact frame.
PointContact(PointContact &&) noexcept=default
Default move constructor.
PointContact(const PointContact &)=default
Default copy constructor.
void computeContactVelocityDerivatives(const pinocchio::Model &model, pinocchio::Data &data, const Eigen::MatrixBase< MatrixType1 > &velocity_partial_dq, const Eigen::MatrixBase< MatrixType2 > &velocity_partial_dv)
Computes the partial derivatives of the contact velocity constraints. Before calling this function,...
Definition: point_contact.hxx:101
PointContact(const pinocchio::Model &model, const ContactModelInfo &contact_model_info)
Construct a point contact model.
const ContactModelInfo & contactModelInfo() const
Gets the contact model info.
void computeContactPositionResidual(const pinocchio::Model &model, const pinocchio::Data &data, const Eigen::MatrixBase< VectorType1 > &desired_contact_position, const Eigen::MatrixBase< VectorType2 > &position_residual) const
Computes the residual of the contact position constraints. Before calling this function,...
Definition: point_contact.hxx:121
void computeJointForceFromContactForce(const Eigen::Vector3d &contact_force, pinocchio::container::aligned_vector< pinocchio::Force > &joint_forces) const
Converts the 3D contact forces in the local coordinate to the corresponding joint spatial forces.
void computeBaumgarteDerivatives(const pinocchio::Model &model, pinocchio::Data &data, const Eigen::MatrixBase< MatrixType1 > &baumgarte_partial_dq, const Eigen::MatrixBase< MatrixType2 > &baumgarte_partial_dv, const Eigen::MatrixBase< MatrixType3 > &baumgarte_partial_da)
Computes the partial derivatives of the contact constraints considered by the Baumgarte's stabilizati...
Definition: point_contact.hxx:37
~PointContact()=default
Default destructor.
void computeContactVelocityResidual(const pinocchio::Model &model, const pinocchio::Data &data, const Eigen::MatrixBase< VectorType > &velocity_residual) const
Computes the residual of the contact velocity constraints. Before calling this function,...
Definition: point_contact.hxx:90
void computeContactPositionDerivative(const pinocchio::Model &model, pinocchio::Data &data, const Eigen::MatrixBase< MatrixType > &position_partial_dq)
Computes the partial derivative of the contact position constraint with respect to the configuratio...
Definition: point_contact.hxx:133
int parentJointId() const
Returns parent joint id, i.e., the index of the parent joint of the contact frame.
PointContact & operator=(const PointContact &)=default
Default copy assign operator.
PointContact()
Default constructor.
void disp(std::ostream &os) const
Displays the point contact onto a ostream.
void setBaumgarteGains(const double baumgarte_position_gain, const double baumgarte_velocity_gain)
Sets the gain parameters of the Baumgarte's stabilization method.
const Eigen::Vector3d & contactPosition(const pinocchio::Data &data) const
Returns the contact position at the current kinematics of the robot. Before calling this function,...
Eigen::Matrix< double, 6, Eigen::Dynamic > Matrix6xd
Definition: point_contact.hpp:24
Definition: constraint_component_base.hpp:17
std::vector< T, Eigen::aligned_allocator< T > > aligned_vector
std vector with Eigen::aligned_allocator.
Definition: aligned_vector.hpp:14
pinocchio::SE3 SE3
Using pinocchio::SE3 without its namespace.
Definition: se3.hpp:15
Info of a contact model.
Definition: contact_model_info.hpp:12