DGtal 2.1.0
Loading...
Searching...
No Matches
DGtal::QuickHull< TKernel >::Facet Struct Reference

#include <DGtal/geometry/tools/QuickHull.h>

Public Member Functions

 Facet ()=default
 
 Facet (const Facet &)=default
 
 Facet (Facet &&)=default
 
Facetoperator= (Facet &&)=default
 
Facetoperator= (const Facet &)=default
 
 Facet (const HalfSpace &aH, Index b)
 
void clear ()
 
void addPointOn (Index p)
 
void display (std::ostream &out) const
 
void addNeighbor (Index n)
 
void subNeighbor (Index n)
 
void swap (Facet &other)
 
Size variableMemory () const
 

Data Fields

HalfSpace H
 the facet geometry
 
IndexRange neighbors
 neighbor facets
 
IndexRange outside_set
 outside set, i.e. points above this facet
 
IndexRange on_set
 on set, i.e. points on this facet, sorted
 
Index below
 index of point that is below this facet
 

Detailed Description

template<typename TKernel>
struct DGtal::QuickHull< TKernel >::Facet

A facet is d-1 dimensional convex cell lying on the boundary of a full dimensional convex set. Its supporting hyperplane defines an half-space touching and enclosing the convex set.

Definition at line 161 of file QuickHull.h.

Constructor & Destructor Documentation

◆ Facet() [1/4]

template<typename TKernel >
DGtal::QuickHull< TKernel >::Facet::Facet ( )
default

◆ Facet() [2/4]

template<typename TKernel >
DGtal::QuickHull< TKernel >::Facet::Facet ( const Facet )
default

◆ Facet() [3/4]

template<typename TKernel >
DGtal::QuickHull< TKernel >::Facet::Facet ( Facet &&  )
default

◆ Facet() [4/4]

template<typename TKernel >
DGtal::QuickHull< TKernel >::Facet::Facet ( const HalfSpace aH,
Index  b 
)
inline

Definition at line 173 of file QuickHull.h.

174 : H( aH ), below( b ) {}
Index below
index of point that is below this facet
Definition QuickHull.h:166
HalfSpace H
the facet geometry
Definition QuickHull.h:162

Member Function Documentation

◆ addNeighbor()

template<typename TKernel >
void DGtal::QuickHull< TKernel >::Facet::addNeighbor ( Index  n)
inline

Definition at line 202 of file QuickHull.h.

203 {
204 const auto it = std::find( neighbors.cbegin(), neighbors.cend(), n );
205 if ( it == neighbors.cend() ) neighbors.push_back( n );
206 }
IndexRange neighbors
neighbor facets
Definition QuickHull.h:163

References DGtal::QuickHull< TKernel >::Facet::neighbors.

◆ addPointOn()

template<typename TKernel >
void DGtal::QuickHull< TKernel >::Facet::addPointOn ( Index  p)
inline

Definition at line 184 of file QuickHull.h.

185 {
186 const auto it = std::find( on_set.cbegin(), on_set.cend(), p );
187 if ( it == on_set.cend() ) on_set.push_back( p );
188 }
IndexRange on_set
on set, i.e. points on this facet, sorted
Definition QuickHull.h:165

References DGtal::QuickHull< TKernel >::Facet::on_set.

◆ clear()

template<typename TKernel >
void DGtal::QuickHull< TKernel >::Facet::clear ( )
inline

◆ display()

template<typename TKernel >
void DGtal::QuickHull< TKernel >::Facet::display ( std::ostream &  out) const
inline

Definition at line 189 of file QuickHull.h.

190 {
191 const auto N = H.internalNormal();
192 out << "[Facet iN=(" << N[0];
193 for ( Dimension i = 1; i < N.dimension; i++ ) out << "," << N[ i ];
194 out << ") c=" << H.internalIntercept() << " b=" << below << " n={";
195 for ( auto&& n : neighbors ) out << " " << n;
196 out << " } #out=" << outside_set.size();
197 out << " on={";
198 for ( auto&& n : on_set ) out << " " << n;
199 out << " }]" << std::endl;
200 }
DGtal::uint32_t Dimension
Definition Common.h:119

References DGtal::QuickHull< TKernel >::Facet::below, DGtal::QuickHull< TKernel >::Facet::H, DGtal::QuickHull< TKernel >::Facet::neighbors, DGtal::QuickHull< TKernel >::Facet::on_set, and DGtal::QuickHull< TKernel >::Facet::outside_set.

Referenced by DGtal::QuickHull< TKernel >::processFacet().

◆ operator=() [1/2]

template<typename TKernel >
Facet & DGtal::QuickHull< TKernel >::Facet::operator= ( const Facet )
default

◆ operator=() [2/2]

template<typename TKernel >
Facet & DGtal::QuickHull< TKernel >::Facet::operator= ( Facet &&  )
default

◆ subNeighbor()

template<typename TKernel >
void DGtal::QuickHull< TKernel >::Facet::subNeighbor ( Index  n)
inline

Definition at line 207 of file QuickHull.h.

208 {
209 auto it = std::find( neighbors.begin(), neighbors.end(), n );
210 if ( it != neighbors.end() ) {
211 std::swap( *it, neighbors.back() );
212 neighbors.pop_back();
213 }
214 }

References DGtal::QuickHull< TKernel >::Facet::neighbors.

◆ swap()

template<typename TKernel >
void DGtal::QuickHull< TKernel >::Facet::swap ( Facet other)
inline

Definition at line 215 of file QuickHull.h.

216 {
217 if ( this != &other ) {
218 std::swap( H, other.H );
219 neighbors.swap ( other.neighbors );
220 outside_set.swap( other.outside_set );
221 on_set.swap ( other.on_set );
222 std::swap( below, other.below );
223 }
224 }

References DGtal::QuickHull< TKernel >::Facet::below, DGtal::QuickHull< TKernel >::Facet::H, DGtal::QuickHull< TKernel >::Facet::neighbors, DGtal::QuickHull< TKernel >::Facet::on_set, and DGtal::QuickHull< TKernel >::Facet::outside_set.

◆ variableMemory()

template<typename TKernel >
Size DGtal::QuickHull< TKernel >::Facet::variableMemory ( ) const
inline

Definition at line 225 of file QuickHull.h.

226 {
227 Size M;
228 M += neighbors.capacity() * sizeof( Index );
229 M += outside_set.capacity() * sizeof( Index );
230 M += on_set.capacity() * sizeof( Index );
231 return M;
232 }
std::size_t Index
Definition QuickHull.h:147
HalfEdgeDataStructure::Size Size

References DGtal::QuickHull< TKernel >::Facet::neighbors, DGtal::QuickHull< TKernel >::Facet::on_set, and DGtal::QuickHull< TKernel >::Facet::outside_set.

Field Documentation

◆ below

template<typename TKernel >
Index DGtal::QuickHull< TKernel >::Facet::below

◆ H

◆ neighbors

◆ on_set

◆ outside_set


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