Element< ELEMENT_DIM, SPACE_DIM > Class Template Reference

#include <Element.hpp>

Inheritance diagram for Element< ELEMENT_DIM, SPACE_DIM >:

Inheritance graph
[legend]
Collaboration diagram for Element< ELEMENT_DIM, SPACE_DIM >:

Collaboration graph
[legend]

List of all members.

Public Member Functions

 Element (unsigned index, std::vector< Node< SPACE_DIM > * > nodes)
 Element (const Element &rElement, const unsigned index)
void RegisterWithNodes ()
void MarkAsDeleted ()
void UpdateNode (const unsigned &rIndex, Node< SPACE_DIM > *pNode)
void ResetIndex (unsigned index)
c_vector< double, SPACE_DIM+1 > CalculateCircumsphere (c_matrix< double, SPACE_DIM, ELEMENT_DIM > &rJacobian, c_matrix< double, ELEMENT_DIM, SPACE_DIM > &rInverseJacobian)
double CalculateCircumsphereVolume ()
double CalculateQuality ()
c_vector< double, SPACE_DIM+1 > CalculateInterpolationWeights (ChastePoint< SPACE_DIM > testPoint)
c_vector< double, SPACE_DIM+1 > CalculateInterpolationWeightsWithProjection (ChastePoint< SPACE_DIM > testPoint)
c_vector< double, SPACE_DIM > CalculatePsi (ChastePoint< SPACE_DIM > testPoint)
bool IncludesPoint (ChastePoint< SPACE_DIM > testPoint, bool strict=false)


Detailed Description

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
class Element< ELEMENT_DIM, SPACE_DIM >

A concrete element class which inherits from AbstractTetrahedralElement.

Definition at line 46 of file Element.hpp.


Constructor & Destructor Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
Element< ELEMENT_DIM, SPACE_DIM >::Element ( unsigned  index,
std::vector< Node< SPACE_DIM > * >  nodes 
) [inline]

Constructor which takes in a vector of nodes.

Parameters:
index the index of the element in the mesh
nodes the nodes owned by the element
Todo:
make rNodes like in AbstractTetrahedralElement? (#991)

Definition at line 42 of file Element.cpp.

References Element< ELEMENT_DIM, SPACE_DIM >::RegisterWithNodes().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
Element< ELEMENT_DIM, SPACE_DIM >::Element ( const Element< ELEMENT_DIM, SPACE_DIM > &  rElement,
const unsigned  index 
) [inline]

Copy constructor which allows a new index to be specified.

Todo:
this is rather dubious; a factory method might be better.
Parameters:
rElement an element to copy
index the index of the new element

Definition at line 49 of file Element.cpp.

References AbstractElement< ELEMENT_DIM, SPACE_DIM >::mIndex, and Element< ELEMENT_DIM, SPACE_DIM >::RegisterWithNodes().


Member Function Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void Element< ELEMENT_DIM, SPACE_DIM >::RegisterWithNodes (  )  [inline, virtual]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void Element< ELEMENT_DIM, SPACE_DIM >::MarkAsDeleted (  )  [inline, virtual]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void Element< ELEMENT_DIM, SPACE_DIM >::UpdateNode ( const unsigned &  rIndex,
Node< SPACE_DIM > *  pNode 
) [inline, virtual]

Update node at the given index.

Parameters:
rIndex is an local index to which node to change
pNode is a pointer to the replacement node
Update node at the given index
Parameters:
rIndex is an local index to which node to change
pNode is a pointer to the replacement node

Implements AbstractElement< ELEMENT_DIM, SPACE_DIM >.

Definition at line 82 of file Element.cpp.

References AbstractElement< ELEMENT_DIM, SPACE_DIM >::mIndex, and AbstractElement< ELEMENT_DIM, SPACE_DIM >::mNodes.

Referenced by MutableMesh< ELEMENT_DIM, SPACE_DIM >::RefineElement().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void Element< ELEMENT_DIM, SPACE_DIM >::ResetIndex ( unsigned  index  )  [inline]

Reset the index of this boundary element in the mesh.

Parameters:
index the new index of the boundary element

Definition at line 97 of file Element.cpp.

References AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::mIndex, AbstractElement< ELEMENT_DIM, SPACE_DIM >::mNodes, and Element< ELEMENT_DIM, SPACE_DIM >::RegisterWithNodes().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
c_vector< double, SPACE_DIM+1 > Element< ELEMENT_DIM, SPACE_DIM >::CalculateCircumsphere ( c_matrix< double, SPACE_DIM, ELEMENT_DIM > &  rJacobian,
c_matrix< double, ELEMENT_DIM, SPACE_DIM > &  rInverseJacobian 
) [inline]

Calculate the circumsphere/circumcircle of this element.

Returns:
a vector containing x_centre, y_centre,...,radius^2 Calculate the circumsphere/circumcircle of this element.
After reconstructing a cylindrical 2d mesh, the jacobian data of the periodic elements is not valid anymore. We want to use the jacobians computed before swapping the nodes.

Returns:
a vector containing x_centre, y_centre,...,radius^2
Parameters:
rJacobian the Jacobian matrix
rInverseJacobian the inverse Jacobian matrix

Definition at line 111 of file Element.cpp.

References AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeLocation().

Referenced by Element< ELEMENT_DIM, SPACE_DIM >::CalculateQuality(), and MutableMesh< ELEMENT_DIM, SPACE_DIM >::CheckVoronoi().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
double Element< ELEMENT_DIM, SPACE_DIM >::CalculateCircumsphereVolume (  ) 

Get the volume of the circumsphere, or area of the circumcircle, of this element.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
double Element< ELEMENT_DIM, SPACE_DIM >::CalculateQuality (  )  [inline]

The quality of a triangle/tetrahedron is the ratio between the volume of the shape and the volume of its circumsphere. This is normalised by dividing through by the Platonic ratio.

Definition at line 156 of file Element.cpp.

References Element< ELEMENT_DIM, SPACE_DIM >::CalculateCircumsphere(), and AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::CalculateInverseJacobian().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
c_vector< double, SPACE_DIM+1 > Element< ELEMENT_DIM, SPACE_DIM >::CalculateInterpolationWeights ( ChastePoint< SPACE_DIM >  testPoint  )  [inline]

Calculate the interpolation weights at a given point.

Parameters:
testPoint the point

Definition at line 196 of file Element.cpp.

References Element< ELEMENT_DIM, SPACE_DIM >::CalculatePsi().

Referenced by Element< ELEMENT_DIM, SPACE_DIM >::CalculateInterpolationWeightsWithProjection(), and Element< ELEMENT_DIM, SPACE_DIM >::IncludesPoint().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
c_vector< double, SPACE_DIM+1 > Element< ELEMENT_DIM, SPACE_DIM >::CalculateInterpolationWeightsWithProjection ( ChastePoint< SPACE_DIM >  testPoint  )  [inline]

Calculate the interpolation weights, but if we are not within the element (one or more negative weights), we project onto the element, rather than extrapolating from it.

Parameters:
testPoint the point
Calculate the interpolation weights, but if we are not within the element (one or more negative weights), we project onto the element, rather than extrapolating from it.

Definition at line 221 of file Element.cpp.

References Element< ELEMENT_DIM, SPACE_DIM >::CalculateInterpolationWeights().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
c_vector< double, SPACE_DIM > Element< ELEMENT_DIM, SPACE_DIM >::CalculatePsi ( ChastePoint< SPACE_DIM >  testPoint  )  [inline]

Calculate psi at a given point.

Parameters:
testPoint the point

Todo:
This method shouldn't need a new Jacobian inverse for every Psi

Definition at line 263 of file Element.cpp.

References AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::CalculateInverseJacobian(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeLocation(), and ChastePoint< DIM >::rGetLocation().

Referenced by Element< ELEMENT_DIM, SPACE_DIM >::CalculateInterpolationWeights().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
bool Element< ELEMENT_DIM, SPACE_DIM >::IncludesPoint ( ChastePoint< SPACE_DIM >  testPoint,
bool  strict = false 
) [inline]

Get whether a given point lies inside this element.

Parameters:
testPoint the point
strict whether the point must not be too close to an edge/face (defaults to false)

Definition at line 284 of file Element.cpp.

References Element< ELEMENT_DIM, SPACE_DIM >::CalculateInterpolationWeights().

Referenced by MutableMesh< ELEMENT_DIM, SPACE_DIM >::RefineElement().


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

Generated on Tue Aug 4 16:11:13 2009 for Chaste by  doxygen 1.5.5