Class SolverIntro

Inheritance Relationships

Base Type

Class Documentation

class SolverIntro : public crocoddyl::SolverFDDP

Public Functions

explicit EIGEN_MAKE_ALIGNED_OPERATOR_NEW SolverIntro(std::shared_ptr<ShootingProblem> problem)

Initialize the INTRO solver.

Parameters:
  • problem[in] Shooting problem

  • reduced[in] Use the reduced Schur-complement approach (default true)

virtual ~SolverIntro()
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 init_reg = NAN)
virtual double tryStep(const double step_length = 1)
virtual double stoppingCriteria()
virtual void resizeData()
virtual double calcDiff()
virtual void computeValueFunction(const std::size_t t, const std::shared_ptr<ActionModelAbstract> &model)
virtual void computeGains(const std::size_t t)
EqualitySolverType get_equality_solver() const

Return the type of solver used for handling the equality constraints.

double get_th_feas() const

Return the threshold for switching to feasibility.

double get_rho() const

Return the rho parameter used in the merit function.

double get_upsilon() const

Return the estimated penalty parameter that balances relative contribution of the cost function and equality constraints.

const std::vector<Eigen::MatrixXd> &get_YZ() const

Return the rank of control-equality constraints \(\mathbf{H_u}\f */ const std::vector<std::size_t>& get_Hu_rank() const; /** @brief Return the span and kernel of control-equality constraints \)\mathbf{H_u}\f.

const std::vector<Eigen::MatrixXd> &get_Qzz() const

Return Hessian of the reduced Hamiltonian \(\mathbf{Q_{zz}}\).

const std::vector<Eigen::MatrixXd> &get_Qxz() const

Return Hessian of the reduced Hamiltonian \(\mathbf{Q_{xz}}\).

const std::vector<Eigen::MatrixXd> &get_Quz() const

Return Hessian of the reduced Hamiltonian \(\mathbf{Q_{uz}}\).

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

Return Jacobian of the reduced Hamiltonian \(\mathbf{Q_{z}}\).

const std::vector<Eigen::MatrixXd> &get_Hy() const

Return span-projected Jacobian of the equality-constraint with respect to the control.

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

Return feedforward term related to the nullspace of \(\mathbf{H_u}\).

const std::vector<Eigen::MatrixXd> &get_Kz() const

Return feedback gain related to the nullspace of \(\mathbf{H_u}\).

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

Return feedforward term related to the equality constraints.

const std::vector<Eigen::MatrixXd> &get_Ks() const

Return feedback gain related to the equality constraints.

bool get_zero_upsilon() const

Return the zero-upsilon label.

True if we set the estimated penalty parameter (upsilon) to zero when solve is called.

void set_equality_solver(const EqualitySolverType type)

Modify the type of solver used for handling the equality constraints.

Note that the default solver is nullspace LU. When we enable parallelization, this strategy is generally faster than others for medium to large systems.

void set_th_feas(const double th_feas)

Modify the threshold for switching to feasibility.

void set_rho(const double rho)

Modify the rho parameter used in the merit function.

void set_zero_upsilon(const bool zero_upsilon)

Modify the zero-upsilon label.

Parameters:

zero_upsilon – True if we set estimated penalty parameter (upsilon) to zero when solve is called.

Protected Attributes

enum EqualitySolverType eq_solver_

Strategy used for handling the equality constraints.

double th_feas_

Threshold for switching to feasibility.

double rho_

Parameter used in the merit function to predict the expected reduction

double upsilon_

Estimated penalty parameter that balances relative contribution of the cost function and equality constraints

bool zero_upsilon_

True if we wish to set estimated penalty parameter (upsilon) to zero when solve is called.

std::vector<std::size_t> Hu_rank_

Rank of the control Jacobian of the equality constraints.

std::vector<Eigen::MatrixXd> KQuu_tmp_
std::vector<Eigen::MatrixXd> YZ_

Span \(\mathbf{Y}\in\mathbb{R}^{rank}\) and kernel \(\mathbf{Z}\in\mathbb{R}^{nullity}\) of the control-equality constraints \(\mathbf{H_u}\)

std::vector<Eigen::MatrixXd> Hy_

Span-projected Jacobian of the equality-constraint with respect to the control

std::vector<Eigen::VectorXd> Qz_

Jacobian of the reduced Hamiltonian \(\mathbf{Q_{z}}\).

std::vector<Eigen::MatrixXd> Qzz_

Hessian of the reduced Hamiltonian \(\mathbf{Q_{zz}}\).

std::vector<Eigen::MatrixXd> Qxz_

Hessian of the reduced Hamiltonian \(\mathbf{Q_{xz}}\).

std::vector<Eigen::MatrixXd> Quz_

Hessian of the reduced Hamiltonian \(\mathbf{Q_{uz}}\).

std::vector<Eigen::VectorXd> kz_

Feedforward term in the nullspace of \(\mathbf{H_u}\).

std::vector<Eigen::MatrixXd> Kz_

Feedback gain in the nullspace of \(\mathbf{H_u}\).

std::vector<Eigen::VectorXd> ks_

Feedforward term related to the equality constraints.

std::vector<Eigen::MatrixXd> Ks_

Feedback gain related to the equality constraints.

std::vector<Eigen::MatrixXd> QuuinvHuT_
std::vector<Eigen::LLT<Eigen::MatrixXd>> Qzz_llt_

Cholesky LLT solver.

std::vector<Eigen::FullPivLU<Eigen::MatrixXd>> Hu_lu_

Full-pivot LU solvers used for computing the span and nullspace matrices

std::vector<Eigen::ColPivHouseholderQR<Eigen::MatrixXd>> Hu_qr_

Column-pivot QR solvers used for computing the span and nullspace matrices

std::vector<Eigen::PartialPivLU<Eigen::MatrixXd>> Hy_lu_

Partial-pivot LU solvers used for computing the feedforward and feedback gain related to the equality constraint