TVM  0.9.4
QuadprogLeastSquareSolver.h
Go to the documentation of this file.
1 /* Copyright 2017-2020 CNRS-AIST JRL and CNRS-UM LIRMM */
2 
3 #pragma once
4 
6 
7 #include <Eigen/QR>
8 
9 #include <eigen-quadprog/QuadProg.h>
10 
11 namespace tvm
12 {
13 
14 namespace solver
15 {
16 class QuadprogLSSolverFactory;
17 
20 {
21  TVM_ADD_NON_DEFAULT_OPTION(big_number, constant::big_number)
22  TVM_ADD_NON_DEFAULT_OPTION(cholesky, false)
23  TVM_ADD_NON_DEFAULT_OPTION(choleskyDamping, 1e-6)
24  TVM_ADD_NON_DEFAULT_OPTION(damping, 1e-12)
25  TVM_ADD_NON_DEFAULT_OPTION(verbose, false)
26 
27 public:
29 };
30 
33 {
34 public:
36 
37 protected:
38  void initializeBuild_(int nObj, int nEq, int nIneq, bool useBounds) override;
39  ImpactFromChanges resize_(int nObj, int nEq, int nIneq, bool useBounds) override;
40  void addBound_(LinearConstraintPtr bound, RangePtr range, bool first) override;
43  void addObjective_(LinearConstraintPtr cstr, SolvingRequirementsPtr req, double additionalWeight) override;
44  void setMinimumNorm_() override;
45  void resetBounds_() override;
46  void preAssignmentProcess_() override;
47  void postAssignmentProcess_() override;
48  bool solve_() override;
49  virtual const Eigen::VectorXd & result_() const override;
50  bool handleDoubleSidedConstraint_() const override { return false; }
54  void removeBounds_(const Range & range) override;
59 
60  void applyImpactLogic(ImpactFromChanges & impact) override;
61 
62  void printProblemData_() const override;
63  void printDiagnostic_() const override;
64 
65 private:
66  using VectorXdSeg = decltype(Eigen::VectorXd().segment(0, 1));
67  using MatrixXdRows = decltype(Eigen::MatrixXd().middleRows(0, 1));
68 
69  Eigen::MatrixXd D_; // We have Q = D^T D
70  Eigen::VectorXd e_; // We have c = D^t e
71  Eigen::MatrixXd Q_;
72  Eigen::VectorXd c_;
73  Eigen::MatrixXd A_;
74  Eigen::VectorXd b_;
75 
76  MatrixXdRows Aineq_; // part of A_ corresponding to inequality constraints
77  VectorXdSeg bineq_; // part of b_ corresponding to inequality constraints
78  VectorXdSeg xl_; // part of b_ corresponding to lower bound constraints
79  VectorXdSeg xu_; // part of b_ corresponding to upper bound constraints
80 
81  Eigen::QuadProgDense qpd_;
82  Eigen::HouseholderQR<Eigen::MatrixXd> qr_; // TODO add option for ColPiv variant
83 
84  bool autoMinNorm_;
85  bool underspecifiedObj_; // true when nObj<n
86  int nIneqInclBounds_; // number of inequality constraints including bounds.
87 
88  // options
89  double big_number_;
90  double damping_; // value added to the diagonal of Q for regularization (non Cholesky case)
91  bool cholesky_; // compute the Cholesky decomposition before calling the solver.
92  double choleskyDamping_; // if nObj<n, the cholesky factor R is trapezoidal. A multiple of
93  // the identity is used to make it triangular using this value.
94 };
95 
100 {
101 public:
102  std::unique_ptr<abstract::LSSolverFactory> clone() const override;
103 
106 
107  std::unique_ptr<abstract::LeastSquareSolver> createSolver() const override;
108 
109 private:
110  QuadprogLSSolverOptions options_;
111 };
112 
113 } // namespace solver
114 
115 } // namespace tvm
#define TVM_ADD_NON_DEFAULT_OPTION(optionName, defaultValue)
Definition: Option.h:30
#define TVM_DLLAPI
Definition: api.h:35
Definition: Range.h:19
Definition: LinearConstraint.h:56
Definition: AssignmentTarget.h:37
Definition: QuadprogLeastSquareSolver.h:100
QuadprogLSSolverFactory(const QuadprogLSSolverOptions &options={})
std::unique_ptr< abstract::LSSolverFactory > clone() const override
std::unique_ptr< abstract::LeastSquareSolver > createSolver() const override
Definition: QuadprogLeastSquareSolver.h:20
Definition: QuadprogLeastSquareSolver.h:33
void updateObjectiveTargetData(scheme::internal::AssignmentTarget &target) override
void initializeBuild_(int nObj, int nEq, int nIneq, bool useBounds) override
void applyImpactLogic(ImpactFromChanges &impact) override
void removeBounds_(const Range &range) override
Range nextInequalityConstraintRange_(const constraint::abstract::LinearConstraint &cstr) const override
Range nextObjectiveRange_(const constraint::abstract::LinearConstraint &cstr) const override
void addEqualityConstraint_(LinearConstraintPtr cstr) override
void updateBoundTargetData(scheme::internal::AssignmentTarget &target) override
void printProblemData_() const override
void addBound_(LinearConstraintPtr bound, RangePtr range, bool first) override
QuadprogLeastSquareSolver(const QuadprogLSSolverOptions &options={})
ImpactFromChanges resize_(int nObj, int nEq, int nIneq, bool useBounds) override
void addObjective_(LinearConstraintPtr cstr, SolvingRequirementsPtr req, double additionalWeight) override
void updateInequalityTargetData(scheme::internal::AssignmentTarget &target) override
virtual const Eigen::VectorXd & result_() const override
void updateEqualityTargetData(scheme::internal::AssignmentTarget &target) override
void addIneqalityConstraint_(LinearConstraintPtr cstr) override
bool handleDoubleSidedConstraint_() const override
Definition: QuadprogLeastSquareSolver.h:50
Range nextEqualityConstraintRange_(const constraint::abstract::LinearConstraint &cstr) const override
void printDiagnostic_() const override
Definition: LeastSquareSolver.h:220
Definition: LeastSquareSolver.h:37
Definition: Clock.h:12
std::shared_ptr< constraint::abstract::LinearConstraint > LinearConstraintPtr
Definition: defs.h:59
std::shared_ptr< Range > RangePtr
Definition: defs.h:61
std::shared_ptr< requirements::SolvingRequirementsWithCallbacks > SolvingRequirementsPtr
Definition: defs.h:63