proxsuite 0.6.7
The Advanced Proximal Optimization Toolbox
Loading...
Searching...
No Matches
proxsuite::proxqp::dense Namespace Reference

Namespaces

namespace  detail
 
namespace  linesearch
 
namespace  nb
 
namespace  preconditioner
 

Classes

struct  BackwardData
 This class stores the jacobians of PROXQP solvers with dense backends at a solutions wrt model parameters. More...
 
struct  BatchQP
 
struct  EigenAllowAlloc
 
struct  Model
 This class stores the model of the QP problem. More...
 
struct  QP
 
struct  QpView
 
struct  QpViewBox
 
struct  QpViewBoxMut
 
struct  QpViewMut
 
struct  Workspace
 This class defines the workspace of the dense solver. More...
 

Typedefs

template<typename T >
using SparseMat = Eigen::SparseMatrix<T, 1>
 
template<typename T >
using Vec = Eigen::Matrix<T, DYN, 1>
 
template<typename T >
using VecRef = Eigen::Ref<Vec<T> const>
 
template<typename T >
using VecRefMut = Eigen::Ref<Vec<T>>
 
template<typename T , int l = layout>
using Mat = Eigen::Matrix<T, DYN, DYN, l>
 
template<typename T , int l = layout>
using MatRef = Eigen::Ref<Mat<T, l> const>
 
template<typename T >
using VecMap = Eigen::Map<Vec<T> const>
 
template<typename T >
using VecMapMut = Eigen::Map<Vec<T>>
 
template<typename T , int l = layout>
using MatMap = Eigen::Map<Mat<T, l> const>
 
template<typename T , int l = layout>
using MatMapMut = Eigen::Map<Mat<T, l>>
 
using VecMapISize = Eigen::Map<Eigen::Matrix<isize, DYN, 1> const>
 
using VecISize = Eigen::Matrix<isize, DYN, 1>
 
using VecMapBool = Eigen::Map<Eigen::Matrix<bool, DYN, 1> const>
 
using VecBool = Eigen::Matrix<bool, DYN, 1>
 

Enumerations

enum  { layout = Eigen::RowMajor }
 

Functions

template<typename T >
void compute_backward (dense::QP< T > &solved_qp, VecRef< T > loss_derivative, T eps=1.E-4, T rho_new=1.E-6, T mu_new=1.E-6)
 
template<typename T >
void compute_backward_loss_ESG (dense::QP< T > &solved_qp, VecRef< T > loss_derivative)
 
template<typename T , typename MatIn , typename VecIn1 , typename VecIn2 , typename VecIn3 >
power_iteration (const Eigen::MatrixBase< MatIn > &H, const Eigen::MatrixBase< VecIn1 > &dw, const Eigen::MatrixBase< VecIn2 > &rhs, const Eigen::MatrixBase< VecIn3 > &err_v, T power_iteration_accuracy, isize nb_power_iteration)
 
template<typename T , typename MatIn , typename VecIn1 , typename VecIn2 , typename VecIn3 >
min_eigen_value_via_modified_power_iteration (const Eigen::MatrixBase< MatIn > &H, const Eigen::MatrixBase< VecIn1 > &dw, const Eigen::MatrixBase< VecIn2 > &rhs, const Eigen::MatrixBase< VecIn3 > &err_v, T max_eigen_value, T power_iteration_accuracy, isize nb_power_iteration)
 
template<typename T , typename MatIn >
estimate_minimal_eigen_value_of_symmetric_matrix (const Eigen::MatrixBase< MatIn > &H, EigenValueEstimateMethodOption estimate_method_option, T power_iteration_accuracy, isize nb_power_iteration)
 
template<typename T >
void update_default_rho_with_minimal_Hessian_eigen_value (optional< T > manual_minimal_H_eigenvalue, Results< T > &results, Settings< T > &settings)
 
template<typename T >
void compute_equality_constrained_initial_guess (Workspace< T > &qpwork, const Settings< T > &qpsettings, const Model< T > &qpmodel, const isize n_constraints, const DenseBackend &dense_backend, const HessianType &hessian_type, Results< T > &qpresults)
 
template<typename T >
void setup_factorization (Workspace< T > &qpwork, const Model< T > &qpmodel, Results< T > &qpresults, const DenseBackend &dense_backend, const HessianType &hessian_type)
 
template<typename T >
void setup_equilibration (Workspace< T > &qpwork, const Settings< T > &qpsettings, const bool box_constraints, const HessianType hessian_type, preconditioner::RuizEquilibration< T > &ruiz, bool execute_preconditioner)
 
template<typename T >
void initial_guess (Workspace< T > &qpwork, Settings< T > &qpsettings, Model< T > &qpmodel, Results< T > &qpresults)
 
template<typename T >
void update (optional< MatRef< T > > H, optional< VecRef< T > > g, optional< MatRef< T > > A, optional< VecRef< T > > b, optional< MatRef< T > > C, optional< VecRef< T > > l, optional< VecRef< T > > u, optional< VecRef< T > > l_box, optional< VecRef< T > > u_box, Model< T > &model, Workspace< T > &work, const bool box_constraints)
 
template<typename T >
void setup (optional< MatRef< T > > H, optional< VecRef< T > > g, optional< MatRef< T > > A, optional< VecRef< T > > b, optional< MatRef< T > > C, optional< VecRef< T > > l, optional< VecRef< T > > u, optional< VecRef< T > > l_box, optional< VecRef< T > > u_box, Settings< T > &qpsettings, Model< T > &qpmodel, Workspace< T > &qpwork, Results< T > &qpresults, const bool box_constraints, preconditioner::RuizEquilibration< T > &ruiz, PreconditionerStatus preconditioner_status, const HessianType hessian_type)
 
template<typename T >
void update_proximal_parameters (Settings< T > &settings, Results< T > &results, Workspace< T > &work, optional< T > rho_new, optional< T > mu_eq_new, optional< T > mu_in_new)
 
template<typename T >
void warm_start (optional< VecRef< T > > x_wm, optional< VecRef< T > > y_wm, optional< VecRef< T > > z_wm, Results< T > &results, Settings< T > &settings, Model< T > &model)
 
template<typename T >
bool operator== (const Model< T > &model1, const Model< T > &model2)
 
template<typename T >
bool operator!= (const Model< T > &model1, const Model< T > &model2)
 
template<typename T >
void refactorize (const Model< T > &qpmodel, Results< T > &qpresults, Workspace< T > &qpwork, const isize n_constraints, const DenseBackend &dense_backend, T rho_new)
 
template<typename T >
void mu_update (const Model< T > &qpmodel, Results< T > &qpresults, Workspace< T > &qpwork, isize n_constraints, const DenseBackend &dense_backend, T mu_eq_new, T mu_in_new)
 
template<typename T >
void iterative_residual (const Model< T > &qpmodel, Results< T > &qpresults, Workspace< T > &qpwork, const isize n_constraints, isize inner_pb_dim, const HessianType &hessian_type)
 
template<typename T >
void solve_linear_system (proxsuite::proxqp::dense::Vec< T > &dw, const Model< T > &qpmodel, Results< T > &qpresults, Workspace< T > &qpwork, const isize n_constraints, const DenseBackend &dense_backend, isize inner_pb_dim, proxsuite::linalg::veg::dynstack::DynStackMut &stack)
 
template<typename T >
void iterative_solve_with_permut_fact (const Settings< T > &qpsettings, const Model< T > &qpmodel, Results< T > &qpresults, Workspace< T > &qpwork, const isize n_constraints, const DenseBackend &dense_backend, const HessianType &hessian_type, T eps, isize inner_pb_dim)
 
template<typename T >
void bcl_update (const Settings< T > &qpsettings, Results< T > &qpresults, Workspace< T > &qpwork, T &primal_feasibility_lhs_new, T &bcl_eta_ext, T &bcl_eta_in, T bcl_eta_ext_init, T eps_in_min, T &new_bcl_mu_in, T &new_bcl_mu_eq, T &new_bcl_mu_in_inv, T &new_bcl_mu_eq_inv)
 
template<typename T >
void Martinez_update (const Settings< T > &qpsettings, Results< T > &qpresults, T &primal_feasibility_lhs_new, T &primal_feasibility_lhs_old, T &bcl_eta_in, T eps_in_min, T &new_bcl_mu_in, T &new_bcl_mu_eq, T &new_bcl_mu_in_inv, T &new_bcl_mu_eq_inv)
 
template<typename T >
auto compute_inner_loop_saddle_point (const Model< T > &qpmodel, Results< T > &qpresults, Workspace< T > &qpwork, const Settings< T > &qpsettings) -> T
 
template<typename T >
void primal_dual_semi_smooth_newton_step (const Settings< T > &qpsettings, const Model< T > &qpmodel, Results< T > &qpresults, Workspace< T > &qpwork, const bool box_constraints, const isize n_constraints, const DenseBackend &dense_backend, const HessianType &hessian_type, T eps)
 
template<typename T >
void primal_dual_newton_semi_smooth (const Settings< T > &qpsettings, const Model< T > &qpmodel, Results< T > &qpresults, Workspace< T > &qpwork, const bool box_constraints, const isize n_constraints, preconditioner::RuizEquilibration< T > &ruiz, const DenseBackend &dense_backend, const HessianType &hessian_type, T eps_int)
 
template<typename T >
void qp_solve (const Settings< T > &qpsettings, const Model< T > &qpmodel, Results< T > &qpresults, Workspace< T > &qpwork, const bool box_constraints, const DenseBackend &dense_backend, const HessianType &hessian_type, preconditioner::RuizEquilibration< T > &ruiz)
 
template<typename T >
void print_setup_header (const Settings< T > &settings, const Results< T > &results, const Model< T > &model, const bool box_constraints, const DenseBackend &dense_backend, const HessianType &hessian_type)
 
template<typename Derived >
void save_data (const std::string &filename, const ::Eigen::MatrixBase< Derived > &mat)
 
template<typename T >
void global_primal_residual (const Model< T > &qpmodel, Results< T > &qpresults, const Settings< T > &qpsettings, Workspace< T > &qpwork, const preconditioner::RuizEquilibration< T > &ruiz, const bool box_constraints, T &primal_feasibility_lhs, T &primal_feasibility_eq_rhs_0, T &primal_feasibility_in_rhs_0, T &primal_feasibility_eq_lhs, T &primal_feasibility_in_lhs)
 
template<typename T >
bool global_primal_residual_infeasibility (VectorViewMut< T > ATdy, VectorViewMut< T > CTdz, VectorViewMut< T > dy, VectorViewMut< T > dz, Workspace< T > &qpwork, const Model< T > &qpmodel, const Settings< T > &qpsettings, const bool box_constraints, const preconditioner::RuizEquilibration< T > &ruiz)
 
template<typename T >
bool global_dual_residual_infeasibility (VectorViewMut< T > Adx, VectorViewMut< T > Cdx, VectorViewMut< T > Hdx, VectorViewMut< T > dx, Workspace< T > &qpwork, const Settings< T > &qpsettings, const Model< T > &qpmodel, const bool box_constraints, const preconditioner::RuizEquilibration< T > &ruiz)
 
template<typename T >
void global_dual_residual (Results< T > &qpresults, Workspace< T > &qpwork, const Model< T > &qpmodel, const bool box_constraints, const preconditioner::RuizEquilibration< T > &ruiz, T &dual_feasibility_lhs, T &dual_feasibility_rhs_0, T &dual_feasibility_rhs_1, T &dual_feasibility_rhs_3, T &rhs_duality_gap, T &duality_gap, const HessianType &hessian_type)
 
 VEG_NIEBLOID (fabs)
 
 VEG_NIEBLOID (sqrt)
 
 VEG_NIEBLOID (pow)
 
 VEG_NIEBLOID (infty_norm)
 
template<typename T >
DenseBackend dense_backend_choice (DenseBackend _dense_backend, isize dim, isize n_eq, isize n_in, bool box_constraints)
 This class defines the API of PROXQP solver with dense backend.
 
template<typename T >
proxqp::Results< T > solve (optional< MatRef< T > > H, optional< VecRef< T > > g, optional< MatRef< T > > A, optional< VecRef< T > > b, optional< MatRef< T > > C, optional< VecRef< T > > l, optional< VecRef< T > > u, optional< VecRef< T > > x=nullopt, optional< VecRef< T > > y=nullopt, optional< VecRef< T > > z=nullopt, optional< T > eps_abs=nullopt, optional< T > eps_rel=nullopt, optional< T > rho=nullopt, optional< T > mu_eq=nullopt, optional< T > mu_in=nullopt, optional< bool > verbose=nullopt, bool compute_preconditioner=true, bool compute_timings=false, optional< isize > max_iter=nullopt, proxsuite::proxqp::InitialGuessStatus initial_guess=proxsuite::proxqp::InitialGuessStatus::EQUALITY_CONSTRAINED_INITIAL_GUESS, bool check_duality_gap=false, optional< T > eps_duality_gap_abs=nullopt, optional< T > eps_duality_gap_rel=nullopt, bool primal_infeasibility_solving=false, optional< T > manual_minimal_H_eigenvalue=nullopt)
 
template<typename T >
proxqp::Results< T > solve (optional< MatRef< T > > H, optional< VecRef< T > > g, optional< MatRef< T > > A, optional< VecRef< T > > b, optional< MatRef< T > > C, optional< VecRef< T > > l, optional< VecRef< T > > u, optional< VecRef< T > > l_box, optional< VecRef< T > > u_box, optional< VecRef< T > > x=nullopt, optional< VecRef< T > > y=nullopt, optional< VecRef< T > > z=nullopt, optional< T > eps_abs=nullopt, optional< T > eps_rel=nullopt, optional< T > rho=nullopt, optional< T > mu_eq=nullopt, optional< T > mu_in=nullopt, optional< bool > verbose=nullopt, bool compute_preconditioner=true, bool compute_timings=false, optional< isize > max_iter=nullopt, proxsuite::proxqp::InitialGuessStatus initial_guess=proxsuite::proxqp::InitialGuessStatus::EQUALITY_CONSTRAINED_INITIAL_GUESS, bool check_duality_gap=false, optional< T > eps_duality_gap_abs=nullopt, optional< T > eps_duality_gap_rel=nullopt, bool primal_infeasibility_solving=false, optional< T > manual_minimal_H_eigenvalue=nullopt)
 
template<typename T >
bool operator== (const QP< T > &qp1, const QP< T > &qp2)
 
template<typename T >
bool operator!= (const QP< T > &qp1, const QP< T > &qp2)
 
template<typename T >
void solve_in_parallel (std::vector< proxqp::dense::QP< T > > &qps, const optional< size_t > num_threads=nullopt)
 
template<typename T >
void solve_in_parallel (proxqp::dense::BatchQP< T > &qps, const optional< size_t > num_threads=nullopt)
 
template<typename T >
void qp_solve_in_parallel (optional< const size_t > num_threads, proxqp::dense::BatchQP< T > &qps)
 
template<typename T >
void qp_solve_backward_in_parallel (optional< const size_t > num_threads, std::vector< proxqp::dense::QP< T > > &qps, std::vector< proxqp::dense::Vec< T > > &loss_derivatives, T eps=1.E-4, T rho_new=1.E-6, T mu_new=1.E-6)
 
template<typename T >
void qp_solve_backward_in_parallel (optional< const size_t > num_threads, proxqp::dense::BatchQP< T > &qps, std::vector< proxqp::dense::Vec< T > > &loss_derivatives, T eps=1.E-4, T rho_new=1.E-6, T mu_new=1.E-6)
 

Variables

static constexpr auto DYN = Eigen::Dynamic
 

Typedef Documentation

◆ SparseMat

template<typename T >
using proxsuite::proxqp::dense::SparseMat = Eigen::SparseMatrix<T, 1>

Definition at line 21 of file fwd.hpp.

◆ Vec

template<typename T >
using proxsuite::proxqp::dense::Vec = Eigen::Matrix<T, DYN, 1>

Definition at line 24 of file fwd.hpp.

◆ VecRef

template<typename T >
using proxsuite::proxqp::dense::VecRef = Eigen::Ref<Vec<T> const>

Definition at line 26 of file fwd.hpp.

◆ VecRefMut

template<typename T >
using proxsuite::proxqp::dense::VecRefMut = Eigen::Ref<Vec<T>>

Definition at line 28 of file fwd.hpp.

◆ Mat

template<typename T , int l = layout>
using proxsuite::proxqp::dense::Mat = Eigen::Matrix<T, DYN, DYN, l>

Definition at line 31 of file fwd.hpp.

◆ MatRef

template<typename T , int l = layout>
using proxsuite::proxqp::dense::MatRef = Eigen::Ref<Mat<T, l> const>

Definition at line 33 of file fwd.hpp.

◆ VecMap

template<typename T >
using proxsuite::proxqp::dense::VecMap = Eigen::Map<Vec<T> const>

Definition at line 38 of file fwd.hpp.

◆ VecMapMut

template<typename T >
using proxsuite::proxqp::dense::VecMapMut = Eigen::Map<Vec<T>>

Definition at line 40 of file fwd.hpp.

◆ MatMap

template<typename T , int l = layout>
using proxsuite::proxqp::dense::MatMap = Eigen::Map<Mat<T, l> const>

Definition at line 43 of file fwd.hpp.

◆ MatMapMut

template<typename T , int l = layout>
using proxsuite::proxqp::dense::MatMapMut = Eigen::Map<Mat<T, l>>

Definition at line 45 of file fwd.hpp.

◆ VecMapISize

using proxsuite::proxqp::dense::VecMapISize = Eigen::Map<Eigen::Matrix<isize, DYN, 1> const>

Definition at line 47 of file fwd.hpp.

◆ VecISize

using proxsuite::proxqp::dense::VecISize = Eigen::Matrix<isize, DYN, 1>

Definition at line 48 of file fwd.hpp.

◆ VecMapBool

using proxsuite::proxqp::dense::VecMapBool = Eigen::Map<Eigen::Matrix<bool, DYN, 1> const>

Definition at line 50 of file fwd.hpp.

◆ VecBool

using proxsuite::proxqp::dense::VecBool = Eigen::Matrix<bool, DYN, 1>

Definition at line 51 of file fwd.hpp.

Enumeration Type Documentation

◆ anonymous enum

anonymous enum
Enumerator
layout 

Definition at line 16 of file fwd.hpp.

Function Documentation

◆ compute_backward()

template<typename T >
void proxsuite::proxqp::dense::compute_backward ( dense::QP< T > & solved_qp,
VecRef< T > loss_derivative,
T eps = 1.E-4,
T rho_new = 1.E-6,
T mu_new = 1.E-6 )

Computes the backward pass for the QP solver.

Parameters
solved_qpThe solved quadratic programming (QP) problem.
loss_derivativeThe derivative of the loss function wrt the solution of the QP.
epsThe epsilon value for the iterative refinement, with a default value of 1.E-4.
rho_newThe new primal proximal parameter, with a default value of 1.E-6.
mu_newThe new dual proximal parameter used for both equality and inequality, with a default value of 1.E-6.

derive solution

Definition at line 31 of file compute_ECJ.hpp.

◆ compute_backward_loss_ESG()

template<typename T >
void proxsuite::proxqp::dense::compute_backward_loss_ESG ( dense::QP< T > & solved_qp,
VecRef< T > loss_derivative )

compute backward derivatives

Definition at line 129 of file compute_ECJ.hpp.

◆ power_iteration()

template<typename T , typename MatIn , typename VecIn1 , typename VecIn2 , typename VecIn3 >
T proxsuite::proxqp::dense::power_iteration ( const Eigen::MatrixBase< MatIn > & H,
const Eigen::MatrixBase< VecIn1 > & dw,
const Eigen::MatrixBase< VecIn2 > & rhs,
const Eigen::MatrixBase< VecIn3 > & err_v,
T power_iteration_accuracy,
isize nb_power_iteration )

Definition at line 30 of file helpers.hpp.

◆ min_eigen_value_via_modified_power_iteration()

template<typename T , typename MatIn , typename VecIn1 , typename VecIn2 , typename VecIn3 >
T proxsuite::proxqp::dense::min_eigen_value_via_modified_power_iteration ( const Eigen::MatrixBase< MatIn > & H,
const Eigen::MatrixBase< VecIn1 > & dw,
const Eigen::MatrixBase< VecIn2 > & rhs,
const Eigen::MatrixBase< VecIn3 > & err_v,
T max_eigen_value,
T power_iteration_accuracy,
isize nb_power_iteration )

Definition at line 71 of file helpers.hpp.

◆ estimate_minimal_eigen_value_of_symmetric_matrix()

template<typename T , typename MatIn >
T proxsuite::proxqp::dense::estimate_minimal_eigen_value_of_symmetric_matrix ( const Eigen::MatrixBase< MatIn > & H,
EigenValueEstimateMethodOption estimate_method_option,
T power_iteration_accuracy,
isize nb_power_iteration )

Estimate minimal eigenvalue of a symmetric Matrix

Parameters
Hsymmetric matrix.
EigenValueEstimateMethodOption
power_iteration_accuracypower iteration algorithm accuracy tracked
nb_power_iterationmaximal number of power iteration executed

Definition at line 125 of file helpers.hpp.

◆ update_default_rho_with_minimal_Hessian_eigen_value()

template<typename T >
void proxsuite::proxqp::dense::update_default_rho_with_minimal_Hessian_eigen_value ( optional< T > manual_minimal_H_eigenvalue,
Results< T > & results,
Settings< T > & settings )

Estimate H minimal eigenvalue

Parameters
settingssolver settings
resultssolver results.
manual_minimal_H_eigenvalueminimal H eigenvalue estimate.

Definition at line 176 of file helpers.hpp.

◆ compute_equality_constrained_initial_guess()

template<typename T >
void proxsuite::proxqp::dense::compute_equality_constrained_initial_guess ( Workspace< T > & qpwork,
const Settings< T > & qpsettings,
const Model< T > & qpmodel,
const isize n_constraints,
const DenseBackend & dense_backend,
const HessianType & hessian_type,
Results< T > & qpresults )

Computes the equality constrained initial guess of a QP problem.

Parameters
qpworkworkspace of the solver.
qpsettingssettings of the solver.
qpmodelQP problem as defined by the user (without any scaling performed).
qpresultssolution results.

Definition at line 201 of file helpers.hpp.

◆ setup_factorization()

template<typename T >
void proxsuite::proxqp::dense::setup_factorization ( Workspace< T > & qpwork,
const Model< T > & qpmodel,
Results< T > & qpresults,
const DenseBackend & dense_backend,
const HessianType & hessian_type )

Setups and performs the first factorization of the regularized KKT matrix of the problem.

Parameters
qpworkworkspace of the solver.
qpmodelQP problem model as defined by the user (without any scaling performed).
qpresultssolution results.

Definition at line 241 of file helpers.hpp.

◆ setup_equilibration()

template<typename T >
void proxsuite::proxqp::dense::setup_equilibration ( Workspace< T > & qpwork,
const Settings< T > & qpsettings,
const bool box_constraints,
const HessianType hessian_type,
preconditioner::RuizEquilibration< T > & ruiz,
bool execute_preconditioner )

Performs the equilibration of the QP problem for reducing its ill-conditionness.

Parameters
qpworkworkspace of the solver.
qpsettingssettings of the solver.
ruizruiz preconditioner.
execute_preconditionerboolean variable for executing or not the ruiz preconditioner. If set to False, it uses the previous preconditioning variables (initialized to the identity preconditioner if it is the first scaling performed).

Definition at line 300 of file helpers.hpp.

◆ initial_guess()

template<typename T >
void proxsuite::proxqp::dense::initial_guess ( Workspace< T > & qpwork,
Settings< T > & qpsettings,
Model< T > & qpmodel,
Results< T > & qpresults )

Setups the solver initial guess.

Parameters
qpworksolver workspace.
qpsettingssolver settings.
qpmodelQP problem model as defined by the user (without any scaling performed).
qpresultssolver results.

Definition at line 342 of file helpers.hpp.

◆ update()

template<typename T >
void proxsuite::proxqp::dense::update ( optional< MatRef< T > > H,
optional< VecRef< T > > g,
optional< MatRef< T > > A,
optional< VecRef< T > > b,
optional< MatRef< T > > C,
optional< VecRef< T > > l,
optional< VecRef< T > > u,
optional< VecRef< T > > l_box,
optional< VecRef< T > > u_box,
Model< T > & model,
Workspace< T > & work,
const bool box_constraints )

Updates the QP solver model.

Parameters
Hquadratic cost input defining the QP model.
glinear cost input defining the QP model.
Aequality constraint matrix input defining the QP model.
bequality constraint vector input defining the QP model.
Cinequality constraint matrix input defining the QP model.
llower inequality constraint vector input defining the QP model.
uupper inequality constraint vector input defining the QP model.
qpworksolver workspace.
qpsettingssolver settings.
qpmodelsolver model.
qpresultssolver result.

Definition at line 374 of file helpers.hpp.

◆ setup()

template<typename T >
void proxsuite::proxqp::dense::setup ( optional< MatRef< T > > H,
optional< VecRef< T > > g,
optional< MatRef< T > > A,
optional< VecRef< T > > b,
optional< MatRef< T > > C,
optional< VecRef< T > > l,
optional< VecRef< T > > u,
optional< VecRef< T > > l_box,
optional< VecRef< T > > u_box,
Settings< T > & qpsettings,
Model< T > & qpmodel,
Workspace< T > & qpwork,
Results< T > & qpresults,
const bool box_constraints,
preconditioner::RuizEquilibration< T > & ruiz,
PreconditionerStatus preconditioner_status,
const HessianType hessian_type )

Setups the QP solver model.

Parameters
Hquadratic cost input defining the QP model.
glinear cost input defining the QP model.
Aequality constraint matrix input defining the QP model.
bequality constraint vector input defining the QP model.
Cinequality constraint matrix input defining the QP model.
llower inequality constraint vector input defining the QP model.
uupper inequality constraint vector input defining the QP model.
qpworksolver workspace.
qpsettingssolver settings.
qpmodelsolver model.
qpresultssolver result.
ruizruiz preconditioner.
preconditioner_statusbool variable for deciding whether executing the preconditioning algorithm, or keeping previous preconditioning variables, or using the identity preconditioner (i.e., no preconditioner).

Definition at line 502 of file helpers.hpp.

◆ update_proximal_parameters()

template<typename T >
void proxsuite::proxqp::dense::update_proximal_parameters ( Settings< T > & settings,
Results< T > & results,
Workspace< T > & work,
optional< T > rho_new,
optional< T > mu_eq_new,
optional< T > mu_in_new )

Update the proximal parameters of the results object.

Parameters
rho_newprimal proximal parameter.
mu_eq_newdual equality proximal parameter.
mu_in_newdual inequality proximal parameter.
resultssolver results.

Definition at line 680 of file helpers.hpp.

◆ warm_start()

template<typename T >
void proxsuite::proxqp::dense::warm_start ( optional< VecRef< T > > x_wm,
optional< VecRef< T > > y_wm,
optional< VecRef< T > > z_wm,
Results< T > & results,
Settings< T > & settings,
Model< T > & model )

Warm start the primal and dual variables.

Parameters
x_wmprimal warm start.
y_wmdual equality warm start.
z_wmdual inequality warm start.
resultssolver result.
settingssolver settings.

Definition at line 717 of file helpers.hpp.

◆ operator==() [1/2]

template<typename T >
bool proxsuite::proxqp::dense::operator== ( const Model< T > & model1,
const Model< T > & model2 )

Definition at line 153 of file model.hpp.

◆ operator!=() [1/2]

template<typename T >
bool proxsuite::proxqp::dense::operator!= ( const Model< T > & model1,
const Model< T > & model2 )

Definition at line 167 of file model.hpp.

◆ refactorize()

template<typename T >
void proxsuite::proxqp::dense::refactorize ( const Model< T > & qpmodel,
Results< T > & qpresults,
Workspace< T > & qpwork,
const isize n_constraints,
const DenseBackend & dense_backend,
T rho_new )

Performs a refactorization of the KKT matrix used by the solver.

Parameters
qpworksolver workspace.
qpmodelQP problem model as defined by the user (without any scaling performed).
qpresultssolver results.
rho_newnew primal proximal parameter used for the refactorization.

Definition at line 40 of file solver.hpp.

◆ mu_update()

template<typename T >
void proxsuite::proxqp::dense::mu_update ( const Model< T > & qpmodel,
Results< T > & qpresults,
Workspace< T > & qpwork,
isize n_constraints,
const DenseBackend & dense_backend,
T mu_eq_new,
T mu_in_new )

Updates the dual proximal parameters of the solver (i.e., penalization parameters of the primal-dual merit function).

Parameters
qpworksolver workspace.
qpmodelQP problem model as defined by the user (without any scaling performed).
qpresultssolver results.
mu_eq_newnew dual equality constrained proximal parameter.
mu_in_newnew dual inequality constrained proximal parameter.

Definition at line 130 of file solver.hpp.

◆ iterative_residual()

template<typename T >
void proxsuite::proxqp::dense::iterative_residual ( const Model< T > & qpmodel,
Results< T > & qpresults,
Workspace< T > & qpwork,
const isize n_constraints,
isize inner_pb_dim,
const HessianType & hessian_type )

Derives the residual of the iterative refinement algorithm used for solving associated linear systems of PROXQP algorithm.

Parameters
qpworksolver workspace.
qpmodelQP problem model as defined by the user (without any scaling performed).
qpresultssolver results.
inner_pb_dimdimension of the linear system.

Definition at line 245 of file solver.hpp.

◆ solve_linear_system()

template<typename T >
void proxsuite::proxqp::dense::solve_linear_system ( proxsuite::proxqp::dense::Vec< T > & dw,
const Model< T > & qpmodel,
Results< T > & qpresults,
Workspace< T > & qpwork,
const isize n_constraints,
const DenseBackend & dense_backend,
isize inner_pb_dim,
proxsuite::linalg::veg::dynstack::DynStackMut & stack )

Definition at line 322 of file solver.hpp.

◆ iterative_solve_with_permut_fact()

template<typename T >
void proxsuite::proxqp::dense::iterative_solve_with_permut_fact ( const Settings< T > & qpsettings,
const Model< T > & qpmodel,
Results< T > & qpresults,
Workspace< T > & qpwork,
const isize n_constraints,
const DenseBackend & dense_backend,
const HessianType & hessian_type,
T eps,
isize inner_pb_dim )

Performs iterative refinement for solving associated linear systems of PROXQP algorithm.

Parameters
qpworksolver workspace.
qpmodelQP problem model as defined by the user (without any scaling performed).
qpsettingssolver settings.
qpresultssolver results.
epsaccuracy required for pursuing or not the iterative refinement.
inner_pb_dimdimension of the linear system.

Definition at line 408 of file solver.hpp.

◆ bcl_update()

template<typename T >
void proxsuite::proxqp::dense::bcl_update ( const Settings< T > & qpsettings,
Results< T > & qpresults,
Workspace< T > & qpwork,
T & primal_feasibility_lhs_new,
T & bcl_eta_ext,
T & bcl_eta_in,
T bcl_eta_ext_init,
T eps_in_min,
T & new_bcl_mu_in,
T & new_bcl_mu_eq,
T & new_bcl_mu_in_inv,
T & new_bcl_mu_eq_inv )

BCL rule for updating penalization parameters and accuracy variables.

Parameters
qpworksolver workspace.
qpsettingssolver settings.
qpresultssolver results.
primal_feasibility_lhs_newprimal infeasibility.
bcl_eta_extBCL variable measuring whether the precisely infeasibility is too large or not.
bcl_eta_inBCL variable setting the accuracy required for solving an associated subproblem.
bcl_eta_ext_initinitial BCL bcl_eta_ext variable value.
eps_in_minminimal possible value for bcl_eta_in.
new_bcl_mu_innew value of the inequality constrained penalization parameter.
new_bcl_mu_eqnew value of the equality constrained penalization parameter.
new_bcl_mu_in_invnew value of the inequality constrained penalization parameter (inverse form).
new_bcl_mu_eq_invnew value of the equality constrained penalization parameter (inverse form).

Definition at line 566 of file solver.hpp.

◆ Martinez_update()

template<typename T >
void proxsuite::proxqp::dense::Martinez_update ( const Settings< T > & qpsettings,
Results< T > & qpresults,
T & primal_feasibility_lhs_new,
T & primal_feasibility_lhs_old,
T & bcl_eta_in,
T eps_in_min,
T & new_bcl_mu_in,
T & new_bcl_mu_eq,
T & new_bcl_mu_in_inv,
T & new_bcl_mu_eq_inv )

Martinez rule for updating penalization parameters and accuracy variables.

Parameters
qpworksolver workspace.
qpsettingssolver settings.
qpresultssolver results.
primal_feasibility_lhs_newprimal infeasibility.
bcl_eta_extBCL variable measuring whether the precisely infeasibility is too large or not.
bcl_eta_inBCL variable setting the accuracy required for solving an associated subproblem.
bcl_eta_ext_initinitial BCL bcl_eta_ext variable value.
eps_in_minminimal possible value for bcl_eta_in.
new_bcl_mu_innew value of the inequality constrained penalization parameter.
new_bcl_mu_eqnew value of the equality constrained penalization parameter.
new_bcl_mu_in_invnew value of the inequality constrained penalization parameter (inverse form).
new_bcl_mu_eq_invnew value of the equality constrained penalization parameter (inverse form).

Definition at line 639 of file solver.hpp.

◆ compute_inner_loop_saddle_point()

template<typename T >
auto proxsuite::proxqp::dense::compute_inner_loop_saddle_point ( const Model< T > & qpmodel,
Results< T > & qpresults,
Workspace< T > & qpwork,
const Settings< T > & qpsettings ) -> T

Derives the stopping criterion value used by the Newton semismooth algorithm to minimize the primal-dual augmented Lagrangian function.

Parameters
qpworksolver workspace.
qpmodelQP problem model as defined by the user (without any scaling performed).
qpresultssolver results.

Definition at line 689 of file solver.hpp.

◆ primal_dual_semi_smooth_newton_step()

template<typename T >
void proxsuite::proxqp::dense::primal_dual_semi_smooth_newton_step ( const Settings< T > & qpsettings,
const Model< T > & qpmodel,
Results< T > & qpresults,
Workspace< T > & qpwork,
const bool box_constraints,
const isize n_constraints,
const DenseBackend & dense_backend,
const HessianType & hessian_type,
T eps )

Derives the Newton semismooth step.

Parameters
qpworksolver workspace.
qpmodelQP problem model as defined by the user (without any scaling performed).
qpsettingssolver settings.
qpresultssolver results.
epsaccuracy required for solving the subproblem.

Definition at line 756 of file solver.hpp.

◆ primal_dual_newton_semi_smooth()

template<typename T >
void proxsuite::proxqp::dense::primal_dual_newton_semi_smooth ( const Settings< T > & qpsettings,
const Model< T > & qpmodel,
Results< T > & qpresults,
Workspace< T > & qpwork,
const bool box_constraints,
const isize n_constraints,
preconditioner::RuizEquilibration< T > & ruiz,
const DenseBackend & dense_backend,
const HessianType & hessian_type,
T eps_int )

Performs the Newton semismooth algorithm to minimize the primal-dual augmented Lagrangian function used by PROXQP algorithm.

Parameters
qpworksolver workspace.
qpmodelQP problem model as defined by the user (without any scaling performed).
qpsettingssolver settings.
qpresultssolver results.
ruizruiz preconditioner.
eps_intaccuracy required for solving the subproblem.

Definition at line 884 of file solver.hpp.

◆ qp_solve()

template<typename T >
void proxsuite::proxqp::dense::qp_solve ( const Settings< T > & qpsettings,
const Model< T > & qpmodel,
Results< T > & qpresults,
Workspace< T > & qpwork,
const bool box_constraints,
const DenseBackend & dense_backend,
const HessianType & hessian_type,
preconditioner::RuizEquilibration< T > & ruiz )

Executes the PROXQP algorithm.

Parameters
qpworksolver workspace.
qpmodelQP problem model as defined by the user (without any scaling performed).
qpsettingssolver settings.
qpresultssolver results.
ruizruiz preconditioner.

\ TODO in a quicker way

\ TODO in a quicker way

\ TODO in a quicker way

\ TODO in a quicker way

effective mu upddate

Definition at line 1090 of file solver.hpp.

◆ print_setup_header()

template<typename T >
void proxsuite::proxqp::dense::print_setup_header ( const Settings< T > & settings,
const Results< T > & results,
const Model< T > & model,
const bool box_constraints,
const DenseBackend & dense_backend,
const HessianType & hessian_type )

Definition at line 33 of file utils.hpp.

◆ save_data()

template<typename Derived >
void proxsuite::proxqp::dense::save_data ( const std::string & filename,
const ::Eigen::MatrixBase< Derived > & mat )

Save a matrix into a CSV format. Used for debug purposes.

Parameters
filenamefilename name for the CSV.
matmatrix to save into CSV format.

Definition at line 133 of file utils.hpp.

◆ global_primal_residual()

template<typename T >
void proxsuite::proxqp::dense::global_primal_residual ( const Model< T > & qpmodel,
Results< T > & qpresults,
const Settings< T > & qpsettings,
Workspace< T > & qpwork,
const preconditioner::RuizEquilibration< T > & ruiz,
const bool box_constraints,
T & primal_feasibility_lhs,
T & primal_feasibility_eq_rhs_0,
T & primal_feasibility_in_rhs_0,
T & primal_feasibility_eq_lhs,
T & primal_feasibility_in_lhs )

Derives the global primal residual of the QP problem.

Parameters
qpworksolver workspace.
qpmodelQP problem model as defined by the user (without any scaling performed).
qpresultssolver results.
ruizruiz preconditioner.
primal_feasibility_lhsprimal infeasibility.
primal_feasibility_eq_rhs_0scalar variable used when using a relative stopping criterion.
primal_feasibility_in_rhs_0scalar variable used when using a relative stopping criterion.
primal_feasibility_eq_lhsscalar variable used when using a relative stopping criterion.
primal_feasibility_in_lhsscalar variable used when using a relative stopping criterion.

Definition at line 166 of file utils.hpp.

◆ global_primal_residual_infeasibility()

template<typename T >
bool proxsuite::proxqp::dense::global_primal_residual_infeasibility ( VectorViewMut< T > ATdy,
VectorViewMut< T > CTdz,
VectorViewMut< T > dy,
VectorViewMut< T > dz,
Workspace< T > & qpwork,
const Model< T > & qpmodel,
const Settings< T > & qpsettings,
const bool box_constraints,
const preconditioner::RuizEquilibration< T > & ruiz )

Check whether the global primal infeasibility criterion is satisfied.

Parameters
qpworksolver workspace.
qpsettingssolver settings.
ruizruiz preconditioner.
ATdyvariable used for testing global primal infeasibility criterion is satisfied.
CTdzvariable used for testing global primal infeasibility criterion is satisfied.
dyvariable used for testing global primal infeasibility criterion is satisfied.
dzvariable used for testing global primal infeasibility criterion is satisfied.

Definition at line 271 of file utils.hpp.

◆ global_dual_residual_infeasibility()

template<typename T >
bool proxsuite::proxqp::dense::global_dual_residual_infeasibility ( VectorViewMut< T > Adx,
VectorViewMut< T > Cdx,
VectorViewMut< T > Hdx,
VectorViewMut< T > dx,
Workspace< T > & qpwork,
const Settings< T > & qpsettings,
const Model< T > & qpmodel,
const bool box_constraints,
const preconditioner::RuizEquilibration< T > & ruiz )

Check whether the global dual infeasibility criterion is satisfied.

Parameters
qpworksolver workspace.
qpsettingssolver settings.
qpmodelQP problem model as defined by the user (without any scaling performed).
ruizruiz preconditioner.
Adxvariable used for testing global dual infeasibility criterion is satisfied.
Cdxvariable used for testing global dual infeasibility criterion is satisfied.
Hdxvariable used for testing global dual infeasibility criterion is satisfied.
dxvariable used for testing global dual infeasibility criterion is satisfied.

Definition at line 345 of file utils.hpp.

◆ global_dual_residual()

template<typename T >
void proxsuite::proxqp::dense::global_dual_residual ( Results< T > & qpresults,
Workspace< T > & qpwork,
const Model< T > & qpmodel,
const bool box_constraints,
const preconditioner::RuizEquilibration< T > & ruiz,
T & dual_feasibility_lhs,
T & dual_feasibility_rhs_0,
T & dual_feasibility_rhs_1,
T & dual_feasibility_rhs_3,
T & rhs_duality_gap,
T & duality_gap,
const HessianType & hessian_type )

Derives the global dual residual of the QP problem.

Parameters
qpworksolver workspace.
qpresultssolver results.
ruizruiz preconditioner.
dual_feasibility_lhsprimal infeasibility.
primal_feasibility_eq_rhs_0scalar variable used when using a relative stopping criterion.
dual_feasibility_rhs_0scalar variable used when using a relative stopping criterion.
dual_feasibility_rhs_1scalar variable used when using a relative stopping criterion.
dual_feasibility_rhs_3scalar variable used when using a relative stopping criterion.

Definition at line 439 of file utils.hpp.

◆ VEG_NIEBLOID() [1/4]

proxsuite::proxqp::dense::VEG_NIEBLOID ( fabs )

◆ VEG_NIEBLOID() [2/4]

proxsuite::proxqp::dense::VEG_NIEBLOID ( sqrt )

◆ VEG_NIEBLOID() [3/4]

proxsuite::proxqp::dense::VEG_NIEBLOID ( pow )

◆ VEG_NIEBLOID() [4/4]

proxsuite::proxqp::dense::VEG_NIEBLOID ( infty_norm )

◆ dense_backend_choice()

template<typename T >
DenseBackend proxsuite::proxqp::dense::dense_backend_choice ( DenseBackend _dense_backend,
isize dim,
isize n_eq,
isize n_in,
bool box_constraints )

This class defines the API of PROXQP solver with dense backend.

Wrapper class for using proxsuite API with dense backend for solving linearly constrained convex QP problem using ProxQp algorithm.

Example usage:

#include <Eigen/Core>
#include <Eigen/Cholesky>
#include <util.hpp>
using T = double;
auto main() -> int {
// Generate a random QP problem with primal variable dimension of size
dim; n_eq equality constraints and n_in inequality constraints
::proxsuite::proxqp::test::rand::set_seed(1);
proxqp::isize dim = 10;
proxqp::isize n_eq(dim / 4);
proxqp::isize n_in(dim / 4);
T strong_convexity_factor(1.e-2);
T sparsity_factor = 0.15; // controls the sparsity of each matrix of the
problem generated T eps_abs = T(1e-9); Qp<T> qp{
random_with_dim_and_neq_and_n_in,
dim,
n_eq,
n_in,
sparsity_factor,
strong_convexity_factor};
// Solve the problem
proxqp::dense::QP<T> Qp{dim, n_eq, n_in}; // creating QP object
Qp.settings.eps_abs = eps_abs; // choose accuracy needed
Qp.init(qp.H, qp.g, qp.A, qp.b, qp.C, qp.u, qp.l); // setup the QP
object Qp.solve(); // solve the problem
// Verify solution accuracy
T pri_res = std::max(
(qp.A * Qp.results.x - qp.b).lpNorm<Eigen::Infinity>(),
(helpers::positive_part(qp.C * Qp.results.x -
qp.u) + helpers::negative_part(qp.C * Qp.results.x - qp.l))
.lpNorm<Eigen::Infinity>());
T dua_res = (qp.H * Qp.results.x + qp.g + qp.A.transpose() *
Qp.results.y + qp.C.transpose() * Qp.results.z) .lpNorm<Eigen::Infinity>();
VEG_ASSERT(pri_res <= eps_abs);
VEG_ASSERT(dua_res <= eps_abs);
// Some solver statistics
std::cout << "------solving qp with dim: " << dim
<< " neq: " << n_eq << " nin: "
<< n_in << std::endl; std::cout << "primal residual: " << pri_res << std::endl;
std::cout << "dual residual: " << dua_res << std::endl;
std::cout << "total number of iteration: " << Qp.results.info.iter
<< std::endl;
}
#define VEG_ASSERT(...)
auto positive_part(T const &expr) PROXSUITE_DEDUCE_RET((expr.array() > 0).select(expr
Returns the positive part of an expression.
auto negative_part(T const &expr) PROXSUITE_DEDUCE_RET((expr.array()< 0).select(expr
Returns the negative part of an expression.

Definition at line 83 of file wrapper.hpp.

◆ solve() [1/2]

template<typename T >
proxqp::Results< T > proxsuite::proxqp::dense::solve ( optional< MatRef< T > > H,
optional< VecRef< T > > g,
optional< MatRef< T > > A,
optional< VecRef< T > > b,
optional< MatRef< T > > C,
optional< VecRef< T > > l,
optional< VecRef< T > > u,
optional< VecRef< T > > x = nullopt,
optional< VecRef< T > > y = nullopt,
optional< VecRef< T > > z = nullopt,
optional< T > eps_abs = nullopt,
optional< T > eps_rel = nullopt,
optional< T > rho = nullopt,
optional< T > mu_eq = nullopt,
optional< T > mu_in = nullopt,
optional< bool > verbose = nullopt,
bool compute_preconditioner = true,
bool compute_timings = false,
optional< isize > max_iter = nullopt,
proxsuite::proxqp::InitialGuessStatus initial_guess = proxsuite::proxqp::InitialGuessStatus::EQUALITY_CONSTRAINED_INITIAL_GUESS,
bool check_duality_gap = false,
optional< T > eps_duality_gap_abs = nullopt,
optional< T > eps_duality_gap_rel = nullopt,
bool primal_infeasibility_solving = false,
optional< T > manual_minimal_H_eigenvalue = nullopt )

Solves the QP problem using PROXQP algorithm without the need to define a QP object, with matrices defined by Dense Eigen matrices. It is possible to set up some of the solver parameters (warm start, initial guess option, proximal step sizes, absolute and relative accuracies, maximum number of iterations, preconditioner execution). There are no box constraints in the model.

Parameters
Hquadratic cost input defining the QP model.
glinear cost input defining the QP model.
Aequality constraint matrix input defining the QP model.
bequality constraint vector input defining the QP model.
Cinequality constraint matrix input defining the QP model.
llower inequality constraint vector input defining the QP model.
uupper inequality constraint vector input defining the QP model.
xprimal warm start.
ydual equality constraint warm start.
zdual inequality constraint warm start.
verboseif set to true, the solver prints more information about each iteration.
compute_preconditionerbool parameter for executing or not the preconditioner.
compute_timingsboolean parameter for computing the solver timings.
rhoproximal step size wrt primal variable.
mu_eqproximal step size wrt equality constrained multiplier.
mu_inproximal step size wrt inequality constrained multiplier.
eps_absabsolute accuracy threshold.
eps_relrelative accuracy threshold.
max_itermaximum number of iteration.
initial_guessinitial guess option for warm starting or not the initial iterate values.
check_duality_gapIf set to true, include the duality gap in absolute and relative stopping criteria.
eps_duality_gap_absabsolute accuracy threshold for the duality-gap criterion.
eps_duality_gap_relrelative accuracy threshold for the duality-gap criterion.

Definition at line 1002 of file wrapper.hpp.

◆ solve() [2/2]

template<typename T >
proxqp::Results< T > proxsuite::proxqp::dense::solve ( optional< MatRef< T > > H,
optional< VecRef< T > > g,
optional< MatRef< T > > A,
optional< VecRef< T > > b,
optional< MatRef< T > > C,
optional< VecRef< T > > l,
optional< VecRef< T > > u,
optional< VecRef< T > > l_box,
optional< VecRef< T > > u_box,
optional< VecRef< T > > x = nullopt,
optional< VecRef< T > > y = nullopt,
optional< VecRef< T > > z = nullopt,
optional< T > eps_abs = nullopt,
optional< T > eps_rel = nullopt,
optional< T > rho = nullopt,
optional< T > mu_eq = nullopt,
optional< T > mu_in = nullopt,
optional< bool > verbose = nullopt,
bool compute_preconditioner = true,
bool compute_timings = false,
optional< isize > max_iter = nullopt,
proxsuite::proxqp::InitialGuessStatus initial_guess = proxsuite::proxqp::InitialGuessStatus::EQUALITY_CONSTRAINED_INITIAL_GUESS,
bool check_duality_gap = false,
optional< T > eps_duality_gap_abs = nullopt,
optional< T > eps_duality_gap_rel = nullopt,
bool primal_infeasibility_solving = false,
optional< T > manual_minimal_H_eigenvalue = nullopt )

Solves the QP problem using PROXQP algorithm without the need to define a QP object, with matrices defined by Dense Eigen matrices. It is possible to set up some of the solver parameters (warm start, initial guess option, proximal step sizes, absolute and relative accuracies, maximum number of iterations, preconditioner execution).

Parameters
Hquadratic cost input defining the QP model.
glinear cost input defining the QP model.
Aequality constraint matrix input defining the QP model.
bequality constraint vector input defining the QP model.
Cinequality constraint matrix input defining the QP model.
llower inequality constraint vector input defining the QP model.
uupper inequality constraint vector input defining the QP model.
l_boxlower box inequality constraint vector input defining the QP model.
u_boxupper box inequality constraint vector input defining the QP model.
xprimal warm start.
ydual equality constraint warm start.
zdual inequality constraint warm start. The upper part must contain a warm start for inequality constraints wrt C matrix, whereas the latter wrt the box inequalities.
verboseif set to true, the solver prints more information about each iteration.
compute_preconditionerbool parameter for executing or not the preconditioner.
compute_timingsboolean parameter for computing the solver timings.
rhoproximal step size wrt primal variable.
mu_eqproximal step size wrt equality constrained multiplier.
mu_inproximal step size wrt inequality constrained multiplier.
eps_absabsolute accuracy threshold.
eps_relrelative accuracy threshold.
max_itermaximum number of iteration.
initial_guessinitial guess option for warm starting or not the initial iterate values.
check_duality_gapIf set to true, include the duality gap in absolute and relative stopping criteria.
eps_duality_gap_absabsolute accuracy threshold for the duality-gap criterion.
eps_duality_gap_relrelative accuracy threshold for the duality-gap criterion.

Definition at line 1132 of file wrapper.hpp.

◆ operator==() [2/2]

template<typename T >
bool proxsuite::proxqp::dense::operator== ( const QP< T > & qp1,
const QP< T > & qp2 )

Definition at line 1237 of file wrapper.hpp.

◆ operator!=() [2/2]

template<typename T >
bool proxsuite::proxqp::dense::operator!= ( const QP< T > & qp1,
const QP< T > & qp2 )

Definition at line 1247 of file wrapper.hpp.

◆ solve_in_parallel() [1/2]

template<typename T >
void proxsuite::proxqp::dense::solve_in_parallel ( std::vector< proxqp::dense::QP< T > > & qps,
const optional< size_t > num_threads = nullopt )

Definition at line 19 of file qp_solve.hpp.

◆ solve_in_parallel() [2/2]

template<typename T >
void proxsuite::proxqp::dense::solve_in_parallel ( proxqp::dense::BatchQP< T > & qps,
const optional< size_t > num_threads = nullopt )

Definition at line 42 of file qp_solve.hpp.

◆ qp_solve_in_parallel()

template<typename T >
void proxsuite::proxqp::dense::qp_solve_in_parallel ( optional< const size_t > num_threads,
proxqp::dense::BatchQP< T > & qps )

Definition at line 64 of file qp_solve.hpp.

◆ qp_solve_backward_in_parallel() [1/2]

template<typename T >
void proxsuite::proxqp::dense::qp_solve_backward_in_parallel ( optional< const size_t > num_threads,
std::vector< proxqp::dense::QP< T > > & qps,
std::vector< proxqp::dense::Vec< T > > & loss_derivatives,
T eps = 1.E-4,
T rho_new = 1.E-6,
T mu_new = 1.E-6 )

Definition at line 86 of file qp_solve.hpp.

◆ qp_solve_backward_in_parallel() [2/2]

template<typename T >
void proxsuite::proxqp::dense::qp_solve_backward_in_parallel ( optional< const size_t > num_threads,
proxqp::dense::BatchQP< T > & qps,
std::vector< proxqp::dense::Vec< T > > & loss_derivatives,
T eps = 1.E-4,
T rho_new = 1.E-6,
T mu_new = 1.E-6 )

Definition at line 114 of file qp_solve.hpp.

Variable Documentation

◆ DYN

constexpr auto proxsuite::proxqp::dense::DYN = Eigen::Dynamic
staticconstexpr

Definition at line 15 of file fwd.hpp.