TVM  0.9.4
HierarchicalLeastSquareSolver.h
Go to the documentation of this file.
1 /* Copyright 2017-2020 CNRS-AIST JRL and CNRS-UM LIRMM */
2 
3 #pragma once
4 
5 #include <tvm/api.h>
6 #include <tvm/defs.h>
7 
12 #include <tvm/utils/internal/map.h>
13 
14 namespace tvm
15 {
16 
17 namespace solver
18 {
19 
20 namespace abstract
21 {
35 {
36 public:
37  HierarchicalLeastSquareSolver(bool verbose = false);
40  virtual ~HierarchicalLeastSquareSolver() = default;
54  void startBuild(const VariableVector & x,
55  const std::vector<int> & nEq,
56  const std::vector<int> & nIneq,
57  bool useBounds,
58  const hint::internal::Substitutions * const subs = nullptr);
60  void finalizeBuild();
61 
66 
69 
76 
80  bool solve();
81 
83  const Eigen::VectorXd & result() const;
84 
90 
96  void process(const internal::SolverEvents & se);
97 
99  int numberOfLevels() const { return useBounds_ ? static_cast<int>(nEq_.size()) + 1 : static_cast<int>(nEq_.size()); }
100 
101 protected:
103  {
104  ImpactFromChanges(int nLvl);
105  ImpactFromChanges(const std::vector<bool> & eq, std::vector<bool> & ineq, bool bounds);
106 
107  std::vector<bool> equalityConstraints_;
108  std::vector<bool> inequalityConstraints_;
109  bool bounds_ = false;
110  int newLevels_ = 0;
111 
112  template<typename Derived>
113  static bool willReallocate(const Eigen::DenseBase<Derived> & M, int rows, int cols = 1);
114 
119 
120  bool any() const
121  {
122  return std::any_of(equalityConstraints_.begin(), equalityConstraints_.end(), [](bool b) { return b; })
123  || std::any_of(inequalityConstraints_.begin(), inequalityConstraints_.end(), [](bool b) { return b; })
124  || bounds_;
125  }
126  };
127 
128  virtual void initializeBuild_(const std::vector<int> & nEq, const std::vector<int> & nIneq, bool useBounds) = 0;
129  virtual ImpactFromChanges resize_(const std::vector<int> & nEq, const std::vector<int> & nIneq, bool useBounds) = 0;
130  virtual void addBound_(LinearConstraintPtr bound, RangePtr range, bool first) = 0;
133  virtual void setMinimumNorm_() = 0;
134  virtual void resetBounds_() = 0;
135  virtual void preAssignmentProcess_() {}
136  virtual void postAssignmentProcess_() {}
137  virtual bool solve_() = 0;
138  virtual const Eigen::VectorXd & result_() const = 0;
139  virtual bool handleDoubleSidedConstraint_() const = 0;
145  virtual void removeBounds_(const Range & range) = 0;
149 
153  virtual void applyImpactLogic(ImpactFromChanges & impact);
154 
155  virtual void printProblemData_() const = 0;
156  virtual void printDiagnostic_() const = 0;
157 
158  const VariableVector & variables() const { return *variables_; }
159  const hint::internal::Substitutions * substitutions() const { return subs_; }
160 
161  template<typename... Args>
162  void addAssignement(Args &&... args);
163 
164 private:
165  void updateWeights(const internal::SolverEvents & se);
166  bool updateVariables(const internal::SolverEvents & se);
167  ImpactFromChanges processRemovedConstraints(const internal::SolverEvents & se);
168  ImpactFromChanges previewAddedConstraints(const internal::SolverEvents & se);
169  void processAddedConstraints(const internal::SolverEvents & se);
170 
171 public:
173  {
174  template<typename... Args>
175  MarkedAssignment(Args &&... args) : assignment(std::forward<Args>(args)...), markedForRemoval(false)
176  {}
179  };
180  template<typename K, typename T>
182  using AssignmentVector = std::vector<std::unique_ptr<MarkedAssignment>>;
183  using AssignmentPtrVector = std::vector<MarkedAssignment *>;
185 
186 protected:
187  bool useBounds_ = false;
188  std::vector<int> nEq_;
189  std::vector<int> nIneq_;
190  std::vector<int> eqSize_;
191  std::vector<int> ineqSize_;
192 
193 private:
194  bool buildInProgress_;
195  bool verbose_;
196  VariableVector const * variables_;
198  AssignmentVector assignments_;
203  std::vector<MapToAssignment> equalityConstraintToAssignments_;
204  std::vector<MapToAssignment> inequalityConstraintToAssignments_;
205  MapToAssignment boundToAssignments_;
206  const hint::internal::Substitutions * subs_;
207 };
208 
215 {
216 protected:
217  HLSSolverFactory(const std::string & solverName) : solverName_(solverName) {}
218 
219 public:
220  virtual ~HLSSolverFactory() = default;
221 
222  virtual std::unique_ptr<HLSSolverFactory> clone() const = 0;
223  virtual std::unique_ptr<HierarchicalLeastSquareSolver> createSolver() const = 0;
224 
225 private:
226  std::string solverName_;
227 };
228 
229 template<typename... Args>
231 { assignments_.emplace_back(new MarkedAssignment(std::forward<Args>(args)...)); }
232 
233 template<typename Derived>
234 inline bool HierarchicalLeastSquareSolver::ImpactFromChanges::willReallocate(const Eigen::DenseBase<Derived> & M,
235  int rows,
236  int cols)
237 { return M.rows() * M.cols() != rows * cols; }
238 
239 } // namespace abstract
240 
241 } // namespace solver
242 
243 } // namespace tvm
#define TVM_DLLAPI
Definition: api.h:35
Definition: Range.h:19
Definition: VariableVector.h:41
Definition: LinearConstraint.h:56
Definition: Substitutions.h:24
Definition: AssignmentTarget.h:37
Definition: Assignment.h:32
Definition: HierarchicalLeastSquareSolver.h:215
virtual std::unique_ptr< HLSSolverFactory > clone() const =0
HLSSolverFactory(const std::string &solverName)
Definition: HierarchicalLeastSquareSolver.h:217
virtual std::unique_ptr< HierarchicalLeastSquareSolver > createSolver() const =0
Definition: HierarchicalLeastSquareSolver.h:35
virtual void addEqualityConstraint_(LinearConstraintPtr cstr, SolvingRequirementsPtr req)=0
void process(const internal::SolverEvents &se)
virtual void updateEqualityTargetData(int lvl, scheme::internal::AssignmentTarget &target)=0
virtual void addIneqalityConstraint_(LinearConstraintPtr cstr, SolvingRequirementsPtr req)=0
virtual ImpactFromChanges resize_(const std::vector< int > &nEq, const std::vector< int > &nIneq, bool useBounds)=0
HierarchicalLeastSquareSolver & operator=(const HierarchicalLeastSquareSolver &)=delete
std::vector< MarkedAssignment * > AssignmentPtrVector
Definition: HierarchicalLeastSquareSolver.h:183
std::vector< int > nIneq_
Definition: HierarchicalLeastSquareSolver.h:189
std::vector< int > eqSize_
Definition: HierarchicalLeastSquareSolver.h:190
const hint::internal::Substitutions * substitutions() const
Definition: HierarchicalLeastSquareSolver.h:159
virtual void updateBoundTargetData(scheme::internal::AssignmentTarget &target)=0
virtual void addBound_(LinearConstraintPtr bound, RangePtr range, bool first)=0
virtual void initializeBuild_(const std::vector< int > &nEq, const std::vector< int > &nIneq, bool useBounds)=0
int numberOfLevels() const
Definition: HierarchicalLeastSquareSolver.h:99
std::vector< std::unique_ptr< MarkedAssignment > > AssignmentVector
Definition: HierarchicalLeastSquareSolver.h:182
void startBuild(const VariableVector &x, const std::vector< int > &nEq, const std::vector< int > &nIneq, bool useBounds, const hint::internal::Substitutions *const subs=nullptr)
virtual Range nextEqualityConstraintRange_(int lvl, const constraint::abstract::LinearConstraint &cstr) const =0
virtual void postAssignmentProcess_()
Definition: HierarchicalLeastSquareSolver.h:136
void addConstraint(LinearConstraintPtr cstr, SolvingRequirementsPtr req)
int constraintSize(const constraint::abstract::LinearConstraint &c) const
virtual const Eigen::VectorXd & result_() const =0
std::vector< int > ineqSize_
Definition: HierarchicalLeastSquareSolver.h:191
virtual Range nextInequalityConstraintRange_(int lvl, const constraint::abstract::LinearConstraint &cstr) const =0
virtual void applyImpactLogic(ImpactFromChanges &impact)
virtual void preAssignmentProcess_()
Definition: HierarchicalLeastSquareSolver.h:135
virtual void removeBounds_(const Range &range)=0
map< constraint::abstract::LinearConstraint *, AssignmentPtrVector > MapToAssignment
Definition: HierarchicalLeastSquareSolver.h:184
HierarchicalLeastSquareSolver(const HierarchicalLeastSquareSolver &)=delete
virtual void updateInequalityTargetData(int lvl, scheme::internal::AssignmentTarget &target)=0
void addAssignement(Args &&... args)
Definition: HierarchicalLeastSquareSolver.h:230
utils::internal::map< K, T > map
Definition: HierarchicalLeastSquareSolver.h:181
const VariableVector & variables() const
Definition: HierarchicalLeastSquareSolver.h:158
std::vector< int > nEq_
Definition: HierarchicalLeastSquareSolver.h:188
Definition: SolverEvents.h:14
Definition: probe.h:44
bool eq(const ObjType &l, const ObjType &r, MemberType ObjType::*member, Args &&... args)
Definition: Log.h:39
std::map< KeyWithId, Value, IdLess< KeyWithId >, Allocator > map
Definition: map.h:41
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
ImpactFromChanges & orAssign(const ImpactFromChanges &other)
bool any() const
Definition: HierarchicalLeastSquareSolver.h:120
ImpactFromChanges(const std::vector< bool > &eq, std::vector< bool > &ineq, bool bounds)
ImpactFromChanges operator||(const ImpactFromChanges &other)
std::vector< bool > inequalityConstraints_
Definition: HierarchicalLeastSquareSolver.h:108
std::vector< bool > equalityConstraints_
Definition: HierarchicalLeastSquareSolver.h:107
static bool willReallocate(const Eigen::DenseBase< Derived > &M, int rows, int cols=1)
Definition: HierarchicalLeastSquareSolver.h:234
Definition: HierarchicalLeastSquareSolver.h:173
bool markedForRemoval
Definition: HierarchicalLeastSquareSolver.h:178
scheme::internal::Assignment assignment
Definition: HierarchicalLeastSquareSolver.h:177
MarkedAssignment(Args &&... args)
Definition: HierarchicalLeastSquareSolver.h:175