robotoc
robotoc - efficient ROBOT Optimal Control solvers
Loading...
Searching...
No Matches
impact_constraint_component_base.hpp
Go to the documentation of this file.
1#ifndef ROBOTOC_IMPACT_CONSTRAINT_COMPONENT_BASE_HPP_
2#define ROBOTOC_IMPACT_CONSTRAINT_COMPONENT_BASE_HPP_
3
4#include <memory>
5
6#include "Eigen/Core"
7
16
17
18namespace robotoc {
19
24class ImpactConstraintComponentBase : public std::enable_shared_from_this<ImpactConstraintComponentBase> {
25public:
34 ImpactConstraintComponentBase(const double barrier=1.0e-03,
35 const double fraction_to_boundary_rule=0.995);
36
41
46 const ImpactConstraintComponentBase&) = default;
47
52 const ImpactConstraintComponentBase&) = default;
53
58 ImpactConstraintComponentBase&&) noexcept = default;
59
64 ImpactConstraintComponentBase&&) noexcept = default;
65
70 virtual KinematicsLevel kinematicsLevel() const = 0;
71
76 virtual void allocateExtraData(ConstraintComponentData& data) const = 0;
77
86 virtual bool isFeasible(Robot& robot, const ImpactStatus& impact_status,
88 const SplitSolution& s) const = 0;
89
97 virtual void setSlack(Robot& robot, const ImpactStatus& impact_status,
99 const SplitSolution& s) const = 0;
100
109 virtual void evalConstraint(Robot& robot, const ImpactStatus& impact_status,
111 const SplitSolution& s) const = 0;
112
124 virtual void evalDerivatives(Robot& robot, const ImpactStatus& impact_status,
126 const SplitSolution& s,
127 SplitKKTResidual& kkt_residual) const = 0;
128
141 const ImpactStatus& impact_status, ConstraintComponentData& data,
142 SplitKKTMatrix& kkt_matrix,
143 SplitKKTResidual& kkt_residual) const = 0;
144
152 virtual void expandSlackAndDual(const ImpactStatus& impact_status,
154 const SplitDirection& d) const = 0;
155
160 virtual int dimc() const = 0;
161
167
174 double maxSlackStepSize(const ConstraintComponentData& data) const;
175
182 double maxDualStepSize(const ConstraintComponentData& data) const;
183
189 static void updateSlack(ConstraintComponentData& data, const double step_size);
190
196 static void updateDual(ConstraintComponentData& data, const double step_size);
197
201 virtual double getBarrierParam() const final;
202
206 virtual double getFractionToBoundaryRule() const final;
207
212 virtual void setBarrierParam(const double barrier_param) final;
213
220 const double fraction_to_boundary_rule) final;
221
228 template <typename Derived>
229 std::shared_ptr<Derived> as_shared_ptr();
230
231protected:
238
247 const int start,
248 const int size) const;
249
257 template <int Size>
259 const int start) const;
260
268 double computeComplementarySlackness(const double slack,
269 const double dual) const;
270
276
284 const int start, const int size);
285
292 template <int Size>
294 const int start);
295
305 static double computeCondensingCoeffcient(const double slack,
306 const double dual,
307 const double residual,
308 const double cmpl);
309
316
325 const int start, const int size);
326
334 template <int Size>
336 const int start);
337
347 static double computeDualDirection(const double slack, const double dual,
348 const double dslack, const double cmpl);
349
355 template <typename VectorType>
356 double logBarrier(const Eigen::MatrixBase<VectorType>& slack) const;
357
358private:
359 double barrier_, fraction_to_boundary_rule_;
360
361};
362
363} // namespace robotoc
364
365#include "robotoc/constraints/impact_constraint_component_base.hxx"
366
367#endif // ROBOTOC_IMPACT_CONSTRAINT_COMPONENT_BASE_HPP_
Data used in constraint components. Composed by slack, dual (Lagrange multiplier),...
Definition: constraint_component_data.hpp:17
Base class for impact constraint components.
Definition: impact_constraint_component_base.hpp:24
static void computeDualDirection(ConstraintComponentData &data)
Computes the direction of the dual variable from slack, primal residual, complementarity slackness,...
Definition: impact_constraint_component_base.hxx:100
virtual void condenseSlackAndDual(const ImpactStatus &impact_status, ConstraintComponentData &data, SplitKKTMatrix &kkt_matrix, SplitKKTResidual &kkt_residual) const =0
Condenses the slack and dual variables, i.e., factorizes the condensed Hessians and KKT residuals....
virtual bool isFeasible(Robot &robot, const ImpactStatus &impact_status, ConstraintComponentData &data, const SplitSolution &s) const =0
Checks whether the current solution s is feasible or not.
static void computeCondensingCoeffcient(ConstraintComponentData &data)
Computes the coefficient of the condensing.
Definition: impact_constraint_component_base.hxx:74
virtual void evalConstraint(Robot &robot, const ImpactStatus &impact_status, ConstraintComponentData &data, const SplitSolution &s) const =0
Computes the primal residual, residual in the complementary slackness, and the log-barrier function o...
virtual void setSlack(Robot &robot, const ImpactStatus &impact_status, ConstraintComponentData &data, const SplitSolution &s) const =0
Sets the slack variables of each constraint components.
virtual double getBarrierParam() const final
Returns the barrier parameter.
virtual double getFractionToBoundaryRule() const final
Returns the parameter of the fraction-to-boundary-rule.
virtual void setFractionToBoundaryRule(const double fraction_to_boundary_rule) final
Sets the parameter of the fraction-to-boundary-rule.
double maxSlackStepSize(const ConstraintComponentData &data) const
Computes and returns the maximum step size by applying fraction-to-boundary-rule to the direction of ...
Definition: impact_constraint_component_base.hxx:12
ImpactConstraintComponentBase & operator=(const ImpactConstraintComponentBase &)=default
Default copy operator.
void setSlackAndDualPositive(ConstraintComponentData &data) const
Sets the slack and dual variables positive.
void computeComplementarySlackness(ConstraintComponentData &data) const
Computes the residual in the complementarity slackness between the slack and dual variables.
Definition: impact_constraint_component_base.hxx:49
ImpactConstraintComponentBase(ImpactConstraintComponentBase &&) noexcept=default
Default move constructor.
double logBarrier(const Eigen::MatrixBase< VectorType > &slack) const
Computes the log barrier function of the slack variable.
Definition: impact_constraint_component_base.hxx:127
ImpactConstraintComponentBase(const ImpactConstraintComponentBase &)=default
Default copy constructor.
double maxDualStepSize(const ConstraintComponentData &data) const
Computes and returns the maximum step size by applying fraction-to-boundary-rule to the direction of ...
Definition: impact_constraint_component_base.hxx:18
virtual void expandSlackAndDual(const ImpactStatus &impact_status, ConstraintComponentData &data, const SplitDirection &d) const =0
Expands the slack and dual, i.e., computes the directions of the slack and dual variables from the di...
virtual int dimc() const =0
Returns the size of the constraint.
ImpactConstraintComponentBase(const double barrier=1.0e-03, const double fraction_to_boundary_rule=0.995)
Constructor.
virtual void allocateExtraData(ConstraintComponentData &data) const =0
Allocates extra data in ConstraintComponentData.
static void updateSlack(ConstraintComponentData &data, const double step_size)
Updates the slack variable according to the step size.
Definition: impact_constraint_component_base.hxx:24
static void updateDual(ConstraintComponentData &data, const double step_size)
Updates the dual variable according to the step size.
Definition: impact_constraint_component_base.hxx:31
virtual void setBarrierParam(const double barrier_param) final
Sets the barrier parameter.
virtual KinematicsLevel kinematicsLevel() const =0
Checks the kinematics level of the constraint component.
virtual void evalDerivatives(Robot &robot, const ImpactStatus &impact_status, ConstraintComponentData &data, const SplitSolution &s, SplitKKTResidual &kkt_residual) const =0
Computes the derivatives of the priaml residual, i.e., the Jacobian of the inequality constraint,...
std::shared_ptr< Derived > as_shared_ptr()
Gets the shared ptr of this object as the specified type. If this fails in dynamic casting,...
Definition: impact_constraint_component_base.hxx:39
virtual ~ImpactConstraintComponentBase()
Destructor.
Definition: impact_constraint_component_base.hpp:40
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
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