Chaste  Release::2018.1
AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM > Class Template Reference

#include <AbstractTetrahedralElement.hpp>

+ Inheritance diagram for AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >:
+ Collaboration diagram for AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >:

Public Member Functions

 AbstractTetrahedralElement (unsigned index, const std::vector< Node< SPACE_DIM > * > &rNodes)
 
 AbstractTetrahedralElement (unsigned index=INDEX_IS_NOT_USED)
 
virtual ~AbstractTetrahedralElement ()
 
c_vector< double, SPACE_DIM > CalculateCentroid () const
 
void CalculateJacobian (c_matrix< double, SPACE_DIM, ELEMENT_DIM > &rJacobian, double &rJacobianDeterminant)
 
void CalculateWeightedDirection (c_vector< double, SPACE_DIM > &rWeightedDirection, double &rJacobianDeterminant)
 
c_vector< double, SPACE_DIM > CalculateNormal ()
 
void CalculateInverseJacobian (c_matrix< double, SPACE_DIM, ELEMENT_DIM > &rJacobian, double &rJacobianDeterminant, c_matrix< double, ELEMENT_DIM, SPACE_DIM > &rInverseJacobian)
 
double GetVolume (double determinant) const
 
void GetStiffnessMatrixGlobalIndices (unsigned problemDim, unsigned *pIndices) const
 
- Public Member Functions inherited from AbstractElement< ELEMENT_DIM, SPACE_DIM >
 AbstractElement (unsigned index, const std::vector< Node< SPACE_DIM > * > &rNodes)
 
 AbstractElement (unsigned index=INDEX_IS_NOT_USED)
 
virtual ~AbstractElement ()
 
virtual void UpdateNode (const unsigned &rIndex, Node< SPACE_DIM > *pNode)=0
 
void ReplaceNode (Node< SPACE_DIM > *pOldNode, Node< SPACE_DIM > *pNewNode)
 
virtual void MarkAsDeleted ()=0
 
virtual void RegisterWithNodes ()=0
 
double GetNodeLocation (unsigned localIndex, unsigned dimension) const
 
c_vector< double, SPACE_DIM > GetNodeLocation (unsigned localIndex) const
 
unsigned GetNodeGlobalIndex (unsigned localIndex) const
 
Node< SPACE_DIM > * GetNode (unsigned localIndex) const
 
unsigned GetNumNodes () const
 
void AddNode (Node< SPACE_DIM > *pNode)
 
bool IsDeleted () const
 
unsigned GetIndex () const
 
void SetIndex (unsigned index)
 
bool GetOwnership () const
 
void SetOwnership (bool ownership)
 
void SetAttribute (double attribute)
 
double GetAttribute ()
 
unsigned GetUnsignedAttribute ()
 
void AddElementAttribute (double attribute)
 
std::vector< double > & rGetElementAttributes ()
 
unsigned GetNumElementAttributes ()
 

Protected Member Functions

void RefreshJacobian (c_matrix< double, SPACE_DIM, ELEMENT_DIM > &rJacobian)
 
- Protected Member Functions inherited from AbstractElement< ELEMENT_DIM, SPACE_DIM >
void ConstructElementAttributes ()
 

Additional Inherited Members

- Protected Attributes inherited from AbstractElement< ELEMENT_DIM, SPACE_DIM >
std::vector< Node< SPACE_DIM > * > mNodes
 
unsigned mIndex
 
bool mIsDeleted
 
bool mOwnership
 
ElementAttributes< ELEMENT_DIM,
SPACE_DIM > * 
mpElementAttributes
 

Detailed Description

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

This abstract class defines a tetrahedral element for use in the Finite Element Method.

Definition at line 52 of file AbstractTetrahedralElement.hpp.

Constructor & Destructor Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::AbstractTetrahedralElement ( unsigned  index,
const std::vector< Node< SPACE_DIM > * > &  rNodes 
)

Constructor which takes in a vector of nodes.

Parameters
indexthe index of the element in the mesh
rNodesthe nodes owned by the element

Definition at line 63 of file AbstractTetrahedralElement.cpp.

References AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::CalculateJacobian(), AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::CalculateWeightedDirection(), and AbstractElement< ELEMENT_DIM, SPACE_DIM >::mNodes.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::AbstractTetrahedralElement ( unsigned  index = INDEX_IS_NOT_USED)

Default constructor, which doesn't fill in any nodes. The nodes must be added later.

Parameters
indexthe index of the element in the mesh (defaults to INDEX_IS_NOT_USED)

Definition at line 104 of file AbstractTetrahedralElement.cpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
virtual AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::~AbstractTetrahedralElement ( )
inlinevirtual

Virtual destructor, since this class has virtual methods.

Definition at line 84 of file AbstractTetrahedralElement.hpp.

Member Function Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
c_vector< double, SPACE_DIM > AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::CalculateCentroid ( ) const
Returns
the location of the centroid of the element.

Definition at line 218 of file AbstractTetrahedralElement.cpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::CalculateInverseJacobian ( c_matrix< double, SPACE_DIM, ELEMENT_DIM > &  rJacobian,
double rJacobianDeterminant,
c_matrix< double, ELEMENT_DIM, SPACE_DIM > &  rInverseJacobian 
)

Compute the Jacobian and the inverse Jacobian for this element.

Parameters
rJacobianthe Jacobian matrix (output)
rJacobianDeterminantthe determinant of the Jacobian (output)
rInverseJacobianthe inverse Jacobian matrix (output)

Definition at line 229 of file AbstractTetrahedralElement.cpp.

References Inverse().

Referenced by AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::AssembleOnCableElement(), and AbstractFunctionalCalculator< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::CalculateOnElement().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::CalculateJacobian ( c_matrix< double, SPACE_DIM, ELEMENT_DIM > &  rJacobian,
double rJacobianDeterminant 
)

Compute the Jacobian for this element.

Parameters
rJacobianthe Jacobian matrix
rJacobianDeterminantthe determinant of the Jacobian

Definition at line 109 of file AbstractTetrahedralElement.cpp.

References Determinant(), and EXCEPTION.

Referenced by AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::AbstractTetrahedralElement(), and MeshBasedCellPopulation< ELEMENT_DIM, SPACE_DIM >::GetVolumeOfCell().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
c_vector< double, SPACE_DIM > AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::CalculateNormal ( )
Returns
computed a unit vector normal to this element, if possible.

Definition at line 197 of file AbstractTetrahedralElement.cpp.

References EXCEPTION.

Referenced by AbstractNonlinearElasticitySolver< DIM >::AssembleOnBoundaryElementForPressureOnDeformedBc().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::CalculateWeightedDirection ( c_vector< double, SPACE_DIM > &  rWeightedDirection,
double rJacobianDeterminant 
)

Compute the weighted direction for this element.

Parameters
rWeightedDirectionthe weighted direction vector
rJacobianDeterminantthe determinant of the Jacobian

Definition at line 149 of file AbstractTetrahedralElement.cpp.

References EXCEPTION, NEVER_REACHED, and SubDeterminant().

Referenced by AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::AbstractTetrahedralElement(), and AbstractTetrahedralElement< 0, SPACE_DIM >::AbstractTetrahedralElement().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::GetStiffnessMatrixGlobalIndices ( unsigned  problemDim,
unsigned pIndices 
) const

Place in the pIndices array, the global indices (within the stiffness matrix) of the degrees of freedom associated with this element.

Parameters
problemDimthe problem dimension e.g. 2 for Bidomain.
pIndiceswhere to store results: an unsigned array with ELEMENT_DIM+1 entries.

Definition at line 262 of file AbstractTetrahedralElement.cpp.

Referenced by AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::DoAssemble(), AbstractFeSurfaceIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::DoAssemble(), and AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::DoAssemble().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
double AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::GetVolume ( double  determinant) const
Returns
the volume of an element (or area in 2d, or length in 1d)
Parameters
determinanta pre-calculated Jacobian determinant for this element.
Returns
volume (which is simply the determinant weighted by the SPACE_DIM)

Definition at line 240 of file AbstractTetrahedralElement.cpp.

Referenced by MeshBasedCellPopulation< ELEMENT_DIM, SPACE_DIM >::GetVolumeOfCell().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::RefreshJacobian ( c_matrix< double, SPACE_DIM, ELEMENT_DIM > &  rJacobian)
protected

Refresh the Jacobian for this element.

Parameters
rJacobianthe Jacobian matrix

Definition at line 47 of file AbstractTetrahedralElement.cpp.

References EXCEPTION.


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