OGS
|
#include <Polygon.h>
Public Member Functions | |
Polygon (const Polyline &ply, bool init=true) | |
Polygon (Polygon const &other) | |
Polygon & | operator= (Polygon const &rhs)=delete |
~Polygon () override | |
bool | initialise () |
bool | isPntInPolygon (const MathLib::Point3d &pnt) const |
bool | containsSegment (GeoLib::LineSegment const &segment) const |
bool | isPolylineInPolygon (const Polyline &ply) const |
bool | isPartOfPolylineInPolygon (const Polyline &ply) const |
bool | getNextIntersectionPointPolygonLine (GeoLib::LineSegment const &seg, GeoLib::Point &intersection_pnt, std::size_t &seg_num) const |
std::list< Polygon * > const & | computeListOfSimplePolygons () |
Public Member Functions inherited from GeoLib::Polyline | |
Polyline (const std::vector< Point * > &pnt_vec) | |
Polyline (const Polyline &ply) | |
Polyline & | operator= (Polyline const &other)=delete |
~Polyline () override=default | |
GEOTYPE | getGeoType () const override |
return a geometry type More... | |
virtual bool | addPoint (std::size_t pnt_id) |
virtual bool | insertPoint (std::size_t pos, std::size_t pnt_id) |
virtual void | removePoint (std::size_t pos) |
void | closePolyline () |
std::size_t | getNumberOfPoints () const |
std::size_t | getNumberOfSegments () const |
bool | isClosed () const |
bool | isCoplanar () const |
bool | isPointIDInPolyline (std::size_t pnt_id) const |
std::size_t | getPointID (std::size_t i) const |
void | setPointID (std::size_t idx, std::size_t id) |
const Point * | getPoint (std::size_t i) const |
returns the i-th point contained in the polyline More... | |
SegmentIterator | begin () const |
SegmentIterator | end () const |
std::vector< Point * > const & | getPointsVec () const |
double | getDistanceAlongPolyline (const MathLib::Point3d &pnt, const double epsilon_radius) const |
Public Member Functions inherited from GeoLib::GeoObject | |
virtual | ~GeoObject ()=default |
Private Member Functions | |
std::vector< GeoLib::Point > | getAllIntersectionPoints (GeoLib::LineSegment const &segment) const |
EdgeType | getEdgeType (std::size_t k, MathLib::Point3d const &pnt) const |
void | ensureCCWOrientation () |
void | splitPolygonAtIntersection (const std::list< Polygon * >::const_iterator &polygon_it) |
void | splitPolygonAtPoint (const std::list< Polygon * >::iterator &polygon_it) |
Private Attributes | |
std::list< Polygon * > | _simple_polygon_list |
AABB | _aabb |
Friends | |
bool | operator== (Polygon const &lhs, Polygon const &rhs) |
Additional Inherited Members | |
Static Public Member Functions inherited from GeoLib::Polyline | |
static Polyline * | constructPolylineFromSegments (const std::vector< Polyline * > &ply_vec, double prox=0.0) |
Protected Member Functions inherited from GeoLib::Polyline | |
Location | getLocationOfPoint (std::size_t k, MathLib::Point3d const &pnt) const |
Protected Attributes inherited from GeoLib::Polyline | |
const std::vector< Point * > & | _ply_pnts |
std::vector< std::size_t > | _ply_pnt_ids |
|
explicit |
constructor checks if the given polyline is closed, and assures that the orientation is clock wise.
ply | closed Polyline |
init | if true, check if polyline is closed, calculate bounding box |
Definition at line 24 of file Polygon.cpp.
References _simple_polygon_list, and initialise().
GeoLib::Polygon::Polygon | ( | Polygon const & | other | ) |
Definition at line 34 of file Polygon.cpp.
References _simple_polygon_list, and GeoLib::Polyline::Polygon.
|
override |
std::list< Polygon * > const & GeoLib::Polygon::computeListOfSimplePolygons | ( | ) |
Subdivides a self-intersecting polygon into a list of non-intersecting shapes.
Definition at line 618 of file Polygon.cpp.
References _simple_polygon_list, splitPolygonAtIntersection(), and splitPolygonAtPoint().
bool GeoLib::Polygon::containsSegment | ( | GeoLib::LineSegment const & | segment | ) | const |
Checks if the straight line segment is contained within the polygon.
segment | the straight line segment that is checked with |
Definition at line 144 of file Polygon.cpp.
References getAllIntersectionPoints(), GeoLib::LineSegment::getBeginPoint(), GeoLib::LineSegment::getEndPoint(), isPntInPolygon(), and MathLib::sqrDist().
|
private |
Definition at line 317 of file Polygon.cpp.
References GeoLib::Polyline::_ply_pnt_ids, GeoLib::CCW, GeoLib::Polyline::getNumberOfPoints(), GeoLib::getOrientation(), GeoLib::Polyline::getPoint(), and GeoLib::rotatePointsToXY().
Referenced by initialise().
|
private |
Computes all intersections of the straight line segment and the polyline boundary
segment | the line segment that will be processed |
Definition at line 128 of file Polygon.cpp.
References GeoLib::lineSegmentIntersect().
Referenced by containsSegment().
|
private |
from book: Computational Geometry and Computer Graphics in C++, page 119 get the type of edge with respect to the given point (2d method!)
k | number of line segment |
pnt | point that is edge type computed for |
Definition at line 282 of file Polygon.cpp.
References GeoLib::BETWEEN, GeoLib::CROSSING, GeoLib::DESTINATION, GeoLib::Polyline::getLocationOfPoint(), GeoLib::Polyline::getPoint(), GeoLib::INESSENTIAL, GeoLib::LEFT, GeoLib::RIGHT, GeoLib::SOURCE, and GeoLib::TOUCHING.
Referenced by isPntInPolygon().
bool GeoLib::Polygon::getNextIntersectionPointPolygonLine | ( | GeoLib::LineSegment const & | seg, |
GeoLib::Point & | intersection_pnt, | ||
std::size_t & | seg_num | ||
) | const |
Calculates the next intersection point between the line segment seg
and the polygon starting with segment seg_num
.
seg | (input) Line segment to compute intersection. |
intersection_pnt | (output) next intersection point |
seg_num | (in/out) the number of the polygon segment that is intersecting |
intersection_pnt
and seg_num
contains new valid values Definition at line 248 of file Polygon.cpp.
References _simple_polygon_list, GeoLib::Polyline::begin(), GeoLib::Polyline::end(), and GeoLib::lineSegmentIntersect().
bool GeoLib::Polygon::initialise | ( | ) |
Definition at line 60 of file Polygon.cpp.
References ensureCCWOrientation(), GeoLib::Polyline::isClosed(), and WARN().
Referenced by Polygon().
bool GeoLib::Polygon::isPartOfPolylineInPolygon | ( | const Polyline & | ply | ) | const |
Method checks first if at least one (end!) point of a line segment of the polyline is inside of the polygon. If this test fails each line segment of the polyline will be tested against each polygon segment for intersection.
ply | the polyline that should be checked |
Definition at line 222 of file Polygon.cpp.
References GeoLib::Polyline::begin(), GeoLib::Polyline::end(), GeoLib::Polyline::getNumberOfPoints(), GeoLib::Polyline::getPoint(), and isPntInPolygon().
bool GeoLib::Polygon::isPntInPolygon | ( | const MathLib::Point3d & | pnt | ) | const |
Method checks if the given point is inside the polygon. The method requires that the polygon has clock wise orientation.
pnt | the Point |
Definition at line 71 of file Polygon.cpp.
References _aabb, _simple_polygon_list, GeoLib::CROSSING, getEdgeType(), GeoLib::AABB::getMaxPoint(), GeoLib::AABB::getMinPoint(), GeoLib::Polyline::getNumberOfPoints(), GeoLib::Polyline::getPoint(), GeoLib::INESSENTIAL, and GeoLib::TOUCHING.
Referenced by containsSegment(), FileIO::GMSH::GMSHPolygonTree::insertPolyline(), isPartOfPolylineInPolygon(), MeshGeoToolsLib::markNodesOutSideOfPolygon(), and FileIO::FEFLOWMeshInterface::setMaterialIDs().
bool GeoLib::Polygon::isPolylineInPolygon | ( | const Polyline & | ply | ) | const |
Method checks if all points of the polyline ply are inside of the polygon.
ply | the polyline that should be checked |
Definition at line 215 of file Polygon.cpp.
References GeoLib::Polyline::begin(), and GeoLib::Polyline::end().
Referenced by GeoLib::SimplePolygonTree::isPolygonInside().
|
private |
Definition at line 398 of file Polygon.cpp.
References GeoLib::Polyline::_ply_pnts, _simple_polygon_list, GeoLib::Polyline::getPointsVec(), GeoLib::Polyline::SegmentIterator::getSegmentNumber(), and GeoLib::lineSegmentsIntersect().
Referenced by computeListOfSimplePolygons().
|
private |
Definition at line 457 of file Polygon.cpp.
References _simple_polygon_list, GeoLib::Polyline::addPoint(), GeoLib::Polyline::Polygon, and BaseLib::quicksort().
Referenced by computeListOfSimplePolygons().
comparison operator for polygons
lhs | the first polygon |
rhs | the second polygon |
Definition at line 521 of file Polygon.cpp.
|
private |
Definition at line 133 of file Polygon.h.
Referenced by isPntInPolygon().
|
private |
Definition at line 132 of file Polygon.h.
Referenced by Polygon(), ~Polygon(), computeListOfSimplePolygons(), getNextIntersectionPointPolygonLine(), isPntInPolygon(), splitPolygonAtIntersection(), and splitPolygonAtPoint().