robotoc
robotoc - efficient ROBOT Optimal Control solvers
Loading...
Searching...
No Matches
robotoc::ImpactWrenchCone Class Referencefinal

Constraint on the wrench firction cone for surface contacts. More...

#include <impact_wrench_cone.hpp>

Inheritance diagram for robotoc::ImpactWrenchCone:
Collaboration diagram for robotoc::ImpactWrenchCone:

Public Member Functions

 ImpactWrenchCone (const Robot &robot, const double X, const double Y)
 Constructor. More...
 
 ImpactWrenchCone ()
 Default constructor. More...
 
 ~ImpactWrenchCone ()
 Destructor. More...
 
 ImpactWrenchCone (const ImpactWrenchCone &)=default
 Default copy constructor. More...
 
ImpactWrenchConeoperator= (const ImpactWrenchCone &)=default
 Default copy operator. More...
 
 ImpactWrenchCone (ImpactWrenchCone &&) noexcept=default
 Default move constructor. More...
 
ImpactWrenchConeoperator= (ImpactWrenchCone &&) noexcept=default
 Default move assign operator. More...
 
void setRectangular (const double X, const double Y)
 
KinematicsLevel kinematicsLevel () const override
 Checks the kinematics level of the constraint component. More...
 
void allocateExtraData (ConstraintComponentData &data) const override
 Allocates extra data in ConstraintComponentData. More...
 
bool isFeasible (Robot &robot, const ImpactStatus &impact_status, ConstraintComponentData &data, const SplitSolution &s) const override
 Checks whether the current solution s is feasible or not. More...
 
void setSlack (Robot &robot, const ImpactStatus &impact_status, ConstraintComponentData &data, const SplitSolution &s) const override
 Sets the slack variables of each constraint components. More...
 
void evalConstraint (Robot &robot, const ImpactStatus &impact_status, ConstraintComponentData &data, const SplitSolution &s) const override
 Computes the primal residual, residual in the complementary slackness, and the log-barrier function of the slack varible. More...
 
void evalDerivatives (Robot &robot, const ImpactStatus &impact_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, and add the product of the Jacobian and the dual variable to the KKT residual. This function is always called just after evalConstraint(). More...
 
void condenseSlackAndDual (const ImpactStatus &impact_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. This function is always called just after evalDerivatives(). More...
 
void expandSlackAndDual (const ImpactStatus &impact_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 directions of the primal variables. More...
 
int dimc () const override
 Returns the size of the constraint. More...
 
- Public Member Functions inherited from robotoc::ImpactConstraintComponentBase
 ImpactConstraintComponentBase (const double barrier=1.0e-03, const double fraction_to_boundary_rule=0.995)
 Constructor. More...
 
virtual ~ImpactConstraintComponentBase ()
 Destructor. More...
 
 ImpactConstraintComponentBase (const ImpactConstraintComponentBase &)=default
 Default copy constructor. More...
 
ImpactConstraintComponentBaseoperator= (const ImpactConstraintComponentBase &)=default
 Default copy operator. More...
 
 ImpactConstraintComponentBase (ImpactConstraintComponentBase &&) noexcept=default
 Default move constructor. More...
 
ImpactConstraintComponentBaseoperator= (ImpactConstraintComponentBase &&) noexcept=default
 Default move assign operator. More...
 
virtual KinematicsLevel kinematicsLevel () const =0
 Checks the kinematics level of the constraint component. More...
 
virtual void allocateExtraData (ConstraintComponentData &data) const =0
 Allocates extra data in ConstraintComponentData. More...
 
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. More...
 
virtual void setSlack (Robot &robot, const ImpactStatus &impact_status, ConstraintComponentData &data, const SplitSolution &s) const =0
 Sets the slack variables of each constraint components. More...
 
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 of the slack varible. More...
 
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, and add the product of the Jacobian and the dual variable to the KKT residual. This function is always called just after evalConstraint(). More...
 
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. This function is always called just after evalDerivatives(). More...
 
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 directions of the primal variables. More...
 
virtual int dimc () const =0
 Returns the size of the constraint. More...
 
void setSlackAndDualPositive (ConstraintComponentData &data) const
 Sets the slack and dual variables positive. More...
 
double maxSlackStepSize (const ConstraintComponentData &data) const
 Computes and returns the maximum step size by applying fraction-to-boundary-rule to the direction of the slack variable. More...
 
double maxDualStepSize (const ConstraintComponentData &data) const
 Computes and returns the maximum step size by applying fraction-to-boundary-rule to the direction of the dual variable. More...
 
virtual double getBarrierParam () const final
 Returns the barrier parameter. More...
 
virtual double getFractionToBoundaryRule () const final
 Returns the parameter of the fraction-to-boundary-rule. More...
 
virtual void setBarrierParam (const double barrier_param) final
 Sets the barrier parameter. More...
 
virtual void setFractionToBoundaryRule (const double fraction_to_boundary_rule) final
 Sets the parameter of the fraction-to-boundary-rule. More...
 
template<typename Derived >
std::shared_ptr< Derived > as_shared_ptr ()
 Gets the shared ptr of this object as the specified type. If this fails in dynamic casting, throws an exception. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from robotoc::ImpactConstraintComponentBase
static void updateSlack (ConstraintComponentData &data, const double step_size)
 Updates the slack variable according to the step size. More...
 
static void updateDual (ConstraintComponentData &data, const double step_size)
 Updates the dual variable according to the step size. More...
 
- Protected Member Functions inherited from robotoc::ImpactConstraintComponentBase
void computeComplementarySlackness (ConstraintComponentData &data) const
 Computes the residual in the complementarity slackness between
the slack and dual variables. More...
 
void computeComplementarySlackness (ConstraintComponentData &data, const int start, const int size) const
 Computes the residual in the complementarity slackness between
the slack and dual variables. More...
 
template<int Size>
void computeComplementarySlackness (ConstraintComponentData &data, const int start) const
 Computes the residual in the complementarity slackness between
the slack and dual variables. More...
 
double computeComplementarySlackness (const double slack, const double dual) const
 Computes the residual in the complementarity slackness between
the slack and dual variables. More...
 
template<typename VectorType >
double logBarrier (const Eigen::MatrixBase< VectorType > &slack) const
 Computes the log barrier function of the slack variable. More...
 
- Static Protected Member Functions inherited from robotoc::ImpactConstraintComponentBase
static void computeCondensingCoeffcient (ConstraintComponentData &data)
 Computes the coefficient of the condensing. More...
 
static void computeCondensingCoeffcient (ConstraintComponentData &data, const int start, const int size)
 Computes the coefficient of the condensing. More...
 
template<int Size>
static void computeCondensingCoeffcient (ConstraintComponentData &data, const int start)
 Computes the coefficient of the condensing. More...
 
static double computeCondensingCoeffcient (const double slack, const double dual, const double residual, const double cmpl)
 Computes the residual in the complementarity slackness between
the slack and dual variables. More...
 
static void computeDualDirection (ConstraintComponentData &data)
 Computes the direction of the dual variable from slack, primal residual, complementarity slackness, and the direction of the slack. More...
 
static void computeDualDirection (ConstraintComponentData &data, const int start, const int size)
 Computes the direction of the dual variable from slack, primal residual, complementarity slackness, and the direction of the slack. More...
 
template<int Size>
static void computeDualDirection (ConstraintComponentData &data, const int start)
 Computes the direction of the dual variable from slack, primal residual, complementary slackness, and the direction of the slack. More...
 
static double computeDualDirection (const double slack, const double dual, const double dslack, const double cmpl)
 Computes the direction of the dual variable from slack, primal residual, complementary slackness, and the direction of the slack. More...
 

Detailed Description

Constraint on the wrench firction cone for surface contacts.

Constructor & Destructor Documentation

◆ ImpactWrenchCone() [1/4]

robotoc::ImpactWrenchCone::ImpactWrenchCone ( const Robot robot,
const double  X,
const double  Y 
)

Constructor.

Parameters
[in]robotRobot model.
[in]XA length of the rectangular. Must be positive.
[in]YA length of the rectangular. Must be positive.

◆ ImpactWrenchCone() [2/4]

robotoc::ImpactWrenchCone::ImpactWrenchCone ( )

Default constructor.

◆ ~ImpactWrenchCone()

robotoc::ImpactWrenchCone::~ImpactWrenchCone ( )

Destructor.

◆ ImpactWrenchCone() [3/4]

robotoc::ImpactWrenchCone::ImpactWrenchCone ( const ImpactWrenchCone )
default

Default copy constructor.

◆ ImpactWrenchCone() [4/4]

robotoc::ImpactWrenchCone::ImpactWrenchCone ( ImpactWrenchCone &&  )
defaultnoexcept

Default move constructor.

Member Function Documentation

◆ allocateExtraData()

void robotoc::ImpactWrenchCone::allocateExtraData ( ConstraintComponentData data) const
overridevirtual

Allocates extra data in ConstraintComponentData.

Parameters
[in]dataConstraint component data.

Implements robotoc::ImpactConstraintComponentBase.

◆ condenseSlackAndDual()

void robotoc::ImpactWrenchCone::condenseSlackAndDual ( const ImpactStatus impact_status,
ConstraintComponentData data,
SplitKKTMatrix kkt_matrix,
SplitKKTResidual kkt_residual 
) const
overridevirtual

Condenses the slack and dual variables, i.e., factorizes the
condensed Hessians and KKT residuals. This function is always called just after evalDerivatives().

Parameters
[in]impact_statusImpact status.
[in]dataConstraints data.
[out]kkt_matrixImpact split KKT matrix. The condensed Hessians
are added to this data.
[out]kkt_residualImpact split KKT residual. The condensed KKT residual are added to this data.

Implements robotoc::ImpactConstraintComponentBase.

◆ dimc()

int robotoc::ImpactWrenchCone::dimc ( ) const
overridevirtual

Returns the size of the constraint.

Returns
Size of the constraints.

Implements robotoc::ImpactConstraintComponentBase.

◆ evalConstraint()

void robotoc::ImpactWrenchCone::evalConstraint ( Robot robot,
const ImpactStatus impact_status,
ConstraintComponentData data,
const SplitSolution s 
) const
overridevirtual

Computes the primal residual, residual in the complementary slackness, and the log-barrier function of the slack varible.

Parameters
[in]robotRobot model.
[in]impact_statusImpact status.
[in]dataConstraints data.
[in]sImpact split solution.

Implements robotoc::ImpactConstraintComponentBase.

◆ evalDerivatives()

void robotoc::ImpactWrenchCone::evalDerivatives ( Robot robot,
const ImpactStatus impact_status,
ConstraintComponentData data,
const SplitSolution s,
SplitKKTResidual kkt_residual 
) const
overridevirtual

Computes the derivatives of the priaml residual, i.e., the Jacobian of the inequality constraint, and add the product of the Jacobian and the dual variable to the KKT residual. This function is always called just after evalConstraint().

Parameters
[in]robotRobot model.
[in]impact_statusImpact status.
[in]dataConstraint data.
[in]sImpact split solution.
[out]kkt_residualImpact split KKT residual.

Implements robotoc::ImpactConstraintComponentBase.

◆ expandSlackAndDual()

void robotoc::ImpactWrenchCone::expandSlackAndDual ( const ImpactStatus impact_status,
ConstraintComponentData data,
const SplitDirection d 
) const
overridevirtual

Expands the slack and dual, i.e., computes the directions of the slack and dual variables from the directions of the primal variables.

Parameters
[in]impact_statusImpact status.
[in,out]dataConstraints data.
[in]dImpact split direction.

Implements robotoc::ImpactConstraintComponentBase.

◆ isFeasible()

bool robotoc::ImpactWrenchCone::isFeasible ( Robot robot,
const ImpactStatus impact_status,
ConstraintComponentData data,
const SplitSolution s 
) const
overridevirtual

Checks whether the current solution s is feasible or not.

Parameters
[in]robotRobot model.
[in]impact_statusImpact status.
[in]dataConstraint data.
[in]sImpact split solution.
Returns
true if s is feasible. false if not.

Implements robotoc::ImpactConstraintComponentBase.

◆ kinematicsLevel()

KinematicsLevel robotoc::ImpactWrenchCone::kinematicsLevel ( ) const
overridevirtual

Checks the kinematics level of the constraint component.

Returns
Kinematics level of the constraint component.

Implements robotoc::ImpactConstraintComponentBase.

◆ operator=() [1/2]

ImpactWrenchCone & robotoc::ImpactWrenchCone::operator= ( const ImpactWrenchCone )
default

Default copy operator.

◆ operator=() [2/2]

ImpactWrenchCone & robotoc::ImpactWrenchCone::operator= ( ImpactWrenchCone &&  )
defaultnoexcept

Default move assign operator.

◆ setRectangular()

void robotoc::ImpactWrenchCone::setRectangular ( const double  X,
const double  Y 
)
Parameters
[in]XA length of the rectangular. Must be positive.
[in]YA length of the rectangular. Must be positive.

◆ setSlack()

void robotoc::ImpactWrenchCone::setSlack ( Robot robot,
const ImpactStatus impact_status,
ConstraintComponentData data,
const SplitSolution s 
) const
overridevirtual

Sets the slack variables of each constraint components.

Parameters
[in]robotRobot model.
[in]impact_statusImpact status.
[out]dataConstraint data.
[in]sImpact split solution.

Implements robotoc::ImpactConstraintComponentBase.


The documentation for this class was generated from the following file: