irafhy
Interval arithmetic based Reachability Analysis Framework for Hybrid Automaton
intervalHull.h
Go to the documentation of this file.
1 #ifndef REPRESENTATION_GEOMETRIC_INTERVAL_HULL_H
2 #define REPRESENTATION_GEOMETRIC_INTERVAL_HULL_H
3 
4 #include <capd/capdlib.h>
5 #include <vector>
6 #include <iostream>
9 
10 namespace irafhy
11 {
12  class IntervalHull : public Geometry<IntervalHull>
13  {
14  private:
18  std::vector<capd::interval> constraints_;
19 
20  public:
26  static IntervalHull Empty(std::size_t dimension = 0);
27 
28  public:
32  IntervalHull();
37  IntervalHull(const IntervalHull& intervalHull) = default;
42  IntervalHull(IntervalHull&& intervalHull) noexcept = default;
47  explicit IntervalHull(const std::vector<capd::interval>& constraints);
52  explicit IntervalHull(const std::vector<Point>& points);
58  IntervalHull(const Point& center, double radius);
63  explicit IntervalHull(const capd::C0Rect2Set& set);
69  IntervalHull(const Eigen::VectorXd& center, double radius);
73  ~IntervalHull() override = default;
78  [[nodiscard]] std::vector<capd::interval> constraints() const;
83  [[nodiscard]] bool empty() const override;
88  [[nodiscard]] int dimension() const override;
93  [[nodiscard]] std::vector<Point> extremeVertices() const;
98  [[nodiscard]] Point randInnerPoint() const;
103  [[nodiscard]] Point centroid() const;
110  bool intersect(const IntervalHull& rhs, IntervalHull& result) const override;
116  [[nodiscard]] IntervalHull unite(const IntervalHull& rhs) const override;
122  [[nodiscard]] bool contains(const Point& point) const override;
128  [[nodiscard]] bool contains(const Eigen::VectorXd& coordinate) const override;
134  [[nodiscard]] bool contains(const IntervalHull& rhs) const;
140  capd::interval& operator[](std::size_t index);
146  const capd::interval& operator[](std::size_t index) const;
152  bool operator==(const IntervalHull& rhs) const;
158  bool operator!=(const IntervalHull& rhs) const;
164  bool operator<(const IntervalHull& rhs) const;
170  bool operator<=(const IntervalHull& rhs) const;
176  bool operator>(const IntervalHull& rhs) const;
182  bool operator>=(const IntervalHull& rhs) const;
188  IntervalHull& operator=(const IntervalHull& rhs) = default;
194  IntervalHull& operator=(IntervalHull&& rhs) = default;
195  };
196 
203  std::ostream& operator<<(std::ostream& out, const IntervalHull& rhs);
204 } // namespace irafhy
205 #endif //REPRESENTATION_GEOMETRIC_INTERVAL_HULL_H
IntervalHull()
default constructor
Definition: intervalHull.cpp:7
bool operator>=(const IntervalHull &rhs) const
relational operator
Definition: intervalHull.cpp:248
bool operator<(const IntervalHull &rhs) const
relational operator
Definition: intervalHull.cpp:226
std::ostream & operator<<(std::ostream &out, const Condition &rhs)
output the right hand side condition to the standard out stream
Definition: condition.cpp:246
IntervalHull & operator=(const IntervalHull &rhs)=default
assignment operator
std::vector< Point > extremeVertices() const
get the extreme vertices of the interval hull
Definition: intervalHull.cpp:99
Definition: intervalHull.h:12
std::vector< capd::interval > constraints_
interval constraints of each dimension related entries
Definition: intervalHull.h:18
bool operator>(const IntervalHull &rhs) const
relational operator
Definition: intervalHull.cpp:246
IntervalHull unite(const IntervalHull &rhs) const override
get the union of the current interval hull and right hand side interval hull
Definition: intervalHull.cpp:162
Definition: point.h:10
bool operator!=(const IntervalHull &rhs) const
relational operator
Definition: intervalHull.cpp:224
bool empty() const override
check if the interval hull is empty or not
Definition: intervalHull.cpp:80
~IntervalHull() override=default
destructor
bool operator<=(const IntervalHull &rhs) const
relational operator
Definition: intervalHull.cpp:236
Definition: condition.cpp:3
std::vector< capd::interval > constraints() const
get the constraints of the interval hull
Definition: intervalHull.cpp:78
bool operator==(const IntervalHull &rhs) const
relational operator
Definition: intervalHull.cpp:217
int dimension() const override
get the dimension of the space which the interval hull in
Definition: intervalHull.cpp:97
Definition: geometry.h:11
capd::interval & operator[](std::size_t index)
get the constraint of specified dimension
Definition: intervalHull.cpp:205
static IntervalHull Empty(std::size_t dimension=0)
static constructor of Empty interval hull in specified space
Definition: intervalHull.cpp:5
bool contains(const Point &point) const override
check if the given point inside the domain defined by the interval hull or not
Definition: intervalHull.cpp:175
bool intersect(const IntervalHull &rhs, IntervalHull &result) const override
check if the current interval hull intersects with the right hand side interval hull or not ...
Definition: intervalHull.cpp:140
Point randInnerPoint() const
randomly sample a point from the domain defined by the interval hull
Definition: intervalHull.cpp:124
Point centroid() const
get the centroid of the interval hull
Definition: intervalHull.cpp:130