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

‍** 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 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.
 
ConstraintSegment2getNearbyBoundaryConstraint (Point2 &p, double tolerance)
 Locates the nearest boundary ConstraintSegment2 of the zone within a specified distance. More...
 
size_t getNumberOfTriangles () const
 Get the number of triangles. More...
 
void getOffsetBoundary (double offset, std::vector< Segment2 > &vOffsetBoundary, double mergeAngleDeg=10.0, double angleStepDeg=20.0)
 Get the boundary of an offset shape. 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...
 
bool hasOn (const Point2 &p)
 Checks if the given point lies on the zone or its boundary. More...
 
bool hasOnBoundary (const Point2 &p)
 Checks if the given point lies on the boundary of the zone. More...
 
Triangle2locate (const Point2 &p)
 Locates the triangle containing the given point within the zone. 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...
 
bool shiftToZone (const Point2 &from, double tolerance, Point2 &result)
 Finds a point close to the input point that lies inside the zone. More...
 
void show (const char *filename, bool bShowFull, bool bWithConstraints) const
 Postscript- and PDF-visualization. More...
 
void show (Visualizer2 *pVisualizer, bool bShowFull, bool bWithConstraints) const
 Postscript- and PDF-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 showVtk (const char *filename, VtkColor zoneColor, VtkColor nonZoneColor=VTK_TRANSPARENT, VtkColor constraintColor=VTK_TRANSPARENT) const
 VTK visualization. More...
 
void showVtk (VtkWriter *pVtk, VtkColor zoneColor, VtkColor nonZoneColor=VTK_TRANSPARENT, VtkColor constraintColor=VTK_TRANSPARENT) const
 VTK 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 smoothing2 (int numIterations, bool bWithXY, bool bWithConstraintZ)
 Smooths the vertices of the zone. 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)
 Computes the difference between two Zone2 objects. More...
 
Zone2zoneIntersection (Zone2 *pZone0, Zone2 *pZone1)
 Computes the intersection of two Zone2 objects. More...
 
Zone2zoneSymmetricDifference (Zone2 *pZone0, Zone2 *pZone1)
 Computes the symmetric difference between two Zone2 objects. More...
 
Zone2zoneUnion (Zone2 *pZone0, Zone2 *pZone1)
 Computes the union of two Zone2 objects. More...
 

Detailed Description

‍**

Connected component with boundary- and hole polygons

Zone2 is a certain defined area of a triangulation

See also
http://www.geom.at/example4-zones-defined-areas-in-triangulations/
http://www.geom.at/boolean-operations/
createZone in the Fade2D class

Member Function Documentation

◆ convertToBoundedZone()

Zone2* GEOM_FADE25D::Zone2::convertToBoundedZone ( )

Convert a zone to a bounded zone.

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

Export triangles from a zone.

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

Get 2.5D Area.

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

Get 2D Area.

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

Get border edges.

Returns
: the CCW oriented border edges of the zone

◆ getBoundarySegments()

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

Compute the boundary segments.

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

Get connected components and their boundary polygons.

This method subdivides the zone into connected components. For each connected component 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 the triangle belonging to the Zone2. Thus the outer boundary polygon is ccw-oriented while the polygons of inner holes are cw-oriented.

◆ getConstraintGraph()

ConstraintGraph2* GEOM_FADE25D::Zone2::getConstraintGraph ( ) const

Get the associated constraint.

Returns
a pointer to the ConstraintGraph2 object which defines the zone.
or NULL for ZL_RESULT-, ZL_GROW and ZL_GLOBAL_-zones.

◆ getNearbyBoundaryConstraint()

ConstraintSegment2* GEOM_FADE25D::Zone2::getNearbyBoundaryConstraint ( Point2 p,
double  tolerance 
)

Locates the nearest boundary ConstraintSegment2 of the zone within a specified distance.

Parameters
pThe point from which to search.
toleranceThe maximum allowed 2D distance between the point and the ConstraintSegment2.
Returns
A pointer to the closest alive ConstraintSegment2 within the tolerance distance, where p has an orthogonal projection, or NULL if no such segment is found.

The first call to Zone2::hasOn(), Zone2::hasOnBoundary(), Zone2::getNearbyBoundaryConstraint(), Zone2::shiftToZone() or Zone2::locate() initializes a search structure. For zones not of type ZL_INSIDE, ZL_BOUNDED, or ZL_GLOBAL, constraint edges are inserted around the zone.

◆ getNumberOfTriangles()

size_t GEOM_FADE25D::Zone2::getNumberOfTriangles ( ) const

Get the number of triangles.

Warning
This method is fast but O(n), so don't call it frequently in a loop.

◆ getOffsetBoundary()

void GEOM_FADE25D::Zone2::getOffsetBoundary ( double  offset,
std::vector< Segment2 > &  vOffsetBoundary,
double  mergeAngleDeg = 10.0,
double  angleStepDeg = 20.0 
)

Get the boundary of an offset shape.

This function computes a zone offset by a specified distance. The result is returned as a vector of oriented Segment2's.

Parameters
offsetThe positive or negative offset distance from the present zone.
[out]vOffsetBoundaryA vector that will contain the output segments, in arbitrary order, oriented counterclockwise around the resulting shape.
mergeAngleDegIf the angle between two offset points of an original vertex is less than this value (in degrees), the points will be merged. Default: 10.0, valid range: greater than 0 and up to 135.
angleStepDegSpecifies the angle interval (in degrees) at which circular arcs are sampled using line segments. Default: 20.0, valid range: greater than 0 and up to 135.

◆ getTriangles()

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

Get the triangles of the zone.

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

Get the zone location.

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

◆ hasOn()

bool GEOM_FADE25D::Zone2::hasOn ( const Point2 p)

Checks if the given point lies on the zone or its boundary.

Parameters
pThe point to check.
Returns
true if the point is on the zone or its boundary, false otherwise.

The first call to Zone2::hasOn(), Zone2::hasOnBoundary(), Zone2::getNearbyBoundaryConstraint(), Zone2::shiftToZone() or Zone2::locate() initializes a search structure. For zones not of type ZL_INSIDE, ZL_BOUNDED, or ZL_GLOBAL, constraint edges are inserted around the zone.

◆ hasOnBoundary()

bool GEOM_FADE25D::Zone2::hasOnBoundary ( const Point2 p)

Checks if the given point lies on the boundary of the zone.

Parameters
pThe point to check.
Returns
true if the point is on the boundary, false otherwise.

The first call to Zone2::hasOn(), Zone2::hasOnBoundary(), Zone2::getNearbyBoundaryConstraint(), Zone2::shiftToZone() or Zone2::locate() initializes a search structure. For zones not of type ZL_INSIDE, ZL_BOUNDED, or ZL_GLOBAL, constraint edges are inserted around the zone.

◆ locate()

Triangle2* GEOM_FADE25D::Zone2::locate ( const Point2 p)

Locates the triangle containing the given point within the zone.

Parameters
pThe point to locate.
Returns
A pointer to the Triangle2 if the point is inside the zone, otherwise NULL.

The first call to Zone2::hasOn(), Zone2::hasOnBoundary(), Zone2::getNearbyBoundaryConstraint(), Zone2::shiftToZone() or Zone2::locate() initializes a search structure. For zones not of type ZL_INSIDE, ZL_BOUNDED, or ZL_GLOBAL, constraint edges are inserted around the zone.

◆ numberOfConstraintGraphs()

size_t GEOM_FADE25D::Zone2::numberOfConstraintGraphs ( ) const

Get a the number of ConstraintGraph2 objects.

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

◆ save() [1/2]

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

Save the zone.

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)

Save the zone.

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.

◆ shiftToZone()

bool GEOM_FADE25D::Zone2::shiftToZone ( const Point2 from,
double  tolerance,
Point2 result 
)

Finds a point close to the input point that lies inside the zone.

Parameters
fromA point outside the zone.
toleranceThe maximum allowed 2D distance to find a point within the zone.
resultA reference to store the point found within the zone, close to from.
Returns
true if a point within the zone is successfully found or created; false otherwise.

The first call to Zone2::hasOn(), Zone2::hasOnBoundary(), Zone2::getNearbyBoundaryConstraint(), Zone2::shiftToZone() or Zone2::locate() initializes a search structure. For zones not of type ZL_INSIDE, ZL_BOUNDED, or ZL_GLOBAL, constraint edges are inserted around the zone.

◆ show() [1/2]

void GEOM_FADE25D::Zone2::show ( const char *  filename,
bool  bShowFull,
bool  bWithConstraints 
) const

Postscript- and PDF-visualization.

Parameters
filenameis 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

Postscript- and PDF-visualization.

Parameters
pVisualizeris a pointer to an existing Visualizer2 object drawing a .ps or .pdf file
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

Geomview visualization.

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

Geomview visualization.

Parameters
pVispoints to a Visualizer3 object
coloris a string ("red green blue alpha"), e.g., "1.0 0.0 0.0 1.0"*

◆ showVtk() [1/2]

void GEOM_FADE25D::Zone2::showVtk ( const char *  filename,
VtkColor  zoneColor,
VtkColor  nonZoneColor = VTK_TRANSPARENT,
VtkColor  constraintColor = VTK_TRANSPARENT 
) const

VTK visualization.

Parameters
filenameThe name of the output file.
zoneColorThe color for the zone's triangles
nonZoneColorThe color of the non-zone-triangles. Use VTK_TRANSPARENT to prevent them from being drawn.
constraintColorThe color of the constraint edges. Use VTK_TRANSPARENT to prevent them from being drawn.

◆ showVtk() [2/2]

void GEOM_FADE25D::Zone2::showVtk ( VtkWriter pVtk,
VtkColor  zoneColor,
VtkColor  nonZoneColor = VTK_TRANSPARENT,
VtkColor  constraintColor = VTK_TRANSPARENT 
) const

VTK visualization.

Parameters
pVtkA VtkWriter object that may already contain other geometric objects
zoneColorThe color for the zone's triangles
nonZoneColorThe color of the non-zone-triangles. Use VTK_TRANSPARENT to prevent them from being drawn.
constraintColorThe color of the constraint edges. Use VTK_TRANSPARENT to prevent them from being drawn.
Note
pVtk->writeFile() finally writes the .vtk file.

◆ slopeValleyRidgeOptimization()

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

Optimize Slopes, Valleys and Ridges.

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 
)

Smoothing.

Deprecated method: This method is deprecated and retained for backwards compatibility. It is recommended to use the new method smoothing2() instead for better results.

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.

◆ smoothing2()

void GEOM_FADE25D::Zone2::smoothing2 ( int  numIterations,
bool  bWithXY,
bool  bWithConstraintZ 
)

Smooths the vertices of the zone.

This function performs a weighted Laplacian smoothing on the vertices of the Zone2. We distinguish between two types of vertices:

  • Static vertices: Belong to constraint-edges or border-edges of the Zone2.
  • Dynamic vertices: All other vertices.
Parameters
numIterationsNumber of smoothing iterations to perform. Higher values yield smoother results; 2-3 passes are recommended.
bWithXYBoolean flag to determine if the x- and y-coordinates of dynamic vertices should also be adjusted.
bConstraintsWithZBoolean flag to control whether the z-coordinates of static vertices are also updated by the smoothing process.
Note
  • The x- and y-coordinates of static vertices are never changed.
  • The z-coordinates of dynamic vertices can always change, regardless of the value of bConstraintsWithZ.

◆ 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 
)

Register a message receiver.

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 
)

Unregister a message receiver.

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

Write the zone to *.obj Writes the triangles of the present Zone2 to an *.obj file (The *.obj format represents a 3D scene).

Parameters
outFilenameis the output filename

◆ writePly() [1/2]

bool GEOM_FADE25D::Zone2::writePly ( const char *  filename,
bool  bASCII = false 
) const

Write the zone to a *.ply file.

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

Write the zone to a *.ply file.

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

Peel off border triangles.

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 if no triangles remain.

◆ peelOffIf [2/2]

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

Peel off border triangles (deprecated)

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 if no triangles remain.

◆ zoneDifference

Zone2* zoneDifference ( Zone2 pZone0,
Zone2 pZone1 
)
friend

Computes the difference between two Zone2 objects.

This function returns a pointer to a new Zone2 object representing the difference between the two input zones. The difference is the area covered by the first zone that is not covered by the second zone. The input zones must belong to the same Fade_2D object.

Parameters
pZone0Pointer to the first Zone2 object.
pZone1Pointer to the second Zone2 object.
Returns
Zone2* A pointer to the resulting Zone2 object representing the difference. If the result is empty, the function returns NULL.

◆ zoneIntersection

Zone2* zoneIntersection ( Zone2 pZone0,
Zone2 pZone1 
)
friend

Computes the intersection of two Zone2 objects.

This function returns a pointer to a new Zone2 object representing the intersection of the two input zones. The intersection is the area covered by both of the input zones. The input zones must belong to the same Fade_2D object.

Parameters
pZone0Pointer to the first Zone2 object.
pZone1Pointer to the second Zone2 object.
Returns
Zone2* A pointer to the resulting Zone2 object representing the intersection. If the result is empty, the function returns NULL.

◆ zoneSymmetricDifference

Zone2* zoneSymmetricDifference ( Zone2 pZone0,
Zone2 pZone1 
)
friend

Computes the symmetric difference between two Zone2 objects.

This function returns a pointer to a new Zone2 object representing the symmetric difference between the two input zones. The symmetric difference is the area covered by either of the input zones but not by both. The input zones must belong to the same Fade_2D object.

Parameters
pZone0Pointer to the first Zone2 object.
pZone1Pointer to the second Zone2 object.
Returns
Zone2* A pointer to the resulting Zone2 object representing the symmetric difference. If the result is empty, the function returns NULL.

◆ zoneUnion

Zone2* zoneUnion ( Zone2 pZone0,
Zone2 pZone1 
)
friend

Computes the union of two Zone2 objects.

This function returns a pointer to a new Zone2 object representing the union of the two input zones. The union is the area covered by either of the input zones. The input zones must belong to the same Fade_2D object.

Parameters
pZone0Pointer to the first Zone2 object.
pZone1Pointer to the second Zone2 object.
Returns
Zone2* A pointer to the resulting Zone2 object representing the union.

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