Fade2.5D Documentation pages v2.14
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 "VtkWriter.h"
36 #include "Zone2.h"
37 #include "ConstraintGraph2.h"
38 #include "Performance.h"
39 #include "MeshGenParams.h"
40 #include "MsgBase.h"
41 #include "SegmentChecker.h"
42 #include "testDataGenerators.h"
43 #include "FadeExport.h"
44 #include "Voronoi2.h"
45 
46 #if GEOM_PSEUDO3D==GEOM_TRUE
47  #include "IsoContours.h"
48  #include "EfficientModel.h"
49  #include "CutAndFill.h"
50  #include "CloudPrepare.h"
51 #endif
52 
53 
54 
55 #if GEOM_PSEUDO3D==GEOM_TRUE
56  namespace GEOM_FADE25D {
57 #elif GEOM_PSEUDO3D==GEOM_FALSE
58  namespace GEOM_FADE2D {
59 #else
60  #error GEOM_PSEUDO3D is not defined
61 #endif
62 
63 
64 
65 class Dt2; // Forward
66 class Visualizer3; // Forward
67 
68 
73 class CLASS_DECLSPEC Fade_2D
74 {
75 public:
83  explicit Fade_2D(unsigned numExpectedVertices=3);
84 
89  Fade_2D(Fade_2D& other,std::vector<Zone2*> vZoneIn=std::vector<Zone2*>(),std::vector<Zone2*> vZoneOut=std::vector<Zone2*>());
90 
93 
94 
112  bool saveTriangulation(const char* filename,std::vector<Zone2*>& vSaveZones);
113 
131  bool saveTriangulation(std::ostream& stream,std::vector<Zone2*>& vSaveZones);
132 
133 
156  bool saveZones(const char* filename,std::vector<Zone2*>& vSaveZones);
157 
158 
181  bool saveZones(std::ostream& stream,std::vector<Zone2*>& vSaveZones);
182 
183 
195  bool load(const char* filename,std::vector<Zone2*>& vZones);
196 
197 
209  bool load(std::istream& stream,std::vector<Zone2*>& vZones);
210 
211 
221  void exportTriangulation(FadeExport& fadeExport,bool bWithCustomIndices,bool bClear);
222 
223 
236  bool checkValidity(bool bCheckEmptyCircleProperty,const char* msg) const;
237 
244  int setNumCPU(int numCPU);
245 
251  Voronoi2* getVoronoiDiagram();
252 
267  void setFastMode(bool bFast);
268 
269 
274  void statistics(const char* s) const;
275 
286  void show(const char* filename,bool bWithConstraints=true) const;
287 
288 
301  void show(Visualizer2* pVis,bool bWithConstraints=true) const;
302 
314  void showVtk(const char* filename,VtkColor color,VtkColor constraintColor=VTK_TRANSPARENT) const;
315 
316 
317 
326  void showVtk(VtkWriter* pVtkWriter,VtkColor color,VtkColor constraintColor=VTK_TRANSPARENT) const;
327 
328 
329 #if GEOM_PSEUDO3D==GEOM_TRUE
337  void showGeomview(const char* filename, const char* color="1 1 1 0.5") const;
338 
346  void showGeomview(Visualizer3* pVis, const char* color="1 1 1 0.5") const;
347 #endif
348 
357  void remove(Point2* pVertex);
358 
366  void remove(std::vector<Point2*>& vPoints);
367 
368 
374  void reset();
375 
387  void getConvexHull(bool bAllVertices,std::vector<Point2*>& vConvexHullPointsOut);
388 
389 #if GEOM_PSEUDO3D==GEOM_TRUE
407  void insert(CloudPrepare* pCloudPrepare,bool bClear=true);
408 #endif
409 
423  Point2* insert(const Point2& p);
424 
433  void insert(const std::vector<Point2>& vInputPoints);
434 
449  void insert(const std::vector<Point2>& vInputPoints,std::vector<Point2*>& vHandles);
450 
451 #if GEOM_PSEUDO3D==GEOM_TRUE
460 #else
469 #endif
470  void insert(int numPoints,double * aCoordinates,Point2 ** aHandles);
471 
472 
489  double measureTriangulationTime(const std::vector<Point2>& vPoints);
490 
491 
492 
507  Triangle2* locate(const Point2& p);
508 
509 
519 
520 #if GEOM_PSEUDO3D==GEOM_TRUE
546  bool getHeight(double x,double y,double& heightOut,Triangle2* pApproxT=NULL,double tolerance=0) const;
547 #endif
548 
564  void refine( Zone2* pZone,
565  double minAngleDegree,
566  double minEdgeLength,
567  double maxEdgeLength,
568  bool bAllowConstraintSplitting
569  );
570 
606  void refineAdvanced(MeshGenParams* pParameters);
607 
608 
609 
610 
611 
612 
622  size_t numberOfPoints() const;
623 
630  size_t numberOfTriangles() const;
631 
632 
644  void getTrianglePointers(std::vector<Triangle2*>& vAllTriangles) const;
645 
656 void getVertexPointers(std::vector<Point2*>& vAllPoints) const;
657 
658 
665 
666 
677 
694  bool hasArea() const;
695 
701  bool is2D() const;
702 
703 
759 #if GEOM_PSEUDO3D==GEOM_TRUE
761  std::vector<Segment2>& vSegments,
763  bool bOrientedSegments=false,
764  bool bUseHeightOfLatest=false
765  );
766 #else
767  ConstraintGraph2* createConstraint(
768  std::vector<Segment2>& vSegments,
770  bool bOrientedSegments=false
771  );
772 #endif
773 
774 
775 
776 
777 
800  Zone2* createZone(ConstraintGraph2* pConstraintGraph,ZoneLocation zoneLoc,bool bVerbose=true);
801 
802 #if GEOM_PSEUDO3D==GEOM_TRUE
819  Zone2* getPeeledZone(double angleThreshold,bool bAvoidSplit);
820 #endif
821 
833  Zone2* createZone(const std::vector<ConstraintGraph2*>& vConstraintGraphs,ZoneLocation zoneLoc,const Point2& startPoint,bool bVerbose=true);
834 
847  Zone2* createZone(ConstraintGraph2* pConstraintGraph,ZoneLocation zoneLoc,const Point2& startPoint,bool bVerbose=true);
848 
859  Zone2* createZone( std::vector<Triangle2*>& vTriangles,bool bVerbose=true );
860 
861 
862 
875  void deleteZone(Zone2* pZone);
876 
877 
898 
899 
907 
915  bool isConstraint(Triangle2* pT,int ith) const;
916 
919  void getAliveConstraintSegments(std::vector<ConstraintSegment2*>& vAliveConstraintSegments) const;
920 
923  void getAliveAndDeadConstraintSegments(std::vector<ConstraintSegment2*>& vAllConstraintSegments) const;
924 
933 
934 
939  void getIncidentTriangles(Point2* pVtx,std::vector<Triangle2*>& vIncidentT) const;
940 
945  void getIncidentVertices(Point2* pVtx,std::vector<Point2*>& vIncidentVertices) const;
946 
947 
953  void writeObj(const char* filename) const;
954 
955 #ifndef SKIPTHREADS
956 #if GEOM_PSEUDO3D==GEOM_TRUE
965  bool writePly(const char* filename,bool bASCII=false) const;
973  bool writePly(std::ostream& os,bool bASCII=false) const;
974 #endif
975 #endif
976 
983  void writeObj(const char* filename,Zone2* pZone) const;
984 
985 
991  void writeWebScene(const char* path) const;
992 
998  void writeWebScene(const char* path,Zone2* pZone) const;
999 
1005 void subscribe(MsgType msgType,MsgBase* pMsg);
1006 
1012 void unsubscribe(MsgType msgType,MsgBase* pMsg);
1013 
1014 
1021  bool isConstraint(Point2* p0,Point2* p1) const;
1022 
1028  bool isConstraint(Point2* pVtx) const;
1029 
1031  void printLicense() const;
1032 
1034  bool checkZoneQuality(Zone2* pZone,double minAngle,const char* name,const AcceptExperimentalFeature& accept);
1036  void setName(const char* s);
1038  std::string getName() const;
1039 
1040 
1061  Zone2* importTriangles( std::vector<Point2>& vPoints,
1062  bool bReorientIfNeeded,
1063  bool bCreateExtendedBoundingBox
1064  );
1065 
1087  Zone2* importTriangles_robust(std::vector<Point2>& vPoints);
1088 
1089 #ifndef SKIPTHREADS
1090 #if GEOM_PSEUDO3D==GEOM_TRUE
1103  Zone2* importTrianglesFromPly(const char* filename);
1104 
1117  Zone2* importTrianglesFromPly(std::istream& is);
1118 #endif
1119 #endif
1120 
1125  Orientation2 getOrientation(const Point2& p0,const Point2& p1,const Point2& p2);
1126 
1157  void cutTriangles( const Point2& knifeStart,
1158  const Point2& knifeEnd,
1159  bool bTurnEdgesIntoConstraints);
1160 
1173  void cutTriangles( std::vector<Segment2>& vSegments,
1174  bool bTurnEdgesIntoConstraints);
1175 
1177  Zone2* createZone_cookieCutter(std::vector<Segment2>& vSegments,bool bProtectEdges);
1178 
1179 
1202  Zone2* createZone_cookieCutter( std::vector<Segment2>& vSegments,
1203  bool bProtectEdges,
1204  ConstraintGraph2*& pProtectedEdgesCG,
1205  ConstraintGraph2*& pBoundaryCG);
1206 
1207 
1257 #if GEOM_PSEUDO3D==GEOM_TRUE
1258  bool drape( std::vector<Segment2>& vSegmentsIn,
1259  std::vector<Segment2>& vSegmentsOut,
1260  double zTolerance,
1261  double offTolerance=0.0) const;
1262 #else
1263  bool drape( std::vector<Segment2>& vSegmentsIn,
1264  std::vector<Segment2>& vSegmentsOut) const;
1265 #endif
1266 
1267 
1268 
1269 
1275 void getDirectedEdges(std::vector<Edge2>& vDirectedEdgesOut) const;
1276 
1282 void getUndirectedEdges(std::vector<Edge2>& vUndirectedEdgesOut) const;
1283 
1285  /* Deprecated: Use setNumCPU() instead. This method is kept for
1286  * compatibility with existing applications. Internally it calls
1287  * setNumCPU(0) to automatically determine and use the number of
1288  * available CPU cores.
1289  */
1290  void enableMultithreading();
1291  // Development functions, not for public use
1293  void internal(int au,int fu,const char* s="");
1295  Dt2* getImpl();
1297  void setDev(const char* s,int ival,double dval);
1298 
1299 protected:
1301  void initFade(unsigned numExpectedVertices);
1303  Fade_2D& operator=(const Fade_2D&); // No assignment allowed
1305  Dt2* pImpl;
1306 }; // end of class
1307 
1308 
1316 inline void copyFade(Fade_2D& dt0,std::vector<Zone2*>& vZones0,Fade_2D& dt1,std::vector<Zone2*>& vZones1)
1317 {
1318  std::stringstream ss;
1319  dt0.saveTriangulation(ss,vZones0);
1320  dt1.load(ss,vZones1);
1321 }
1322 
1323 
1324 
1325 
1326 } // (namespace)
1327 
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:1316
VtkColor
Enumeration of colors used by the VTKWriter.
Definition: VtkWriter.h:39
Bbox2 is an axis aligned 2D bounding box.
Definition: Bbox2.h:37
CloudPrepare simplifies overdense point clouds and helps to avoid memory-usage-peaks during data tran...
Definition: CloudPrepare.h:87
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:74
void getTrianglePointers(std::vector< Triangle2 * > &vAllTriangles) const
Get pointers to all triangles.
void writeWebScene(const char *path, Zone2 *pZone) const
Write a zone to an *.obj file.
void showVtk(VtkWriter *pVtkWriter, VtkColor color, VtkColor constraintColor=VTK_TRANSPARENT) const
Draws the triangulation using an existing VtkWriter.
bool getHeight(double x, double y, double &heightOut, Triangle2 *pApproxT=NULL, double tolerance=0) const
Compute the height of a certain point.
bool load(const char *filename, std::vector< Zone2 * > &vZones)
Load a triangulation.
void getVertexPointers(std::vector< Point2 * > &vAllPoints) const
Get pointers to all vertices.
void statistics(const char *s) const
Statistics.
void cutTriangles(const Point2 &knifeStart, const Point2 &knifeEnd, bool bTurnEdgesIntoConstraints)
Cut through a triangulation.
void getConvexHull(bool bAllVertices, std::vector< Point2 * > &vConvexHullPointsOut)
Compute the convex hull.
void getAliveConstraintSegments(std::vector< ConstraintSegment2 * > &vAliveConstraintSegments) const
Get active (alive) constraint segments.
void insert(const std::vector< Point2 > &vInputPoints)
Insert a vector of points.
void showGeomview(Visualizer3 *pVis, const char *color="1 1 1 0.5") const
Draws the triangulation in 3D.
void writeWebScene(const char *path) const
Write the current triangulation to an *.obj file.
Voronoi2 * getVoronoiDiagram()
Get the Voronoi diagram.
void unsubscribe(MsgType msgType, MsgBase *pMsg)
Unregister a message receiver.
Point2 * getVertexPointer(const Point2 &p)
Get pointer.
void setFastMode(bool bFast)
Set fast mode.
void remove(std::vector< Point2 * > &vPoints)
Remove vertices.
void insert(CloudPrepare *pCloudPrepare, bool bClear=true)
Insert points from a CloudPrepare object.
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 insert(int numPoints, double *aCoordinates, Point2 **aHandles)
Insert points from an array.
void refineAdvanced(MeshGenParams *pParameters)
Delaunay refinement and grid meshing.
size_t numberOfTriangles() const
Number of triangles.
void writeObj(const char *filename, Zone2 *pZone) const
Write a zone to an *.obj file.
void insert(const std::vector< Point2 > &vInputPoints, std::vector< Point2 * > &vHandles)
Insert points from a std::vector and store pointers in vHandles.
ConstraintGraph2 * createConstraint(std::vector< Segment2 > &vSegments, ConstraintInsertionStrategy cis, bool bOrientedSegments=false, bool bUseHeightOfLatest=false)
Add constraint edges (edges, polyline, polygon) .
ConstraintSegment2 * getConstraintSegment(Point2 *p0, Point2 *p1) const
Retrieve a ConstraintSegment2.
bool drape(std::vector< Segment2 > &vSegmentsIn, std::vector< Segment2 > &vSegmentsOut, double zTolerance, double offTolerance=0.0) const
Drape segments along a surface.
bool saveZones(std::ostream &stream, std::vector< Zone2 * > &vSaveZones)
Save zones.
Zone2 * getPeeledZone(double angleThreshold, bool bAvoidSplit)
Get a peeled zone The getPeeledZone() function returns a zone that includes all triangles of dt excep...
bool isConstraint(Point2 *p0, Point2 *p1) const
Check if an edge is a constraint edge.
Point2 * insert(const Point2 &p)
Insert a single point.
void applyConstraintsAndZones()
Apply conforming constraints and zones (deprecated!)
bool hasArea() const
Check if the triangulation contains triangles (which is the case if at least 3 non-collinear points e...
void refine(Zone2 *pZone, double minAngleDegree, double minEdgeLength, double maxEdgeLength, bool bAllowConstraintSplitting)
Delaunay refinement.
void remove(Point2 *pVertex)
Remove a single vertex.
Zone2 * createZone(ConstraintGraph2 *pConstraintGraph, ZoneLocation zoneLoc, bool bVerbose=true)
Create a zone.
bool saveTriangulation(const char *filename, std::vector< Zone2 * > &vSaveZones)
Save a triangulation.
void showGeomview(const char *filename, const char *color="1 1 1 0.5") const
Draws the triangulation in 3D.
bool writePly(const char *filename, bool bASCII=false) const
Write the current triangulation to a *.ply file.
Point2 * getNearestNeighbor(const Point2 &p)
Get nearest neighbor.
Zone2 * importTrianglesFromPly(const char *filename)
Import Triangles from *.PLY.
void deleteZone(Zone2 *pZone)
Delete a Zone2 object.
bool writePly(std::ostream &os, bool bASCII=false) const
Write the current triangulation to a *.ply file.
size_t numberOfPoints() const
Number of points.
bool saveZones(const char *filename, std::vector< Zone2 * > &vSaveZones)
Save zones.
void getIncidentTriangles(Point2 *pVtx, std::vector< Triangle2 * > &vIncidentT) const
Get incident triangles.
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 show(Visualizer2 *pVis, bool bWithConstraints=true) const
Draws the triangulation using an existing Visualizer2 object.
Zone2 * importTriangles(std::vector< Point2 > &vPoints, bool bReorientIfNeeded, bool bCreateExtendedBoundingBox)
Import Triangles.
Zone2 * importTrianglesFromPly(std::istream &is)
Import Triangles from *.PLY.
Triangle2 * getAdjacentTriangle(Point2 *p0, Point2 *p1) const
Get adjacent triangle.
Triangle2 * locate(const Point2 &p)
Locate a triangle which contains p.
bool isConstraint(Point2 *pVtx) const
Check if a vertex is a constraint vertex.
void showVtk(const char *filename, VtkColor color, VtkColor constraintColor=VTK_TRANSPARENT) const
Draws the triangulation using the VTK file format.
void exportTriangulation(FadeExport &fadeExport, bool bWithCustomIndices, bool bClear)
Export triangulation data from Fade.
Fade_2D(unsigned numExpectedVertices=3)
Constructor of the triangulation class.
void getUndirectedEdges(std::vector< Edge2 > &vUndirectedEdgesOut) const
Get undirected edges.
Zone2 * createZone(std::vector< Triangle2 * > &vTriangles, bool bVerbose=true)
Create a zone defined by a vector of triangles.
void show(const char *filename, bool bWithConstraints=true) const
Draws the triangulation as postscript or pdf file.
Zone2 * importTriangles_robust(std::vector< Point2 > &vPoints)
Import Triangles - Robust (Error-tolerant version)
void getDirectedEdges(std::vector< Edge2 > &vDirectedEdgesOut) const
Get directed edges Edges are counterclockwise oriented around their triangle. The present method retu...
bool checkValidity(bool bCheckEmptyCircleProperty, const char *msg) const
Checks if a triangulation is valid.
bool isConstraint(Triangle2 *pT, int ith) const
Check if an edge is a constraint edge.
void printLicense() const
Prints license information.
Fade_2D(Fade_2D &other, std::vector< Zone2 * > vZoneIn=std::vector< Zone2 * >(), std::vector< Zone2 * > vZoneOut=std::vector< Zone2 * >())
Copy constructor, copy also zones.
int setNumCPU(int numCPU)
Set the number CPU cores for multithreading (deprecated)
void subscribe(MsgType msgType, MsgBase *pMsg)
Register a message receiver.
void writeObj(const char *filename) const
Write the current triangulation to an *.obj file.
bool saveTriangulation(std::ostream &stream, std::vector< Zone2 * > &vSaveZones)
Save a triangulation.
void cutTriangles(std::vector< Segment2 > &vSegments, bool bTurnEdgesIntoConstraints)
Cut through a triangulation.
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.
bool load(std::istream &stream, std::vector< Zone2 * > &vZones)
Load a triangulation.
Orientation2 getOrientation(const Point2 &p0, const Point2 &p1, const Point2 &p2)
Compute the orientation of 3 points.
void getIncidentVertices(Point2 *pVtx, std::vector< Point2 * > &vIncidentVertices) const
Get incident vertices.
double measureTriangulationTime(const std::vector< Point2 > &vPoints)
Measure the Delaunay triangulation time.
~Fade_2D()
Destructor.
void getAliveAndDeadConstraintSegments(std::vector< ConstraintSegment2 * > &vAllConstraintSegments) const
Get all (alive and dead) constraint segments.
Bbox2 computeBoundingBox() const
Compute the axis-aligned bounding box of the points.
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:53
Triangle.
Definition: Triangle2.h:60
Visualizer2 is a PDF- and Postscript writer.
Definition: Visualizer2.h:54
Visualizer3 is a 3D scene writer for the Geomview viewer.
Definition: Visualizer3.h:37
A writer for the VTK file format.
Definition: VtkWriter.h:62
Zone2 is a certain defined area of a triangulation.
Definition: Zone2.h:96
FadeExport is a simple struct to export triangulation data.
Definition: FadeExport.h:43