GranOO  3.0
A robust and versatile workbench to build 3D dynamic simulations based on the Discrete Element Method
Public Types | Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
PolyhedronCore Class Reference

#include <PolyhedronCore.hpp>

Public Types

typedef std::pair< unsigned int, unsigned int > VertexPrimitivePair
 

Public Member Functions

 PolyhedronCore ()
 instanciates a new PolyhedronCore object More...
 
 PolyhedronCore (std::string filePath)
 instanciates a PolyhedronCore object with given OFF file More...
 
 PolyhedronCore (std::stringstream &stringStream)
 instanciates a PolyhedronCore object with given OFF file contents, as std::stringstream. More...
 
 ~PolyhedronCore ()
 Default PolyhedronCore destructor. More...
 
void read_off_file (std::string filePath)
 Reads an OFF file given a file path. More...
 
void read_off_file (std::stringstream &stringStream)
 Reads an OFF file from a stream. More...
 
void to_off_format (std::stringstream &offString, bool basefacet=true) const
 Creates an std::stringstream with the polyhedron contents, in OFF file format. More...
 
void to_vtk_ugrid (vtkUnstructuredGrid &uGrid)
 Creates a VTK unstructured grid with the polyhedron contents. More...
 
void clear ()
 Clears the geometric entities of the PolyhedronCore object. More...
 
void clear_metric ()
 Clears the all computations performed on the polyhedron. More...
 
void build_data_structure ()
 Builds the data structure which can be used to perform queries, e.g. adjacency. More...
 
void compute_metric ()
 Computes geometric parameters related to the polyhedron. More...
 
void build_aabb_tree ()
 Builds a static AABB tree based on the polyhedron primitives. More...
 
unsigned int vertex_number () const
 
std::vector< double > & get_hitPoint ()
 
unsigned int basefacet_number () const
 
unsigned int facet_number () const
 
unsigned int edge_number () const
 
bool is_watertight () const
 
const double * min_corner_ptr () const
 
const double * max_corner_ptr () const
 
bool is_convex () const
 
double volume () const
 
const double * initial_centroid_ptr () const
 
const double * centroid_ptr () const
 
void inertia_tensor (double(&tensor)[3][3]) const
 update the tensor to the inertia operator of the polyhedron. More...
 
unsigned int vertex_number_of_basefacet (int basefacetIndex) const
 
const unsigned int vertex_index_of_basefacet (int basefacetIndex, int vertexIndex) const
 Allows to get a the global vertex index for a given base facet, with a given local vertex index. More...
 
const unsigned int vertex_index_of_facet (int facetIndex, int vertexIndex) const
 Allows to get a the global vertex index for a given facet, with a given local vertex index. More...
 
const double * vertex_ptr_of_facet (unsigned int facetIndex, unsigned int vectorIndex) const
 Allows to get the in-plane vectors of a given facet. More...
 
const double * normal_ptr_of_facet (int facetIndex) const
 Allows to get the normal of a given facet. More...
 
const double * vertex_ptr (int vertexIndex) const
 Allows to get a vertex coordinates array from its index in the _vertices vector. More...
 
const double * normal_ptr_of_vertex (int vertexIndex) const
 Allows to get a refernce to a vertex coordinates array from its index in the _vertices vector. More...
 
unsigned int basefacet_index_for_facet (unsigned int facetIndex) const
 Returns the base facet index for a triangle facet. More...
 
double basefacet_area (unsigned int facetIndex) const
 Computes the area of a base facet. More...
 
double facet_area (unsigned int facetIndex) const
 Computes the area of a triangle facet with index facetIndex. More...
 
void apply_offset (double offsetX, double offsetY, double offsetZ)
 Applies an offset to all vertices coordinates of the polyhedron. More...
 
void translate_to_centroid ()
 translates the polyhedron to its center of mass. More...
 
void apply_scale_factor (double factor, double yFactor, double zFactor)
 Applies a scale factor to all vertices coordinates of the polyhedron. More...
 
bool contain (double aPoint[3]) const
 Checks if the polyhedron contains a 3D point. More...
 
std::vector< unsigned int > facet_hit_by_ray (Ray aRay) const
 get_s the list of facets hit by a ray. More...
 
const double * closest_vertex_ptr (double aPoint[3]) const
 Finds the vertex in the polyhedron which is the closest to the query point. More...
 
VertexPrimitivePair closest_vertex_and_primitive (double aPoint[3]) const
 Finds the vertex and the primitive it is associated to in the polyhedron which are the closest to the query point. More...
 
void closest_point (double aPoint[3], double(&closestPoint)[3]) const
 Finds the point in the polyhedron which is the closest to the query point. More...
 
void furthest_vertex_along_axis (const double *anAxis, unsigned int &furthestVertex) const
 Finds the vertex in the polyhedron which is the furthest in the given direction. More...
 
void closest_vertex_callback (unsigned int aFacet, const double *aPoint, double &distance) const
 Finds the vertex in a polyhedron facets which is the closest to the requested point. More...
 
bool raycast_aabb_tree_callback (unsigned int aFacet, Ray aRay) const
 Checks if a primitive is hit by a ray. More...
 
void closest_vertex_and_primitive_callback (unsigned int aFacet, const double *aPoint, double &distance) const
 Finds the vertex in a polyhedron facet and the associated facet which realizes the shortest distance to the requested point. More...
 
void closest_point_callback (unsigned int aFacet, const double *aPoint, double &distance) const
 Finds the point in a polyhedron facets which is the closest to the requested point. More...
 
void furthest_vertex_callback (unsigned int aFacet, const double *aPoint, const double *anAxis, double &distance) const
 Finds the vertex in a polyhedron facets which is the furthest, in the requested direction. More...
 
bool is_facet_crossed_by_ray (unsigned int aFacet, Ray aRay, double(&hitPoint)[3]) const
 Checks if a facet is crossed by a given ray. More...
 
void closest_point_in_facet (unsigned int aFacet, const double *aPoint, double(&closestPoint)[3]) const
 Finds the closest point lying on a facet closest to aPoint. More...
 
const std::vector< unsigned int > & facet_to_vertex () const
 
std::vector< unsigned int > & facet_to_vertex ()
 

Private Member Functions

unsigned int find_local_edge_index_of_facet (unsigned int aFacet, unsigned int anEdge)
 Finds the local index of an edge of a given global facet. More...
 
double closest_point_to_edge (const double *point0, const double *point1, const double *, double(&closestPoint)[3]) const
 

Private Attributes

std::vector< double > _vertices
 The list of vertices coordinates, given as [x0,y0,z0,...,xn,yn,zn]. More...
 
std::vector< unsigned int > _basefacetToVertex
 The relationship between base facets (or raw facets, directly from the geometry file that is, which may have any number of vertices). This vector must be used along with the _basefacetToVertexOffset vector. More...
 
std::vector< unsigned int > _basefacetToVertexOffset
 stores the offset required to access the first vertex of a base facet in _basefacetToVertex vector. More...
 
std::vector< unsigned int > _verticesPerBaseFacet
 The number of vertices that form a base facet. More...
 
std::vector< unsigned int > _facetToVertex
 The relationship between triangle facets and vertices. More...
 
std::vector< unsigned int > _basefacetToFacetOffset
 stores the offset required to access the first triangle facet of a base facet in _facetToBaseFacet vector. More...
 
std::vector< unsigned int > _facetToBaseFacet
 The relationship between a triangle facet to its base facet. More...
 
std::vector< unsigned int > _facetToFacet
 The facet to facets connectivity list. More...
 
std::vector< unsigned int > _facetToEdge
 The facet to edges connectivity list. More...
 
std::vector< unsigned int > _edgeToVertex
 The edge to vertices connectivity list. More...
 
std::vector< unsigned int > _edgeToFacet
 The edge to facets connectivity list. More...
 
std::vector< unsigned int > _nonManifoldEdges
 The list of the non-manifold edges. More...
 
std::vector< double > _facetsVectors
 The in-plane vectors of each facet in the polyhedron. More...
 
std::vector< double > _facetsNormals
 The normals of each facet in the polyhedron. More...
 
std::vector< double > _verticesNormals
 The normals of each vertex in the polyhedron. More...
 
bool _isWatertight
 A flag that say if the polyhedron is watertight or not. More...
 
double _minCorner [3]
 The min corner of the polyhedron bounding box. More...
 
double _maxCorner [3]
 The max corner of the polyhedron bounding box. More...
 
double _greatestSquaredDiagonal
 Greatest squared diagonal used to normalize errors. More...
 
StaticAABBTree _aabbTree
 The Axis Aligned Bounding Box tree built from the facets of the polyhedron. More...
 
bool _isConvex
 The convexity of the polyhedron. More...
 
double _volume
 The volume of the tetrahedron. More...
 
double _initialCentroid [3]
 The initial centroid of the polyhedron, before any transformation. More...
 
double _centroid [3]
 The centroid of the polyhedron. More...
 
double _iMatrix [3][3]
 The inertia matrix of the polyhedron. More...
 
unsigned int _closestVertexIndex
 The index of the closest vertex found during a closest point query on the AABB tree, and set by the closest_point_callback() method. More...
 
unsigned int _primitivesHitByRayCount
 The number of primitives hit by a ray during a intersections count query on the AABB tree, and set by the raycast_aabb_tree_callback() method. More...
 
double _closestPoint [3]
 contains the closest point found during a call to closest_point() method. More...
 
unsigned int _furthestVertex
 contains the furthest vertex found during a call to furthest_vertex_along_axis() method. More...
 
std::vector< double > _hitPoints
 The list of points coordinates on primitives hit by a ray during raycast on the AABB tree, and set by the raycast_aabb_tree_callback() method. More...
 
std::vector< unsigned int > _hitPrimitives
 The list primitives hit by a ray during a raycast on the AABB tree, set by the raycast_aabb_tree_callback() callback method. More...
 
VertexPrimitivePair _closestVertexAndPrimitiveIndices
 The index of the closest primitive associated to the closest vertex found during a closest and primitive point query on the AABB tree, and set by the closest_pointAndPrimitiveCallback() callback method. More...
 
std::vector< unsigned int > _candidateFacets
 A vector containing candidate facets during the closest_pointAndPrimitive() procedure. It is set by the closest_pointAndPrimitiveCallback() callback method, and post-processed in the closest_pointAndPrimitive() method. More...
 

Static Private Attributes

static double SquaredDistanceErrorThreshold = 1e-16
 Small value used to detect squared distance mismatch between two points. More...
 
static unsigned char FacetEdgeVertex [3][2]
 contains the local connectivity of a facet vertices by edge. More...
 
static unsigned char FacetEdgeDirection [3][3]
 contains flags that allow to know if an edge is oriented towards a facet. or its incident counterpart. More...
 
static const unsigned int NotInList = (unsigned int)(-1)
 A flag to determine if an index is in a list. More...
 

Member Typedef Documentation

◆ VertexPrimitivePair

typedef std::pair<unsigned int, unsigned int> PolyhedronCore::VertexPrimitivePair

Constructor & Destructor Documentation

◆ PolyhedronCore() [1/3]

PolyhedronCore::PolyhedronCore ( )

instanciates a new PolyhedronCore object

instanciates a new empty PolyhedronCore object. By using this constuctor, the created object is empty and must be built manually using read_off_file().

◆ PolyhedronCore() [2/3]

PolyhedronCore::PolyhedronCore ( std::string  filePath)

instanciates a PolyhedronCore object with given OFF file

instanciates a new PolyhedronCore object with the OFF file specified by filePath. This constructor also prepares the data structure, checks for mesh sanity and build the AABB tree for fast queries.

◆ PolyhedronCore() [3/3]

PolyhedronCore::PolyhedronCore ( std::stringstream &  stringStream)

instanciates a PolyhedronCore object with given OFF file contents, as std::stringstream.

instanciates a new PolyhedronCore object with the OFF file contained in stringStream. This constructor also prepares the data structure, checks for mesh sanity and build the AABB tree for fast queries.

◆ ~PolyhedronCore()

PolyhedronCore::~PolyhedronCore ( )

Default PolyhedronCore destructor.

Member Function Documentation

◆ apply_offset()

void PolyhedronCore::apply_offset ( double  offsetX,
double  offsetY,
double  offsetZ 
)

Applies an offset to all vertices coordinates of the polyhedron.

After applying the offsets, all metrics and the AABB tree are recomputed.

◆ apply_scale_factor()

void PolyhedronCore::apply_scale_factor ( double  factor,
double  yFactor,
double  zFactor 
)

Applies a scale factor to all vertices coordinates of the polyhedron.

The scale factors are applied with respect to origin, that is point {0, 0, 0}. After applying the scale factors, all metrics and the AABB tree are recomputed.

◆ basefacet_area()

double PolyhedronCore::basefacet_area ( unsigned int  facetIndex) const

Computes the area of a base facet.

Computes the area of the base facet with index facetIndex. In case of gerometry file where facets are not necessarily triangles, this method computes the sum of areas of triangular facets composing the base facet.

◆ basefacet_index_for_facet()

unsigned int PolyhedronCore::basefacet_index_for_facet ( unsigned int  facetIndex) const
inline

Returns the base facet index for a triangle facet.

Returns the base facet index for a triangle facet with index facetIndex. In case of gerometry file where facets are not necessarily triangles, this method returns the facet index as seen in the geometry file. In case of triangle facets based geometry, the returned index is the same as the input.

◆ basefacet_number()

unsigned int PolyhedronCore::basefacet_number ( ) const
inline
Returns
The number of base facets in the polyhedron.

◆ build_aabb_tree()

void PolyhedronCore::build_aabb_tree ( )

Builds a static AABB tree based on the polyhedron primitives.

Building the AABB tree of the polyhedron primitives enables queries in the polyhedron, such as closest point search or ray casting (e.g. for point containment check).

◆ build_data_structure()

void PolyhedronCore::build_data_structure ( )

Builds the data structure which can be used to perform queries, e.g. adjacency.

The variables defined in this method are the _facetToFacet, _facetToEdge, _edgeToVertex, _edgeToFacet and _nonManifoldEdges vectors, that contain relationships between facets, edges and vertices.

◆ centroid_ptr()

const double* PolyhedronCore::centroid_ptr ( ) const
inline
Returns
The centroid of the polyhedron.

◆ clear()

void PolyhedronCore::clear ( )

Clears the geometric entities of the PolyhedronCore object.

Removes all facets, edges and vertices, and non-manifold edges, and invokes the clear_metric() method.

◆ clear_metric()

void PolyhedronCore::clear_metric ( )

Clears the all computations performed on the polyhedron.

Removes all normals and vectors, and reset the centroid, inertia matrix, volume, watertightness and convexity flags.

◆ closest_point()

void PolyhedronCore::closest_point ( double  aPoint[3],
double(&)  closestPoint[3] 
) const

Finds the point in the polyhedron which is the closest to the query point.

This method takes advantage of the AABB tree to find which point is closest to aPoint. Since the AABB tree does not know anything about the primitives it stores (except their Axis Aligned Bounding Box), the tree sends back a request to the Polyhedron object, through the closestPointCallback() callback method, to find which vertex is closest to the requested point.

◆ closest_point_callback()

void PolyhedronCore::closest_point_callback ( unsigned int  aFacet,
const double *  aPoint,
double &  distance 
) const

Finds the point in a polyhedron facets which is the closest to the requested point.

This method is intented to be invoked by the AABB tree during a closest point query. Given a global facet index aFacet and a comparison point Point, the method updates the distance parameter if the distance between aPoint and one of the facet vertices is smaller than the current value. Additionally, this method updates the index of the closest vertex in the facet, _closestVertexIndex.

◆ closest_point_in_facet()

void PolyhedronCore::closest_point_in_facet ( unsigned int  aFacet,
const double *  aPoint,
double(&)  closestPoint[3] 
) const

Finds the closest point lying on a facet closest to aPoint.

Finds the closest point lying in facet aFacet, which is closest to point aPoint, and sets closestPoint.

Note that Points lying exactly on a facet edge are not yet taken into account. If such a case appears, an assertion is raised.

◆ closest_point_to_edge()

double PolyhedronCore::closest_point_to_edge ( const double *  point0,
const double *  point1,
const double *  aPoint,
double(&)  closestPoint[3] 
) const
private

◆ closest_vertex_and_primitive()

PolyhedronCore::VertexPrimitivePair PolyhedronCore::closest_vertex_and_primitive ( double  aPoint[3]) const

Finds the vertex and the primitive it is associated to in the polyhedron which are the closest to the query point.

This method takes advantage of the AABB tree to find which vertex and the primitive it is associated to are closest to aPoint. Since the AABB tree does not know anything about the primitives it stores (except their Axis Aligned Bounding Box), the tree sends back a request to the Polyhedron object, through the closestVertexAndPrimitiveCallback() callback method, to find which vertex is closest to the requested point and the facet which realizes the shortest distance between the query point and its center of mass.

◆ closest_vertex_and_primitive_callback()

void PolyhedronCore::closest_vertex_and_primitive_callback ( unsigned int  aFacet,
const double *  aPoint,
double &  distance 
) const

Finds the vertex in a polyhedron facet and the associated facet which realizes the shortest distance to the requested point.

This method is intented to be invoked by the AABB tree during a closest point query. Given a global facet index aFacet and a comparison point Point, the method updates the distance parameter if the distance between aPoint and one of the facet vertices is smaller than the current value. Additionally, this method updates the index of the closest vertex and the closest facet, _closestVertexAnPrimitiveIndices. The facet is chosen as the closest to the point query with respect to its center of mass.

◆ closest_vertex_callback()

void PolyhedronCore::closest_vertex_callback ( unsigned int  aFacet,
const double *  aPoint,
double &  distance 
) const

Finds the vertex in a polyhedron facets which is the closest to the requested point.

This method is intented to be invoked by the AABB tree during a closest point query. Given a global facet index aFacet and a comparison point Point, the method updates the distance parameter if the distance between aPoint and one of the facet vertices is smaller than the current value. Additionally, this method updates the index of the closest vertex in the facet, _closestVertexIndex.

◆ closest_vertex_ptr()

const double * PolyhedronCore::closest_vertex_ptr ( double  aPoint[3]) const

Finds the vertex in the polyhedron which is the closest to the query point.

This method takes advantage of the AABB tree to find which vertex is closest to aPoint. Since the AABB tree does not know anything about the primitives it stores (except their Axis Aligned Bounding Box), the tree sends back a request to the Polyhedron object, through the closestVertexCallback() method, to find which vertex is closest to the query point.

◆ compute_metric()

void PolyhedronCore::compute_metric ( )

Computes geometric parameters related to the polyhedron.

Computes informations such as min and max corners, volume, inertia matrix, etc.

◆ contain()

bool PolyhedronCore::contain ( double  aPoint[3]) const

Checks if the polyhedron contains a 3D point.

This method takes advantage of the AABB tree to perform a raycasting on the tree AABBs. The ray is created from the origin 3D point aPoint, towards the +X direction to avoid useless AABB faces intersection computations. Each time an AABB is hit, the facet it contains is tested wether it is also hit by the ray and, if yes, the intersection point and facet index are stored in _hitPoints and _hitPrimitives vectors, and the intersections counter _primitivesHitByRayCount is incremented. The point is considered inside the polyhedron if _primitivesHitByRayCount is greater than 0, and if intersections modulo 2 equals 0.

◆ edge_number()

unsigned int PolyhedronCore::edge_number ( ) const
inline
Returns
The number of edges in the polyhedron.

◆ facet_area()

double PolyhedronCore::facet_area ( unsigned int  facetIndex) const

Computes the area of a triangle facet with index facetIndex.

◆ facet_hit_by_ray()

std::vector< unsigned int > PolyhedronCore::facet_hit_by_ray ( Ray  aRay) const

get_s the list of facets hit by a ray.

This method takes advantage of the AABB tree to perform a raycasting on the tree AABBs. Each time an AABB is hit, the facet it contains is tested wether it is also hit by the ray and, if yes, the intersection point and facet index are stored in _hitPoints and _hitPrimitives vectors, and the intersections counter _primitivesHitByRayCount is incremented.

◆ facet_number()

unsigned int PolyhedronCore::facet_number ( ) const
inline
Returns
The number of facets in the polyhedron.

◆ facet_to_vertex() [1/2]

std::vector<unsigned int>& PolyhedronCore::facet_to_vertex ( )
inline

◆ facet_to_vertex() [2/2]

const std::vector<unsigned int>& PolyhedronCore::facet_to_vertex ( ) const
inline

◆ find_local_edge_index_of_facet()

unsigned int PolyhedronCore::find_local_edge_index_of_facet ( unsigned int  aFacet,
unsigned int  anEdge 
)
private

Finds the local index of an edge of a given global facet.

Finds the local index of an edge anEdge in the given facet with global index aFacet.

Returns
The local index (in 0,1,2) of the given edge.

◆ furthest_vertex_along_axis()

void PolyhedronCore::furthest_vertex_along_axis ( const double *  anAxis,
unsigned int &  furthestVertex 
) const

Finds the vertex in the polyhedron which is the furthest in the given direction.

This method takes advantage of the AABB tree to find which vertex in the polyhedron is the furthest in the direction given by axis. Since the AABB tree does not know anything about the primitives it stores (except their Axis Aligned Bounding Box), the tree sends back a request to the Polyhedron object, through the furthestVertexCallback() callback method, to find which vertex is the furthest in the requested direction, for a given facet.

◆ furthest_vertex_callback()

void PolyhedronCore::furthest_vertex_callback ( unsigned int  aFacet,
const double *  aPoint,
const double *  anAxis,
double &  distance 
) const

Finds the vertex in a polyhedron facets which is the furthest, in the requested direction.

This method is intented to be invoked by the AABB tree during a furthest vertex query. Given a global facet index aFacet, a comparison point aPoint and a direction anAxis, the method updates the distance parameter if the distance between aPoint and one of the facet vertices is greater than the current value. Additionally, this method updates the index of the furthest vertex in the facet, _furthestVertexIndex.

◆ get_hitPoint()

std::vector<double>& PolyhedronCore::get_hitPoint ( )
inline
Returns
The list of hit point (see raycasting).

◆ inertia_tensor()

void PolyhedronCore::inertia_tensor ( double(&)  tensor[3][3]) const
inline

update the tensor to the inertia operator of the polyhedron.

◆ initial_centroid_ptr()

const double* PolyhedronCore::initial_centroid_ptr ( ) const
inline
Returns
The initial centroid of the polyhedron, before any transformation.

◆ is_convex()

bool PolyhedronCore::is_convex ( ) const
inline
Returns
The convexity of the polyhedron.

◆ is_facet_crossed_by_ray()

bool PolyhedronCore::is_facet_crossed_by_ray ( unsigned int  aFacet,
Ray  aRay,
double(&)  hitPoint[3] 
) const

Checks if a facet is crossed by a given ray.

Checks if facet with index aFacet is crossed by the ray aRay.

Returns
true or false whether the facet is crossed by the ray or not.

◆ is_watertight()

bool PolyhedronCore::is_watertight ( ) const
inline
Returns
The watertightness of the polyhedron.

◆ max_corner_ptr()

const double* PolyhedronCore::max_corner_ptr ( ) const
inline
Returns
The max corner of the polyhedron (axis aligned) bounding box.

◆ min_corner_ptr()

const double* PolyhedronCore::min_corner_ptr ( ) const
inline
Returns
The min corner of the polyhedron (axis aligned) bounding box.

◆ normal_ptr_of_facet()

const double* PolyhedronCore::normal_ptr_of_facet ( int  facetIndex) const
inline

Allows to get the normal of a given facet.

Returns
A double array pointer to the components of the facet normal with index facetIndex in the _facetsNormals vector.

◆ normal_ptr_of_vertex()

const double* PolyhedronCore::normal_ptr_of_vertex ( int  vertexIndex) const
inline

Allows to get a refernce to a vertex coordinates array from its index in the _vertices vector.

Returns
A reference double array pointer to the coordinates of the vertex with index vertexIndex in the _vertices vector.

Allows to get a vertex normal array from its index in the _verticesNormals vector.

Returns
A double array pointer to the components of the vertex normal with index vertexIndex in the _verticesNormals vector.

◆ raycast_aabb_tree_callback()

bool PolyhedronCore::raycast_aabb_tree_callback ( unsigned int  aFacet,
Ray  aRay 
) const

Checks if a primitive is hit by a ray.

This method is intented to be invoked by the AABB tree during counting primtives hit by a ray. Given a global facet index aFacet and a ray aRay, the method checks if the ray hits the primitive. Additionally, this method updates the global number of primitives hit by the ray, _primitivesHitByRayCount.

◆ read_off_file() [1/2]

void PolyhedronCore::read_off_file ( std::string  filePath)

Reads an OFF file given a file path.

This method reads an Object File Format (OFF) mesh file. It only reads the vertices and facets (not the edges).

◆ read_off_file() [2/2]

void PolyhedronCore::read_off_file ( std::stringstream &  stringStream)

Reads an OFF file from a stream.

This method reads an Object File Format (OFF) mesh file. It only reads the vertices and facets (not the edges).

◆ to_off_format()

void PolyhedronCore::to_off_format ( std::stringstream &  offString,
bool  basefacet = true 
) const

Creates an std::stringstream with the polyhedron contents, in OFF file format.

This method creates a standard string stream with the polyhedron contents, in Object File Format (OFF) mesh file. It only writes the vertices and facets (not the edges).

◆ to_vtk_ugrid()

void PolyhedronCore::to_vtk_ugrid ( vtkUnstructuredGrid &  uGrid)

Creates a VTK unstructured grid with the polyhedron contents.

This method creates a VTK unstructured grid with the contents of the polyhedron. It allows to further save it in .vtu file.

◆ translate_to_centroid()

void PolyhedronCore::translate_to_centroid ( )

translates the polyhedron to its center of mass.

After applying the translation, all metrics and the AABB tree are recomputed.

◆ vertex_index_of_basefacet()

const unsigned int PolyhedronCore::vertex_index_of_basefacet ( int  basefacetIndex,
int  vertexIndex 
) const
inline

Allows to get a the global vertex index for a given base facet, with a given local vertex index.

Returns
An unsigned int corresponding to the global index of local vertex vertexIndex, for the base facet with index facetIndex.

◆ vertex_index_of_facet()

const unsigned int PolyhedronCore::vertex_index_of_facet ( int  facetIndex,
int  vertexIndex 
) const
inline

Allows to get a the global vertex index for a given facet, with a given local vertex index.

Returns
An unsigned int corresponding to the global index of local vertex vertexIndex, for the facet with index facetIndex.

◆ vertex_number()

unsigned int PolyhedronCore::vertex_number ( ) const
inline
Returns
The number of vertices in the polyhedron.

◆ vertex_number_of_basefacet()

unsigned int PolyhedronCore::vertex_number_of_basefacet ( int  basefacetIndex) const
inline
Returns
The number of vertices in the base facet defined by basefacet.

◆ vertex_ptr()

const double* PolyhedronCore::vertex_ptr ( int  vertexIndex) const
inline

Allows to get a vertex coordinates array from its index in the _vertices vector.

Returns
A double array pointer to the coordinates of the vertex with index vertexIndex in the _vertices vector.

◆ vertex_ptr_of_facet()

const double* PolyhedronCore::vertex_ptr_of_facet ( unsigned int  facetIndex,
unsigned int  vectorIndex 
) const
inline

Allows to get the in-plane vectors of a given facet.

Returns
A double array pointer to the components of the vector with index vectorIndex of facet facetIndex in the _facetsVectors vector.

◆ volume()

double PolyhedronCore::volume ( ) const
inline
Returns
The volume of the polyhedron.

Member Data Documentation

◆ _aabbTree

StaticAABBTree PolyhedronCore::_aabbTree
private

The Axis Aligned Bounding Box tree built from the facets of the polyhedron.

It is used to perform queries, such as closest point, containment checks, etc.

◆ _basefacetToFacetOffset

std::vector<unsigned int> PolyhedronCore::_basefacetToFacetOffset
private

stores the offset required to access the first triangle facet of a base facet in _facetToBaseFacet vector.

◆ _basefacetToVertex

std::vector<unsigned int> PolyhedronCore::_basefacetToVertex
private

The relationship between base facets (or raw facets, directly from the geometry file that is, which may have any number of vertices). This vector must be used along with the _basefacetToVertexOffset vector.

◆ _basefacetToVertexOffset

std::vector<unsigned int> PolyhedronCore::_basefacetToVertexOffset
private

stores the offset required to access the first vertex of a base facet in _basefacetToVertex vector.

◆ _candidateFacets

std::vector<unsigned int> PolyhedronCore::_candidateFacets
mutableprivate

A vector containing candidate facets during the closest_pointAndPrimitive() procedure. It is set by the closest_pointAndPrimitiveCallback() callback method, and post-processed in the closest_pointAndPrimitive() method.

◆ _centroid

double PolyhedronCore::_centroid[3]
private

The centroid of the polyhedron.

◆ _closestPoint

double PolyhedronCore::_closestPoint[3]
mutableprivate

contains the closest point found during a call to closest_point() method.

◆ _closestVertexAndPrimitiveIndices

VertexPrimitivePair PolyhedronCore::_closestVertexAndPrimitiveIndices
mutableprivate

The index of the closest primitive associated to the closest vertex found during a closest and primitive point query on the AABB tree, and set by the closest_pointAndPrimitiveCallback() callback method.

The closest primitive realizes the sortest distance between its center of mass and the point query.

◆ _closestVertexIndex

unsigned int PolyhedronCore::_closestVertexIndex
mutableprivate

The index of the closest vertex found during a closest point query on the AABB tree, and set by the closest_point_callback() method.

◆ _edgeToFacet

std::vector<unsigned int> PolyhedronCore::_edgeToFacet
private

The edge to facets connectivity list.

In a sane mesh, an edge may be connected to one or two facets, which are given in this array as global indices. When populating this conectivity list, if an edge is found to be adjacent to more than two facets, this edge is added to the _nonManifoldEdges vector for further sanity check.

◆ _edgeToVertex

std::vector<unsigned int> PolyhedronCore::_edgeToVertex
private

The edge to vertices connectivity list.

An edge is defined by two vertices, which are given in this array as global indices.

◆ _facetsNormals

std::vector<double> PolyhedronCore::_facetsNormals
private

The normals of each facet in the polyhedron.

The facets normals can be used to check for intersections or for 3D (flat) rendering.

◆ _facetsVectors

std::vector<double> PolyhedronCore::_facetsVectors
private

The in-plane vectors of each facet in the polyhedron.

The facets vectors can be used to check for intersections.

◆ _facetToBaseFacet

std::vector<unsigned int> PolyhedronCore::_facetToBaseFacet
private

The relationship between a triangle facet to its base facet.

◆ _facetToEdge

std::vector<unsigned int> PolyhedronCore::_facetToEdge
private

The facet to edges connectivity list.

A facet has three edges, which are given in this array as global indices.

◆ _facetToFacet

std::vector<unsigned int> PolyhedronCore::_facetToFacet
private

The facet to facets connectivity list.

A facet may have up to three adjacent facets through its three edges, that are listed in this vector. If a facet has no adjacent facet through a given local edge, it is marked as NotInList.

◆ _facetToVertex

std::vector<unsigned int> PolyhedronCore::_facetToVertex
private

The relationship between triangle facets and vertices.

◆ _furthestVertex

unsigned int PolyhedronCore::_furthestVertex
mutableprivate

contains the furthest vertex found during a call to furthest_vertex_along_axis() method.

◆ _greatestSquaredDiagonal

double PolyhedronCore::_greatestSquaredDiagonal
private

Greatest squared diagonal used to normalize errors.

◆ _hitPoints

std::vector<double> PolyhedronCore::_hitPoints
mutableprivate

The list of points coordinates on primitives hit by a ray during raycast on the AABB tree, and set by the raycast_aabb_tree_callback() method.

◆ _hitPrimitives

std::vector<unsigned int> PolyhedronCore::_hitPrimitives
mutableprivate

The list primitives hit by a ray during a raycast on the AABB tree, set by the raycast_aabb_tree_callback() callback method.

◆ _iMatrix

double PolyhedronCore::_iMatrix[3][3]
private

The inertia matrix of the polyhedron.

◆ _initialCentroid

double PolyhedronCore::_initialCentroid[3]
private

The initial centroid of the polyhedron, before any transformation.

◆ _isConvex

bool PolyhedronCore::_isConvex
private

The convexity of the polyhedron.

◆ _isWatertight

bool PolyhedronCore::_isWatertight
private

A flag that say if the polyhedron is watertight or not.

A polyhedron is assumed watertight if it does not contain free edges (such as holes in the mesh), and non-manifold edges (edges connected to 3 or more facets).

◆ _maxCorner

double PolyhedronCore::_maxCorner[3]
private

The max corner of the polyhedron bounding box.

◆ _minCorner

double PolyhedronCore::_minCorner[3]
private

The min corner of the polyhedron bounding box.

◆ _nonManifoldEdges

std::vector<unsigned int> PolyhedronCore::_nonManifoldEdges
private

The list of the non-manifold edges.

An edge is considered as non-manifold when is adjacent two more than two facets.

◆ _primitivesHitByRayCount

unsigned int PolyhedronCore::_primitivesHitByRayCount
mutableprivate

The number of primitives hit by a ray during a intersections count query on the AABB tree, and set by the raycast_aabb_tree_callback() method.

◆ _vertices

std::vector<double> PolyhedronCore::_vertices
private

The list of vertices coordinates, given as [x0,y0,z0,...,xn,yn,zn].

◆ _verticesNormals

std::vector<double> PolyhedronCore::_verticesNormals
private

The normals of each vertex in the polyhedron.

The vertices normals can be used for 3D (smooth) rendering.

◆ _verticesPerBaseFacet

std::vector<unsigned int> PolyhedronCore::_verticesPerBaseFacet
private

The number of vertices that form a base facet.

◆ _volume

double PolyhedronCore::_volume
private

The volume of the tetrahedron.

◆ FacetEdgeDirection

unsigned char PolyhedronCore::FacetEdgeDirection
staticprivate
Initial value:
= {
{ 0, 0, 1},
{ 1, 0, 0},
{ 0, 1, 0}
}

contains flags that allow to know if an edge is oriented towards a facet. or its incident counterpart.

Given two local vertices index of a facet, this table gives 0 or 1 whether the edge formed by these two vertices is oriented towards the facet or the adjacent one. Invoking FacetEdgeDirection[0][1] gives 0 (towards the facet), while FacetEdgeDirection[1][0] gives 1 (towards the adjacent facet).

◆ FacetEdgeVertex

unsigned char PolyhedronCore::FacetEdgeVertex
staticprivate
Initial value:
= {
{0, 1},
{1, 2},
{2, 0}
}

contains the local connectivity of a facet vertices by edge.

A facet is defined by three vertices labeled 0,1,2. It contains 3 edges (in 0,1,2), and an edge is defined by two vertices (in 0,1). Invoking FacetEdgeVertex[0][1] gives the second vertex of edge 0 in the facet, corresponding to the local facet vertex 1.

◆ NotInList

const unsigned int PolyhedronCore::NotInList = (unsigned int)(-1)
staticprivate

A flag to determine if an index is in a list.

◆ SquaredDistanceErrorThreshold

double PolyhedronCore::SquaredDistanceErrorThreshold = 1e-16
staticprivate

Small value used to detect squared distance mismatch between two points.


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