hpp-constraints  6.0.0
Definition of basic geometric constraints for motion planning
hpp::constraints::solver::HierarchicalIterative Class Reference

More...

#include <hpp/constraints/solver/hierarchical-iterative.hh>

Inheritance diagram for hpp::constraints::solver::HierarchicalIterative:
Collaboration diagram for hpp::constraints::solver::HierarchicalIterative:

Classes

struct  Data
 

Public Types

enum  Status { ERROR_INCREASED , MAX_ITERATION_REACHED , INFEASIBLE , SUCCESS }
 
typedef Eigen::RowBlockIndices Indices_t
 
typedef lineSearch::FixedSequence DefaultLineSearch
 
typedef shared_ptr< saturation::BaseSaturation_t
 

Public Member Functions

 HierarchicalIterative (const LiegroupSpacePtr_t &configSpace)
 
 HierarchicalIterative (const HierarchicalIterative &other)
 
virtual ~HierarchicalIterative ()
 
value_type residualError () const
 Returns the squared norm of the error vector. More...
 
void residualError (vectorOut_t error) const
 Returns the error vector. More...
 
bool definesSubmanifoldOf (const HierarchicalIterative &solver) const
 
virtual std::ostream & print (std::ostream &os) const
 
template<typename LineSearchType >
solver::HierarchicalIterative::Status solve (vectorOut_t arg, LineSearchType lineSearch) const
 
Problem definition
const LiegroupSpacePtr_tconfigSpace () const
 Get configuration space on which constraints are defined. More...
 
virtual bool contains (const ImplicitPtr_t &numericalConstraint) const
 
virtual bool add (const ImplicitPtr_t &constraint, const std::size_t &priority)
 
virtual void merge (const HierarchicalIterative &other)
 
void saturation (const Saturation_t &saturate)
 Set the saturation function. More...
 
const Saturation_tsaturation () const
 Get the saturation function. More...
 
Problem resolution
template<typename LineSearchType >
Status solve (vectorOut_t arg, LineSearchType ls=LineSearchType()) const
 
Status solve (vectorOut_t arg) const
 
bool isSatisfied (vectorIn_t arg) const
 
bool isSatisfied (vectorIn_t arg, value_type errorThreshold) const
 
bool isConstraintSatisfied (const ImplicitPtr_t &constraint, vectorIn_t arg, vectorOut_t error, bool &constraintFound) const
 
const value_typesigma () const
 
Parameters
void freeVariables (const segments_t intervals)
 
void freeVariables (const Indices_t &indices)
 
const Indices_tfreeVariables () const
 Get free velocity variables. More...
 
void maxIterations (size_type iterations)
 Set maximal number of iterations. More...
 
size_type maxIterations () const
 Get maximal number of iterations in config projector. More...
 
void errorThreshold (const value_type &threshold)
 Set error threshold. More...
 
value_type errorThreshold () const
 Get error threshold. More...
 
value_type squaredErrorThreshold () const
 Get error threshold. More...
 
value_type inequalityThreshold () const
 Get the inequality threshold. More...
 
void inequalityThreshold (const value_type &it)
 set the inequality threshold More...
 
void lastIsOptional (bool optional)
 
bool lastIsOptional () const
 
void solveLevelByLevel (bool solveLevelByLevel)
 
bool solveLevelByLevel () const
 
Stack
const ImplicitConstraintSetconstraints (const std::size_t priority)
 Get set of constraints for a give priority level. More...
 
const NumericalConstraints_tconstraints () const
 Get constraints (implicit and explicit) More...
 
std::size_t numberStacks () const
 
const size_typedimension () const
 
const size_typereducedDimension () const
 
ArrayXb activeParameters () const
 Configuration parameters involved in the constraint resolution. More...
 
ArrayXb activeDerivativeParameters () const
 Velocity parameters involved in the constraint resolution. More...
 
Right hand side accessors
vector_t rightHandSideFromConfig (ConfigurationIn_t config)
 
virtual bool rightHandSideFromConfig (const ImplicitPtr_t &constraint, ConfigurationIn_t config)
 
virtual bool rightHandSide (const ImplicitPtr_t &constraint, vectorIn_t rhs)
 
virtual bool getRightHandSide (const ImplicitPtr_t &constraint, vectorOut_t rhs) const
 Get right hand side of a constraints. More...
 
virtual void rightHandSide (vectorIn_t rhs)
 
void rightHandSideAt (const value_type &s)
 
vector_t rightHandSide () const
 
size_type rightHandSideSize () const
 
Access to internal data

You should know what you do when you call these functions

template<bool ComputeJac>
void computeValue (vectorIn_t arg) const
 Compute the value of each level, and the jacobian if ComputeJac is true. More...
 
void computeSaturation (vectorIn_t arg) const
 
void getValue (vectorOut_t v) const
 
void getReducedJacobian (matrixOut_t J) const
 
void computeError () const
 
const vector_tlastStep () const
 Accessor to the last step done. More...
 
virtual bool integrate (vectorIn_t from, vectorIn_t velocity, vectorOut_t result) const
 

Protected Types

typedef Eigen::JacobiSVD< matrix_tSVD_t
 

Protected Member Functions

void update ()
 
virtual void computeActiveRowsOfJ (std::size_t iStack)
 
void computeDescentDirection () const
 
void expandDqSmall () const
 
void saturate (vectorOut_t arg) const
 
 HierarchicalIterative ()
 

Protected Attributes

value_type squaredErrorThreshold_
 
value_type inequalityThreshold_
 
size_type maxIterations_
 
std::vector< ImplicitConstraintSetstacks_
 
LiegroupSpacePtr_t configSpace_
 
size_type dimension_
 
size_type reducedDimension_
 
bool lastIsOptional_
 
bool solveLevelByLevel_
 
Indices_t freeVariables_
 Unknown of the set of implicit constraints. More...
 
Saturation_t saturate_
 
NumericalConstraints_t constraints_
 Members moved from core::ConfigProjector. More...
 
std::map< DifferentiableFunctionPtr_t, size_typeiq_
 Value rank of constraint in its priority level. More...
 
std::map< DifferentiableFunctionPtr_t, size_typeiv_
 Derivative rank of constraint in its priority level. More...
 
std::map< DifferentiableFunctionPtr_t, std::size_t > priority_
 Priority level of constraint. More...
 
value_type sigma_
 The smallest non-zero singular value. More...
 
vector_t dq_
 
vector_t dqSmall_
 
matrix_t reducedJ_
 
Eigen::VectorXi saturation_
 
Eigen::VectorXi reducedSaturation_
 
Configuration_t qSat_
 
ArrayXb tmpSat_
 
value_type squaredNorm_
 
std::vector< Datadatas_
 
SVD_t svd_
 
vector_t OM_
 
vector_t OP_
 

Friends

struct lineSearch::Backtracking
 

Detailed Description

Solve a system of non-linear equations on a robot configuration

The non-linear system of equations is built by adding equations with method HierarchicalIterative::add.

Note that a hierarchy between the equations can be provided. In this case, the solver will try to solve the highest priority equations first, and then to solve the lower priority equations. Note that priorities are in decreasing order: 0 has higher priority than 1.

The algorithm used is a Newton-Raphson like algorithm that works as follows: for a single level of priority, let \(f (\mathbf{q}) = 0\) be the system of equations where \(f\) is a \(C^1\) mapping from the robot configuration space to a Lie group space \(\mathcal{L}\).

Starting from initial guess \(\mathbf{q}_0\), the method HierarchicalIterative::solve builds a sequence of configurations \(\mathbf{q}_i\) as follows:

\begin{eqnarray*} \mathbf{q}_{i+1} = \mathbf{q}_i - \alpha_i \frac{\partial f}{\partial \mathbf{q}}(\mathbf{q}_i)^{+} f (\mathbf{q}_i) \end{eqnarray*}

where

The computation of the direction of descent in the case of multiple level of hierarchy is described in this paper.

The error threshold can be accessed by methods HierarchicalIterative::errorThreshold. The maximal number of iterations can be accessed by methods HierarchicalIterative::maxIterations.

Note
Solving equations one after the other

For some applications, it can be more efficient to solve a set of equations one after the other. In other words, an equation is ignored until the previous one is solved (norm below the threshold). To do so, introduce the equations using method HierarchicalIterative::add with increasing value of priority, and call method HierarchicalIterative::solveLevelByLevel(true).

Note
Lie group

The unknowns \(\mathbf{q}\) may take values in a more general set than the configuration space of a robot. This set should be a Cartesian product of Lie groups: hpp::pinocchio::LiegroupSpace.

Note
Saturation
Right hand side and comparison types

Instead of \(f(\mathbf{q}) = 0\), other constraints can be defined. Several comparison types are available:

  • Equality: \(f(\mathbf{q}) = rhs\), where \(rhs\) is a parameterizable right hand side,
  • EqualToZero: \(f(\mathbf{q}) = 0\),
  • Superior: \(f(\mathbf{q}) > 0\)
  • Inferior: \(f(\mathbf{q}) < 0\) If several constraint are of type equality, the right hand side of the system of equations can be modified by methods HierarchicalIterative::rightHandSideFromInput, HierarchicalIterative::rightHandSide.
Note
Free variables

Some variables can be locked, or computed explicitely. In this case, the iterative resolution will only change the other variables called free variables.

See also
methods

Member Typedef Documentation

◆ DefaultLineSearch

◆ Indices_t

◆ Saturation_t

◆ SVD_t

Member Enumeration Documentation

◆ Status

Enumerator
ERROR_INCREASED 
MAX_ITERATION_REACHED 
INFEASIBLE 
SUCCESS 

Constructor & Destructor Documentation

◆ HierarchicalIterative() [1/3]

hpp::constraints::solver::HierarchicalIterative::HierarchicalIterative ( const LiegroupSpacePtr_t configSpace)

◆ HierarchicalIterative() [2/3]

hpp::constraints::solver::HierarchicalIterative::HierarchicalIterative ( const HierarchicalIterative other)

◆ ~HierarchicalIterative()

virtual hpp::constraints::solver::HierarchicalIterative::~HierarchicalIterative ( )
inlinevirtual

◆ HierarchicalIterative() [3/3]

hpp::constraints::solver::HierarchicalIterative::HierarchicalIterative ( )
inlineprotected

Member Function Documentation

◆ activeDerivativeParameters()

ArrayXb hpp::constraints::solver::HierarchicalIterative::activeDerivativeParameters ( ) const

Velocity parameters involved in the constraint resolution.

◆ activeParameters()

ArrayXb hpp::constraints::solver::HierarchicalIterative::activeParameters ( ) const

Configuration parameters involved in the constraint resolution.

◆ add()

virtual bool hpp::constraints::solver::HierarchicalIterative::add ( const ImplicitPtr_t constraint,
const std::size_t &  priority 
)
virtual

Add an implicit constraint

Parameters
constraintimplicit constraint
prioritylevel of priority of the constraint: priority are in decreasing order: 0 is the highest priority level,

Reimplemented in hpp::constraints::solver::BySubstitution.

◆ computeActiveRowsOfJ()

virtual void hpp::constraints::solver::HierarchicalIterative::computeActiveRowsOfJ ( std::size_t  iStack)
protectedvirtual

Compute which rows of the jacobian of stack_[iStack] are not zero, using the activeDerivativeParameters of the functions. The result is stored in datas_[i].activeRowsOfJ

Reimplemented in hpp::constraints::solver::BySubstitution.

◆ computeDescentDirection()

void hpp::constraints::solver::HierarchicalIterative::computeDescentDirection ( ) const
protected

Compute a SVD decomposition of each level and find the best descent direction at the first order. Linearization of the system of equations rhs - v_{i} = J (q_i) (dq_{i+1} - q_{i}) q_{i+1} - q_{i} = J(q_i)^{+} ( rhs - v_{i} ) dq = J(q_i)^{+} ( rhs - v_{i} )

Warning
computeValue<true> must have been called first.

◆ computeError()

void hpp::constraints::solver::HierarchicalIterative::computeError ( ) const

If lastIsOptional() is true, then the last level is ignored.

Warning
computeValue must have been called first.

◆ computeSaturation()

void hpp::constraints::solver::HierarchicalIterative::computeSaturation ( vectorIn_t  arg) const

◆ computeValue()

template<bool ComputeJac>
void hpp::constraints::solver::HierarchicalIterative::computeValue ( vectorIn_t  arg) const

Compute the value of each level, and the jacobian if ComputeJac is true.

◆ configSpace()

const LiegroupSpacePtr_t& hpp::constraints::solver::HierarchicalIterative::configSpace ( ) const
inline

Get configuration space on which constraints are defined.

◆ constraints() [1/2]

const NumericalConstraints_t& hpp::constraints::solver::HierarchicalIterative::constraints ( ) const
inline

Get constraints (implicit and explicit)

◆ constraints() [2/2]

const ImplicitConstraintSet& hpp::constraints::solver::HierarchicalIterative::constraints ( const std::size_t  priority)
inline

Get set of constraints for a give priority level.

◆ contains()

virtual bool hpp::constraints::solver::HierarchicalIterative::contains ( const ImplicitPtr_t numericalConstraint) const
virtual

Check whether a numerical constraint has been added

Parameters
numericalConstraintnumerical constraint
Returns
true if numerical constraint is already in the solver

◆ definesSubmanifoldOf()

bool hpp::constraints::solver::HierarchicalIterative::definesSubmanifoldOf ( const HierarchicalIterative solver) const

Inclusion of manifolds

Parameters
solveranother solver

This function returns true if the solution manifold of the solver given as input is a submanifold of the solution manifold of this solver. The function tests whether constraints of the input solver are also in this solver.

◆ dimension()

const size_type& hpp::constraints::solver::HierarchicalIterative::dimension ( ) const
inline

◆ errorThreshold() [1/2]

value_type hpp::constraints::solver::HierarchicalIterative::errorThreshold ( ) const
inline

Get error threshold.

◆ errorThreshold() [2/2]

void hpp::constraints::solver::HierarchicalIterative::errorThreshold ( const value_type threshold)
inline

Set error threshold.

◆ expandDqSmall()

void hpp::constraints::solver::HierarchicalIterative::expandDqSmall ( ) const
protected

◆ freeVariables() [1/3]

const Indices_t& hpp::constraints::solver::HierarchicalIterative::freeVariables ( ) const
inline

Get free velocity variables.

◆ freeVariables() [2/3]

void hpp::constraints::solver::HierarchicalIterative::freeVariables ( const Indices_t indices)
inline

Set free velocity variables

The other variables will be left unchanged by the iterative resolution.

◆ freeVariables() [3/3]

void hpp::constraints::solver::HierarchicalIterative::freeVariables ( const segments_t  intervals)
inline

Set free velocity variables

The other variables will be left unchanged by the iterative resolution.

Parameters
intervalsset of index intervals

◆ getReducedJacobian()

void hpp::constraints::solver::HierarchicalIterative::getReducedJacobian ( matrixOut_t  J) const

◆ getRightHandSide()

virtual bool hpp::constraints::solver::HierarchicalIterative::getRightHandSide ( const ImplicitPtr_t constraint,
vectorOut_t  rhs 
) const
virtual

Get right hand side of a constraints.

Reimplemented in hpp::constraints::solver::BySubstitution.

◆ getValue()

void hpp::constraints::solver::HierarchicalIterative::getValue ( vectorOut_t  v) const

◆ inequalityThreshold() [1/2]

value_type hpp::constraints::solver::HierarchicalIterative::inequalityThreshold ( ) const
inline

Get the inequality threshold.

◆ inequalityThreshold() [2/2]

void hpp::constraints::solver::HierarchicalIterative::inequalityThreshold ( const value_type it)
inline

set the inequality threshold

◆ integrate()

virtual bool hpp::constraints::solver::HierarchicalIterative::integrate ( vectorIn_t  from,
vectorIn_t  velocity,
vectorOut_t  result 
) const
virtual

◆ isConstraintSatisfied()

bool hpp::constraints::solver::HierarchicalIterative::isConstraintSatisfied ( const ImplicitPtr_t constraint,
vectorIn_t  arg,
vectorOut_t  error,
bool &  constraintFound 
) const

Whether a constraint is satisfied for an input vector

Parameters
constraint,theconstraint in the solver,
argthe input vector
Return values
errorthe error of the constraint.
constraintFoundwhether the constraint belongs to the solver,
Returns
true if constraint belongs to solver and error is below the threshold, false otherwise.

◆ isSatisfied() [1/2]

bool hpp::constraints::solver::HierarchicalIterative::isSatisfied ( vectorIn_t  arg) const
inline

Whether input vector satisfies the constraints of the solver

Parameters
arginput vector Compares to internal error threshold.

◆ isSatisfied() [2/2]

bool hpp::constraints::solver::HierarchicalIterative::isSatisfied ( vectorIn_t  arg,
value_type  errorThreshold 
) const
inline

Whether input vector satisfies the constraints of the solver

Parameters
arginput vector
errorThresholdthreshold to use instead of the value stored in the solver.

◆ lastIsOptional() [1/2]

bool hpp::constraints::solver::HierarchicalIterative::lastIsOptional ( ) const
inline

◆ lastIsOptional() [2/2]

void hpp::constraints::solver::HierarchicalIterative::lastIsOptional ( bool  optional)
inline

◆ lastStep()

const vector_t& hpp::constraints::solver::HierarchicalIterative::lastStep ( ) const
inline

Accessor to the last step done.

◆ maxIterations() [1/2]

size_type hpp::constraints::solver::HierarchicalIterative::maxIterations ( ) const
inline

Get maximal number of iterations in config projector.

◆ maxIterations() [2/2]

void hpp::constraints::solver::HierarchicalIterative::maxIterations ( size_type  iterations)
inline

Set maximal number of iterations.

◆ merge()

virtual void hpp::constraints::solver::HierarchicalIterative::merge ( const HierarchicalIterative other)
virtual

add constraints of another solver

Parameters
otherother solver

Add constraints of other to this solver.

Note
right hand side of other is not copied.

◆ numberStacks()

std::size_t hpp::constraints::solver::HierarchicalIterative::numberStacks ( ) const
inline

◆ print()

virtual std::ostream& hpp::constraints::solver::HierarchicalIterative::print ( std::ostream &  os) const
virtual

◆ reducedDimension()

const size_type& hpp::constraints::solver::HierarchicalIterative::reducedDimension ( ) const
inline

Dimension of the problem after removing the rows of the jacobian which do not influence the error (only zeros along those lines).

◆ residualError() [1/2]

value_type hpp::constraints::solver::HierarchicalIterative::residualError ( ) const
inline

Returns the squared norm of the error vector.

◆ residualError() [2/2]

void hpp::constraints::solver::HierarchicalIterative::residualError ( vectorOut_t  error) const

Returns the error vector.

◆ rightHandSide() [1/3]

vector_t hpp::constraints::solver::HierarchicalIterative::rightHandSide ( ) const

Get the right hand side

Returns
the right hand side
Note
size of result is equal to total dimension of parameterizable constraints (type Equality).

◆ rightHandSide() [2/3]

virtual bool hpp::constraints::solver::HierarchicalIterative::rightHandSide ( const ImplicitPtr_t constraint,
vectorIn_t  rhs 
)
virtual

Set right hand side of a constraints

Parameters
constraintthe constraint,
rhsright hand side.
Note
Size of rhs should be equal to the total dimension of parameterizable constraints (type Equality) .

Reimplemented in hpp::constraints::solver::BySubstitution.

◆ rightHandSide() [3/3]

virtual void hpp::constraints::solver::HierarchicalIterative::rightHandSide ( vectorIn_t  rhs)
virtual

Set the right hand side

Parameters
rhsthe right hand side
Note
Size of rhs should be equal to the total dimension of parameterizable constraints (type Equality).

Reimplemented in hpp::constraints::solver::BySubstitution.

◆ rightHandSideAt()

void hpp::constraints::solver::HierarchicalIterative::rightHandSideAt ( const value_type s)

Set the right hand side at a given parameter.

Parameters
sparameter passed to Implicit::rightHandSideAt

◆ rightHandSideFromConfig() [1/2]

vector_t hpp::constraints::solver::HierarchicalIterative::rightHandSideFromConfig ( ConfigurationIn_t  config)

Compute right hand side of equality constraints from a configuration

Parameters
configa configuration.

for each constraint of type Equality, set right hand side as \(rhs = f(\mathbf{q})\).

Note
Only parameterizable constraints (type Equality) are set

◆ rightHandSideFromConfig() [2/2]

virtual bool hpp::constraints::solver::HierarchicalIterative::rightHandSideFromConfig ( const ImplicitPtr_t constraint,
ConfigurationIn_t  config 
)
virtual

Compute right hand side of a constraint from a configuration

Parameters
constraintthe constraint,
configa configuration.

Set right hand side as \(rhs = f(\mathbf{q})\).

Note
Only parameterizable constraints (type Equality) are set

Reimplemented in hpp::constraints::solver::BySubstitution.

◆ rightHandSideSize()

size_type hpp::constraints::solver::HierarchicalIterative::rightHandSideSize ( ) const

Get size of the right hand side

Returns
sum of dimensions of parameterizable constraints (type Equality)

◆ saturate()

void hpp::constraints::solver::HierarchicalIterative::saturate ( vectorOut_t  arg) const
protected

◆ saturation() [1/2]

const Saturation_t& hpp::constraints::solver::HierarchicalIterative::saturation ( ) const
inline

Get the saturation function.

◆ saturation() [2/2]

void hpp::constraints::solver::HierarchicalIterative::saturation ( const Saturation_t saturate)
inline

Set the saturation function.

◆ sigma()

const value_type& hpp::constraints::solver::HierarchicalIterative::sigma ( ) const
inline

Returns the lowest singular value. If the jacobian has maximum rank r, then it corresponds to r-th greatest singular value. This value is zero when the jacobian is singular.

◆ solve() [1/3]

Status hpp::constraints::solver::HierarchicalIterative::solve ( vectorOut_t  arg) const
inline

Solve the system of non linear equations

Parameters
arginitial guess,

Use Newton Rhapson like iterative method until the error is below the threshold, or until the maximal number of iterations has been reached. Use the default line search method (fixed sequence of \(\alpha_i\)).

Note
Explicit constraints are expressed in their implicit form: \(\mathbf{q}_2 = f (\mathbf{q}_1)\) is replaced by \(\mathbf{q}_2 - f (\mathbf{q}_1) = 0\), where \(\mathbf{q}_1\) and \(\mathbf{q}_2\) are vectors composed of the components of \(\mathbf{q}\).

◆ solve() [2/3]

template<typename LineSearchType >
solver::HierarchicalIterative::Status hpp::constraints::solver::HierarchicalIterative::solve ( vectorOut_t  arg,
LineSearchType  lineSearch 
) const
inline

◆ solve() [3/3]

template<typename LineSearchType >
Status hpp::constraints::solver::HierarchicalIterative::solve ( vectorOut_t  arg,
LineSearchType  ls = LineSearchType() 
) const

Solve the system of non linear equations

Parameters
arginitial guess,
lsline search method used.

Use Newton Rhapson like iterative method until the error is below the threshold, or until the maximal number of iterations has been reached.

Note
Explicit constraints are expressed in their implicit form: \(\mathbf{q}_2 = f (\mathbf{q}_1)\) is replaced by \(\mathbf{q}_2 - f (\mathbf{q}_1) = 0\), where \(\mathbf{q}_1\) and \(\mathbf{q}_2\) are vectors composed of the components of \(\mathbf{q}\).

◆ solveLevelByLevel() [1/2]

bool hpp::constraints::solver::HierarchicalIterative::solveLevelByLevel ( ) const
inline

◆ solveLevelByLevel() [2/2]

void hpp::constraints::solver::HierarchicalIterative::solveLevelByLevel ( bool  solveLevelByLevel)
inline

◆ squaredErrorThreshold()

value_type hpp::constraints::solver::HierarchicalIterative::squaredErrorThreshold ( ) const
inline

Get error threshold.

◆ update()

void hpp::constraints::solver::HierarchicalIterative::update ( )
protected

Allocate datas and update sizes of the problem Should be called whenever the stack is modified.

Friends And Related Function Documentation

◆ lineSearch::Backtracking

friend struct lineSearch::Backtracking
friend

Member Data Documentation

◆ configSpace_

LiegroupSpacePtr_t hpp::constraints::solver::HierarchicalIterative::configSpace_
protected

◆ constraints_

NumericalConstraints_t hpp::constraints::solver::HierarchicalIterative::constraints_
protected

Members moved from core::ConfigProjector.

◆ datas_

std::vector<Data> hpp::constraints::solver::HierarchicalIterative::datas_
mutableprotected

◆ dimension_

size_type hpp::constraints::solver::HierarchicalIterative::dimension_
protected

◆ dq_

vector_t hpp::constraints::solver::HierarchicalIterative::dq_
mutableprotected

◆ dqSmall_

vector_t hpp::constraints::solver::HierarchicalIterative::dqSmall_
protected

◆ freeVariables_

Indices_t hpp::constraints::solver::HierarchicalIterative::freeVariables_
protected

Unknown of the set of implicit constraints.

◆ inequalityThreshold_

value_type hpp::constraints::solver::HierarchicalIterative::inequalityThreshold_
protected

◆ iq_

std::map<DifferentiableFunctionPtr_t, size_type> hpp::constraints::solver::HierarchicalIterative::iq_
protected

Value rank of constraint in its priority level.

◆ iv_

std::map<DifferentiableFunctionPtr_t, size_type> hpp::constraints::solver::HierarchicalIterative::iv_
protected

Derivative rank of constraint in its priority level.

◆ lastIsOptional_

bool hpp::constraints::solver::HierarchicalIterative::lastIsOptional_
protected

◆ maxIterations_

size_type hpp::constraints::solver::HierarchicalIterative::maxIterations_
protected

◆ OM_

vector_t hpp::constraints::solver::HierarchicalIterative::OM_
mutableprotected

◆ OP_

vector_t hpp::constraints::solver::HierarchicalIterative::OP_
mutableprotected

◆ priority_

std::map<DifferentiableFunctionPtr_t, std::size_t> hpp::constraints::solver::HierarchicalIterative::priority_
protected

Priority level of constraint.

◆ qSat_

Configuration_t hpp::constraints::solver::HierarchicalIterative::qSat_
mutableprotected

◆ reducedDimension_

size_type hpp::constraints::solver::HierarchicalIterative::reducedDimension_
protected

◆ reducedJ_

matrix_t hpp::constraints::solver::HierarchicalIterative::reducedJ_
mutableprotected

◆ reducedSaturation_

Eigen::VectorXi hpp::constraints::solver::HierarchicalIterative::reducedSaturation_
protected

◆ saturate_

Saturation_t hpp::constraints::solver::HierarchicalIterative::saturate_
protected

◆ saturation_

Eigen::VectorXi hpp::constraints::solver::HierarchicalIterative::saturation_
mutableprotected

◆ sigma_

value_type hpp::constraints::solver::HierarchicalIterative::sigma_
mutableprotected

The smallest non-zero singular value.

◆ solveLevelByLevel_

bool hpp::constraints::solver::HierarchicalIterative::solveLevelByLevel_
protected

◆ squaredErrorThreshold_

value_type hpp::constraints::solver::HierarchicalIterative::squaredErrorThreshold_
protected

◆ squaredNorm_

value_type hpp::constraints::solver::HierarchicalIterative::squaredNorm_
mutableprotected

◆ stacks_

std::vector<ImplicitConstraintSet> hpp::constraints::solver::HierarchicalIterative::stacks_
protected

◆ svd_

SVD_t hpp::constraints::solver::HierarchicalIterative::svd_
mutableprotected

◆ tmpSat_

ArrayXb hpp::constraints::solver::HierarchicalIterative::tmpSat_
mutableprotected

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