Chaste  Build::
Facet< DIM > Class Template Reference

#include <Facet.hpp>

+ Collaboration diagram for Facet< DIM >:

Public Member Functions

 Facet (std::vector< boost::shared_ptr< Polygon< DIM > > > polygons)
 
 Facet (boost::shared_ptr< Polygon< DIM > > pPolygon)
 
 ~Facet ()
 
void AddPolygons (std::vector< boost::shared_ptr< Polygon< DIM > > > polygons)
 
void AddPolygon (boost::shared_ptr< Polygon< DIM > > pPolygon)
 
bool ContainsPoint (const DimensionalChastePoint< DIM > &location)
 
std::vector< units::quantity< unit::length > > GetBoundingBox ()
 
DimensionalChastePoint< DIM > GetCentroid ()
 
units::quantity< unit::length > GetDistance (const DimensionalChastePoint< DIM > &rLocation)
 
std::string GetLabel ()
 
vtkSmartPointer< vtkPlane > GetPlane ()
 
c_vector< double, DIM > GetNormal ()
 
std::vector< boost::shared_ptr< Polygon< DIM > > > GetPolygons ()
 
std::vector< boost::shared_ptr< DimensionalChastePoint< DIM > > > GetVertices ()
 
std::pair< vtkSmartPointer< vtkPoints >, vtkSmartPointer< vtkIdTypeArray > > GetVtkVertices ()
 
void RotateAboutAxis (c_vector< double, 3 > axis, double angle)
 
void SetLabel (const std::string &label)
 
void Translate (DimensionalChastePoint< DIM > translationVector)
 
void UpdateVertices ()
 

Static Public Member Functions

static boost::shared_ptr< Facet< DIM > > Create (std::vector< boost::shared_ptr< Polygon< DIM > > > polygons)
 
static boost::shared_ptr< Facet< DIM > > Create (boost::shared_ptr< Polygon< DIM > > pPolygon)
 

Private Attributes

std::vector< boost::shared_ptr< Polygon< DIM > > > mPolygons
 
std::vector< boost::shared_ptr< DimensionalChastePoint< DIM > > > mVertices
 
bool mVerticesUpToDate
 
std::map< std::string, double > mData
 
std::string mLabel
 
units::quantity< unit::length > mReferenceLength
 

Detailed Description

template<unsigned DIM>
class Facet< DIM >

A collection of planar polygons

Definition at line 54 of file Facet.hpp.

Constructor & Destructor Documentation

template<unsigned DIM>
Facet< DIM >::Facet ( std::vector< boost::shared_ptr< Polygon< DIM > > >  polygons)

Constructor

Parameters
polygonsa facet is made from these polygons

Definition at line 46 of file Facet.cpp.

template<unsigned DIM>
Facet< DIM >::Facet ( boost::shared_ptr< Polygon< DIM > >  pPolygon)

Constructor

Parameters
pPolygona single polygon for the facet

Definition at line 58 of file Facet.cpp.

References Facet< DIM >::mPolygons.

template<unsigned DIM>
Facet< DIM >::~Facet ( )

Desctructor

Definition at line 84 of file Facet.cpp.

Member Function Documentation

template<unsigned DIM>
void Facet< DIM >::AddPolygon ( boost::shared_ptr< Polygon< DIM > >  pPolygon)

Add polygon

Parameters
pPolygona polygon

Definition at line 96 of file Facet.cpp.

References Facet< DIM >::mPolygons, and Facet< DIM >::mVerticesUpToDate.

template<unsigned DIM>
void Facet< DIM >::AddPolygons ( std::vector< boost::shared_ptr< Polygon< DIM > > >  polygons)

Add polygons

Parameters
polygonsplanar polygons

Definition at line 89 of file Facet.cpp.

References Facet< DIM >::mPolygons, and Facet< DIM >::mVerticesUpToDate.

template<unsigned DIM>
bool Facet< DIM >::ContainsPoint ( const DimensionalChastePoint< DIM > &  location)

Return true if the specified location is in the facet

Parameters
locationthe location to be tested
Returns
true if the location is in the facet

Definition at line 103 of file Facet.cpp.

References Facet< DIM >::mPolygons.

template<unsigned DIM>
boost::shared_ptr< Facet< DIM > > Facet< DIM >::Create ( std::vector< boost::shared_ptr< Polygon< DIM > > >  polygons)
static

Factory constructor method

Parameters
polygonsplanar polygons
Returns
a shared pointer to a new facet

Definition at line 70 of file Facet.cpp.

template<unsigned DIM>
boost::shared_ptr< Facet< DIM > > Facet< DIM >::Create ( boost::shared_ptr< Polygon< DIM > >  pPolygon)
static

Factory constructor method

Parameters
pPolygona polygon
Returns
a smart pointer to a new facet

Definition at line 77 of file Facet.cpp.

template<unsigned DIM>
std::vector< units::quantity< unit::length > > Facet< DIM >::GetBoundingBox ( )

Return the bounding box of the facet

Returns
the bounding box (xmin, xmax, ymin, ymax, zmin, zmax)

Definition at line 118 of file Facet.cpp.

References Facet< DIM >::GetVertices(), and Facet< DIM >::mReferenceLength.

template<unsigned DIM>
DimensionalChastePoint< DIM > Facet< DIM >::GetCentroid ( )

Return the centroid of the facet

Returns
the centroid of the facet

Definition at line 157 of file Facet.cpp.

References Facet< DIM >::GetVtkVertices(), and Facet< DIM >::mReferenceLength.

Referenced by Facet< DIM >::GetPlane().

template<unsigned DIM>
units::quantity< unit::length > Facet< DIM >::GetDistance ( const DimensionalChastePoint< DIM > &  rLocation)

Return the distance to the facet

Parameters
rLocationreference to the location of the point for distance calculation
Returns
the distance to the facet

Definition at line 184 of file Facet.cpp.

References DimensionalChastePoint< DIM >::GetLocation(), Facet< DIM >::GetPlane(), and Facet< DIM >::mReferenceLength.

template<unsigned DIM>
std::string Facet< DIM >::GetLabel ( )

Get the label for boundary conditions

Returns
label the boundary condition label

Definition at line 178 of file Facet.cpp.

References Facet< DIM >::mLabel.

template<unsigned DIM>
c_vector< double, DIM > Facet< DIM >::GetNormal ( )

Return the normal to the facet

Returns
the normal to the facet

Definition at line 202 of file Facet.cpp.

References Facet< DIM >::GetVertices(), and Facet< DIM >::GetVtkVertices().

Referenced by Facet< DIM >::GetPlane().

template<unsigned DIM>
vtkSmartPointer< vtkPlane > Facet< DIM >::GetPlane ( )

Return the facet's plane

Returns
a vtk plane on the facet's plane

Definition at line 233 of file Facet.cpp.

References Facet< DIM >::GetCentroid(), Facet< DIM >::GetNormal(), and Facet< DIM >::mReferenceLength.

Referenced by Facet< DIM >::GetDistance().

template<unsigned DIM>
std::vector< boost::shared_ptr< Polygon< DIM > > > Facet< DIM >::GetPolygons ( )

Return the polygons

Returns
the polygons making up the facet

Definition at line 252 of file Facet.cpp.

References Facet< DIM >::mPolygons.

template<unsigned DIM>
std::vector< boost::shared_ptr< DimensionalChastePoint< DIM > > > Facet< DIM >::GetVertices ( )
template<unsigned DIM>
std::pair< vtkSmartPointer< vtkPoints >, vtkSmartPointer< vtkIdTypeArray > > Facet< DIM >::GetVtkVertices ( )

Return the facet vertices as a set of VtkPoints.

Returns
the facet vertices as a set of VtkPoints.

Definition at line 268 of file Facet.cpp.

References Facet< DIM >::GetVertices(), and Facet< DIM >::mReferenceLength.

Referenced by Facet< DIM >::GetCentroid(), and Facet< DIM >::GetNormal().

template<unsigned DIM>
void Facet< DIM >::RotateAboutAxis ( c_vector< double, 3 >  axis,
double  angle 
)

Rotate about the specified axis by the specified angle

Parameters
axisthe rotation axis
anglethe rotation angle

Definition at line 292 of file Facet.cpp.

References Facet< DIM >::GetVertices().

template<unsigned DIM>
void Facet< DIM >::SetLabel ( const std::string &  label)

Set the label for boundary conditions

Parameters
labelthe boundary condition label

Definition at line 302 of file Facet.cpp.

References Facet< DIM >::mLabel.

template<unsigned DIM>
void Facet< DIM >::Translate ( DimensionalChastePoint< DIM >  translationVector)

Move the facet along the translation vector

Parameters
translationVectorthe new location is the original + the translationVector

Definition at line 308 of file Facet.cpp.

References Facet< DIM >::GetVertices().

template<unsigned DIM>
void Facet< DIM >::UpdateVertices ( )

Update the mVertices member

Definition at line 318 of file Facet.cpp.

References Facet< DIM >::mPolygons, Facet< DIM >::mVertices, and Facet< DIM >::mVerticesUpToDate.

Referenced by Facet< DIM >::GetVertices().

Member Data Documentation

template<unsigned DIM>
std::map<std::string, double> Facet< DIM >::mData
private

Data container, useful for specifying boundary conditions on facets

Definition at line 76 of file Facet.hpp.

template<unsigned DIM>
std::string Facet< DIM >::mLabel
private

A label for the application of boundary conditions

Definition at line 81 of file Facet.hpp.

Referenced by Facet< DIM >::GetLabel(), and Facet< DIM >::SetLabel().

template<unsigned DIM>
std::vector<boost::shared_ptr<Polygon<DIM> > > Facet< DIM >::mPolygons
private
template<unsigned DIM>
units::quantity<unit::length> Facet< DIM >::mReferenceLength
private
template<unsigned DIM>
std::vector<boost::shared_ptr<DimensionalChastePoint<DIM> > > Facet< DIM >::mVertices
private

Unique vertices in the facet. This is not always up-to-date. Use GetVertices() to ensure up-to-date vertices are used.

Definition at line 65 of file Facet.hpp.

Referenced by Facet< DIM >::GetVertices(), and Facet< DIM >::UpdateVertices().

template<unsigned DIM>
bool Facet< DIM >::mVerticesUpToDate
private

Whether mVertices is up-to-date. This should be set false when new polygons are added.

Definition at line 71 of file Facet.hpp.

Referenced by Facet< DIM >::AddPolygon(), Facet< DIM >::AddPolygons(), Facet< DIM >::GetVertices(), and Facet< DIM >::UpdateVertices().


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