Fade2D Documentation pages v2.12
Delaunay Features
Fade_2D.h
Go to the documentation of this file.
1 // Copyright (C) Geom Software e.U, Bernhard Kornberger, Graz/Austria
2 //
3 // This file is part of the Fade2D library. The student license is free
4 // of charge and covers personal non-commercial research. Licensees
5 // holding a commercial license may use this file in accordance with
6 // the Commercial License Agreement.
7 //
8 // This software is provided AS IS with NO WARRANTY OF ANY KIND,
9 // INCLUDING THE WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS
10 // FOR A PARTICULAR PURPOSE.
11 //
12 // Please contact the author if any conditions of this licensing are
13 // not clear to you.
14 //
15 // Author: Bernhard Kornberger, bkorn (at) geom.at
16 // http://www.geom.at
17 
18 
24 #pragma once
25 
26 #include "common.h"
27 #include "freeFunctions.h"
28 #ifndef FADE2D_EXPORT
29  #include "License.h"
30 #endif
31 #include "Point2.h"
32 #include "Triangle2.h"
34 #include "Visualizer2.h"
35 #include "Zone2.h"
36 #include "ConstraintGraph2.h"
37 #include "Performance.h"
38 #include "MeshGenParams.h"
39 #include "MsgBase.h"
40 #include "SegmentChecker.h"
41 #include "testDataGenerators.h"
42 #include "FadeExport.h"
43 #include "Voronoi2.h"
44 
45 #if GEOM_PSEUDO3D==GEOM_TRUE
46  #include "IsoContours.h"
47  #include "EfficientModel.h"
48  #include "CutAndFill.h"
49  #include "CloudPrepare.h"
50 #endif
51 
52 
53 
54 #if GEOM_PSEUDO3D==GEOM_TRUE
55  namespace GEOM_FADE25D {
56 #elif GEOM_PSEUDO3D==GEOM_FALSE
57  namespace GEOM_FADE2D {
58 #else
59  #error GEOM_PSEUDO3D is not defined
60 #endif
61 
62 
63 
64 class Dt2; // Forward
65 class Visualizer3; // Forward
66 
67 
72 class CLASS_DECLSPEC Fade_2D
73 {
74 public:
82  explicit Fade_2D(unsigned numExpectedVertices=3);
83 
88  Fade_2D(Fade_2D& other,std::vector<Zone2*> vZoneIn=std::vector<Zone2*>(),std::vector<Zone2*> vZoneOut=std::vector<Zone2*>());
89 
92 
93 
111  bool saveTriangulation(const char* filename,std::vector<Zone2*>& vSaveZones);
112 
130  bool saveTriangulation(std::ostream& stream,std::vector<Zone2*>& vSaveZones);
131 
132 
155  bool saveZones(const char* filename,std::vector<Zone2*>& vSaveZones);
156 
157 
180  bool saveZones(std::ostream& stream,std::vector<Zone2*>& vSaveZones);
181 
182 
194  bool load(const char* filename,std::vector<Zone2*>& vZones);
195 
196 
208  bool load(std::istream& stream,std::vector<Zone2*>& vZones);
209 
210 
220  void exportTriangulation(FadeExport& fadeExport,bool bWithCustomIndices,bool bClear);
221 
222 
235  bool checkValidity(bool bCheckEmptyCircleProperty,const char* msg) const;
236 
243  int setNumCPU(int numCPU);
244 
251 
266  void setFastMode(bool bFast);
267 
268 
273  void statistics(const char* s) const;
274 
285  void show(const char* postscriptFilename,bool bWithConstraints=true) const;
286 
287 
300  void show(Visualizer2* pVis,bool bWithConstraints=true) const;
301 
302 #if GEOM_PSEUDO3D==GEOM_TRUE
310  void showGeomview(const char* filename, const char* color="1 1 1 0.5") const;
311 
319  void showGeomview(Visualizer3* pVis, const char* color="1 1 1 0.5") const;
320 
321 
322 #endif
323 
332  void remove(Point2* pVertex);
333 
341  void remove(std::vector<Point2*>& vPoints);
342 
343 
349  void reset();
350 
362  void getConvexHull(bool bAllVertices,std::vector<Point2*>& vConvexHullPointsOut);
363 
364 #if GEOM_PSEUDO3D==GEOM_TRUE
382  void insert(CloudPrepare* pCloudPrepare,bool bClear=true);
383 #endif
384 
398  Point2* insert(const Point2& p);
399 
408  void insert(const std::vector<Point2>& vInputPoints);
409 
424  void insert(const std::vector<Point2>& vInputPoints,std::vector<Point2*>& vHandles);
425 
426 #if GEOM_PSEUDO3D==GEOM_TRUE
435 #else
444 #endif
445  void insert(int numPoints,double * aCoordinates,Point2 ** aHandles);
446 
447 
464  double measureTriangulationTime(const std::vector<Point2>& vPoints);
465 
466 
467 
482  Triangle2* locate(const Point2& p);
483 
484 
494 
495 #if GEOM_PSEUDO3D==GEOM_TRUE
521  bool getHeight(double x,double y,double& heightOut,Triangle2* pApproxT=NULL,double tolerance=0) const;
522 #endif
523 
539  void refine( Zone2* pZone,
540  double minAngleDegree,
541  double minEdgeLength,
542  double maxEdgeLength,
543  bool bAllowConstraintSplitting
544  );
545 
581  void refineAdvanced(MeshGenParams* pParameters);
582 
583 
584 
585 
586 
587 
597  size_t numberOfPoints() const;
598 
605  size_t numberOfTriangles() const;
606 
607 
619  void getTrianglePointers(std::vector<Triangle2*>& vAllTriangles) const;
620 
631 void getVertexPointers(std::vector<Point2*>& vAllPoints) const;
632 
633 
640 
641 
652 
669  bool hasArea() const;
670 
676  bool is2D() const;
677 
678 
734 #if GEOM_PSEUDO3D==GEOM_TRUE
735  ConstraintGraph2* createConstraint(
736  std::vector<Segment2>& vSegments,
738  bool bOrientedSegments=false,
739  bool bUseHeightOfLatest=false
740  );
741 #else
743  std::vector<Segment2>& vSegments,
745  bool bOrientedSegments=false
746  );
747 #endif
748 
749 
750 
751 
752 
775  Zone2* createZone(ConstraintGraph2* pConstraintGraph,ZoneLocation zoneLoc,bool bVerbose=true);
776 
777 #if GEOM_PSEUDO3D==GEOM_TRUE
794  Zone2* getPeeledZone(double angleThreshold,bool bAvoidSplit);
795 #endif
796 
808  Zone2* createZone(const std::vector<ConstraintGraph2*>& vConstraintGraphs,ZoneLocation zoneLoc,const Point2& startPoint,bool bVerbose=true);
809 
822  Zone2* createZone(ConstraintGraph2* pConstraintGraph,ZoneLocation zoneLoc,const Point2& startPoint,bool bVerbose=true);
823 
834  Zone2* createZone( std::vector<Triangle2*>& vTriangles,bool bVerbose=true );
835 
836 
837 
850  void deleteZone(Zone2* pZone);
851 
852 
873 
874 
882 
890  bool isConstraint(Triangle2* pT,int ith) const;
891 
894  void getAliveConstraintSegments(std::vector<ConstraintSegment2*>& vAliveConstraintSegments) const;
895 
898  void getAliveAndDeadConstraintSegments(std::vector<ConstraintSegment2*>& vAllConstraintSegments) const;
899 
908 
909 
914  void getIncidentTriangles(Point2* pVtx,std::vector<Triangle2*>& vIncidentT) const;
915 
920  void getIncidentVertices(Point2* pVtx,std::vector<Point2*>& vIncidentVertices) const;
921 
922 
928  void writeObj(const char* filename) const;
929 
930 #ifndef SKIPTHREADS
931 #if GEOM_PSEUDO3D==GEOM_TRUE
940  bool writePly(const char* filename,bool bASCII=false) const;
948  bool writePly(std::ostream& os,bool bASCII=false) const;
949 #endif
950 #endif
951 
958  void writeObj(const char* filename,Zone2* pZone) const;
959 
960 
966  void writeWebScene(const char* path) const;
967 
973  void writeWebScene(const char* path,Zone2* pZone) const;
974 
980 void subscribe(MsgType msgType,MsgBase* pMsg);
981 
987 void unsubscribe(MsgType msgType,MsgBase* pMsg);
988 
989 
996  bool isConstraint(Point2* p0,Point2* p1) const;
997 
1003  bool isConstraint(Point2* pVtx) const;
1004 
1006  void printLicense() const;
1007 
1009  bool checkZoneQuality(Zone2* pZone,double minAngle,const char* name,const AcceptExperimentalFeature& accept);
1011  void setName(const char* s);
1013  std::string getName() const;
1014 
1015 
1036  Zone2* importTriangles( std::vector<Point2>& vPoints,
1037  bool bReorientIfNeeded,
1038  bool bCreateExtendedBoundingBox
1039  );
1040 
1057  Zone2* importTriangles_robust(std::vector<Point2>& vPoints);
1058 
1059 #ifndef SKIPTHREADS
1060 #if GEOM_PSEUDO3D==GEOM_TRUE
1073  Zone2* importTrianglesFromPly(const char* filename);
1074 
1087  Zone2* importTrianglesFromPly(std::istream& is);
1088 #endif
1089 #endif
1090 
1095  Orientation2 getOrientation(const Point2& p0,const Point2& p1,const Point2& p2);
1096 
1127  void cutTriangles( const Point2& knifeStart,
1128  const Point2& knifeEnd,
1129  bool bTurnEdgesIntoConstraints);
1130 
1143  void cutTriangles( std::vector<Segment2>& vSegments,
1144  bool bTurnEdgesIntoConstraints);
1145 
1147  Zone2* createZone_cookieCutter(std::vector<Segment2>& vSegments,bool bProtectEdges);
1148 
1149 
1172  Zone2* createZone_cookieCutter( std::vector<Segment2>& vSegments,
1173  bool bProtectEdges,
1174  ConstraintGraph2*& pProtectedEdgesCG,
1175  ConstraintGraph2*& pBoundaryCG);
1176 
1177 
1227 #if GEOM_PSEUDO3D==GEOM_TRUE
1228  bool drape( std::vector<Segment2>& vSegmentsIn,
1229  std::vector<Segment2>& vSegmentsOut,
1230  double zTolerance,
1231  double offTolerance=0.0) const;
1232 #else
1233  bool drape( std::vector<Segment2>& vSegmentsIn,
1234  std::vector<Segment2>& vSegmentsOut) const;
1235 #endif
1236 
1237 
1238 
1239 
1245 void getDirectedEdges(std::vector<Edge2>& vDirectedEdgesOut) const;
1246 
1252 void getUndirectedEdges(std::vector<Edge2>& vUndirectedEdgesOut) const;
1253 
1255  /* Deprecated: Use setNumCPU() instead. This method is kept for
1256  * compatibility with existing applications. Internally it calls
1257  * setNumCPU(0) to automatically determine and use the number of
1258  * available CPU cores.
1259  */
1260  void enableMultithreading();
1261  // Development functions, not for public use
1263  void internal(int au,int fu,const char* s="");
1265  Dt2* getImpl();
1267  void setDev(const char* s,int ival,double dval);
1268 
1269 protected:
1271  void initFade(unsigned numExpectedVertices);
1273  Fade_2D& operator=(const Fade_2D&); // No assignment allowed
1275  Dt2* pImpl;
1276 }; // end of class
1277 
1278 
1286 inline void copyFade(Fade_2D& dt0,std::vector<Zone2*>& vZones0,Fade_2D& dt1,std::vector<Zone2*>& vZones1)
1287 {
1288  std::stringstream ss;
1289  dt0.saveTriangulation(ss,vZones0);
1290  dt1.load(ss,vZones1);
1291 }
1292 
1293 
1294 
1295 
1296 } // (namespace)
1297 
ConstraintInsertionStrategy
Constraint Insertion Strategy determines how a constraint edge shall be inserted:
Definition: ConstraintSegment2.h:52
void copyFade(Fade_2D &dt0, std::vector< Zone2 * > &vZones0, Fade_2D &dt1, std::vector< Zone2 * > &vZones1)
Copy a Fade_2D object and selected Zone2 objects.
Definition: Fade_2D.h:1286
Bbox2 is an axis aligned 2D bounding box.
Definition: Bbox2.h:37
ConstraintGraph2 is a set of Constraint Edges (ConstraintSegment2)
Definition: ConstraintGraph2.h:52
A ConstraintSegment2 represents a Constraint Edge.
Definition: ConstraintSegment2.h:69
Fade_2D is the Delaunay triangulation class.
Definition: Fade_2D.h:73
Triangle2 * getAdjacentTriangle(Point2 *p0, Point2 *p1) const
Get adjacent triangle.
void insert(const std::vector< Point2 > &vInputPoints, std::vector< Point2 * > &vHandles)
Insert points from a std::vector and store pointers in vHandles.
void unsubscribe(MsgType msgType, MsgBase *pMsg)
Unregister a message receiver.
int setNumCPU(int numCPU)
Set the number CPU cores for multithreading (deprecated)
void setFastMode(bool bFast)
Set fast mode.
bool drape(std::vector< Segment2 > &vSegmentsIn, std::vector< Segment2 > &vSegmentsOut) const
Drape segments along a surface.
void getConvexHull(bool bAllVertices, std::vector< Point2 * > &vConvexHullPointsOut)
Compute the convex hull.
void subscribe(MsgType msgType, MsgBase *pMsg)
Register a message receiver.
bool hasArea() const
Check if the triangulation contains triangles (which is the case if at least 3 non-collinear points e...
Zone2 * importTriangles_robust(std::vector< Point2 > &vPoints)
Import Triangles - Robust (Error-tolerant version)
Orientation2 getOrientation(const Point2 &p0, const Point2 &p1, const Point2 &p2)
Compute the orientation of 3 points.
void show(Visualizer2 *pVis, bool bWithConstraints=true) const
Draws the triangulation as postscript file using an existing Visualizer2 object.
void getUndirectedEdges(std::vector< Edge2 > &vUndirectedEdgesOut) const
Get undirected edges.
ConstraintGraph2 * createConstraint(std::vector< Segment2 > &vSegments, ConstraintInsertionStrategy cis, bool bOrientedSegments=false)
Add constraint edges (edges, polyline, polygon) .
void insert(const std::vector< Point2 > &vInputPoints)
Insert a vector of points.
double measureTriangulationTime(const std::vector< Point2 > &vPoints)
Measure the Delaunay triangulation time.
Zone2 * createZone(ConstraintGraph2 *pConstraintGraph, ZoneLocation zoneLoc, const Point2 &startPoint, bool bVerbose=true)
Create a zone limited by a ConstraintGraph by growing from a start point.
void writeWebScene(const char *path, Zone2 *pZone) const
Write a zone to an *.obj file.
void writeObj(const char *filename, Zone2 *pZone) const
Write a zone to an *.obj file.
bool load(std::istream &stream, std::vector< Zone2 * > &vZones)
Load a triangulation.
void getVertexPointers(std::vector< Point2 * > &vAllPoints) const
Get pointers to all vertices.
void writeObj(const char *filename) const
Write the current triangulation to an *.obj file.
Bbox2 computeBoundingBox() const
Compute the axis-aligned bounding box of the points.
Zone2 * createZone(const std::vector< ConstraintGraph2 * > &vConstraintGraphs, ZoneLocation zoneLoc, const Point2 &startPoint, bool bVerbose=true)
Create a zone limited by multiple ConstraintGraph2 objects by growing from a start point.
void getTrianglePointers(std::vector< Triangle2 * > &vAllTriangles) const
Get pointers to all triangles.
Fade_2D(Fade_2D &other, std::vector< Zone2 * > vZoneIn=std::vector< Zone2 * >(), std::vector< Zone2 * > vZoneOut=std::vector< Zone2 * >())
Copy constructor, copy also zones.
Point2 * insert(const Point2 &p)
Insert a single point.
void printLicense() const
Prints license information.
Fade_2D(unsigned numExpectedVertices=3)
Constructor of the triangulation class.
ConstraintSegment2 * getConstraintSegment(Point2 *p0, Point2 *p1) const
Retrieve a ConstraintSegment2.
void remove(Point2 *pVertex)
Remove a single vertex.
void getIncidentTriangles(Point2 *pVtx, std::vector< Triangle2 * > &vIncidentT) const
Get incident triangles.
bool saveTriangulation(const char *filename, std::vector< Zone2 * > &vSaveZones)
Save a triangulation.
Zone2 * createZone(ConstraintGraph2 *pConstraintGraph, ZoneLocation zoneLoc, bool bVerbose=true)
Create a zone.
Voronoi2 * getVoronoiDiagram()
Get the Voronoi diagram.
Point2 * getNearestNeighbor(const Point2 &p)
Get nearest neighbor.
bool saveZones(std::ostream &stream, std::vector< Zone2 * > &vSaveZones)
Save zones.
bool saveTriangulation(std::ostream &stream, std::vector< Zone2 * > &vSaveZones)
Save a triangulation.
void applyConstraintsAndZones()
Apply conforming constraints and zones (deprecated!)
void reset()
Reset.
size_t numberOfPoints() const
Number of points.
~Fade_2D()
Destructor.
Zone2 * createZone_cookieCutter(std::vector< Segment2 > &vSegments, bool bProtectEdges, ConstraintGraph2 *&pProtectedEdgesCG, ConstraintGraph2 *&pBoundaryCG)
Cookie Cutter The Cookie Cutter cuts out a part of a triangulation and returns it as a Zone2 object.
void getAliveAndDeadConstraintSegments(std::vector< ConstraintSegment2 * > &vAllConstraintSegments) const
Get all (alive and dead) constraint segments.
void getDirectedEdges(std::vector< Edge2 > &vDirectedEdgesOut) const
Get directed edges Edges are counterclockwise oriented around their triangle. The present method retu...
void refineAdvanced(MeshGenParams *pParameters)
Delaunay refinement and grid meshing.
bool isConstraint(Point2 *pVtx) const
Check if a vertex is a constraint vertex.
bool load(const char *filename, std::vector< Zone2 * > &vZones)
Load a triangulation.
void cutTriangles(const Point2 &knifeStart, const Point2 &knifeEnd, bool bTurnEdgesIntoConstraints)
Cut through a triangulation.
void exportTriangulation(FadeExport &fadeExport, bool bWithCustomIndices, bool bClear)
Export triangulation data from Fade.
void statistics(const char *s) const
Statistics.
Point2 * getVertexPointer(const Point2 &p)
Get pointer.
void getAliveConstraintSegments(std::vector< ConstraintSegment2 * > &vAliveConstraintSegments) const
Get active (alive) constraint segments.
void refine(Zone2 *pZone, double minAngleDegree, double minEdgeLength, double maxEdgeLength, bool bAllowConstraintSplitting)
Delaunay refinement.
void cutTriangles(std::vector< Segment2 > &vSegments, bool bTurnEdgesIntoConstraints)
Cut through a triangulation.
bool checkValidity(bool bCheckEmptyCircleProperty, const char *msg) const
Checks if a triangulation is valid.
bool saveZones(const char *filename, std::vector< Zone2 * > &vSaveZones)
Save zones.
bool isConstraint(Triangle2 *pT, int ith) const
Check if an edge is a constraint edge.
bool isConstraint(Point2 *p0, Point2 *p1) const
Check if an edge is a constraint edge.
void writeWebScene(const char *path) const
Write the current triangulation to an *.obj file.
void getIncidentVertices(Point2 *pVtx, std::vector< Point2 * > &vIncidentVertices) const
Get incident vertices.
void remove(std::vector< Point2 * > &vPoints)
Remove vertices.
void insert(int numPoints, double *aCoordinates, Point2 **aHandles)
Insert points from an array.
void show(const char *postscriptFilename, bool bWithConstraints=true) const
Draws the triangulation as postscript file.
Zone2 * importTriangles(std::vector< Point2 > &vPoints, bool bReorientIfNeeded, bool bCreateExtendedBoundingBox)
Import Triangles.
size_t numberOfTriangles() const
Number of triangles.
void deleteZone(Zone2 *pZone)
Delete a Zone2 object.
Zone2 * createZone(std::vector< Triangle2 * > &vTriangles, bool bVerbose=true)
Create a zone defined by a vector of triangles.
Triangle2 * locate(const Point2 &p)
Locate a triangle which contains p.
Parameters for the mesh generator.
Definition: MeshGenParams.h:59
MsgBase, a base class for message subscriber classes.
Definition: MsgBase.h:47
Point.
Definition: Point2.h:52
Triangle.
Definition: Triangle2.h:60
Visualizer2 is a general Postscript writer. It draws the objects Point2, Segment2,...
Definition: Visualizer2.h:57
Voronoi diagram.
Definition: Voronoi2.h:56
Zone2 is a certain defined area of a triangulation.
Definition: Zone2.h:95
FadeExport is a simple struct to export triangulation data.
Definition: FadeExport.h:43