robotoc
robotoc - efficient ROBOT Optimal Control solvers
Loading...
Searching...
No Matches
joint_acceleration_upper_limit.hpp
Go to the documentation of this file.
1#ifndef ROBOTOC_JOINT_ACCELERATION_UPPER_LIMIT_HPP_
2#define ROBOTOC_JOINT_ACCELERATION_UPPER_LIMIT_HPP_
3
4#include "Eigen/Core"
5
14
15
16namespace robotoc {
17
23public:
29 JointAccelerationUpperLimit(const Robot& robot, const Eigen::VectorXd& amax);
30
35
40
45
50 const JointAccelerationUpperLimit&) = default;
51
56
61 JointAccelerationUpperLimit&&) noexcept = default;
62
64
65 void allocateExtraData(ConstraintComponentData& data) const override {}
66
67 bool isFeasible(Robot& robot, const ContactStatus& contact_status,
69 const SplitSolution& s) const override;
70
71 void setSlack(Robot& robot, const ContactStatus& contact_status,
73 const SplitSolution& s) const override;
74
75 void evalConstraint(Robot& robot, const ContactStatus& contact_status,
77 const SplitSolution& s) const override;
78
79 void evalDerivatives(Robot& robot, const ContactStatus& contact_status,
81 SplitKKTResidual& kkt_residual) const override;
82
83 void condenseSlackAndDual(const ContactStatus& contact_status,
85 SplitKKTMatrix& kkt_matrix,
86 SplitKKTResidual& kkt_residual) const override;
87
88 void expandSlackAndDual(const ContactStatus& contact_status,
90 const SplitDirection& d) const override;
91
92 int dimc() const override;
93
94private:
95 int dimc_, dim_passive_;
96 Eigen::VectorXd amax_;
97
98};
99
100} // namespace robotoc
101
102#endif // ROBOTOC_JOINT_ACCELERATION_UPPER_LIMIT_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
Contact status of robot model.
Definition: contact_status.hpp:32
Constraint on the upper limits of the joint acceleration.
Definition: joint_acceleration_upper_limit.hpp:22
JointAccelerationUpperLimit(JointAccelerationUpperLimit &&) noexcept=default
Default move constructor.
JointAccelerationUpperLimit(const Robot &robot, const Eigen::VectorXd &amax)
Constructor.
JointAccelerationUpperLimit()
Default constructor.
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,...
KinematicsLevel kinematicsLevel() const override
Checks the kinematics level of the constraint component.
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 allocateExtraData(ConstraintComponentData &data) const override
Allocates extra data in ConstraintComponentData.
Definition: joint_acceleration_upper_limit.hpp:65
JointAccelerationUpperLimit & operator=(const JointAccelerationUpperLimit &)=default
Default copy operator.
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.
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...
int dimc() const override
Returns the size of the constraint.
JointAccelerationUpperLimit(const JointAccelerationUpperLimit &)=default
Default copy constructor.
void setSlack(Robot &robot, const ContactStatus &contact_status, ConstraintComponentData &data, const SplitSolution &s) const override
Sets the slack variables of each constraint components.
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