irafhy
Interval arithmetic based Reachability Analysis Framework for Hybrid Automaton
polytope.h
Go to the documentation of this file.
1 #ifndef REPRESENTATION_GEOMETRIC_POLYTOPE_H
2 #define REPRESENTATION_GEOMETRIC_POLYTOPE_H
3 
4 #include <vector>
5 #include <iostream>
6 #include <capd/capdlib.h>
11 
12 namespace irafhy
13 {
14  class Polytope : public Geometry<Polytope>
15  {
16  private:
20  std::vector<Point> pointConstraints_;
24  std::vector<HalfSpace> halfSpaceConstraints_;
28  double volume_ = 0.0;
32  int dimension_ = 0;
33 
34  private:
42  Polytope(const std::vector<Point>& pointConstraints,
43  const std::vector<HalfSpace>& halfSpaceConstraints,
44  double volume,
45  int dimension);
50  void init(const convexConstraints& constraints);
56  std::vector<Point> verticesEnumerationDualPly(const std::vector<HalfSpace>& halfSpaces);
62  std::vector<Point> verticesEnumerationIntersect(const std::vector<HalfSpace>& halfSpaces);
68  std::vector<Point> verticesEnumeration(const std::vector<capd::interval>& constraints);
75  std::vector<Point> verticesEnumeration(const Eigen::VectorXd& center, double radius);
76 
77  public:
83  static Polytope Empty(std::size_t dimension = 0);
84 
85  public:
89  Polytope();
94  Polytope(const Polytope& polytope) = default;
99  Polytope(Polytope&& polytope) noexcept = default;
104  explicit Polytope(const std::vector<Point>& points);
109  explicit Polytope(const std::vector<HalfSpace>& halfSpaces);
114  explicit Polytope(const std::vector<capd::interval>& constraints);
119  explicit Polytope(const capd::C0Rect2Set& set);
125  Polytope(const Point& center, double radius);
131  Polytope(const Eigen::VectorXd& center, double radius);
135  ~Polytope() override = default;
140  [[nodiscard]] std::vector<Point> pointConstraints() const;
145  [[nodiscard]] std::vector<HalfSpace> halfSpaceConstraints() const;
150  [[nodiscard]] Point centroid() const;
155  [[nodiscard]] std::vector<capd::interval> constraints() const;
160  [[nodiscard]] bool empty() const override;
165  [[nodiscard]] int dimension() const override;
170  [[nodiscard]] double volume() const;
177  bool intersect(const Polytope& rhs, Polytope& result) const override;
183  [[nodiscard]] Polytope unite(const Polytope& rhs) const override;
189  [[nodiscard]] bool contains(const Point& point) const override;
195  [[nodiscard]] bool contains(const Eigen::VectorXd& coordinate) const override;
201  Polytope& operator=(const Polytope& rhs) = default;
207  Polytope& operator=(Polytope&& rhs) noexcept = default;
208  };
209 
216  std::ostream& operator<<(std::ostream& out, const Polytope& rhs);
217 } // namespace irafhy
218 #endif //REPRESENTATION_GEOMETRIC_POLYTOPE_H
std::ostream & operator<<(std::ostream &out, const Condition &rhs)
output the right hand side condition to the standard out stream
Definition: condition.cpp:246
bool empty() const override
check if the polytope is empty or not
Definition: polytope.cpp:285
static Polytope Empty(std::size_t dimension=0)
static constructor of Empty polytope
Definition: polytope.cpp:166
Definition: convexHull.hpp:21
std::vector< HalfSpace > halfSpaceConstraints() const
get the half space constraints of the polytope
Definition: polytope.cpp:254
std::vector< capd::interval > constraints() const
get the variance constraints of all dimensions
Definition: polytope.cpp:265
Definition: point.h:10
int dimension_
dimension of the space which current convex hull in
Definition: polytope.h:32
~Polytope() override=default
destructor
void init(const convexConstraints &constraints)
initialize the polytope with given convex hull
Definition: polytope.cpp:17
std::vector< Point > verticesEnumerationIntersect(const std::vector< HalfSpace > &halfSpaces)
enumerate vertices using intersection computation
Definition: polytope.cpp:78
std::vector< Point > verticesEnumeration(const std::vector< capd::interval > &constraints)
enumerate extreme vertices with given constraints in interval form
Definition: polytope.cpp:144
Definition: condition.cpp:3
double volume_
volume of the convex hull
Definition: polytope.h:28
Polytope unite(const Polytope &rhs) const override
get the union of the current polytope and the given right hand side one
Definition: polytope.cpp:386
std::vector< Point > pointConstraints() const
get the extreme vertices of the polytope
Definition: polytope.cpp:252
int dimension() const override
get the dimension of the space which the polytope in
Definition: polytope.cpp:290
bool intersect(const Polytope &rhs, Polytope &result) const override
check if the current polytope intersect with the given right hand side one or not ...
Definition: polytope.cpp:302
Definition: geometry.h:11
double volume() const
get the volume of the polytope
Definition: polytope.cpp:296
bool contains(const Point &point) const override
check if the given point inside the domain which defined by the current polytope or not ...
Definition: polytope.cpp:396
Point centroid() const
get the centroid of the polytope
Definition: polytope.cpp:256
Polytope()
constructor
Definition: polytope.cpp:158
std::vector< Point > pointConstraints_
point constraints which define the convex hull entity
Definition: polytope.h:20
Definition: polytope.h:14
std::vector< HalfSpace > halfSpaceConstraints_
half spaces constraints which define the convex hull entity
Definition: polytope.h:24
Polytope & operator=(const Polytope &rhs)=default
assignment operator
std::vector< Point > verticesEnumerationDualPly(const std::vector< HalfSpace > &halfSpaces)
enumerate vertices using dual polytope method
Definition: polytope.cpp:25