1#ifndef ROBOTOC_ROBOT_HPP_
2#define ROBOTOC_ROBOT_HPP_
10#include "pinocchio/multibody/model.hpp"
11#include "pinocchio/multibody/data.hpp"
12#include "pinocchio/container/aligned-vector.hpp"
13#include "pinocchio/spatial/force.hpp"
80 template <typename TangentVectorType, typename ConfigVectorType>
82 const Eigen::MatrixBase<TangentVectorType>& v,
83 const
double integration_length,
84 const Eigen::MatrixBase<ConfigVectorType>& q) const;
95 template <typename ConfigVectorType1, typename TangentVectorType,
96 typename ConfigVectorType2>
98 const Eigen::MatrixBase<ConfigVectorType1>& q,
99 const Eigen::MatrixBase<TangentVectorType>& v,
100 const
double integration_length,
101 const Eigen::MatrixBase<ConfigVectorType2>& q_integrated) const;
113 template <typename ConfigVectorType, typename TangentVectorType,
114 typename MatrixType1, typename MatrixType2>
116 const Eigen::MatrixBase<ConfigVectorType>& q,
117 const Eigen::MatrixBase<TangentVectorType>& v,
118 const Eigen::MatrixBase<MatrixType1>& Jin,
119 const Eigen::MatrixBase<MatrixType2>& Jout) const;
131 template <typename ConfigVectorType, typename TangentVectorType,
132 typename MatrixType1, typename MatrixType2>
134 const Eigen::MatrixBase<ConfigVectorType>& q,
135 const Eigen::MatrixBase<TangentVectorType>& v,
136 const Eigen::MatrixBase<MatrixType1>& Jin,
137 const Eigen::MatrixBase<MatrixType2>& Jout) const;
148 template <typename ConfigVectorType1, typename ConfigVectorType2,
149 typename TangentVectorType>
151 const Eigen::MatrixBase<ConfigVectorType1>& qf,
152 const Eigen::MatrixBase<ConfigVectorType2>& q0,
153 const Eigen::MatrixBase<TangentVectorType>& qdiff) const;
163 template <typename ConfigVectorType1, typename ConfigVectorType2,
166 const Eigen::MatrixBase<ConfigVectorType1>& qf,
167 const Eigen::MatrixBase<ConfigVectorType2>& q0,
168 const Eigen::MatrixBase<MatrixType>& dqdiff_dqf) const;
178 template <typename ConfigVectorType1, typename ConfigVectorType2,
181 const Eigen::MatrixBase<ConfigVectorType1>& qf,
182 const Eigen::MatrixBase<ConfigVectorType2>& q0,
183 const Eigen::MatrixBase<MatrixType>& dqdiff_dq0) const;
193 template <typename ConfigVectorType1, typename ConfigVectorType2,
194 typename ConfigVectorType3>
196 const Eigen::MatrixBase<ConfigVectorType1>& q1,
197 const Eigen::MatrixBase<ConfigVectorType2>& q2, const
double t,
198 const Eigen::MatrixBase<ConfigVectorType3>& qout) const;
205 template <typename ConfigVectorType, typename MatrixType>
207 const Eigen::MatrixBase<MatrixType>& J) const;
216 template <typename ConfigVectorType, typename TangentVectorType1,
217 typename TangentVectorType2>
219 const Eigen::MatrixBase<TangentVectorType1>& v,
220 const Eigen::MatrixBase<TangentVectorType2>& a);
228 template <typename ConfigVectorType, typename TangentVectorType>
230 const Eigen::MatrixBase<TangentVectorType>& v);
237 template <typename ConfigVectorType>
247 template <typename ConfigVectorType, typename TangentVectorType1,
248 typename TangentVectorType2>
250 const Eigen::MatrixBase<TangentVectorType1>& v,
251 const Eigen::MatrixBase<TangentVectorType2>& a);
259 template <typename ConfigVectorType, typename TangentVectorType>
261 const Eigen::MatrixBase<TangentVectorType>& v);
268 template <typename ConfigVectorType>
323 const Eigen::Vector3d&
CoM() const;
335 const pinocchio::ReferenceFrame reference_frame=pinocchio::LOCAL_WORLD_ALIGNED) const;
346 const std::
string& frame_name,
347 const pinocchio::ReferenceFrame reference_frame=pinocchio::LOCAL_WORLD_ALIGNED) const;
359 const pinocchio::ReferenceFrame reference_frame=pinocchio::LOCAL_WORLD_ALIGNED) const;
370 const std::
string& frame_name,
371 const pinocchio::ReferenceFrame reference_frame=pinocchio::LOCAL_WORLD_ALIGNED) const;
383 const pinocchio::ReferenceFrame reference_frame=pinocchio::LOCAL_WORLD_ALIGNED) const;
394 const std::
string& frame_name,
395 const pinocchio::ReferenceFrame reference_frame=pinocchio::LOCAL_WORLD_ALIGNED) const;
411 template <typename MatrixType>
413 const Eigen::MatrixBase<MatrixType>& J);
420 template <typename MatrixType>
421 void getCoMJacobian(const Eigen::MatrixBase<MatrixType>& J) const;
432 template <typename Vector3dType>
434 const
int frame_id, const Eigen::Vector3d& vec_local,
435 const Eigen::MatrixBase<Vector3dType>& vec_world) const;
447 template <typename Vector3dType, typename MatrixType>
449 const
int frame_id, const Eigen::MatrixBase<Vector3dType>& vec_world,
450 const Eigen::MatrixBase<MatrixType>& J);
460 template <typename VectorType>
463 const Eigen::MatrixBase<VectorType>& baumgarte_residual);
477 template <typename MatrixType1, typename MatrixType2, typename MatrixType3>
480 const Eigen::MatrixBase<MatrixType1>& baumgarte_partial_dq,
481 const Eigen::MatrixBase<MatrixType2>& baumgarte_partial_dv,
482 const Eigen::MatrixBase<MatrixType3>& baumgarte_partial_da);
491 template <typename VectorType>
494 const Eigen::MatrixBase<VectorType>& velocity_residual) const;
505 template <typename MatrixType1, typename MatrixType2>
508 const Eigen::MatrixBase<MatrixType1>& velocity_partial_dq,
509 const Eigen::MatrixBase<MatrixType2>& velocity_partial_dv);
518 template <typename VectorType>
521 const Eigen::MatrixBase<VectorType>& position_residual);
531 template <typename MatrixType>
534 const Eigen::MatrixBase<MatrixType>& position_partial_dq);
565 template <typename ConfigVectorType, typename TangentVectorType1,
566 typename TangentVectorType2, typename TangentVectorType3>
567 void RNEA(const Eigen::MatrixBase<ConfigVectorType>& q,
568 const Eigen::MatrixBase<TangentVectorType1>& v,
569 const Eigen::MatrixBase<TangentVectorType2>& a,
570 const Eigen::MatrixBase<TangentVectorType3>& tau);
590 template <typename ConfigVectorType, typename TangentVectorType1,
591 typename TangentVectorType2, typename MatrixType1,
592 typename MatrixType2, typename MatrixType3>
594 const Eigen::MatrixBase<TangentVectorType1>& v,
595 const Eigen::MatrixBase<TangentVectorType2>& a,
596 const Eigen::MatrixBase<MatrixType1>& dRNEA_partial_dq,
597 const Eigen::MatrixBase<MatrixType2>& dRNEA_partial_dv,
598 const Eigen::MatrixBase<MatrixType3>& dRNEA_partial_da);
609 template <typename ConfigVectorType, typename TangentVectorType1,
610 typename TangentVectorType2>
611 void RNEAImpact(const Eigen::MatrixBase<ConfigVectorType>& q,
612 const Eigen::MatrixBase<TangentVectorType1>& dv,
613 const Eigen::MatrixBase<TangentVectorType2>& res);
628 template <typename ConfigVectorType, typename TangentVectorType,
629 typename MatrixType1, typename MatrixType2>
631 const Eigen::MatrixBase<ConfigVectorType>& q,
632 const Eigen::MatrixBase<TangentVectorType>& dv,
633 const Eigen::MatrixBase<MatrixType1>& dRNEA_partial_dq,
634 const Eigen::MatrixBase<MatrixType2>& dRNEA_partial_ddv);
643 template <typename MatrixType1, typename MatrixType2>
644 void computeMinv(const Eigen::MatrixBase<MatrixType1>& M,
645 const Eigen::MatrixBase<MatrixType2>& Minv);
657 template <typename MatrixType1, typename MatrixType2, typename MatrixType3>
659 const Eigen::MatrixBase<MatrixType2>& J,
660 const Eigen::MatrixBase<MatrixType3>& MJtJinv);
672 template <typename ConfigVectorType>
674 const Eigen::MatrixBase<ConfigVectorType>& q) const;
682 int frameId(const std::
string& frame_name) const;
873 const Eigen::VectorXd& lower_joint_position_limit);
881 const Eigen::VectorXd& upper_joint_position_limit);
905 void disp(std::ostream& os) const;
907 friend std::ostream& operator<<(std::ostream& os, const
Robot& robot);
909 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
915 pinocchio::Model model_, impact_model_;
916 pinocchio::Data data_, impact_data_;
918 Eigen::MatrixXd dimpact_dv_;
923 int dimq_, dimv_, dimu_, dim_passive_, max_dimf_, max_num_contacts_;
925 Eigen::VectorXd joint_effort_limit_, joint_velocity_limit_,
926 lower_joint_position_limit_, upper_joint_position_limit_;
931#include "robotoc/robot/robot.hxx"
Impact status of robot model. Wrapper of ContactStatus to treat impacts.
Definition: impact_status.hpp:21
Dynamics and kinematics model of robots. Wraps pinocchio::Model and pinocchio::Data....
Definition: robot.hpp:32
void RNEA(const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType1 > &v, const Eigen::MatrixBase< TangentVectorType2 > &a, const Eigen::MatrixBase< TangentVectorType3 > &tau)
Computes inverse dynamics, i.e., generalized torques corresponding for given configuration,...
Definition: robot.hxx:524
std::vector< std::string > pointContactFrameNames() const
Retruns the names of the frames involving the contacts.
void integrateConfiguration(const Eigen::MatrixBase< TangentVectorType > &v, const double integration_length, const Eigen::MatrixBase< ConfigVectorType > &q) const
Integrates the generalized velocity, that is, performs q <- q + integration_length * v like computati...
Definition: robot.hxx:23
const RobotModelInfo & robotModelInfo() const
Gets the robot model info.
int dimv() const
Returns the dimensiton of the velocity, i.e, tangent space.
void computeImpactVelocityDerivatives(const ImpactStatus &impact_status, const Eigen::MatrixBase< MatrixType1 > &velocity_partial_dq, const Eigen::MatrixBase< MatrixType2 > &velocity_partial_dv)
Computes the partial derivatives of the impact velocity constraint. Before calling this function,...
Definition: robot.hxx:396
const Eigen::Vector3d & CoM() const
Returns the position of the center of mass. Before calling this function, updateKinematics() or updat...
ImpactStatus createImpactStatus() const
Creates a ImpactStatus for this robot model.
void computeImpactVelocityResidual(const ImpactStatus &impact_status, const Eigen::MatrixBase< VectorType > &velocity_residual) const
Computes the residual of the impact velocity constraint. Before calling this function,...
Definition: robot.hxx:366
const Eigen::Vector3d & CoMVelocity() const
Returns the linear velocity of the center of mass expressed in the world frame. Before calling this f...
void subtractConfiguration(const Eigen::MatrixBase< ConfigVectorType1 > &qf, const Eigen::MatrixBase< ConfigVectorType2 > &q0, const Eigen::MatrixBase< TangentVectorType > &qdiff) const
Computes qf - q0 at the tangent space. The result means that the unit velocity from initial configura...
Definition: robot.hxx:97
const RobotProperties & robotProperties() const
Gets a collectio of the properties of this robot model.
std::vector< int > surfaceContactFrames() const
Retruns the indices of the frames involving the surface contacts.
double totalWeight() const
Returns the total weight of this robot model.
void RNEADerivatives(const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType1 > &v, const Eigen::MatrixBase< TangentVectorType2 > &a, const Eigen::MatrixBase< MatrixType1 > &dRNEA_partial_dq, const Eigen::MatrixBase< MatrixType2 > &dRNEA_partial_dv, const Eigen::MatrixBase< MatrixType3 > &dRNEA_partial_da)
Computes the partial dervatives of the function of inverse dynamics with respect to the configuration...
Definition: robot.hxx:550
Vector6d frameSpatialVelocity(const int frame_id, const pinocchio::ReferenceFrame reference_frame=pinocchio::LOCAL_WORLD_ALIGNED) const
Returns the spatial velocity of the frame. Before calling this function, updateKinematics() or upda...
void getJacobianTransformFromLocalToWorld(const int frame_id, const Eigen::MatrixBase< Vector3dType > &vec_world, const Eigen::MatrixBase< MatrixType > &J)
Jacobian of the transformation mapping of a 3D quantity from local frame to the world frame.
Definition: robot.hxx:275
const Eigen::Vector3d & framePosition(const int frame_id) const
Returns the position of the frame. Before calling this function, updateKinematics() or updateFrameKin...
void setImpactForces(const ImpactStatus &impact_status, const std::vector< Vector6d > &f)
Set impact forces in this robot model for each active impacts.
Definition: robot.hxx:494
Eigen::Vector3d frameAngularVelocity(const int frame_id, const pinocchio::ReferenceFrame reference_frame=pinocchio::LOCAL_WORLD_ALIGNED) const
Returns the angular velocity of the frame. Before calling this function, updateKinematics() or upda...
int dim_passive() const
Returns the dimensiton of the generalized torques corresponding to the passive joints.
void dSubtractConfiguration_dq0(const Eigen::MatrixBase< ConfigVectorType1 > &qf, const Eigen::MatrixBase< ConfigVectorType2 > &q0, const Eigen::MatrixBase< MatrixType > &dqdiff_dq0) const
Computes the partial derivative of the function of subtraction qf - q0 w.r.t. the initial configurati...
Definition: robot.hxx:128
void dIntegrateTransport_dv(const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType > &v, const Eigen::MatrixBase< MatrixType1 > &Jin, const Eigen::MatrixBase< MatrixType2 > &Jout) const
Transport the Jacobian w.r.t. the generalized velocity evaluated at the integrated configuration q + ...
Definition: robot.hxx:78
void setLowerJointPositionLimit(const Eigen::VectorXd &lower_joint_position_limit)
Sets the lower limit of the position of each joints.
int maxNumSurfaceContacts() const
Returns the maximum number of the surface contacts.
void setContactForces(const ContactStatus &contact_status, const std::vector< Vector6d > &f)
Set contact forces in this robot model for each active contacts.
const SE3 & framePlacement(const int frame_id) const
Returns the SE(3) of the frame. Before calling this function, updateKinematics() or updateFrameKinema...
Eigen::VectorXd generateFeasibleConfiguration() const
Generates feasible configuration randomly.
void computeMJtJinv(const Eigen::MatrixBase< MatrixType1 > &M, const Eigen::MatrixBase< MatrixType2 > &J, const Eigen::MatrixBase< MatrixType3 > &MJtJinv)
Computes the inverse of the contact dynamics matrix [[M J^T], [J O]].
Definition: robot.hxx:642
int maxNumContacts() const
Returns the maximum number of the contacts.
void computeMinv(const Eigen::MatrixBase< MatrixType1 > &M, const Eigen::MatrixBase< MatrixType2 > &Minv)
Computes the inverse of the joint inertia matrix M.
Definition: robot.hxx:628
const Eigen::Matrix3d & frameRotation(const int frame_id) const
Returns the rotation matrix of the frame. Before calling this function, updateKinematics() or updat...
void computeContactPositionDerivative(const ImpactStatus &impact_status, const Eigen::MatrixBase< MatrixType > &position_partial_dq)
Computes the partial derivative of the contact position at the impact. Before calling this function,...
Definition: robot.hxx:464
void updateFrameKinematics(const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType1 > &v, const Eigen::MatrixBase< TangentVectorType2 > &a)
Updates the frame kinematics of the robot. The frame placements, velocities, and accelerations are ca...
Definition: robot.hxx:212
std::vector< ContactType > contactTypes() const
Returns the types of the contacts.
void setRobotProperties(const RobotProperties &properties)
Sets a collection of the properties for this robot model.
Eigen::VectorXd lowerJointPositionLimit() const
Returns the lower limit of the position of each joints.
ContactType contactType(const int contact_index) const
Returns the type of the contact.
void dIntegrateTransport_dq(const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType > &v, const Eigen::MatrixBase< MatrixType1 > &Jin, const Eigen::MatrixBase< MatrixType2 > &Jout) const
Transport the Jacobian w.r.t. the configuration evaluated at the integrated configuration q + v to q.
Definition: robot.hxx:59
Eigen::VectorXd jointVelocityLimit() const
Returns the joint velocity limit of each joints.
void setJointVelocityLimit(const Eigen::VectorXd &joint_velocity_limit)
Sets the joint velocity limit of each joints.
std::vector< int > pointContactFrames() const
Retruns the indices of the frames involving the contacts.
std::vector< int > contactFrames() const
Retruns the indices of the frames involving the point or surface contacts.
~Robot()=default
Default destructor.
int maxNumPointContacts() const
Returns the maximum number of the point contacts.
void integrateCoeffWiseJacobian(const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< MatrixType > &J) const
Transforms the Jacobian w.r.t. configuration to tangentspace.
Definition: robot.hxx:159
void computeBaumgarteResidual(const ContactStatus &contact_status, const Eigen::MatrixBase< VectorType > &baumgarte_residual)
Computes the residual of the contact constriants represented by Baumgarte's stabilization method....
Definition: robot.hxx:291
void initializeJointLimits()
Initializes the results of jointEffortLimit(), jointVelocityLimit(), lowerJointPositionLimit(),...
Eigen::VectorXd upperJointPositionLimit() const
Returns the upper limit of the position of each joints.
void transformFromLocalToWorld(const int frame_id, const Eigen::Vector3d &vec_local, const Eigen::MatrixBase< Vector3dType > &vec_world) const
Transforms 3D quantity from local frame to the world frame.
Definition: robot.hxx:265
void computeContactPositionResidual(const ImpactStatus &impact_status, const Eigen::MatrixBase< VectorType > &position_residual)
Computes the residual of the contact position constraint at the impact. Before calling this function,...
Definition: robot.hxx:434
std::vector< std::string > surfaceContactFrameNames() const
Retruns the names of the frames involving the surface contacts.
Robot & operator=(const Robot &)=default
Default copy assign operator.
void disp(std::ostream &os) const
Displays the robot model onto a ostream.
Eigen::Vector3d frameLinearVelocity(const int frame_id, const pinocchio::ReferenceFrame reference_frame=pinocchio::LOCAL_WORLD_ALIGNED) const
Returns the linear (translational) velocity of the frame. Before calling this function,...
int max_dimf() const
Returns the maximum dimension of the contacts.
Robot()
Default constructor.
int frameId(const std::string &frame_name) const
Gets the id of the specified frame.
void RNEAImpactDerivatives(const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType > &dv, const Eigen::MatrixBase< MatrixType1 > &dRNEA_partial_dq, const Eigen::MatrixBase< MatrixType2 > &dRNEA_partial_ddv)
Computes the partial dervatives of the function of impuse dynamics with respect to configuration and ...
Definition: robot.hxx:604
void getFrameJacobian(const int frame_id, const Eigen::MatrixBase< MatrixType > &J)
Computes the Jacobian of the frame position expressed in the local coordinate. Before calling this fu...
Definition: robot.hxx:247
Eigen::VectorXd jointEffortLimit() const
Returns the effort limit of each joints.
void updateKinematics(const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType1 > &v, const Eigen::MatrixBase< TangentVectorType2 > &a)
Updates the kinematics of the robot. The frame placements, frame velocity, frame acceleration,...
Definition: robot.hxx:172
Robot(Robot &&) noexcept=default
Default move constructor.
bool hasFloatingBase() const
Returns true if the robot has a floating base and false if not.
void getCoMJacobian(const Eigen::MatrixBase< MatrixType > &J) const
Gets the Jacobian of the position of the center of mass. Before calling this function,...
Definition: robot.hxx:257
Eigen::Matrix< double, 6, 1 > Vector6d
Definition: robot.hpp:34
void setUpperJointPositionLimit(const Eigen::VectorXd &upper_joint_position_limit)
Sets the upper limit of the position of each joints.
int dimu() const
Returns the dimensiton of the velocity, i.e, tangent space.
void setJointEffortLimit(const Eigen::VectorXd &joint_effort_limit)
Sets the effort limit of each joints.
std::vector< std::string > contactFrameNames() const
Retruns the names of the contact frames involving the point or.
std::string frameName(const int frame_id) const
Gets the name of the specified frame.
Robot(const RobotModelInfo &robot_model_info)
Constructs a robot model according to the input model info.
int dimq() const
Returns the dimensiton of the configuration.
Robot(const Robot &)=default
Default copy constructor.
void RNEAImpact(const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType1 > &dv, const Eigen::MatrixBase< TangentVectorType2 > &res)
Computes the residual of the impact dynamics for given configuration and impact change in the general...
Definition: robot.hxx:589
double totalMass() const
Returns the total mass of this robot model.
void computeBaumgarteDerivatives(const ContactStatus &contact_status, 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 constriants represented by the Baumgarte's stabilizat...
Definition: robot.hxx:321
void normalizeConfiguration(const Eigen::MatrixBase< ConfigVectorType > &q) const
Normalizes a configuration vector.
Definition: robot.hxx:688
ContactStatus createContactStatus() const
Creates a ContactStatus for this robot model.
void dSubtractConfiguration_dqf(const Eigen::MatrixBase< ConfigVectorType1 > &qf, const Eigen::MatrixBase< ConfigVectorType2 > &q0, const Eigen::MatrixBase< MatrixType > &dqdiff_dqf) const
Computes the partial derivative of the function of subtraction qf - q0 w.r.t. the terminal configurat...
Definition: robot.hxx:112
void interpolateConfiguration(const Eigen::MatrixBase< ConfigVectorType1 > &q1, const Eigen::MatrixBase< ConfigVectorType2 > &q2, const double t, const Eigen::MatrixBase< ConfigVectorType3 > &qout) const
Computes an interpolation of two configurations by (1.0 - t) * s1 + t * s2.
Definition: robot.hxx:144
Definition: constraint_component_base.hpp:17
ContactType
Types of contacts.
Definition: contact_status.hpp:22
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 robot model.
Definition: robot_model_info.hpp:24
Collection of the robot properties, which can change after constructing robot models.
Definition: robot_properties.hpp:30