robotoc
robotoc - efficient ROBOT Optimal Control solvers
Loading...
Searching...
No Matches
friction_cone.hpp
Go to the documentation of this file.
1#ifndef ROBOTOC_FRICTION_CONE_HPP_
2#define ROBOTOC_FRICTION_CONE_HPP_
3
4#include "Eigen/Core"
5
14
15
16namespace robotoc {
17
23public:
24 using Vector5d = Eigen::Matrix<double, 5, 1>;
25
30 FrictionCone(const Robot& robot);
31
36
41
45 FrictionCone(const FrictionCone&) = default;
46
51
55 FrictionCone(FrictionCone&&) noexcept = default;
56
60 FrictionCone& operator=(FrictionCone&&) noexcept = default;
61
63
64 void allocateExtraData(ConstraintComponentData& data) const override;
65
66 bool isFeasible(Robot& robot, const ContactStatus& contact_status,
68 const SplitSolution& s) const override;
69
70 void setSlack(Robot& robot, const ContactStatus& contact_status,
72 const SplitSolution& s) const override;
73
74 void evalConstraint(Robot& robot, const ContactStatus& contact_status,
76 const SplitSolution& s) const override;
77
78 void evalDerivatives(Robot& robot, const ContactStatus& contact_status,
80 SplitKKTResidual& kkt_residual) const override;
81
82 void condenseSlackAndDual(const ContactStatus& contact_status,
84 SplitKKTMatrix& kkt_matrix,
85 SplitKKTResidual& kkt_residual) const override;
86
87 void expandSlackAndDual(const ContactStatus& contact_status,
89 const SplitDirection& d) const override;
90
91 int dimc() const override;
92
101 template <typename VectorType1, typename VectorType2>
102 static void frictionConeResidual(const double mu,
103 const Eigen::MatrixBase<VectorType1>& f_world,
104 const Eigen::Matrix3d& R_surface,
105 const Eigen::MatrixBase<VectorType2>& res) {
106 assert(mu > 0);
107 assert(f_world.size() == 3);
108 assert((R_surface*R_surface.transpose()).isIdentity());
109 assert(res.size() == 5);
110 const Eigen::Vector3d f_local = R_surface.transpose() * f_world;
111 const_cast<Eigen::MatrixBase<VectorType2>&>(res).coeffRef(0) = - f_local.coeff(2);
112 const_cast<Eigen::MatrixBase<VectorType2>&>(res).coeffRef(1)
113 = f_local.coeff(0) - mu * f_local.coeff(2) / std::sqrt(2);
114 const_cast<Eigen::MatrixBase<VectorType2>&>(res).coeffRef(2)
115 = - f_local.coeff(0) - mu * f_local.coeff(2) / std::sqrt(2);
116 const_cast<Eigen::MatrixBase<VectorType2>&>(res).coeffRef(3)
117 = f_local.coeff(1) - mu * f_local.coeff(2) / std::sqrt(2);
118 const_cast<Eigen::MatrixBase<VectorType2>&>(res).coeffRef(4)
119 = - f_local.coeff(1) - mu * f_local.coeff(2) / std::sqrt(2);
120 }
121
122 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
123
124private:
125 int dimv_, dimc_, max_num_contacts_;
126 std::vector<int> contact_frame_;
127 std::vector<ContactType> contact_types_;
128
129 static Eigen::VectorXd& fW(ConstraintComponentData& data,
130 const int contact_idx) {
131 return data.r[contact_idx];
132 }
133
134 Eigen::VectorXd& r(ConstraintComponentData& data,
135 const int contact_idx) const {
136 return data.r[max_num_contacts_+contact_idx];
137 }
138
139 static Eigen::MatrixXd& dg_dq(ConstraintComponentData& data,
140 const int contact_idx) {
141 return data.J[contact_idx];
142 }
143
144 Eigen::MatrixXd& dg_df(ConstraintComponentData& data,
145 const int contact_idx) const {
146 return data.J[max_num_contacts_+contact_idx];
147 }
148
149 Eigen::MatrixXd& dfW_dq(ConstraintComponentData& data,
150 const int contact_idx) const {
151 return data.J[2*max_num_contacts_+contact_idx];
152 }
153
154 Eigen::MatrixXd& r_dg_df(ConstraintComponentData& data,
155 const int contact_idx) const {
156 return data.J[3*max_num_contacts_+contact_idx];
157 }
158
159 Eigen::MatrixXd& cone_local(ConstraintComponentData& data,
160 const int contact_idx) const {
161 return data.J[4*max_num_contacts_+contact_idx];
162 }
163
164 Eigen::MatrixXd& cone_world(ConstraintComponentData& data,
165 const int contact_idx) const {
166 return data.J[5*max_num_contacts_+contact_idx];
167 }
168
169};
170
171} // namespace robotoc
172
173#endif // ROBOTOC_FRICTION_CONE_HPP_
Base class for constraint components.
Definition: constraint_component_base.hpp:34
Data used in constraint components. Composed by slack, dual (Lagrange multiplier),...
Definition: constraint_component_data.hpp:17
std::vector< Eigen::VectorXd > r
std vector of Eigen::VectorXd used to store residual temporaly. Only be allocated in ConstraintCompon...
Definition: constraint_component_data.hpp:108
Contact status of robot model.
Definition: contact_status.hpp:32
Constraint on the inner-approximated firction cone.
Definition: friction_cone.hpp:22
void expandSlackAndDual(const ContactStatus &contact_status, ConstraintComponentData &data, const SplitDirection &d) const override
Expands the slack and dual, i.e., computes the directions of the slack and dual variables from the di...
void evalConstraint(Robot &robot, const ContactStatus &contact_status, ConstraintComponentData &data, const SplitSolution &s) const override
Computes the primal residual, residual in the complementary slackness, and the log-barrier function o...
FrictionCone(const Robot &robot)
Constructor.
int dimc() const override
Returns the size of the constraint.
bool isFeasible(Robot &robot, const ContactStatus &contact_status, ConstraintComponentData &data, const SplitSolution &s) const override
Checks whether the current solution s is feasible or not.
KinematicsLevel kinematicsLevel() const override
Checks the kinematics level of the constraint component.
FrictionCone()
Default constructor.
Eigen::Matrix< double, 5, 1 > Vector5d
Definition: friction_cone.hpp:24
static void frictionConeResidual(const double mu, const Eigen::MatrixBase< VectorType1 > &f_world, const Eigen::Matrix3d &R_surface, const Eigen::MatrixBase< VectorType2 > &res)
Computes the friction cone residual.
Definition: friction_cone.hpp:102
void allocateExtraData(ConstraintComponentData &data) const override
Allocates extra data in ConstraintComponentData.
void condenseSlackAndDual(const ContactStatus &contact_status, ConstraintComponentData &data, SplitKKTMatrix &kkt_matrix, SplitKKTResidual &kkt_residual) const override
Condenses the slack and dual variables, i.e., factorizes the condensed Hessians and KKT residuals....
void evalDerivatives(Robot &robot, const ContactStatus &contact_status, ConstraintComponentData &data, const SplitSolution &s, SplitKKTResidual &kkt_residual) const override
Computes the derivatives of the priaml residual, i.e., the Jacobian of the inequality constraint,...
~FrictionCone()
Destructor.
void setSlack(Robot &robot, const ContactStatus &contact_status, ConstraintComponentData &data, const SplitSolution &s) const override
Sets the slack variables of each constraint components.
FrictionCone(FrictionCone &&) noexcept=default
Default move constructor.
FrictionCone & operator=(const FrictionCone &)=default
Default copy operator.
FrictionCone(const FrictionCone &)=default
Default copy constructor.
Dynamics and kinematics model of robots. Wraps pinocchio::Model and pinocchio::Data....
Definition: robot.hpp:32
Newton direction of the solution to the optimal control problem split into a time stage.
Definition: split_direction.hpp:20
The KKT matrix split into a time stage.
Definition: split_kkt_matrix.hpp:18
KKT residual split into each time stage.
Definition: split_kkt_residual.hpp:18
Solution to the optimal control problem split into a time stage.
Definition: split_solution.hpp:20
Definition: constraint_component_base.hpp:17
KinematicsLevel
Kinematics level of the constraint component used in ConstraintComponentBase.
Definition: constraint_component_base.hpp:24