crocoddyl 1.9.0
Contact RObot COntrol by Differential DYnamic programming Library (Crocoddyl)
 
Loading...
Searching...
No Matches
SolverKKT Class Reference
Inheritance diagram for SolverKKT:
Collaboration diagram for SolverKKT:

Public Member Functions

EIGEN_MAKE_ALIGNED_OPERATOR_NEW SolverKKT (boost::shared_ptr< ShootingProblem > problem)
 
virtual void computeDirection (const bool recalc=true)
 Compute the search direction \((\delta\mathbf{x}^k,\delta\mathbf{u}^k)\) for the current guess \((\mathbf{x}^k_s,\mathbf{u}^k_s)\). More...
 
virtual const Eigen::Vector2d & expectedImprovement ()
 Return the expected improvement \(dV_{exp}\) from a given current search direction \((\delta\mathbf{x}^k,\delta\mathbf{u}^k)\). More...
 
const std::vector< Eigen::VectorXd > & get_dus () const
 
const std::vector< Eigen::VectorXd > & get_dxs () const
 
const Eigen::MatrixXd & get_kkt () const
 
const Eigen::VectorXd & get_kktref () const
 
const std::vector< Eigen::VectorXd > & get_lambdas () const
 
std::size_t get_ndx () const
 
std::size_t get_nu () const
 
std::size_t get_nx () const
 
const Eigen::VectorXd & get_primaldual () const
 
virtual bool solve (const std::vector< Eigen::VectorXd > &init_xs=DEFAULT_VECTOR, const std::vector< Eigen::VectorXd > &init_us=DEFAULT_VECTOR, const std::size_t maxiter=100, const bool is_feasible=false, const double regInit=1e-9)
 Compute the optimal trajectory \(\mathbf{x}^*_s,\mathbf{u}^*_s\) as lists of \(T+1\) and \(T\) terms. More...
 
virtual double stoppingCriteria ()
 Return a positive value that quantifies the algorithm termination. More...
 
virtual double tryStep (const double steplength=1)
 Try a predefined step length \(\alpha\) and compute its cost improvement \(dV\). More...
 
- Public Member Functions inherited from SolverAbstract
EIGEN_MAKE_ALIGNED_OPERATOR_NEW SolverAbstract (boost::shared_ptr< ShootingProblem > problem)
 Initialize the solver. More...
 
virtual void computeDirection (const bool recalc)=0
 Compute the search direction \((\delta\mathbf{x}^k,\delta\mathbf{u}^k)\) for the current guess \((\mathbf{x}^k_s,\mathbf{u}^k_s)\). More...
 
double computeDynamicFeasibility ()
 Compute the dynamic feasibility \(\|\mathbf{f}_{\mathbf{s}}\|_{\infty,1}\) for the current guess \((\mathbf{x}^k,\mathbf{u}^k)\). More...
 
virtual const Eigen::Vector2d & expectedImprovement ()=0
 Return the expected improvement \(dV_{exp}\) from a given current search direction \((\delta\mathbf{x}^k,\delta\mathbf{u}^k)\). More...
 
double get_cost () const
 Return the total cost. More...
 
const Eigen::Vector2d & get_d () const
 Return the LQ approximation of the expected improvement. More...
 
double get_dV () const
 Return the cost reduction \(dV\). More...
 
double get_dVexp () const
 Return the expected cost reduction \(dV_{exp}\). More...
 
double get_ffeas () const
 Return the feasibility of the dynamic constraints \(\|\mathbf{f}_{\mathbf{s}}\|_{\infty,1}\) of the current guess. More...
 
const std::vector< Eigen::VectorXd > & get_fs () const
 Return the gaps \(\mathbf{f}_{s}\). More...
 
bool get_inffeas () const
 Return the norm used for the computing the feasibility (true for \(\ell_\infty\), false for \(\ell_1\)) More...
 
bool get_is_feasible () const
 Return the feasibility status of the \((\mathbf{x}_s,\mathbf{u}_s)\) trajectory. More...
 
std::size_t get_iter () const
 Return the number of iterations performed by the solver. More...
 
const boost::shared_ptr< ShootingProblem > & get_problem () const
 Return the shooting problem. More...
 
double get_steplength () const
 Return the step length \(\alpha\). More...
 
double get_stop () const
 Return the value computed by stoppingCriteria() More...
 
double get_th_acceptstep () const
 Return the threshold used for accepting a step. More...
 
double get_th_gaptol () const
 Return the threshold for accepting a gap as non-zero. More...
 
double get_th_stop () const
 Return the tolerance for stopping the algorithm. More...
 
double get_ureg () const
 Return the control regularization value. More...
 
const std::vector< Eigen::VectorXd > & get_us () const
 Return the control trajectory \(\mathbf{u}_s\). More...
 
double get_xreg () const
 Return the state regularization value. More...
 
const std::vector< Eigen::VectorXd > & get_xs () const
 Return the state trajectory \(\mathbf{x}_s\). More...
 
const std::vector< boost::shared_ptr< CallbackAbstract > > & getCallbacks () const
 Return the list of callback functions using for diagnostic. More...
 
virtual void resizeData ()
 Resizing the solver data. More...
 
void set_inffeas (const bool inffeas)
 Modify the current norm used for computed the feasibility. More...
 
void set_th_acceptstep (const double th_acceptstep)
 Modify the threshold used for accepting step. More...
 
void set_th_gaptol (const double th_gaptol)
 Modify the threshold for accepting a gap as non-zero. More...
 
void set_th_stop (const double th_stop)
 Modify the tolerance for stopping the algorithm. More...
 
void set_ureg (const double ureg)
 Modify the control regularization value. More...
 
void set_us (const std::vector< Eigen::VectorXd > &us)
 Modify the control trajectory \(\mathbf{u}_s\). More...
 
void set_xreg (const double xreg)
 Modify the state regularization value. More...
 
void set_xs (const std::vector< Eigen::VectorXd > &xs)
 Modify the state trajectory \(\mathbf{x}_s\). More...
 
void setCallbacks (const std::vector< boost::shared_ptr< CallbackAbstract > > &callbacks)
 Set a list of callback functions using for the solver diagnostic. More...
 
void setCandidate (const std::vector< Eigen::VectorXd > &xs_warm=DEFAULT_VECTOR, const std::vector< Eigen::VectorXd > &us_warm=DEFAULT_VECTOR, const bool is_feasible=false)
 Set the solver candidate trajectories \((\mathbf{x}_s,\mathbf{u}_s)\). More...
 
virtual bool solve (const std::vector< Eigen::VectorXd > &init_xs=DEFAULT_VECTOR, const std::vector< Eigen::VectorXd > &init_us=DEFAULT_VECTOR, const std::size_t maxiter=100, const bool is_feasible=false, const double reg_init=1e-9)=0
 Compute the optimal trajectory \(\mathbf{x}^*_s,\mathbf{u}^*_s\) as lists of \(T+1\) and \(T\) terms. More...
 
virtual double stoppingCriteria ()=0
 Return a positive value that quantifies the algorithm termination. More...
 
virtual double tryStep (const double steplength=1)=0
 Try a predefined step length \(\alpha\) and compute its cost improvement \(dV\). More...
 

Protected Attributes

double cost_try_
 
double reg_decfactor_
 
double reg_incfactor_
 
double reg_max_
 
double reg_min_
 
std::vector< Eigen::VectorXd > us_try_
 
std::vector< Eigen::VectorXd > xs_try_
 
- Protected Attributes inherited from SolverAbstract
std::vector< boost::shared_ptr< CallbackAbstract > > callbacks_
 Callback functions. More...
 
double cost_
 Total cost. More...
 
Eigen::Vector2d d_
 LQ approximation of the expected improvement. More...
 
double dV_
 Cost reduction obtained by tryStep() More...
 
double dVexp_
 Expected cost reduction. More...
 
double ffeas_
 Feasibility of the dynamic constraints. More...
 
std::vector< Eigen::VectorXd > fs_
 Gaps/defects between shooting nodes. More...
 
bool inffeas_
 
bool is_feasible_
 Label that indicates is the iteration is feasible. More...
 
std::size_t iter_
 Number of iteration performed by the solver. More...
 
boost::shared_ptr< ShootingProblemproblem_
 optimal control problem More...
 
double steplength_
 Current applied step-length. More...
 
double stop_
 Value computed by stoppingCriteria() More...
 
double th_acceptstep_
 Threshold used for accepting step. More...
 
double th_gaptol_
 Threshold limit to check non-zero gaps. More...
 
double th_stop_
 Tolerance for stopping the algorithm. More...
 
double tmp_feas_
 Temporal variables used for computed the feasibility. More...
 
double ureg_
 Current control regularization values. More...
 
std::vector< Eigen::VectorXd > us_
 Control trajectory. More...
 
bool was_feasible_
 Label that indicates in the previous iterate was feasible. More...
 
double xreg_
 Current state regularization value. More...
 
std::vector< Eigen::VectorXd > xs_
 State trajectory. More...
 

Detailed Description

Definition at line 20 of file kkt.hpp.

Constructor & Destructor Documentation

◆ SolverKKT()

SolverKKT ( boost::shared_ptr< ShootingProblem problem)
explicit

Definition at line 14 of file kkt.cpp.

◆ ~SolverKKT()

~SolverKKT ( )
virtual

Definition at line 33 of file kkt.cpp.

Member Function Documentation

◆ solve()

bool solve ( const std::vector< Eigen::VectorXd > &  init_xs = DEFAULT_VECTOR,
const std::vector< Eigen::VectorXd > &  init_us = DEFAULT_VECTOR,
const std::size_t  maxiter = 100,
const bool  is_feasible = false,
const double  reg_init = 1e-9 
)
virtual

Compute the optimal trajectory \(\mathbf{x}^*_s,\mathbf{u}^*_s\) as lists of \(T+1\) and \(T\) terms.

From an initial guess init_xs, init_us (feasible or not), iterate over computeDirection() and tryStep() until stoppingCriteria() is below threshold. It also describes the globalization strategy used during the numerical optimization.

Parameters
[in]init_xsinitial guess for state trajectory with \(T+1\) elements (default [])
[in]init_usinitial guess for control trajectory with \(T\) elements (default [])
[in]maxitermaximum allowed number of iterations (default 100)
[in]isFeasibletrue if the init_xs are obtained from integrating the init_us (rollout) (default false)
[in]regInitinitial guess for the regularization value. Very low values are typical used with very good guess points (init_xs, init_us)
Returns
A boolean that describes if convergence was reached.

Implements SolverAbstract.

Definition at line 35 of file kkt.cpp.

◆ computeDirection()

void computeDirection ( const bool  recalc = true)
virtual

Compute the search direction \((\delta\mathbf{x}^k,\delta\mathbf{u}^k)\) for the current guess \((\mathbf{x}^k_s,\mathbf{u}^k_s)\).

You must call setCandidate() first in order to define the current guess. A current guess defines a state and control trajectory \((\mathbf{x}^k_s,\mathbf{u}^k_s)\) of \(T+1\) and \(T\) elements, respectively.

Parameters
[in]recalctrue for recalculating the derivatives at current state and control
Returns
The search direction \((\delta\mathbf{x},\delta\mathbf{u})\) and the dual lambdas as lists of \(T+1\), \(T\) and \(T+1\) lengths, respectively

Implements SolverAbstract.

Definition at line 85 of file kkt.cpp.

◆ tryStep()

double tryStep ( const double  steplength = 1)
virtual

Try a predefined step length \(\alpha\) and compute its cost improvement \(dV\).

It uses the search direction found by computeDirection() to try a determined step length \(\alpha\). Therefore, it assumes that we have run computeDirection() first. Additionally, it returns the cost improvement \(dV\) along the predefined step length \(\alpha\).

Parameters
[in]steplengthapplied step length ( \(0\leq\alpha\leq1\))
Returns
the cost improvement

Implements SolverAbstract.

Definition at line 111 of file kkt.cpp.

◆ stoppingCriteria()

double stoppingCriteria ( )
virtual

Return a positive value that quantifies the algorithm termination.

These values typically represents the gradient norm which tell us that it's been reached the local minima. The stopping criteria strictly speaking depends on the search direction (calculated by computeDirection()) but it could also depend on the chosen step length, tested by tryStep().

Implements SolverAbstract.

Definition at line 129 of file kkt.cpp.

◆ expectedImprovement()

const Eigen::Vector2d & expectedImprovement ( )
virtual

Return the expected improvement \(dV_{exp}\) from a given current search direction \((\delta\mathbf{x}^k,\delta\mathbf{u}^k)\).

For computing the expected improvement, you need to compute the search direction first via computeDirection().

Implements SolverAbstract.

Definition at line 152 of file kkt.cpp.

◆ get_kkt()

const Eigen::MatrixXd & get_kkt ( ) const

Definition at line 162 of file kkt.cpp.

◆ get_kktref()

const Eigen::VectorXd & get_kktref ( ) const

Definition at line 164 of file kkt.cpp.

◆ get_primaldual()

const Eigen::VectorXd & get_primaldual ( ) const

Definition at line 166 of file kkt.cpp.

◆ get_dxs()

const std::vector< Eigen::VectorXd > & get_dxs ( ) const

Definition at line 168 of file kkt.cpp.

◆ get_dus()

const std::vector< Eigen::VectorXd > & get_dus ( ) const

Definition at line 170 of file kkt.cpp.

◆ get_lambdas()

const std::vector< Eigen::VectorXd > & get_lambdas ( ) const

Definition at line 172 of file kkt.cpp.

◆ get_nx()

std::size_t get_nx ( ) const

Definition at line 174 of file kkt.cpp.

◆ get_ndx()

std::size_t get_ndx ( ) const

Definition at line 176 of file kkt.cpp.

◆ get_nu()

std::size_t get_nu ( ) const

Definition at line 178 of file kkt.cpp.

Member Data Documentation

◆ reg_incfactor_

double reg_incfactor_
protected

Definition at line 46 of file kkt.hpp.

◆ reg_decfactor_

double reg_decfactor_
protected

Definition at line 47 of file kkt.hpp.

◆ reg_min_

double reg_min_
protected

Definition at line 48 of file kkt.hpp.

◆ reg_max_

double reg_max_
protected

Definition at line 49 of file kkt.hpp.

◆ cost_try_

double cost_try_
protected

Definition at line 50 of file kkt.hpp.

◆ xs_try_

std::vector<Eigen::VectorXd> xs_try_
protected

Definition at line 51 of file kkt.hpp.

◆ us_try_

std::vector<Eigen::VectorXd> us_try_
protected

Definition at line 52 of file kkt.hpp.


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