Fade2.5D Documentation pages v2.12
Delaunay Features
GEOM_FADE25D::Zone2 Class Reference

Zone2 is a certain defined area of a triangulation. More...

#include <Zone2.h>

Public Member Functions

Zone2convertToBoundedZone ()
 Convert a zone to a bounded zone. More...
 
void debug (const char *name="")
 Development function.
 
void exportZone (FadeExport &fadeExport, bool bWithCustomIndices) const
 Export triangles from a zone. More...
 
double getArea25D () const
 Get 2.5D Area. More...
 
double getArea2D () const
 Get 2D Area. More...
 
void getBorderEdges (std::vector< Edge2 > &vBorderEdgesOut) const
 Get border edges. More...
 
void getBoundaryEdges (std::vector< Edge2 > &vEdges) const
 Compute the boundary edges. More...
 
void getBoundarySegments (std::vector< Segment2 > &vSegments) const
 Compute the boundary segments. More...
 
Bbox2 getBoundingBox () const
 Compute the bounding box.
 
void getComponentPolygons (std::vector< CompPolygon > &vCompPolygons) const
 Get connected components and their boundary polygons. More...
 
ConstraintGraph2getConstraintGraph () const
 Get the associated constraint. More...
 
void getConstraintGraphs (std::vector< ConstraintGraph2 * > &vConstraintGraphs_) const
 Get the associated constraint graphs.
 
size_t getNumberOfTriangles () const
 Get the number of triangles. More...
 
void getTriangles (std::vector< Triangle2 * > &vTriangles_) const
 Get the triangles of the zone. More...
 
void getVertices (std::vector< Point2 * > &vVertices_) const
 Get the vertices of the zone.
 
ZoneLocation getZoneLocation () const
 Get the zone location. More...
 
size_t numberOfConstraintGraphs () const
 Get a the number of ConstraintGraph2 objects. More...
 
bool save (const char *filename)
 Save the zone. More...
 
bool save (std::ostream &stream)
 Save the zone. More...
 
void show (const char *postscriptFilename, bool bShowFull, bool bWithConstraints) const
 Postscript visualization. More...
 
void show (Visualizer2 *pVisualizer, bool bShowFull, bool bWithConstraints) const
 Postscript visualization. More...
 
void showGeomview (const char *filename, const char *color) const
 Geomview visualization. More...
 
void showGeomview (Visualizer3 *pVis, const char *color) const
 Geomview visualization. More...
 
void slopeValleyRidgeOptimization (OptimizationMode om=OPTMODE_BETTER)
 Optimize Slopes, Valleys and Ridges. More...
 
void smoothing (int numIterations=2, bool bWithXY=true)
 Smoothing. More...
 
void statistics (const char *s) const
 
void subscribe (MsgType msgType, MsgBase *pMsg)
 Register a message receiver. More...
 
void unifyGrid (double tolerance)
 
void unsubscribe (MsgType msgType, MsgBase *pMsg)
 Unregister a message receiver. More...
 
void writeObj (const char *outFilename) const
 Write the zone to *.obj Writes the triangles of the present Zone2 to an *.obj file (The *.obj format represents a 3D scene). More...
 
bool writePly (const char *filename, bool bASCII=false) const
 Write the zone to a *.ply file. More...
 
bool writePly (std::ostream &os, bool bASCII=false) const
 Write the zone to a *.ply file. More...
 

Protected Member Functions

Zone2operator= (const Zone2 &)
 

Friends

Zone2peelOffIf (Zone2 *pZone, bool bAvoidSplit, PeelPredicateTS *pPredicate)
 Peel off border triangles. More...
 
Zone2peelOffIf (Zone2 *pZone, UserPredicateT *pPredicate, bool bVerbose)
 Peel off border triangles (deprecated) More...
 
Zone2zoneDifference (Zone2 *pZone0, Zone2 *pZone1)
 Compute the difference of two zones. More...
 
Zone2zoneIntersection (Zone2 *pZone0, Zone2 *pZone1)
 Compute the intersection of two zones. More...
 
Zone2zoneSymmetricDifference (Zone2 *pZone0, Zone2 *pZone1)
 Compute the symmetric difference of two zones. More...
 
Zone2zoneUnion (Zone2 *pZone0, Zone2 *pZone1)
 Compute the union of two zones. More...
 

Detailed Description

Member Function Documentation

◆ convertToBoundedZone()

Zone2* GEOM_FADE25D::Zone2::convertToBoundedZone ( )

The mesh generation algorithms refine() and refineAdvanced() require a zone object that is bounded by constraint segments. This is always the case for zones with zoneLocation ZL_INSIDE but other types of zones may be unbounded. For convenience this method is provided to create a bounded zone from a possibly unbounded one.

Returns
a pointer to a new Zone2 object with zoneLocation ZL_RESULT_BOUNDED or a pointer to the present zone if this->getZoneLocation() is ZL_INSIDE.

◆ exportZone()

void GEOM_FADE25D::Zone2::exportZone ( FadeExport fadeExport,
bool  bWithCustomIndices 
) const
Parameters
fadeExportis a struct that will hold the requested triangulation data
bWithCustomIndicesdetermines whether the custom indices of the points are also stored

◆ getArea25D()

double GEOM_FADE25D::Zone2::getArea25D ( ) const

Returns the 2.5D area of the zone.

Note: The getArea() method is deprecated and replaced by getArea2D() and getArea25D()

◆ getArea2D()

double GEOM_FADE25D::Zone2::getArea2D ( ) const

Returns the 2D area of the zone.

Note: The getArea() method is deprecated and replaced by getArea2D() and getArea25D()

◆ getBorderEdges()

void GEOM_FADE25D::Zone2::getBorderEdges ( std::vector< Edge2 > &  vBorderEdgesOut) const
Returns
: the CCW oriented border edges of the zone

◆ getBoundaryEdges()

void GEOM_FADE25D::Zone2::getBoundaryEdges ( std::vector< Edge2 > &  vEdges) const

Outputs the boundary edges of the zone. Edge2 elements are always ccw-oriented but the edges are not returned in any specific order.

◆ getBoundarySegments()

void GEOM_FADE25D::Zone2::getBoundarySegments ( std::vector< Segment2 > &  vSegments) const

Outputs the boundary segments of the zone. These are ccw-oriented but not returned in any specific order.

◆ getComponentPolygons()

void GEOM_FADE25D::Zone2::getComponentPolygons ( std::vector< CompPolygon > &  vCompPolygons) const

This method subdivides the zone into connected components. For each connected components it then returns a CompPolygon object consisting of the triangles, their outer boundary polygon and the hole polygons. Edges are represented by a triangle and an index and they are always counterclockwise (ccw) around their triangle. Thus the outer boundary polygon is also ccw-oriented while the polygons of inner holes are cw-oriented.

◆ getConstraintGraph()

ConstraintGraph2* GEOM_FADE25D::Zone2::getConstraintGraph ( ) const
Returns
a pointer to the ConstraintGraph2 object which defines the zone.
or NULL for ZL_RESULT-, ZL_GROW and ZL_GLOBAL_-zones.

◆ getNumberOfTriangles()

size_t GEOM_FADE25D::Zone2::getNumberOfTriangles ( ) const
Warning
This method is fast but O(n), so don't call it frequently in a loop.

◆ getTriangles()

void GEOM_FADE25D::Zone2::getTriangles ( std::vector< Triangle2 * > &  vTriangles_) const

This command fetches the existing triangles of the zone.

Note
Fade_2D::void applyConstraintsAndZones() must be called after the last insertion of points and constraints.
that the lifetime of data from the Fade2D datastructures does exceed the lifetime of the Fade2D object.

◆ getZoneLocation()

ZoneLocation GEOM_FADE25D::Zone2::getZoneLocation ( ) const
Returns
ZL_INSIDE if the zone applies to the triangles inside one or more ConstraintGraph2 objects
ZL_OUTSIDE if the zone applies to the outside triangles
ZL_GLOBAL if the zone applies (dynamically) to all triangles
ZL_RESULT if the zone is the result of a set operation
ZL_GROW if the zone is specified by a set of constraint graphs and an inner point
An ouside zone and in inside zone

◆ numberOfConstraintGraphs()

size_t GEOM_FADE25D::Zone2::numberOfConstraintGraphs ( ) const

A Zone2 object might be defined by zero, one or more ConstraintGraph2 objects.

◆ save() [1/2]

bool GEOM_FADE25D::Zone2::save ( const char *  filename)

This command saves the present Zone2 to a binary file. Any constraint edges and custom indices in the domain are retained.

Parameters
[in]filenameis the output filename
Returns
whether the operation was successful
Note
A Delaunay triangulation is convex without holes but this may not hold for the zone to be saved. Thus extra triangles may be saved to fill concavities. These extra-triangles will belong to the Fade_2D instance but not to the Zone2 object when reloaded.
See also
save(std::ostream& stream). Use the similar command Fade_2D::saveZones(const char* file, std::vector<Zone2*>& vZones) to store more than just one zone. Use Fade_2D::saveTriangulation() to store all triangles of the triangulation plus any specified zones. Use Fade_2D::load() to reload the data from such files.

◆ save() [2/2]

bool GEOM_FADE25D::Zone2::save ( std::ostream &  stream)

This command saves the present Zone2 to an ostream. Any constraint edges and custom indices in the domain are retained.

Parameters
streamis the output stream
Returns
whether the operation was successful
Note
A Delaunay triangulation is convex without holes but this may not hold for the zone to be saved. Thus extra triangles may be saved to fill concavities. These extra-triangles will belong to the Fade_2D instance but not to the Zone2 object when reloaded.
See also
Use the similar command Fade_2D::saveZones(const char* file, std::vector<Zone2*>& vZones) to store more than just one zone. Use Fade_2D::saveTriangulation() to store all triangles of the triangulation plus any specified zones. Use Fade_2D::load() to reload the data from such files.

◆ show() [1/2]

void GEOM_FADE25D::Zone2::show ( const char *  postscriptFilename,
bool  bShowFull,
bool  bWithConstraints 
) const
Parameters
postscriptFilenameis the name of the output file.
bShowFullspecifies if only the zone or the full triangulation shall be drawn
bWithConstraintsspecifies if constraint edges shall be drawn

◆ show() [2/2]

void GEOM_FADE25D::Zone2::show ( Visualizer2 pVisualizer,
bool  bShowFull,
bool  bWithConstraints 
) const
Parameters
pVisualizeris a pointer to an existing Visualizer2 object.
Note
You must call pVisualizer->writeFile() before program end
Parameters
bShowFullspecifies if only the zone or the full triangulation shall be drawn
bWithConstraintsspecifies if constraint edges shall be drawn

◆ showGeomview() [1/2]

void GEOM_FADE25D::Zone2::showGeomview ( const char *  filename,
const char *  color 
) const
Parameters
filenameis the name of the output file.
coloris a string ("red green blue alpha"), e.g., "1.0 0.0 0.0 1.0"*

◆ showGeomview() [2/2]

void GEOM_FADE25D::Zone2::showGeomview ( Visualizer3 pVis,
const char *  color 
) const
Parameters
pVispoints to a Visualizer3 object
coloris a string ("red green blue alpha"), e.g., "1.0 0.0 0.0 1.0"*

◆ slopeValleyRidgeOptimization()

void GEOM_FADE25D::Zone2::slopeValleyRidgeOptimization ( OptimizationMode  om = OPTMODE_BETTER)

A pure Delaunay triangulation takes only the x and y coordinates into account. However, for terrain scans, it is important to consider the z coordinate as well, otherwise ridges, valleys and rivers will look unnatural. This method leaves the points constant, but uses edge flips to change the connectivity, making the surface smoother overall.

Parameters
omis the optimization mode: OPTMODE_NORMAL is the fastest. OPTMODE_BETTER provides significantly better results while still taking a moderate amount of time. OPTMODE_BEST delivers the best results, but also has a significantly higher time requirement. This method supports the progress-bar mechanism.
Note
Flipping edges makes the triangulation non-delaunay, i.e. the empty-circle-property is then no longer given. Improving the smoothness of the surface by edge flips also means degrading the interior angles of the triangles (to a certain degree).

◆ smoothing()

void GEOM_FADE25D::Zone2::smoothing ( int  numIterations = 2,
bool  bWithXY = true 
)

Weighted laplacian smoothing for the z-coordinate of all zone vertices. The x and y coordinates can also be optimized but only for vertices not belonging to the boundary of the zone or to constraint edges. This method is very fast but does nevertheless support the progress bar mechanism.

Parameters
numIterationsis the number of smoothing passes.
bWithXYspecifies if the x and y coordinates are also adapted.

◆ statistics()

void GEOM_FADE25D::Zone2::statistics ( const char *  s) const

Statistics

Prints statistics to stdout.

◆ subscribe()

void GEOM_FADE25D::Zone2::subscribe ( MsgType  msgType,
MsgBase pMsg 
)
Parameters
msgTypeis the type of message the subscriber shall receive, e.g. MSG_PROGRESS or MSG_WARNING
pMsgis a pointer to a custom class derived from MsgBase

◆ unifyGrid()

void GEOM_FADE25D::Zone2::unifyGrid ( double  tolerance)

Unify Grid

A Delaunay triangulation is not unique when when 2 or more triangles lie on a common circle. As a consequence the four corners of a rectangle can be triangulated in two different ways: Either the diagonal proceeds from the lower left to the upper right corner or it connects the other two corners. Both solutions are valid and an arbitrary one is applied when points are triangulated. To improve the repeatability and for reasons of visual appearance this method unifies such diagonals to point from the lower left to the upper right corner (or in horizontal direction).

Parameters
toleranceis 0 when only exact cases of more than 3 points on a common circumcircle shall be changed. But in practice input data can be disturbed by noise and tiny rounding errors such that grid points are not exactly on a grid. The numeric error is computed as $error=\frac{abs(diagonalA-diagonalB)}{max(diagonalA,diagonalB)}$. and tolerance is an upper threshold to allow modification despite such tiny inaccuracies. Use with caution, such flips break the empty circle property and this may or may not fit your setting.

◆ unsubscribe()

void GEOM_FADE25D::Zone2::unsubscribe ( MsgType  msgType,
MsgBase pMsg 
)
Parameters
msgTypeis the type of message the subscriber shall not receive anymore
pMsgis a pointer to a custom class derived from MsgBase

◆ writeObj()

void GEOM_FADE25D::Zone2::writeObj ( const char *  outFilename) const
Parameters
outFilenameis the output filename

◆ writePly() [1/2]

bool GEOM_FADE25D::Zone2::writePly ( const char *  filename,
bool  bASCII = false 
) const
Parameters
filenameis the output filename
bASCIIspecifies whether to write the *.ply in ASCII or binary format
Note
Method available for platforms with C++11

◆ writePly() [2/2]

bool GEOM_FADE25D::Zone2::writePly ( std::ostream &  os,
bool  bASCII = false 
) const
Parameters
osis the output file
bASCIIspecifies whether to write the *.ply in ASCII or binary format
Note
Method available for platforms with C++11

Friends And Related Function Documentation

◆ peelOffIf [1/2]

Zone2* peelOffIf ( Zone2 pZone,
bool  bAvoidSplit,
PeelPredicateTS pPredicate 
)
friend
Parameters
pZoneis the input zone
bAvoidSplitif true, then the algorithm removes a triangle only if it does not break the zone into independent components.
pPredicateis a user-defined predicate that decides if a triangle shall be removed.
Returns
a new zone containing a subset of the triangles of pZone or NULL when no triangles remain.
Attention
Check whether NULL is returned!

◆ peelOffIf [2/2]

Zone2* peelOffIf ( Zone2 pZone,
UserPredicateT pPredicate,
bool  bVerbose 
)
friend

This function is DEPRECATED but kept for backwards compatibility. The new and better function is: peelOffIf(Zone2* pZone, bool bAvoidSplit,PeelPredicateTS* pPredicate)

Parameters
pZone
pPredicate
bVerbose
Returns
a new zone containing a subset of the triangles of pZone or NULL when no triangles remain.

◆ zoneDifference

Zone2* zoneDifference ( Zone2 pZone0,
Zone2 pZone1 
)
friend
Returns
a new zone containing the triangles of *pZone0 minus the ones of *pZone1
Note
pZone0 and pZone1 must belong to the same Fade_2D object.

◆ zoneIntersection

Zone2* zoneIntersection ( Zone2 pZone0,
Zone2 pZone1 
)
friend
Returns
a new zone containing the intersection of *pZone0 and *pZone1
Note
pZone0 and pZone1 must belong to the same Fade_2D object.

◆ zoneSymmetricDifference

Zone2* zoneSymmetricDifference ( Zone2 pZone0,
Zone2 pZone1 
)
friend
Returns
a new zone containing the triangles that are present in one of the zones but not in the other one.
Note
pZone0 and pZone1 must belong to the same Fade_2D object.

◆ zoneUnion

Zone2* zoneUnion ( Zone2 pZone0,
Zone2 pZone1 
)
friend
Returns
a new zone containing the union of the triangles of *pZone0 and *pZone1
Note
pZone0 and pZone1 must belong to the same Fade_2D object.

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