Fade2D Documentation pages v2.16.8
Delaunay Features
GEOM_FADE2D::Fade_2D Class Reference

Fade_2D is a class that contains a Delaunay triangulation. More...

#include <Fade_2D.h>

Public Member Functions

 Fade_2D (Fade_2D &other, std::vector< Zone2 * > vZoneIn=std::vector< Zone2 * >(), std::vector< Zone2 * > vZoneOut=std::vector< Zone2 * >())
 Copy constructor. More...
 
 Fade_2D (unsigned numExpectedVertices=3)
 Constructor of the triangulation class. More...
 
 ~Fade_2D ()
 Destructor.
 
void applyConstraintsAndZones ()
 Apply conforming constraints and zones (deprecated!) More...
 
bool checkValidity (bool bCheckEmptyCircleProperty, const char *msg) const
 Checks the validity of the triangulation. More...
 
Bbox2 computeBoundingBox () const
 Computes the axis-aligned bounding box of the triangulation. More...
 
ConstraintGraph2createConstraint (std::vector< Segment2 > &vSegments, ConstraintInsertionStrategy cis, bool bOrientedSegments=false)
 Insert constraint edges (edges, polyline or polygon) . More...
 
Zone2createZone (const std::vector< ConstraintGraph2 * > &vConstraintGraphs, ZoneLocation zoneLoc, const Point2 &startPoint, bool bVerbose=true)
 Create a Zone2 bounded by multiple ConstraintGraph2 objects by growing from a seed point. More...
 
Zone2createZone (ConstraintGraph2 *pConstraintGraph, ZoneLocation zoneLoc, bool bVerbose=true)
 Create a Zone2 . More...
 
Zone2createZone (ConstraintGraph2 *pConstraintGraph, ZoneLocation zoneLoc, const Point2 &startPoint, bool bVerbose=true)
 Create a Zone2 bounded by a ConstraintGraph2 object by growing from a seed point. More...
 
Zone2createZone (std::vector< Triangle2 * > &vTriangles, bool bVerbose=true)
 Create a Zone2 defined by a vector of triangles. More...
 
Zone2createZone_cookieCutter (std::vector< Segment2 > &vSegments, bool bProtectEdges)
 Cookie Cutter (deprecated) More...
 
Zone2createZone_cookieCutter (std::vector< Segment2 > &vSegments, bool bProtectEdges, ConstraintGraph2 *&pProtectedEdgesCG, ConstraintGraph2 *&pBoundaryCG)
 Cookie Cutter (deprecated) More...
 
Zone2createZone_cookieCutter2 (std::vector< Segment2 > &vSegments, bool bProtectEdges)
 Cookie Cutter2. More...
 
Zone2createZone_cookieCutter2 (std::vector< Segment2 > &vSegments, bool bProtectEdges, ConstraintGraph2 *&pProtectedEdgesCG, ConstraintGraph2 *&pBoundaryCG)
 Cookie Cutter2. More...
 
void cutTriangles (const Point2 &knifeStart, const Point2 &knifeEnd, bool bTurnEdgesIntoConstraints)
 Cut through a triangulation. More...
 
void cutTriangles (std::vector< Segment2 > &vSegments, bool bTurnEdgesIntoConstraints)
 Cut through a triangulation. More...
 
void deleteConstraintsAndZones ()
 Deletes all constraints and zones. More...
 
void deleteZone (Zone2 *pZone)
 Delete a Zone2 object. More...
 
bool drape (std::vector< Segment2 > &vSegmentsIn, std::vector< Segment2 > &vSegmentsOut, double offTolerance=0.0, bool bKeepOffSegments=false) const
 Drape segments onto a surface. More...
 
void exportTriangulation (FadeExport &fadeExport, bool bWithCustomIndices, bool bClear)
 Export triangulation data. More...
 
Triangle2getAdjacentTriangle (Point2 *p0, Point2 *p1) const
 Get the adjacent triangle of an oriented edge. More...
 
void getAliveAndDeadConstraintSegments (std::vector< ConstraintSegment2 * > &vAllConstraintSegments) const
 Get all (alive and dead) constraint segments.
 
void getAliveConstraintSegments (std::vector< ConstraintSegment2 * > &vAliveConstraintSegments) const
 Get active (alive) constraint segments. More...
 
ConstraintSegment2getConstraintSegment (Point2 *p0, Point2 *p1) const
 Retrieve a ConstraintSegment2. More...
 
void getConvexHull (bool bAllVertices, std::vector< Point2 * > &vConvexHullPointsOut)
 Compute the convex hull. More...
 
void getDirectedEdges (std::vector< Edge2 > &vDirectedEdgesOut) const
 Get directed edges Edges are counterclockwise oriented around their triangle. The present method returns directed edges. That means each edge(a,b) is returend twice, as edge(a,b) and as edge(b,a).
 
void getIncidentTriangles (Point2 *pVtx, std::vector< Triangle2 * > &vIncidentT) const
 Get incident triangles. More...
 
void getIncidentVertices (Point2 *pVtx, std::vector< Point2 * > &vAdjVtx) const
 Get adjacent vertices. More...
 
double getMedianEdgeLength2D ()
 Get the median 2D edge length. More...
 
Point2getNearestNeighbor (const Point2 &p)
 Get nearest neighbor. More...
 
Orientation2 getOrientation (const Point2 &p0, const Point2 &p1, const Point2 &p2)
 Compute the orientation of 3 points. More...
 
void getTrianglePointers (std::vector< Triangle2 * > &vAllTriangles) const
 Get pointers to all triangles. More...
 
void getUndirectedEdges (std::vector< Edge2 > &vUndirectedEdgesOut) const
 Get undirected edges. More...
 
Point2getVertexPointer (const Point2 &p)
 Get vertex pointer. More...
 
void getVertexPointers (std::vector< Point2 * > &vAllPoints) const
 Get pointers to all vertices. More...
 
Voronoi2getVoronoiDiagram ()
 Get the Voronoi diagram of the triangulation. More...
 
void getZones (std::vector< Zone2 * > &vZones) const
 Retrieves all zones of the triangulation. More...
 
bool hasArea () const
 Check if the triangulation contains any triangles. More...
 
Zone2importTriangles (std::vector< Point2 > &vPoints, bool bReorientIfNeeded, bool bCreateExtendedBoundingBox)
 Import Triangles (deprecated in favor of the robust version) More...
 
Zone2importTriangles_robust (std::vector< Point2 > &vPoints)
 Import Triangles - Robust (Error-tolerant version) More...
 
Point2insert (const Point2 &p)
 Insert a single point. More...
 
void insert (const std::vector< Point2 > &vInputPoints)
 Insert a vector of points. More...
 
void insert (const std::vector< Point2 > &vInputPoints, std::vector< Point2 * > &vHandles)
 Insert points from a std::vector<Point2> and store pointers in vHandles. More...
 
void insert (int numPoints, double *aCoordinates, Point2 **aHandles)
 Insert points from an array. More...
 
bool isConstraint (Point2 *p0, Point2 *p1) const
 Check if an edge is a constraint edge. More...
 
bool isConstraint (Point2 *pVtx) const
 Check if a vertex is a constraint vertex. More...
 
bool isConstraint (Triangle2 *pT, int ith) const
 Check if an edge is a constraint edge. More...
 
bool load (const char *filename, std::vector< Zone2 * > &vZones)
 Load a triangulation from a binary file. More...
 
bool load (std::istream &stream, std::vector< Zone2 * > &vZones)
 Load a triangulation from a binary file. More...
 
Triangle2locate (const Point2 &p)
 Locate a triangle that contains p. More...
 
double measureTriangulationTime (const std::vector< Point2 > &vPoints)
 Measure the Delaunay triangulation time. More...
 
size_t numberOfPoints () const
 Number of points. More...
 
size_t numberOfTriangles () const
 Number of triangles. More...
 
void printLicense () const
 Prints license information.
 
void refine (Zone2 *pZone, double minAngleDegree, double minEdgeLength, double maxEdgeLength, bool bAllowConstraintSplitting)
 Delaunay refinement. More...
 
void refineAdvanced (MeshGenParams *pParameters)
 Delaunay refinement and grid meshing. More...
 
void remove (Point2 *pVertex)
 Remove a single vertex. More...
 
void remove (std::vector< Point2 * > &vPoints)
 Remove vertices. More...
 
Zone2replicateZoneBoundary (Zone2 *pZoneIn)
 Replicate the boundaries of a zone in the current Fade_2D instance. More...
 
void reset ()
 Reset the triangulation. More...
 
bool saveTriangulation (const char *filename, std::vector< Zone2 * > &vSaveZones)
 Save a triangulation to a binary file. More...
 
bool saveTriangulation (std::ostream &stream, std::vector< Zone2 * > &vSaveZones)
 Save a triangulation to a binary file. More...
 
bool saveZones (const char *filename, std::vector< Zone2 * > &vSaveZones)
 Save zones to a binary file. More...
 
bool saveZones (std::ostream &stream, std::vector< Zone2 * > &vSaveZones)
 Save zones to a binary file. More...
 
void setFastMode (bool bFast)
 Set fast mode. More...
 
int setNumCPU (int numCPU)
 Set the number of CPU cores for multithreading (deprecated) More...
 
void show (const char *filename, bool bWithConstraints=true) const
 Draws the triangulation as a PostScript or PDF file. More...
 
void show (Visualizer2 *pVis, bool bWithConstraints=true) const
 Draws the triangulation as a PostScript or PDF file. More...
 
void showVtk (const char *filename, VtkColor color, VtkColor constraintColor=VTK_TRANSPARENT) const
 Draws the triangulation using the VTK file format. More...
 
void showVtk (VtkWriter *pVtkWriter, VtkColor color, VtkColor constraintColor=VTK_TRANSPARENT) const
 Draws the triangulation using the VTK file format. More...
 
void statistics (const char *s) const
 Statistics. More...
 
void subscribe (MsgType msgType, MsgBase *pMsg)
 Registers a message receiver. More...
 
void unsubscribe (MsgType msgType, MsgBase *pMsg)
 Unregisters a message receiver. More...
 
void writeObj (const char *filename) const
 Write the current triangulation to an .obj file. More...
 
void writeObj (const char *filename, Zone2 *pZone) const
 Write a Zone2 an .obj file. More...
 
void writeWebScene (const char *path) const
 Create a scene of a triangulation that can be viewed in a web browser. More...
 
void writeWebScene (const char *path, Zone2 *pZone) const
 Create a scene of a Zone2 that can be viewed in a web browser. More...
 

Detailed Description

The library comes in two versions: 2D and 2.5D. The 2D version computes classic Delaunay triangulations in the xy-plane, while the 2.5D version includes a z-coordinate and provides additional algorithms for lifted triangulations.

A 2D Delaunay triangulation with constraint segments
2D Delaunay triangulation with constraint edges.
A 2.5D Delaunay triangulation or lifted Delaunay triangulation, computed with the 2.5D version of Fade
2.5D Delaunay triangulation (Fade 2.5D)

You are currently viewing the 2D documentation. If you are interested in lifted triangulations, you can switch to the 2.5D documentation.

Constructor & Destructor Documentation

◆ Fade_2D() [1/2]

GEOM_FADE2D::Fade_2D::Fade_2D ( unsigned  numExpectedVertices = 3)
explicit
Parameters
numExpectedVerticesThe number of points that will be inserted. This is an optional parameter.

◆ Fade_2D() [2/2]

GEOM_FADE2D::Fade_2D::Fade_2D ( Fade_2D other,
std::vector< Zone2 * >  vZoneIn = std::vector< Zone2 * >(),
std::vector< Zone2 * >  vZoneOut = std::vector< Zone2 * >() 
)

This constructor creates a copy of the provided Fade_2D instance, including the option to specify zones to copy into the new instance.

Parameters
[in]otherThe Fade_2D instance to copy from.
[in]vZoneInA vector of Zone2* objects to copy.
[out]vZoneOutUsed to return pointers of the copied Zone2* objects in the new triangulation.

Member Function Documentation

◆ applyConstraintsAndZones()

void GEOM_FADE2D::Fade_2D::applyConstraintsAndZones ( )
Deprecated:
The present function applyConstraintsAndZones() as well as the two constraint insertion strategies CIS_CONFORMING_DELAUNAY and CIS_CONFORMING_DELAUNAY_SEGMENT_LEVEL are deprecated. These are only kept for backwards compatibilty. The replacement is CIS_CONSTRAINED_DELAUNAY along with the methods Fade_2D::drape() and/or ConstraintGraph2::makeDelaunay(). See the example code in examples_25D/terrain.cpp

This method establishes conforming constraint segments and zones which depend on them. For technical reasons conforming constraint segments are not immediately established but inserted at the end of the triangulation process. This step must be triggered manually i.e., it is up to the user to call applyConstraintsAndZones() before the resulting triangulation is used. If afterwards the triangulation is changed in any way, applyConstraintsAndZones() must be called again.

◆ checkValidity()

bool GEOM_FADE2D::Fade_2D::checkValidity ( bool  bCheckEmptyCircleProperty,
const char *  msg 
) const

This method checks the validity of the triangulation's data structure. It verifies various properties, and can optionally recheck the empty circle property using slow multiprecision arithmetic for higher accuracy.

Parameters
[in]bCheckEmptyCirclePropertyA boolean flag that specifies whether multiprecision arithmetic should be used to recheck the empty circle property. Enabling this option will slow down the check.
[in]msgA debug message that will be displayed in the terminal output to indicate which call to checkValidity is currently running.

◆ computeBoundingBox()

Bbox2 GEOM_FADE2D::Fade_2D::computeBoundingBox ( ) const

This method calculates the smallest axis-aligned bounding box (AABB) that encloses all points currently present in the triangulation.

Returns
Bbox2 An object representing the computed bounding box.
Note
If no points have been inserted into the triangulation yet, the returned Bbox2 object is invalid. In this case, calling Bbox2::isValid() on the returned object will return false.

◆ createConstraint()

ConstraintGraph2* GEOM_FADE2D::Fade_2D::createConstraint ( std::vector< Segment2 > &  vSegments,
ConstraintInsertionStrategy  cis,
bool  bOrientedSegments = false 
)
Parameters
vSegmentsSegments to be inserted as enforced edges of the triangulation.
cisUse CIS_CONSTRAINED_DELAUNAY, the supported constraint-insertion-strategy.
bOrientedSegmentsSpecifies whether the input segments in vSegments are oriented with a counterclockwise direction around an enclosed area. This property is required if the returned ConstraintGraph2 object shall be used for later Zone2 creation. If bOrientedSegments is false and vSegments represents a simple polygon, the edges are automatically reoriented to establish this property. For complex or uncertain polygon orientations, use the PolygonClipper to repair and orient vSegments upfront.
Returns
A pointer to the created ConstraintGraph2 object

The present method subdivides the segments only if required, which is the case if existing vertices or constraint edges are crossed. If you want to subdivide the created constraint edges subsequently to achieve better triangle shapes (see the third image), use ConstraintGraph2::makeDelaunay() on the returned ConstraintGraph2.

See also
ConstraintGraph2::makeDelaunay() subdivides a ConstraintGraph2 to achieve better triangle shapes.
https://www.geom.at/example4-zones-defined-areas-in-triangulations/ provides examples of Zone2 creation.
https://www.geom.at/polygon-clipper-repairing-polygons/ contains examples of repairing polygons.
Delaunay triangulation without constraints
Constrained Delaunay triangulation
Conforming Delaunay triangulation using ConstraintGraph::makeDelaunay()

◆ createZone() [1/4]

Zone2* GEOM_FADE2D::Fade_2D::createZone ( const std::vector< ConstraintGraph2 * > &  vConstraintGraphs,
ZoneLocation  zoneLoc,
const Point2 startPoint,
bool  bVerbose = true 
)

Creates a Zone2 object bounded by multiple ConstraintGraph2 objects, grown from a seed point.

Zone2, grown from a seed point until the constraint edges limit the growing
Parameters
vConstraintGraphsA vector of ConstraintGraph2 objects.
zoneLocThe ZoneLocation must be ZL_GROW.
startPointThe point from which the area grows until the borders defined in vConstraintGraphs are reached.
bVerboseIs true by default and triggers a warning if NULL is returned.
Returns
A pointer to the created Zone2 object, or NULL if zoneLoc != ZL_GROW or no triangles exist.

◆ createZone() [2/4]

Zone2* GEOM_FADE2D::Fade_2D::createZone ( ConstraintGraph2 pConstraintGraph,
ZoneLocation  zoneLoc,
bool  bVerbose = true 
)

Creates a Zone2 object which represents an area of a triangulation.

Parameters
[in]zoneLoccan be ZL_INSIDE, ZL_OUTSIDE or ZL_GLOBAL.
[in]pConstraintGraphpoints to an existing ConstraintGraph2 object (which must be oriented) or is NULL if zoneLoc is ZL_GLOBAL.
[in]bVerboseis true by default and causes a warning if NULL is returned.
Returns
A pointer to the created Zone2 object, or NULL if no triangles exist, or if pConstraintGraph->isOriented() returns false. Make sure to check for NULL before proceeding.
A zone in a 2D Delaunay triangulation
A Zone in a 2D Delaunay triangulation

◆ createZone() [3/4]

Zone2* GEOM_FADE2D::Fade_2D::createZone ( ConstraintGraph2 pConstraintGraph,
ZoneLocation  zoneLoc,
const Point2 startPoint,
bool  bVerbose = true 
)

Creates a Zone2 object bounded by a ConstraintGraph2 object, grown from a seed point.

Zone2, grown from a seed point until the constraint edges limit the growing
Parameters
pConstraintGraphA pointer to the ConstraintGraph2 that defines the boundary.
zoneLocThe ZoneLocation must be ZL_GROW.
startPointThe point from which the area grows until the borders defined in pConstraintGraph are reached.
bVerboseIs true by default and triggers a warning if NULL is returned.
Returns
A pointer to the created Zone2 object, or NULL if zoneLoc != ZL_GROW or no triangles exist.

◆ createZone() [4/4]

Zone2* GEOM_FADE2D::Fade_2D::createZone ( std::vector< Triangle2 * > &  vTriangles,
bool  bVerbose = true 
)

This method creates a Zone2 defined by the triangles in @p vTriangles.

Parameters
vTrianglesA vector of Triangle2* pointers defining the area.
bVerboseTriggers a warning if NULL is returned. Default: true.
Returns
A pointer to the created Zone2 object, or NULL if vTriangles is empty.
Zone2, defined by triangles

◆ createZone_cookieCutter() [1/2]

Zone2* GEOM_FADE2D::Fade_2D::createZone_cookieCutter ( std::vector< Segment2 > &  vSegments,
bool  bProtectEdges 
)
Deprecated:
This method is deprecated. It works as expected, but does not support zones with holes. To avoid behavioral changes, this method is left unchanged and has been deprecated in favor of the new method Fade_2D::createZone_cookieCutter2(), which allows for zones with holes, is faster, and provides improved numerical accuracy.

◆ createZone_cookieCutter() [2/2]

Zone2* GEOM_FADE2D::Fade_2D::createZone_cookieCutter ( std::vector< Segment2 > &  vSegments,
bool  bProtectEdges,
ConstraintGraph2 *&  pProtectedEdgesCG,
ConstraintGraph2 *&  pBoundaryCG 
)
Deprecated:
This method is deprecated. It works as expected, but does not support zones with holes. To avoid behavioral changes, this method is left unchanged and has been deprecated in favor of the new method Fade_2D::createZone_cookieCutter2(), which allows for zones with holes, is faster, and provides improved numerical accuracy.

◆ createZone_cookieCutter2() [1/2]

Zone2* GEOM_FADE2D::Fade_2D::createZone_cookieCutter2 ( std::vector< Segment2 > &  vSegments,
bool  bProtectEdges 
)

The Cookie Cutter cuts out a part of a 2D or 2.5D triangulation and returns it as a Zone2 object.

Cookie Cutter: Cuts a part of the triangulation, with the red polygon (showing vSegments) defining the shape.
Parameters
[in]vSegmentsspecifies the input polygon.
[in]bProtectEdgesspecifies if existing edges shall be locked by turning them into constraint segments
Returns
a Zone2 object consisting of all triangles inside the polygon or NULL when the operation has failed due to unmet preconditions.

◆ createZone_cookieCutter2() [2/2]

Zone2* GEOM_FADE2D::Fade_2D::createZone_cookieCutter2 ( std::vector< Segment2 > &  vSegments,
bool  bProtectEdges,
ConstraintGraph2 *&  pProtectedEdgesCG,
ConstraintGraph2 *&  pBoundaryCG 
)

The Cookie Cutter cuts out a part of a triangulation and returns it as a Zone2 object.

Cookie Cutter: Cuts a part of the triangulation, with the red polygon (showing vSegments) defining the shape.
Parameters
[in]vSegmentsspecifies the input polygon.
[in]bProtectEdgesspecifies if existing edges shall be locked by turning them into constraint segments
[out]pProtectedEdgesCGis used to return the ConstraintGraph2* of edges that have been protected
[out]pBoundaryCGis used to return the ConstraintGraph2* of the boundary
Returns
a Zone2 object consisting of all triangles inside the polygon or NULL when the operation has failed due to unmet preconditions.

◆ cutTriangles() [1/2]

void GEOM_FADE2D::Fade_2D::cutTriangles ( const Point2 knifeStart,
const Point2 knifeEnd,
bool  bTurnEdgesIntoConstraints 
)
Parameters
knifeStartis one point of the knife segment
knifeEndis the second point of the knife segment
bTurnEdgesIntoConstraintsturns all 3 edges of each intersected triangle into constraint segments.

This method inserts a constraint edge knife(knifeStart,knifeEnd). If existing edges E are intersected by knife, then knife is subdivided at the intersection points P.

In any case knife will exist (in a possibly subdivided form) in the result. But a consequence of the insertion of the points P is that the edges E and even edges which are not intersected by knife may be flipped. Use bTurnEdgesIntoConstraints=true to avoid that.

Note
The intersection point of two line segments may not be exactly representable in double precision floating point arithmetic and thus tiny rounding errors may occur. As a consequence two very close intersection points may be rounded to the same coordinates.
When more than one knife segment is inserted then the method void cutTriangles(std::vector<Segment2>& vSegments,bool bTurnEdgesIntoConstraints) should be used. The reason is that each individual cut operation changes the triangulation and thus iterative calls to the present version of the method can lead to a different result.

◆ cutTriangles() [2/2]

void GEOM_FADE2D::Fade_2D::cutTriangles ( std::vector< Segment2 > &  vSegments,
bool  bTurnEdgesIntoConstraints 
)
Parameters
vSegmentsare the knife segments
bTurnEdgesIntoConstraintsspecifies if intersected edges shall automatically be turned into constraints

Same method as void cutTriangles(const Point2& knifeStart,const Point2& knifeEnd,bool bTurnEdgesIntoConstraints) but it takes a vector of segments instead of a single segment. This is the recommended method to cut through a triangulation when more than one knife segment exists.

◆ deleteConstraintsAndZones()

void GEOM_FADE2D::Fade_2D::deleteConstraintsAndZones ( )

This method deletes all ConstraintGraph2 objects, all ConstraintSegment2 objects, and all Zone2 objects. It does not modify existing triangles; existing edges that have been created as constraint segments before are not flipped back to satisfy the 'empty circle property' of a Delaunay triangulation. However, these edges become unprotected, and any subsequent modifications to the triangulation may cause them to be flipped.

Warning
This method is dangerous. It violates the 'empty circle property', and it invalidates the pointers of the deleted objects. Use it only if you are fully aware of the consequences.

◆ deleteZone()

void GEOM_FADE2D::Fade_2D::deleteZone ( Zone2 pZone)

Zone2 objects are automatically destroyed with their Fade_2D objects. In addition this method provides the possibility to eliminate Zone2 objects earlier.

Note
Zones are designed transparently: When two zones Z1 and Z2 are combined to a new one Z3 (for example through a boolean operation) then Z1, Z2 and Z3 form a tree datastructure such that changes in the leaf nodes Z1 and Z2 can propagate up to the root node Z3. For this reason Z1 and Z2 must be alive as long as Z3 is used.

◆ drape()

bool GEOM_FADE2D::Fade_2D::drape ( std::vector< Segment2 > &  vSegmentsIn,
std::vector< Segment2 > &  vSegmentsOut,
double  offTolerance = 0.0,
bool  bKeepOffSegments = false 
) const

Projects the segments from vSegmentsIn onto the triangulation. Degenerate input segments are ignored. Segments are subdivided at intersections with triangulation edges.

Parameters
[in]vSegmentsInInput segments.
[out]vSegmentsOutOutput segments.
[in]offToleranceUsed to prevent truncation of segments that are slightly outside the convex hull due to tiny rounding errors. If this value is greater than 0, segment endpoints slightly outside the convex hull of the triangulation will be inserted into the triangulation.
[in]bKeepOffSegmentsSpecifies whether segments outside the triangulation should be kept. Default: false.
Returns
true if all input segments are inside the convex hull of the triangulation. Otherwise, false is returned, and the result is still valid but only contains the segment parts inside the convex hull unless bKeepOffSegments is true.

  • Drape: Input segments are draped (red) onto a TIN. They are subdivided (blue points) at intersections with triangulation edges
Note
Draping segments onto a triangulation does not insert them. Use Fade_2D::createConstraint() for that purpose.

◆ exportTriangulation()

void GEOM_FADE2D::Fade_2D::exportTriangulation ( FadeExport fadeExport,
bool  bWithCustomIndices,
bool  bClear 
)

This method exports triangulation data from the current Fade instance to a simple FadeExport structure, which is defined in a header and consists solely of arrays containing basic data types. The method also offers the option to clear the Fade object during the export to prevent memory usage peaks.

Parameters
[out]fadeExportA FadeExport structure that will hold the requested triangulation data after the export operation.
[in]bWithCustomIndicesA boolean flag that determines whether custom indices for the points should be included in the export.
[in]bClearA boolean flag that specifies whether the Fade instance should be cleared during the export to free memory. If true, all memory associated with the Fade instance is released, and all existing pointers to its objects will become invalid.

◆ getAdjacentTriangle()

Triangle2* GEOM_FADE2D::Fade_2D::getAdjacentTriangle ( Point2 p0,
Point2 p1 
) const

A method to retrieve the counterclockwise adjacent triangle of an oriented edge of the triangulation.

The adjacent triangle of edge (p0,p1)
Parameters
[in]p0,p1Pointers to the vertices of the triangulation that form the endpoints of the oriented edge.
Returns
A pointer to the Triangle2 containing the oriented edge (p0, p1), or NULL if no such triangle is present.
Note
Recall the counter-clockwise enumeration of vertices in a Triangle2: While an unoriented edge (p0, p1) can have two adjacent triangles, only one of them contains the oriented edge (p0, p1), while the other contains the oriented edge (p1, p0).

◆ getAliveConstraintSegments()

void GEOM_FADE2D::Fade_2D::getAliveConstraintSegments ( std::vector< ConstraintSegment2 * > &  vAliveConstraintSegments) const

This method retrieves the active (alive) constraint segments and stores them in the provided vector.

Parameters
vAliveConstraintSegmentsA reference to a vector that will be populated with pointers to the active constraint segments.

◆ getConstraintSegment()

ConstraintSegment2* GEOM_FADE2D::Fade_2D::getConstraintSegment ( Point2 p0,
Point2 p1 
) const
Parameters
[in]p0,p1The endpoints of the ConstraintSegment2
Returns
A pointer to the ConstraintSegment2 between p0 and p1, or NULL if the segment is not a constraint edge of the triangulation (or is dead because it has been split).

◆ getConvexHull()

void GEOM_FADE2D::Fade_2D::getConvexHull ( bool  bAllVertices,
std::vector< Point2 * > &  vConvexHullPointsOut 
)

Computes the convex hull of the triangulation.

Parameters
[in]bAllVerticesA boolean that determines whether all convex hull points are returned or if collinear ones should be removed.
[out]vConvexHullPointsOutA vector used to return the convex hull vertices in counterclockwise order.

◆ getIncidentTriangles()

void GEOM_FADE2D::Fade_2D::getIncidentTriangles ( Point2 pVtx,
std::vector< Triangle2 * > &  vIncidentT 
) const

Stores pointers to all triangles around the vertex pVtx into the vector vIncidentT.

Parameters
[in]pVtxA pointer to a vertex of the triangulation.
[out]vIncidentTA vector used to return pointers to the incident triangles of pVtx.
Incident triangles of Point2* pVtx

◆ getIncidentVertices()

void GEOM_FADE2D::Fade_2D::getIncidentVertices ( Point2 pVtx,
std::vector< Point2 * > &  vAdjVtx 
) const

Stores pointers to all vertices around the vertex pVtx into the vector vAdjVtx

Parameters
[in]pVtxA pointer to a vertex of the triangulation.
[out]vAdjVtxA vector used to return pointers to all vertices adjacent to pVtx.
Adjacent vertices of Point2* pVtx

◆ getMedianEdgeLength2D()

double GEOM_FADE2D::Fade_2D::getMedianEdgeLength2D ( )

Computes the median edge length of all edges in the triangulation in 2D.

Returns
The median 2D edge length of all edges in the triangulation.

◆ getNearestNeighbor()

Point2* GEOM_FADE2D::Fade_2D::getNearestNeighbor ( const Point2 p)

This method returns the closest vertex to the point p.

Nearest neighbor query: The query point (red) and the result vertex (blue)
Parameters
[in]pA pointer to the query point.
Returns
A pointer to the closest vertex to p.

◆ getOrientation()

Orientation2 GEOM_FADE2D::Fade_2D::getOrientation ( const Point2 p0,
const Point2 p1,
const Point2 p2 
)

Determines if the ordered points p0, p1, and p2 are collinear, clockwise, or counterclockwise.

Parameters
[in]p0,p1,p2The points whose orientation is computed.
Returns
ORIENTATION2_COLLINEAR, ORIENTATION2_CW, or ORIENTATION2_CCW.

◆ getTrianglePointers()

void GEOM_FADE2D::Fade_2D::getTrianglePointers ( std::vector< Triangle2 * > &  vAllTriangles) const

This method retrieves pointers to all existing triangles in the triangulation.

Parameters
[out]vAllTrianglesA vector used to return the pointers to all triangles.

◆ getUndirectedEdges()

void GEOM_FADE2D::Fade_2D::getUndirectedEdges ( std::vector< Edge2 > &  vUndirectedEdgesOut) const

Retrieves all undirected edges from the triangulation.

Parameters
[out]vUndirectedEdgesOutA vector used to return a unique set of undirected edges.

◆ getVertexPointer()

Point2* GEOM_FADE2D::Fade_2D::getVertexPointer ( const Point2 p)

Checks if the given point is a vertex of the triangulation. If is a vertex, the function returns a pointer to that vertex. Otherwise, it returns NULL.

Parameters
[in]pThe query point.
Returns
A pointer to the vertex if is a vertex of the triangulation, or NULL if is not a vertex.

◆ getVertexPointers()

void GEOM_FADE2D::Fade_2D::getVertexPointers ( std::vector< Point2 * > &  vAllPoints) const

Returns all vertex pointers of the triangulation. The order in which the points are returned is arbitrary.

Parameters
[out]vAllPointsA vector to return the vertex pointers.
Note
For duplicate point insertions, only one vertex pointer exists and it is stored only once. Thus, the number of returned points is smaller than the number of inserted points if duplicate points have been inserted.

◆ getVoronoiDiagram()

Voronoi2* GEOM_FADE2D::Fade_2D::getVoronoiDiagram ( )

This method returns the Voronoi diagram of the existing vertices, which is the dual diagram of their Delaunay triangulation.

Returns
A pointer to a Voronoi2 instance, representing the dual diagram of the triangulation.
Note
The Voronoi diagram is dynamically updated when the triangulation changes.
Warning
While the Voronoi diagram is the dual diagram of a Delaunay triangulation, it is NOT dual to a Constrained Delaunay triangulation, so do not create any constraints in it.
Voronoi diagram

◆ getZones()

void GEOM_FADE2D::Fade_2D::getZones ( std::vector< Zone2 * > &  vZones) const

This method retrieves all the zones in the triangulation and stores their pointers in the provided vector.

Parameters
[out]vZonesA vector used to return pointers to all existing Zone2 objects.

◆ hasArea()

bool GEOM_FADE2D::Fade_2D::hasArea ( ) const

A triangulation only contains triangles if at least three non-collinear points exist. If all inserted points are collinear, no triangles can be formed, regardless of how many points have been inserted. This method checks whether the triangulation currently contains any triangles, i.e., if an area exists.

Collinear points (black), triangles are generated as soon as the first non-collinear point (red) is inserted.
Returns
true if at least one triangle exists, false otherwise

◆ importTriangles()

Zone2* GEOM_FADE2D::Fade_2D::importTriangles ( std::vector< Point2 > &  vPoints,
bool  bReorientIfNeeded,
bool  bCreateExtendedBoundingBox 
)

This method imports triangles into an empty triangulation.

Deprecated:
Using this older version of the method is only recommended if you are certain that your input triangulation is 100% correct. If there is any possibility that the input points have been rounded (e.g., due to conversion from 64-bit to 32-bit precision or simply because they were parsed from an ASCII file), you should use the robust version, importTriangles_robust(), instead.
Parameters
vPointsA vector containing the input triangulation, with 3 subsequent points per triangle.
bReorientIfNeededSpecifies whether the orientations of the point triples should be checked and corrected. If you are confident that the point triples are already oriented in counterclockwise order, you can skip the orientation test.
bCreateExtendedBoundingBoxIf true, the four corners of an extended bounding box will be inserted as dummy points. This can be useful in some cases. If unsure, set this to false.
Returns
A pointer to a Zone2 object, or NULL if the input data is invalid.

◆ importTriangles_robust()

Zone2* GEOM_FADE2D::Fade_2D::importTriangles_robust ( std::vector< Point2 > &  vPoints)

This method imports triangles into an empty Fade object. Unlike the classic importTriangles() method, this robust version handles errors in the input data and automatically corrects overlapping triangles. This robust version is very fast, and due to the additional checks, it is also very safe. Thus, it is highly recommended to use it always rather than the classic importTriangles() method.

Parameters
vPointsA vector containing the input triangulation, with 3 subsequent points per triangle.
Returns
A pointer to a Zone2 object, or NULL if the input data is invalid.
Note
A Delaunay triangulation is always convex, and thus additional fill-triangles are created, if necessary, to complete the convex hull. Consequently, when you later call Fade_2D::getTrianglePointers(), you will also receive the additional triangles. If you are only interested in the triangles that you have imported, use the Zone2* returned by the present command, and call Zone2::getTriangles() instead.

◆ insert() [1/4]

Point2* GEOM_FADE2D::Fade_2D::insert ( const Point2 p)

This method inserts a single point into the triangulation.

Parameters
pThe point to be inserted.
Returns
A pointer to the vertex in the triangulation.

The triangulation keeps a copy of p, and the return value is a pointer to this copy. If duplicate points are inserted, the triangulation does not create new copies but instead returns a pointer to the copy of the first insertion.

Note
This method offers very good performance, but it is still faster to pass all points at once, if possible.

◆ insert() [2/4]

void GEOM_FADE2D::Fade_2D::insert ( const std::vector< Point2 > &  vInputPoints)

This method inserts a vector of points into the triangulation.

Parameters
vInputPointsA vector containing the points to be inserted.
Note
To activate multithreading, use setGlobalNumCPU() before using this method.

◆ insert() [3/4]

void GEOM_FADE2D::Fade_2D::insert ( const std::vector< Point2 > &  vInputPoints,
std::vector< Point2 * > &  vHandles 
)

This method inserts points from vInputPoints and stores the corresponding pointers to the inserted points in vHandles.

Parameters
vInputPointsA vector containing the points to be inserted.
vHandlesAn empty vector that will be populated with Point2 pointers.

Internally, the triangulation keeps copies of the inserted points, and these copies are returned in vHandles (in the same order as in vInputPoints). If duplicate points are found in vInputPoints, only one copy will be created, and a pointer to this unique copy will be stored in vHandles for each occurrence of that point.

Note
To activate multithreading, use setGlobalNumCPU() before using this method.

◆ insert() [4/4]

void GEOM_FADE2D::Fade_2D::insert ( int  numPoints,
double *  aCoordinates,
Point2 **  aHandles 
)

This method inserts points from an array of coordinates and stores the corresponding pointers to the inserted points in aHandles.

Parameters
numPointsThe number of points to be inserted.
aCoordinatesAn array of 2n double values, representing the coordinates of the points in the format: {x0, y0, x1, y1, ..., xn, yn}.
aHandlesAn empty array of size n where pointers to the inserted points will be stored by Fade.
Note
To activate multithreading, use Fade_2D::setNumCPU() before using this method.

◆ isConstraint() [1/3]

bool GEOM_FADE2D::Fade_2D::isConstraint ( Point2 p0,
Point2 p1 
) const

This method checks whether the edge formed by the points p0 and p1 is a constraint edge.

Parameters
p0,p1Pointers to endpoints of the edge
Returns
A boolean value indicating whether the edge (p0, p1) is a constraint edge (true) or not (false).

◆ isConstraint() [2/3]

bool GEOM_FADE2D::Fade_2D::isConstraint ( Point2 pVtx) const

This method checks whether the vertex pVtx belongs to a constraint edge.

Parameters
pVtxPointer to the vertex to check.
Returns
A boolean value indicating whether the vertex pVtx belongs to a constraint edge (true) or not (false).

◆ isConstraint() [3/3]

bool GEOM_FADE2D::Fade_2D::isConstraint ( Triangle2 pT,
int  ith 
) const

This method checks whether the edge opposite to the ith vertex in the triangle pT is a constraint edge.

Parameters
pTThe triangle to check.
ithThe index of the vertex in the triangle.
Returns
A boolean value indicating whether the edge opposite to the ith vertex is a constraint edge (true) or not (false).

◆ load() [1/2]

bool GEOM_FADE2D::Fade_2D::load ( const char *  filename,
std::vector< Zone2 * > &  vZones 
)

This method loads a triangulation from a binary file, including any custom indices, constraint edges, and zones. The file must have been saved using the saveTriangulation(), saveZones(), or Zone2::save() methods.

Parameters
[in]filenameThe name of the input file.
[out]vZonesA vector of Zone2* that will be populated with pointers to the zones (if any) stored in the file. The order of the pointers will match the order at the time of storage.
Returns
true if the operation was successful, false otherwise.

◆ load() [2/2]

bool GEOM_FADE2D::Fade_2D::load ( std::istream &  stream,
std::vector< Zone2 * > &  vZones 
)

This method loads a triangulation from a binary file, including any custom indices, constraint edges, and zones. The file must have been saved using the saveTriangulation(), saveZones(), or Zone2::save() methods.

Parameters
[in]streamis an input stream
[out]vZonesA vector of Zone2* that will be populated with pointers to the zones (if any) stored in the file. The order of the pointers will match the order at the time of storage.
Returns
true if the operation was successful, false otherwise.

◆ locate()

Triangle2* GEOM_FADE2D::Fade_2D::locate ( const Point2 p)
Point location

This method searches for and returns a pointer to the Triangle2 that contains the query point p.

Parameters
[in]pThe query point whose containing Triangle2 is to be found.
Returns
A pointer to a Triangle2 containing p, or NULL if p lies outside the triangulation.

◆ measureTriangulationTime()

double GEOM_FADE2D::Fade_2D::measureTriangulationTime ( const std::vector< Point2 > &  vPoints)

This method measures the time taken to create a Delaunay triangulation from a given set of points.

Benchmark of single- and multithreaded Delaunay triangulation (Intel i9 10980XE, Linux)
Parameters
[in]vPointsA vector of points to be inserted into the triangulation.
Returns
The total wall time for point insertion in seconds.
Note
To enable multithreading, use setGlobalNumCPU().

◆ numberOfPoints()

size_t GEOM_FADE2D::Fade_2D::numberOfPoints ( ) const

Computes the number of points in the triangulation.

Returns
the number of points in the triangulation.
Note
Due to possibly duplicate input points the number of points is not stored somewhere but freshly computed in O(n) time. This is fast but you are adviced to avoid calling this method over-frequently in a loop. Duplicate point insertions count only once.

◆ numberOfTriangles()

size_t GEOM_FADE2D::Fade_2D::numberOfTriangles ( ) const

Computes the number of points in the triangulation

Returns
the number of triangles in the triangulation.

◆ refine()

void GEOM_FADE2D::Fade_2D::refine ( Zone2 pZone,
double  minAngleDegree,
double  minEdgeLength,
double  maxEdgeLength,
bool  bAllowConstraintSplitting 
)

Creates a mesh inside the area given by a Zone2 object.

Parameters
pZoneis the zone whose triangles are refined. Allowed zoneLocation values are ZL_INSIDE and ZL_BOUNDED.
minAngleDegree(up to 30) is the minimum interior triangle angle
minEdgeLengthis a lower threshold on the edge length. Triangles with smaller edges are not refined.
maxEdgeLengthis an upper threshold on the edge length. Triangles with larger edges are always refined.
bAllowConstraintSplittingspecifies if constraint edges may be splitted
Note
You can convert any Zone2 into a bounded zone using Zone2::convertToBoundedZone().
See also
void refineAdvanced() This method allows more control over the refinement process

◆ refineAdvanced()

void GEOM_FADE2D::Fade_2D::refineAdvanced ( MeshGenParams pParameters)

This method calls an advanced Delaunay mesh generator and grid mesher. The parameters are encapsulated in the MeshGenParams class, which provides default parameters that can be used as is. Alternatively, client code can derive from MeshGenParams and override the methods and parameters to gain full control over the mesh generation process.

See also
https://www.geom.at/advanced-mesh-generation/ This article provides examples on advanced mesh generation.

Control over the local mesh density in a Delaunay Mesh

Control over the local mesh density

Delaunay meshing with maximum edge length

Delaunay meshing with maximum edge length

Axis-aligned grid mesh

Grid mesh, axis-aligned

Grid mesh aligned to a certain direction

Grid mesh, aligned to a certain direction

◆ remove() [1/2]

void GEOM_FADE2D::Fade_2D::remove ( Point2 pVertex)
Parameters
pVertexThe vertex to be removed.
Note
pVertex must not belong to a constraint edge. Otherwise the vertex will not be removed, and a warning will be issued.

◆ remove() [2/2]

void GEOM_FADE2D::Fade_2D::remove ( std::vector< Point2 * > &  vPoints)
Parameters
vPointsare the points to be removed
Note
vPoints must not be vertices of any constraint edge. Otherwise they will not be removed.

◆ replicateZoneBoundary()

Zone2* GEOM_FADE2D::Fade_2D::replicateZoneBoundary ( Zone2 pZoneIn)

This method replicates only the boundaries of the given Zone2 pZoneIn in the current Fade_2D instance. It is primarily intended to facilitate boolean operations with zones, as these operations require the zones to belong to the same Fade_2D instance.

Parameters
pZoneInThe input zone whose boundaries are to be replicated.
Returns
A new zone, constrained by the replicated boundaries.

◆ reset()

void GEOM_FADE2D::Fade_2D::reset ( )

Resets the triangulation, making it reusable.

◆ saveTriangulation() [1/2]

bool GEOM_FADE2D::Fade_2D::saveTriangulation ( const char *  filename,
std::vector< Zone2 * > &  vSaveZones 
)

The saveTriangulation() method saves all triangles of the current triangulation to a binary file. It retains constraint edges and custom vertex indices, if available. If Zone2 pointers are provided, these zones will also be saved, and their order will be preserved.

Parameters
[in]filenameThe name of the output file.
[in]vSaveZonesA vector of Zone2 pointers representing zones to be saved.
See also
To store zones alone, use Zone2::save() or Fade_2D::saveZones(). To reload data from saved files, use Fade_2D::load().
Returns
true if the operation was successful, false otherwise.

◆ saveTriangulation() [2/2]

bool GEOM_FADE2D::Fade_2D::saveTriangulation ( std::ostream &  stream,
std::vector< Zone2 * > &  vSaveZones 
)

The saveTriangulation() method saves all triangles of the current triangulation to a binary file. It retains constraint edges and custom vertex indices, if available. If Zone2* pointers are provided, these zones will also be saved, and their order will be preserved.

Parameters
streamis the output stream
[in]vSaveZonesA vector of Zone2 pointers representing zones to be saved.
See also
To store zones alone, use Zone2::save() or Fade_2D::saveZones(). To reload data from saved files, use Fade_2D::load().
Returns
true if the operation was successful, false otherwise.

◆ saveZones() [1/2]

bool GEOM_FADE2D::Fade_2D::saveZones ( const char *  filename,
std::vector< Zone2 * > &  vSaveZones 
)

The saveZones() method saves the triangles of the zones in vSaveZones to a binary file. It retains the order of the zones, and maintains the constraint edges and custom vertex indices within the domain.

Note
A Delaunay triangulation is convex without holes, but this may not be the case for the zones being saved. As a result, additional triangles may be added to fill concavities. These extra triangles will belong to the Fade_2D instance, when reloaded, but will not be associated with any reloaded Zone2.
Parameters
[in]filenameThe name of the output file
[in]vSaveZonesA non-empty vector of Zone2 pointers specifying the zones to be saved.
Returns
true if the operation was successful, false otherwise.
See also
Use the saveTriangulation() method to store all triangles of a triangulation along with any specified zones. The Zone2::save() method can be used to save a single zone. To reload data from a saved file, use Fade_2D::load().

◆ saveZones() [2/2]

bool GEOM_FADE2D::Fade_2D::saveZones ( std::ostream &  stream,
std::vector< Zone2 * > &  vSaveZones 
)

The saveZones() method saves the triangles of the zones in vSaveZones to a binary file. It retains the order of the zones, and maintains the constraint edges and custom vertex indices within the domain.

Note
A Delaunay triangulation is convex without holes, but this may not be the case for the zones being saved. As a result, additional triangles may be added to fill concavities. These extra triangles will belong to the Fade_2D instance, when reloaded, but will not be associated with any reloaded Zone2.
Parameters
[in]streamis the name of output stream
[in]vSaveZonesA non-empty vector of Zone2 pointers specifying the zones to be saved.
Returns
true if the operation was successful, false otherwise.
See also
Use the saveTriangulation() method to store all triangles of a triangulation along with any specified zones. The Zone2::save() method can be used to save a single zone. To reload data from a saved file, use Fade_2D::load().

◆ setFastMode()

void GEOM_FADE2D::Fade_2D::setFastMode ( bool  bFast)

By default, multiple-precision arithmetic may be used (if required) to compute a 100% perfect Delaunay triangulation. However, the quality difference compared to using double precision arithmetic is hardly noticeable. Use this method to skip computationally expensive calculations and improve performance.

Parameters
bFastSet to true to avoid using computationally expensive multiple-precision arithmetic in favor of better performance.

Depending on the position of the input points, enabling setFastMode(true) may have no effect if multiple-precision arithmetic wouldn't be used anyway, or provide considerable acceleration if it would otherwise be needed.

◆ setNumCPU()

int GEOM_FADE2D::Fade_2D::setNumCPU ( int  numCPU)
Deprecated:
Use setGlobalNumCPU() instead. This method is kept for backward compatibility. It internally forwards the call to setGlobalNumCPU().
Parameters
[in]numCPUThe number of CPU cores to use for multithreading.

◆ show() [1/2]

void GEOM_FADE2D::Fade_2D::show ( const char *  filename,
bool  bWithConstraints = true 
) const

show() is a convenience function for quickly generating outputs with a default appearance. It is also possible to use the Visualizer2 class directly to draw arbitrary circles, line segments, vertices, and labels with custom colors.

Parameters
filenameThe output file name, e.g., "file.ps" for PostScript or "file.pdf" for PDF.
bWithConstraintsSpecifies whether constraint segments should be shown (default: true).
See also
https://www.geom.at/example2-traversing/ for examples of how to use the Visualizer2 class.

◆ show() [2/2]

void GEOM_FADE2D::Fade_2D::show ( Visualizer2 pVis,
bool  bWithConstraints = true 
) const

show() is a convenience function for quickly generating outputs with a default appearance. It is also possible to use the Visualizer2 class directly to draw arbitrary circles, line segments, vertices, and labels with custom colors.

Parameters
pVisA pointer to a Visualizer2 object, which may already contain geometric elements or may later be used to draw additional elements.
bWithConstraintsSpecifies whether constraint segments should be shown (default: true).
Note
After adding the last elements, finalize the visualization by calling pVis->writeFile().
See also
https://www.geom.at/example2-traversing/ for examples on how to use the Visualizer2 class.

◆ showVtk() [1/2]

void GEOM_FADE2D::Fade_2D::showVtk ( const char *  filename,
VtkColor  color,
VtkColor  constraintColor = VTK_TRANSPARENT 
) const

This method generates a .VTK file representing the triangulation, which can be viewed with tools like Paraview.

Parameters
filenameThe output filename, e.g., "file.vtk".
colorThe color for the triangles
constraintColorThe color for the constraint edges. VTK_TRANSPARENT hides them.
Note
It is also possible to use the VtkWriter class directly to visualize custom points, line segments, and triangles.

◆ showVtk() [2/2]

void GEOM_FADE2D::Fade_2D::showVtk ( VtkWriter pVtkWriter,
VtkColor  color,
VtkColor  constraintColor = VTK_TRANSPARENT 
) const

This method generates a .VTK file representing the triangulation, which can be viewed with tools like Paraview.

Parameters
pVtkWriterA pointer to a VtkWriter object, which may already contain geometric elements or can be used to add further elements.
colorThe color for the triangles.
constraintColorThe color for the constraint edges. Use VTK_TRANSPARENT to hide them.
Note
Once all elements have been drawn, you need to call pVtkWriter->writeFile() to save the file.
It is also possible to use the VtkWriter class directly to visualize custom points, line segments, and triangles.

◆ statistics()

void GEOM_FADE2D::Fade_2D::statistics ( const char *  s) const

Prints mesh quality statistics.

◆ subscribe()

void GEOM_FADE2D::Fade_2D::subscribe ( MsgType  msgType,
MsgBase pMsg 
)

This method registers a subscriber to receive messages of a specified type.

Parameters
msgTypeThe type of message the subscriber will receive, e.g., MSG_PROGRESS or MSG_WARNING.
pMsgA pointer to a custom class derived from MsgBase that will handle the messages.

◆ unsubscribe()

void GEOM_FADE2D::Fade_2D::unsubscribe ( MsgType  msgType,
MsgBase pMsg 
)

This method removes a subscriber from receiving messages of a specified type.

Parameters
msgTypeThe type of message that the subscriber will no longer receive.
pMsgA pointer to a custom class derived from MsgBase representing the subscriber to be removed.

◆ writeObj() [1/2]

void GEOM_FADE2D::Fade_2D::writeObj ( const char *  filename) const

This method writes the current triangulation data to an .obj file, which can be used for visualization and further processing in tools like Meshlab.

Parameters
filenameThe output filename, e.g., "triangulation.obj".

◆ writeObj() [2/2]

void GEOM_FADE2D::Fade_2D::writeObj ( const char *  filename,
Zone2 pZone 
) const

This method writes a Zone2 to an .obj file, which can be used for visualization and further processing in tools like Meshlab.

Parameters
filenameThe output filename, e.g., "triangulation.obj".
pZoneA pointer to the pZone to be written.

◆ writeWebScene() [1/2]

void GEOM_FADE2D::Fade_2D::writeWebScene ( const char *  path) const

This method generates a scene of the triangulation, which can be viewed in a web browser.

Parameters
pathThe output filename or file path where the scene will be saved.

◆ writeWebScene() [2/2]

void GEOM_FADE2D::Fade_2D::writeWebScene ( const char *  path,
Zone2 pZone 
) const

This method generates a scene of a Zone2 object, which can be viewed in a web browser.

Parameters
pathThe output filename or file path where the scene will be saved.
pZoneA pointer to the Zone2 object to be visualized.

The documentation for this class was generated from the following file: