Chaste
Release::3.4
|
#include <AbstractTetrahedralElement.hpp>
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 |
This abstract class defines a tetrahedral element for use in the Finite Element Method.
Definition at line 52 of file AbstractTetrahedralElement.hpp.
AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::AbstractTetrahedralElement | ( | unsigned | index, |
const std::vector< Node< SPACE_DIM > * > & | rNodes | ||
) |
Constructor which takes in a vector of nodes.
index | the index of the element in the mesh |
rNodes | the 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.
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.
index | the index of the element in the mesh (defaults to INDEX_IS_NOT_USED) |
Definition at line 104 of file AbstractTetrahedralElement.cpp.
|
inlinevirtual |
Virtual destructor, since this class has virtual methods.
Definition at line 84 of file AbstractTetrahedralElement.hpp.
c_vector< double, SPACE_DIM > AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::CalculateCentroid | ( | ) | const |
Definition at line 219 of file AbstractTetrahedralElement.cpp.
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.
rJacobian | the Jacobian matrix (output) |
rJacobianDeterminant | the determinant of the Jacobian (output) |
rInverseJacobian | the inverse Jacobian matrix (output) |
Definition at line 230 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().
void AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::CalculateJacobian | ( | c_matrix< double, SPACE_DIM, ELEMENT_DIM > & | rJacobian, |
double & | rJacobianDeterminant | ||
) |
Compute the Jacobian for this element.
rJacobian | the Jacobian matrix |
rJacobianDeterminant | the 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().
c_vector< double, SPACE_DIM > AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::CalculateNormal | ( | ) |
Definition at line 198 of file AbstractTetrahedralElement.cpp.
References EXCEPTION.
Referenced by AbstractNonlinearElasticitySolver< DIM >::AssembleOnBoundaryElementForPressureOnDeformedBc().
void AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::CalculateWeightedDirection | ( | c_vector< double, SPACE_DIM > & | rWeightedDirection, |
double & | rJacobianDeterminant | ||
) |
Compute the weighted direction for this element.
rWeightedDirection | the weighted direction vector |
rJacobianDeterminant | the 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().
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.
problemDim | the problem dimension e.g. 2 for Bidomain. |
pIndices | where to store results: an unsigned array with ELEMENT_DIM+1 entries. |
Definition at line 263 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().
double AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::GetVolume | ( | double | determinant | ) | const |
determinant | a pre-calculated Jacobian determinant for this element. |
Definition at line 241 of file AbstractTetrahedralElement.cpp.
Referenced by MeshBasedCellPopulation< ELEMENT_DIM, SPACE_DIM >::GetVolumeOfCell().
|
protected |
Refresh the Jacobian for this element.
rJacobian | the Jacobian matrix |
Definition at line 47 of file AbstractTetrahedralElement.cpp.
References EXCEPTION.