hpp::constraints::solver::HierarchicalIterative Class Reference

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

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

Inheritance diagram for hpp::constraints::solver::HierarchicalIterative:
[legend]
Collaboration diagram for hpp::constraints::solver::HierarchicalIterative:
[legend]

Classes

struct  Data
 

Public Types

enum  Status {
  ERROR_INCREASED,
  MAX_ITERATION_REACHED,
  INFEASIBLE,
  SUCCESS
}
 
typedef Eigen::RowBlockIndices Indices_t
 
typedef lineSearch::FixedSequence DefaultLineSearch
 
typedef boost::function< bool(vectorIn_t q, vectorOut_t qSat, Eigen::VectorXi &saturation)> Saturation_t
 This function checks which degrees of freedom are saturated. More...
 

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...
 
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
 Check whether a numerical constraint has been added. More...
 
void add (const DifferentiableFunctionPtr_t &f, const std::size_t &priority) HPP_CONSTRAINTS_DEPRECATED
 Add an implicit equality constraint. More...
 
void add (const DifferentiableFunctionPtr_t &f, const std::size_t &priority, const ComparisonTypes_t &comp) HPP_CONSTRAINTS_DEPRECATED
 Add an implicit constraint. More...
 
void add (const ImplicitPtr_t &constraint, const std::size_t &priority)
 Add an implicit constraint. More...
 
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
 Solve the system of non linear equations. More...
 
Status solve (vectorOut_t arg) const
 Solve the system of non linear equations. More...
 
bool isSatisfied (vectorIn_t arg) const
 
const value_typesigma () const
 Returns the lowest singular value. More...
 
Parameters
void freeVariables (const segments_t intervals)
 Set free velocity variables. More...
 
void freeVariables (const Indices_t &indices)
 Set free velocity variables. More...
 
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
 
Stack
const ImplicitConstraintSetconstraints (const std::size_t priority)
 Get set of constraints for a give priority level. More...
 
std::size_t numberStacks () const
 
const size_typedimension () const
 
const size_typereducedDimension () const
 Dimension of the problem after removing the rows of the jacobian which do not influence the error (only zeros along those lines). More...
 
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)
 Compute right hand side of equality constraints from a configuration. More...
 
virtual bool rightHandSideFromConfig (const ImplicitPtr_t &constraint, ConfigurationIn_t config)
 Compute right hand side of a constraint from a configuration. More...
 
virtual bool rightHandSide (const ImplicitPtr_t &constraint, vectorIn_t rhs)
 Set right hand side of a constraints. More...
 
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)
 Set the right hand side. More...
 
void rightHandSideAt (const value_type &s)
 Set the right hand side at a given parameter. More...
 
vector_t rightHandSide () const
 Get the right hand side. More...
 
size_type rightHandSideSize () const
 Get size of the right hand side. More...
 
Access to internal datas

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
 If lastIsOptional() is true, then the last level is ignored. More...
 
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 ()
 Allocate datas and update sizes of the problem Should be called whenever the stack is modified. More...
 
virtual void computeActiveRowsOfJ (std::size_t iStack)
 Compute which rows of the jacobian of stack_[iStack] are not zero, using the activeDerivativeParameters of the functions. More...
 
void computeDescentDirection () const
 Compute a SVD decomposition of each level and find the best descent direction at the first order. More...
 
void expandDqSmall () const
 
void saturate (vectorOut_t arg) const
 

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_
 
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_
 
mutable::hpp::statistics::SuccessStatistics statistics_
 

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: 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 error threshold can be accessed by methods HierarchicalIterative::errorThreshold. The maximal number of iterations can be accessed by methods HierarchicalIterative::maxIterations.

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

To prevent configuration variables to get out of joint limits during Newton Raphson iterations, the user may provide a method of type HierarchicalIterative::Saturation_t using setter and getter HierarchicalIterative::saturation.

Note
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

typedef boost::function<bool (vectorIn_t q, vectorOut_t qSat, Eigen::VectorXi& saturation)> hpp::constraints::solver::HierarchicalIterative::Saturation_t

This function checks which degrees of freedom are saturated.

Parameters
qa configuration,
Return values
qSatconfiguration after saturing values out of bounds
saturationvector: for each degree of freedom, saturation is set to
  • -1 if the lower bound is reached,
  • 1 if the upper bound is reached,
  • 0 otherwise.
Returns
true if and only if at least one degree of freedom has been saturated

Member Enumeration Documentation

Enumerator
ERROR_INCREASED 
MAX_ITERATION_REACHED 
INFEASIBLE 
SUCCESS 

Constructor & Destructor Documentation

hpp::constraints::solver::HierarchicalIterative::HierarchicalIterative ( const LiegroupSpacePtr_t configSpace)
hpp::constraints::solver::HierarchicalIterative::HierarchicalIterative ( const HierarchicalIterative other)
virtual hpp::constraints::solver::HierarchicalIterative::~HierarchicalIterative ( )
inlinevirtual

Member Function Documentation

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

Velocity parameters involved in the constraint resolution.

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

Configuration parameters involved in the constraint resolution.

void hpp::constraints::solver::HierarchicalIterative::add ( const DifferentiableFunctionPtr_t f,
const std::size_t &  priority 
)

Add an implicit equality constraint.

Parameters
fdifferentiable function from the robot configuration space to a Lie group (See hpp::pinocchio::LiegroupSpace),
prioritylevel of priority of the constraint: priority are in decreasing order: 0 is the highest priority level.

Constraint is defined by \(f (\mathbf{q}) = 0\).

void hpp::constraints::solver::HierarchicalIterative::add ( const DifferentiableFunctionPtr_t f,
const std::size_t &  priority,
const ComparisonTypes_t comp 
)

Add an implicit constraint.

Parameters
fdifferentiable function from the robot configuration space to a Lie group (See hpp::pinocchio::LiegroupSpace),
prioritylevel of priority of the constraint: priority are in decreasing order: 0 is the highest priority level,
compcomparison type. See class documentation for details.
void hpp::constraints::solver::HierarchicalIterative::add ( const ImplicitPtr_t constraint,
const std::size_t &  priority 
)

Add an implicit constraint.

Parameters
constraintimplicit constraint
prioritylevel of priority of the constraint: priority are in decreasing order: 0 is the highest priority level,
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.

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.

Referenced by hpp::constraints::solver::BySubstitution::impl_solve().

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

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

Warning
computeValue must have been called first.

Referenced by hpp::constraints::solver::BySubstitution::impl_solve().

void hpp::constraints::solver::HierarchicalIterative::computeSaturation ( vectorIn_t  arg) const
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.

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

Get configuration space on which constraints are defined.

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

Get set of constraints for a give priority level.

References Eigen::assert().

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 whatever the passive dofs are.

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

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

Set error threshold.

value_type hpp::constraints::solver::HierarchicalIterative::errorThreshold ( ) const
inline
void hpp::constraints::solver::HierarchicalIterative::expandDqSmall ( ) const
protected
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
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.

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

Get free velocity variables.

void hpp::constraints::solver::HierarchicalIterative::getReducedJacobian ( matrixOut_t  J) const
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.

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

Get the inequality threshold.

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

set the inequality threshold

virtual bool hpp::constraints::solver::HierarchicalIterative::integrate ( vectorIn_t  from,
vectorIn_t  velocity,
vectorOut_t  result 
) const
virtual
bool hpp::constraints::solver::HierarchicalIterative::isSatisfied ( vectorIn_t  arg) const
inline
void hpp::constraints::solver::HierarchicalIterative::lastIsOptional ( bool  optional)
inline
bool hpp::constraints::solver::HierarchicalIterative::lastIsOptional ( ) const
inline
const vector_t& hpp::constraints::solver::HierarchicalIterative::lastStep ( ) const
inline

Accessor to the last step done.

References integrate().

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

Set maximal number of iterations.

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

Get maximal number of iterations in config projector.

std::size_t hpp::constraints::solver::HierarchicalIterative::numberStacks ( ) const
inline
virtual std::ostream& hpp::constraints::solver::HierarchicalIterative::print ( std::ostream &  os) const
virtual
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).

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

Returns the squared norm of the error vector.

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

Returns the error vector.

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.

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.

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).
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
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
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.

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

Get size of the right hand side.

Returns
sum of dimensions of parameterizable constraints (type Equality)
void hpp::constraints::solver::HierarchicalIterative::saturate ( vectorOut_t  arg) const
protected
void hpp::constraints::solver::HierarchicalIterative::saturation ( const Saturation_t saturate)
inline

Set the saturation function.

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

Get the saturation function.

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.

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

References Eigen::assert().

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}\).
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}\).
value_type hpp::constraints::solver::HierarchicalIterative::squaredErrorThreshold ( ) const
inline

Get error threshold.

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

friend struct lineSearch::Backtracking
friend

Member Data Documentation

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

Members moved from core::ConfigProjector.

std::vector<Data> hpp::constraints::solver::HierarchicalIterative::datas_
mutableprotected
size_type hpp::constraints::solver::HierarchicalIterative::dimension_
protected
vector_t hpp::constraints::solver::HierarchicalIterative::dq_
mutableprotected
vector_t hpp::constraints::solver::HierarchicalIterative::dqSmall_
mutableprotected
Indices_t hpp::constraints::solver::HierarchicalIterative::freeVariables_
protected

Unknown of the set of implicit constraints.

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

Value rank of constraint in its priority level.

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

Derivative rank of constraint in its priority level.

bool hpp::constraints::solver::HierarchicalIterative::lastIsOptional_
protected
size_type hpp::constraints::solver::HierarchicalIterative::maxIterations_
protected
vector_t hpp::constraints::solver::HierarchicalIterative::OM_
mutableprotected
vector_t hpp::constraints::solver::HierarchicalIterative::OP_
mutableprotected
std::map<DifferentiableFunctionPtr_t, std::size_t> hpp::constraints::solver::HierarchicalIterative::priority_
protected

Priority level of constraint.

Configuration_t hpp::constraints::solver::HierarchicalIterative::qSat_
mutableprotected
size_type hpp::constraints::solver::HierarchicalIterative::reducedDimension_
protected
matrix_t hpp::constraints::solver::HierarchicalIterative::reducedJ_
mutableprotected
Eigen::VectorXi hpp::constraints::solver::HierarchicalIterative::reducedSaturation_
mutableprotected
Saturation_t hpp::constraints::solver::HierarchicalIterative::saturate_
protected
Eigen::VectorXi hpp::constraints::solver::HierarchicalIterative::saturation_
mutableprotected
value_type hpp::constraints::solver::HierarchicalIterative::sigma_
mutableprotected

The smallest non-zero singular value.

value_type hpp::constraints::solver::HierarchicalIterative::squaredErrorThreshold_
protected
value_type hpp::constraints::solver::HierarchicalIterative::squaredNorm_
mutableprotected
std::vector<ImplicitConstraintSet> hpp::constraints::solver::HierarchicalIterative::stacks_
protected
mutable ::hpp::statistics::SuccessStatistics hpp::constraints::solver::HierarchicalIterative::statistics_
protected
SVD_t hpp::constraints::solver::HierarchicalIterative::svd_
mutableprotected
ArrayXb hpp::constraints::solver::HierarchicalIterative::tmpSat_
mutableprotected