Quadratic Programming problem
scs-eigen can be used to solve Quadratic Programming (QP) problems in the for
where is the optimization variable. The objective function is defined by a positive semidefinite matrix and vector . The linear constraints are defined by matrix and vectors and so that and for all .
The following example shows how scs-eigen can be used to solve the QP problem:
First of all you should include ScsEigen
#include <ScsEigen/ScsEigen.h>
You can also define the Hessian, gradient the constraint matrix and vectors.
Eigen::Matrix2d H; H << 3, 2, 2, 4; Eigen::Vector2d gradient; gradient << 3, 1; Eigen::MatrixXd A(3,2); A.setZero(); A(0,0) = 1; A(0,1) = 1; A(1,0) = 1; A(2,1) = 1; Eigen::Vector3d lowerBound; lowerBound << 1, 0, 0; Eigen::Vector3d upperBound; upperBound << 1, 0.7, 0.7;
Once the matrices used to described the QP problem have been defined you can create ScsEigen::
and initialize the number of variables.
ScsEigen::Solver solver; solver.setNumberOfVariabels(2);
Now you can set the constraints and the cost using ScsEigen::
class
solver.mathematicalProgram().addLinearConstraint( std::make_shared<ScsEigen::LinearConstraint>(A, lowerBound, upperBound), "linear constraint"); solver.mathematicalProgram().addQuadraticCost( std::make_shared<ScsEigen::QuadraticCost>(H, gradient), "quadratic cost");
You can finally solve the problem and get the solution
solver.solve(); Eigen::Vector2d solution = solver.solution().solution;
The complete example follows
### CMakeLists.txt project(QP) find_package(ScsEigen REQUIRED) add_executable(QP qp.cpp) target_link_libraries(QP ScsEigen::ScsEigen)
/// qp.cpp #include <ScsEigen/ScsEigen.h> #include <Eigen/Dense> int main() { Eigen::Matrix2d H; H << 3, 2, 2, 4; Eigen::Vector2d gradient; gradient << 3, 1; Eigen::MatrixXd A(3,2); A.setZero(); A(0,0) = 1; A(0,1) = 1; A(1,0) = 1; A(2,1) = 1; Eigen::Vector3d lowerBound; lowerBound << 1, 0, 0; Eigen::Vector3d upperBound; upperBound << 1, 0.7, 0.7; ScsEigen::Solver solver; solver.mathematicalProgram().setNumberOVfariables(2); if (!solver.mathematicalProgram().addQuadraticCost( std::make_shared<ScsEigen::QuadraticCost>(H, gradient), "quadratic cost")) return EXIT_FAILURE; if (!solver.mathematicalProgram().addLinearConstraint( std::make_shared<ScsEigen::LinearConstraint>(A, lowerBound, upperBound), "linear constraint")) return EXIT_FAILURE; if (!solver.solve()) return EXIT_FAILURE; if (!solver.solution().isValid()) return EXIT_FAILURE; if (solver.solution().status != ScsEigen::Solution::Status::solved) return EXIT_FAILURE; Eigen::Vector2d solution = solver.solution().solution; return EXIT_SUCCESS; }