Chaste Commit::f2ff7ee04e70ac9d06c57344df8d017dbb12b97b
AbstractElement< ELEMENT_DIM, SPACE_DIM > Class Template Referenceabstract

#include <AbstractElement.hpp>

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

Public Member Functions

 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 ConstructElementAttributes ()
 

Protected Attributes

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 AbstractElement< ELEMENT_DIM, SPACE_DIM >

An abstract element class for use in finite element meshes.

Definition at line 55 of file AbstractElement.hpp.

Constructor & Destructor Documentation

◆ AbstractElement() [1/2]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
AbstractElement< ELEMENT_DIM, SPACE_DIM >::AbstractElement ( 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 48 of file AbstractElement.cpp.

◆ AbstractElement() [2/2]

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

Default constructor, which doesn't add any nodes: they must be added later.

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

Definition at line 60 of file AbstractElement.cpp.

◆ ~AbstractElement()

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

Virtual destructor, since this class has virtual methods. Deletes added attributes (when they exist)

Definition at line 68 of file AbstractElement.cpp.

Member Function Documentation

◆ AddElementAttribute()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractElement< ELEMENT_DIM, SPACE_DIM >::AddElementAttribute ( double  attribute)

Add an attribute to the list of element attributes.

Parameters
attributethe value of the attribute to be added

Definition at line 203 of file AbstractElement.cpp.

◆ AddNode()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractElement< ELEMENT_DIM, SPACE_DIM >::AddNode ( Node< SPACE_DIM > *  pNode)

Add a node to this element.

Parameters
pNodepointer to the new node

Definition at line 129 of file AbstractElement.cpp.

Referenced by QuadraticMeshHelper< DIM >::AddInternalNodesToElements(), and QuadraticMeshHelper< DIM >::AddNodeToBoundaryElement().

◆ ConstructElementAttributes()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractElement< ELEMENT_DIM, SPACE_DIM >::ConstructElementAttributes ( )
protected

Construct an empty ElementAttributes container.

Definition at line 194 of file AbstractElement.cpp.

◆ GetAttribute()

◆ GetIndex()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetIndex ( ) const
Returns
the index of this element

Definition at line 141 of file AbstractElement.cpp.

Referenced by AbstractCorrectionTermAssembler< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::AbstractCorrectionTermAssembler(), ImmersedBoundaryMesh< ELEMENT_DIM, SPACE_DIM >::ImmersedBoundaryMesh(), MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >::MutableVertexMesh(), PottsMesh< DIM >::PottsMesh(), VertexMesh< ELEMENT_DIM, SPACE_DIM >::VertexMesh(), MutableMesh< ELEMENT_DIM, SPACE_DIM >::AddElement(), PottsMesh< DIM >::AddElement(), MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >::AddElement(), FarhadifarForce< DIM >::AddForceContribution(), NagaiHondaForce< DIM >::AddForceContribution(), QuadraticMeshHelper< DIM >::AddInternalNodesToElements(), AbstractNonlinearElasticitySolver< DIM >::AssembleOnBoundaryElement(), ContinuumMechanicsNeumannBcsAssembler< DIM >::AssembleOnBoundaryElement(), AbstractNonlinearElasticitySolver< DIM >::AssembleOnBoundaryElementForPressureOnDeformedBc(), CompressibleNonlinearElasticitySolver< DIM >::AssembleOnElement(), IncompressibleNonlinearElasticitySolver< DIM >::AssembleOnElement(), AbstractContinuumMechanicsAssembler< DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >::AssembleOnElement(), AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::AssembleOnElement(), AbstractFeSurfaceIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::AssembleOnSurfaceElement(), ImmersedBoundaryLinearMembraneForce< DIM >::CalculateForcesOnElement(), ImmersedBoundaryMorseMembraneForce< DIM >::CalculateForcesOnElement(), MutableMesh< ELEMENT_DIM, SPACE_DIM >::CheckIsVoronoi(), MutableVertexMesh< 2, 2 >::ClearLocationsOfT1Swaps(), ElectrodesStimulusFactory< DIM >::ComputeElectrodeTotalFlux(), AveragedSourceEllipticPde< DIM >::ComputeLinearInUCoeffInSourceTerm(), MonodomainStiffnessMatrixAssembler< ELEMENT_DIM, SPACE_DIM >::ComputeMatrixTerm(), BidomainAssembler< ELEMENT_DIM, SPACE_DIM >::ComputeMatrixTerm(), ExtendedBidomainAssembler< ELEMENT_DIM, SPACE_DIM >::ComputeMatrixTerm(), AveragedSourceParabolicPde< DIM >::ComputeSourceTerm(), MutableMesh< ELEMENT_DIM, SPACE_DIM >::DeleteBoundaryNodeAt(), ImmersedBoundaryMesh< ELEMENT_DIM, SPACE_DIM >::DivideElement(), ImmersedBoundaryMesh< ELEMENT_DIM, SPACE_DIM >::DivideElementAlongGivenAxis(), ImmersedBoundaryMesh< ELEMENT_DIM, SPACE_DIM >::DivideElementAlongShortAxis(), AbstractContinuumMechanicsAssembler< DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >::DoAssemble(), AbstractCorrectionTermAssembler< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::ElementAssemblyCriterion(), AbstractNonlinearElasticitySolver< DIM >::GetElementCentroidStrain(), AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::Initialise(), MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >::PerformIntersectionSwap(), MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >::PerformT2Swap(), MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >::PerformT3Swap(), Cylindrical2dMesh::ReconstructCylindricalMesh(), VertexMeshOperationRecorder< ELEMENT_DIM, SPACE_DIM >::RecordCellDivideOperation(), VertexMeshOperationRecorder< ELEMENT_DIM, SPACE_DIM >::RecordEdgeMergeOperation(), VertexMeshOperationRecorder< ELEMENT_DIM, SPACE_DIM >::RecordEdgeSplitOperation(), VertexMeshOperationRecorder< ELEMENT_DIM, SPACE_DIM >::RecordNewEdgeOperation(), VertexMeshOperationRecorder< ELEMENT_DIM, SPACE_DIM >::RecordNodeMergeOperation(), ImmersedBoundaryMesh< ELEMENT_DIM, SPACE_DIM >::TagBoundaryElements(), StressPerElementWriter< DIM >::Visit(), and StreeterFibreGenerator< SPACE_DIM >::Visit().

◆ GetNode()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
Node< SPACE_DIM > * AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNode ( unsigned  localIndex) const

Get the node with a given local index in this element.

Parameters
localIndexlocal index of the node in this element
Returns
a pointer to the node.

Definition at line 116 of file AbstractElement.cpp.

Referenced by ImmersedBoundaryMesh< ELEMENT_DIM, SPACE_DIM >::ImmersedBoundaryMesh(), MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >::MutableVertexMesh(), PottsMesh< DIM >::PottsMesh(), VertexElement< ELEMENT_DIM, SPACE_DIM >::VertexElement(), VertexElement< ELEMENT_DIM - 1, SPACE_DIM >::VertexElement(), VertexMesh< ELEMENT_DIM, SPACE_DIM >::VertexMesh(), ImmersedBoundaryMesh< DIM, DIM >::~ImmersedBoundaryMesh(), FarhadifarForce< DIM >::AddForceContribution(), NagaiHondaForce< DIM >::AddForceContribution(), BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::AddNeumannBoundaryCondition(), QuadraticMeshHelper< DIM >::AddNodeToBoundaryElement(), BidomainProblem< DIM >::AnalyseMeshForBath(), AbstractNonlinearElasticitySolver< DIM >::AssembleOnBoundaryElementForPressureOnDeformedBc(), AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::AssembleOnCableElement(), AbstractContinuumMechanicsAssembler< DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >::AssembleOnElement(), AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::AssembleOnElement(), AbstractFeSurfaceIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::AssembleOnSurfaceElement(), ImmersedBoundaryLinearMembraneForce< DIM >::CalculateForcesOnElement(), ImmersedBoundaryMorseMembraneForce< DIM >::CalculateForcesOnElement(), AbstractFunctionalCalculator< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::CalculateOnElement(), MutableMesh< ELEMENT_DIM, SPACE_DIM >::CheckIsVoronoi(), ImmersedBoundaryMesh< DIM, DIM >::Clear(), MutableMesh< ELEMENT_DIM, SPACE_DIM >::DeleteBoundaryNodeAt(), MutableVertexMesh< 2, 2 >::DeleteElementPriorToReMesh(), ImmersedBoundaryMesh< ELEMENT_DIM, SPACE_DIM >::DivideElement(), PottsMesh< DIM >::DivideElement(), StreeterFibreGenerator< SPACE_DIM >::GetAveragedThicknessLocalNode(), MixedDimensionMesh< ELEMENT_DIM, ELEMENT_DIM >::GetCablesAtNode(), PottsMesh< DIM >::GetNeighbouringElementIndices(), VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetNeighbouringElementIndices(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::EdgeIterator::GetNodeA(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::EdgeIterator::GetNodeB(), PottsMesh< SPACE_DIM >::GetNumNodes(), VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetRosetteRankOfElement(), ImmersedBoundaryMesh< DIM, DIM >::GetVolumeOfElement(), Cylindrical2dMesh::ReMesh(), Toroidal2dMesh::ReMesh(), ImmersedBoundaryMesh< ELEMENT_DIM, SPACE_DIM >::ReMeshElement(), ImmersedBoundaryMesh< ELEMENT_DIM, SPACE_DIM >::ReMeshLamina(), FineCoarseMeshPair< DIM >::SetUpBoxes(), VertexMesh< 2, 2 >::SolveBoundaryElementMapping(), StreeterFibreGenerator< SPACE_DIM >::Visit(), and DistanceMapCalculator< ELEMENT_DIM, SPACE_DIM >::WorkOnLocalQueue().

◆ GetNodeGlobalIndex()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex ( unsigned  localIndex) const

Given the local index of a node owned by this element, return the global index of the node in the mesh.

Parameters
localIndexthe node's local index in this element
Returns
the global index

Definition at line 109 of file AbstractElement.cpp.

Referenced by VoltageInterpolaterOntoMechanicsMesh< DIM >::VoltageInterpolaterOntoMechanicsMesh(), QuadraticMeshHelper< DIM >::AddExtraBoundaryNodes(), QuadraticMeshHelper< DIM >::AddNodesToBoundaryElements(), AbstractNonlinearElasticitySolver< DIM >::AssembleOnBoundaryElement(), AbstractNonlinearElasticitySolver< DIM >::AssembleOnBoundaryElementForPressureOnDeformedBc(), AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::AssembleOnCableElement(), CompressibleNonlinearElasticitySolver< DIM >::AssembleOnElement(), IncompressibleNonlinearElasticitySolver< DIM >::AssembleOnElement(), AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::AssembleOnElement(), CompressibleNonlinearElasticitySolver< DIM >::AssembleSystem(), IncompressibleNonlinearElasticitySolver< DIM >::AssembleSystem(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::CalculateMaximumNodeConnectivityPerProcess(), AbstractFunctionalCalculator< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::CalculateOnElement(), MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >::CheckForIntersections(), MutableMesh< ELEMENT_DIM, SPACE_DIM >::CheckIsVoronoi(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::CheckOutwardNormals(), AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::ComputeDeformationGradientAndStretchInEachElement(), MonodomainCorrectionTermAssembler< ELEM_DIM, SPACE_DIM >::ComputeVectorTerm(), BidomainCorrectionTermAssembler< ELEM_DIM, SPACE_DIM >::ComputeVectorTerm(), EllipticBoxDomainPdeModifier< DIM >::ConstructBoundaryConditionsContainer(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ConstructFromMesh(), Toroidal2dMesh::CorrectCylindricalNonPeriodicMesh(), Cylindrical2dMesh::CorrectNonPeriodicMesh(), Toroidal2dMesh::CorrectToroidalNonPeriodicMesh(), AbstractContinuumMechanicsAssembler< DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >::DoAssemble(), AbstractCorrectionTermAssembler< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::ElementAssemblyCriterion(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ExportToMesher(), AbstractNonlinearElasticitySolver< DIM >::GetElementCentroidStrain(), MeshBasedCellPopulation< ELEMENT_DIM, SPACE_DIM >::GetNeighbouringNodeIndices(), VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetNeighbouringNodeNotAlsoInElement(), AbstractTetrahedralMeshWriter< ELEMENT_DIM, SPACE_DIM >::GetNextCableElement(), AbstractTetrahedralMeshWriter< ELEMENT_DIM, SPACE_DIM >::GetNextElement(), ImmersedBoundaryMesh< ELEMENT_DIM, SPACE_DIM >::GetSurfaceAreaOfElement(), VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetSurfaceAreaOfElement(), ImmersedBoundaryCellPopulation< DIM >::GetTetrahedralMeshForPdeModifier(), VertexBasedCellPopulation< DIM >::GetTetrahedralMeshForPdeModifier(), ImmersedBoundaryMesh< ELEMENT_DIM, SPACE_DIM >::GetVoronoiSurfaceAreaOfElement(), QuadraticMeshHelper< DIM >::HelperMethod1(), MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >::IdentifySwapType(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::EdgeIterator::operator++(), MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >::PerformIntersectionSwap(), MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >::PerformProtorosetteResolution(), MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >::PerformRosetteRankDecrease(), MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >::PerformT3Swap(), Cylindrical2dMesh::ReconstructCylindricalMesh(), AbstractCorrectionTermAssembler< ELEM_DIM, SPACE_DIM, 2 >::ResetInterpolatedQuantities(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::SetElementOwnerships(), MutableMesh< 2, 2 >::SetNode(), CellwiseDataGradient< DIM >::SetupGradients(), CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::Solve(), MutableMesh< ELEMENT_DIM, SPACE_DIM >::SplitEdge(), AbstractBoxDomainPdeModifier< DIM >::UpdateCellData(), Cylindrical2dMesh::UseTheseElementsToDecideMeshing(), CellPopulationElementWriter< ELEMENT_DIM, SPACE_DIM >::Visit(), CellPopulationElementWriter< ELEMENT_DIM, SPACE_DIM >::Visit(), and CellPopulationElementWriter< ELEMENT_DIM, SPACE_DIM >::Visit().

◆ GetNodeLocation() [1/2]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
c_vector< double, SPACE_DIM > AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeLocation ( unsigned  localIndex) const
Returns
the location in space of one of the nodes in this element.
Parameters
localIndexthe index of the node to query, in [0,N) where N is the number of nodes in this element.

Definition at line 102 of file AbstractElement.cpp.

◆ GetNodeLocation() [2/2]

◆ GetNumElementAttributes()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNumElementAttributes ( )
Returns
the number of node attributes associated with this node.

Definition at line 222 of file AbstractElement.cpp.

Referenced by ImmersedBoundaryMesh< ELEMENT_DIM, SPACE_DIM >::DivideElement().

◆ GetNumNodes()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNumNodes ( ) const
Returns
the number of nodes owned by this element.

Definition at line 123 of file AbstractElement.cpp.

Referenced by ImmersedBoundaryMesh< ELEMENT_DIM, SPACE_DIM >::ImmersedBoundaryMesh(), MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >::MutableVertexMesh(), PottsMesh< DIM >::PottsMesh(), VertexMesh< ELEMENT_DIM, SPACE_DIM >::VertexMesh(), VoltageInterpolaterOntoMechanicsMesh< DIM >::VoltageInterpolaterOntoMechanicsMesh(), FarhadifarForce< DIM >::AddForceContribution(), NagaiHondaForce< DIM >::AddForceContribution(), QuadraticMeshHelper< DIM >::AddInternalNodesToElements(), BidomainProblem< DIM >::AnalyseMeshForBath(), AbstractNonlinearElasticitySolver< DIM >::AssembleOnBoundaryElementForPressureOnDeformedBc(), AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::AssembleOnCableElement(), AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::AssembleOnElement(), AbstractFeSurfaceIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::AssembleOnSurfaceElement(), ImmersedBoundaryLinearMembraneForce< DIM >::CalculateForcesOnElement(), ImmersedBoundaryMorseMembraneForce< DIM >::CalculateForcesOnElement(), ImmersedBoundaryMesh< ELEMENT_DIM, SPACE_DIM >::CalculateMomentsOfElement(), VertexMesh< ELEMENT_DIM, SPACE_DIM >::CalculateMomentsOfElement(), AbstractFunctionalCalculator< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::CalculateOnElement(), MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >::CheckForIntersections(), MutableMesh< ELEMENT_DIM, SPACE_DIM >::CheckIsVoronoi(), MutableVertexMesh< 2, 2 >::ClearLocationsOfIntersectionSwaps(), ElectrodesStimulusFactory< DIM >::ComputeElectrodeTotalFlux(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ConstructFromMesh(), MutableMesh< ELEMENT_DIM, SPACE_DIM >::DeleteBoundaryNodeAt(), MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >::DivideEdge(), ImmersedBoundaryMesh< ELEMENT_DIM, SPACE_DIM >::DivideElement(), PottsMesh< DIM >::DivideElement(), ImmersedBoundaryMesh< ELEMENT_DIM, SPACE_DIM >::DivideElementAlongGivenAxis(), VertexMesh< ELEMENT_DIM, SPACE_DIM >::ElementIncludesPoint(), StreeterFibreGenerator< SPACE_DIM >::GetAveragedThicknessLocalNode(), PottsMesh< DIM >::GetCentroidOfElement(), ImmersedBoundaryMesh< ELEMENT_DIM, SPACE_DIM >::GetCentroidOfElement(), VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetCentroidOfElement(), VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetLocalIndexForElementEdgeClosestToPoint(), PottsMesh< DIM >::GetNeighbouringElementIndices(), VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetNeighbouringElementIndices(), MeshBasedCellPopulation< ELEMENT_DIM, SPACE_DIM >::GetNeighbouringNodeIndices(), VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetNeighbouringNodeNotAlsoInElement(), VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetRosetteRankOfElement(), ImmersedBoundaryMesh< ELEMENT_DIM, SPACE_DIM >::GetSurfaceAreaOfElement(), VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetSurfaceAreaOfElement(), ImmersedBoundaryCellPopulation< DIM >::GetTetrahedralMeshForPdeModifier(), VertexBasedCellPopulation< DIM >::GetTetrahedralMeshForPdeModifier(), PottsMesh< DIM >::GetVolumeOfElement(), ImmersedBoundaryMesh< ELEMENT_DIM, SPACE_DIM >::GetVolumeOfElement(), VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetVolumeOfElement(), ImmersedBoundaryMesh< ELEMENT_DIM, SPACE_DIM >::GetVoronoiSurfaceAreaOfElement(), MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >::IdentifySwapType(), MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >::PerformIntersectionSwap(), MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >::PerformProtorosetteResolution(), MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >::PerformRosetteRankDecrease(), MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >::PerformT2Swap(), MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >::PerformT3Swap(), ImmersedBoundaryMesh< ELEMENT_DIM, SPACE_DIM >::ReMeshElement(), ImmersedBoundaryMesh< ELEMENT_DIM, SPACE_DIM >::ReMeshLamina(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::SetElementOwnerships(), CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::Solve(), VertexMesh< 2, 2 >::SolveBoundaryElementMapping(), CellPopulationElementWriter< ELEMENT_DIM, SPACE_DIM >::Visit(), CellPopulationElementWriter< ELEMENT_DIM, SPACE_DIM >::Visit(), CellPopulationElementWriter< ELEMENT_DIM, SPACE_DIM >::Visit(), and DistanceMapCalculator< ELEMENT_DIM, SPACE_DIM >::WorkOnLocalQueue().

◆ GetOwnership()

◆ GetUnsignedAttribute()

◆ IsDeleted()

◆ MarkAsDeleted()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
virtual void AbstractElement< ELEMENT_DIM, SPACE_DIM >::MarkAsDeleted ( )
pure virtual

◆ RegisterWithNodes()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
virtual void AbstractElement< ELEMENT_DIM, SPACE_DIM >::RegisterWithNodes ( )
pure virtual

◆ ReplaceNode()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractElement< ELEMENT_DIM, SPACE_DIM >::ReplaceNode ( Node< SPACE_DIM > *  pOldNode,
Node< SPACE_DIM > *  pNewNode 
)

Replace one of the nodes in this element with another.

Parameters
pOldNodepointer to the current node
pNewNodepointer to the replacement node

Definition at line 74 of file AbstractElement.cpp.

Referenced by MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >::PerformT2Swap(), Cylindrical2dMesh::ReconstructCylindricalMesh(), and MutableMesh< ELEMENT_DIM, SPACE_DIM >::SplitEdge().

◆ rGetElementAttributes()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
std::vector< double > & AbstractElement< ELEMENT_DIM, SPACE_DIM >::rGetElementAttributes ( )
Returns
reference to a vector containing the element attributes.

Definition at line 211 of file AbstractElement.cpp.

Referenced by ImmersedBoundaryMesh< ELEMENT_DIM, SPACE_DIM >::DivideElement().

◆ SetAttribute()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractElement< ELEMENT_DIM, SPACE_DIM >::SetAttribute ( double  attribute)

Set an attribute (a value associated with the element)

Parameters
attributethe value of an attribute

Definition at line 165 of file AbstractElement.cpp.

Referenced by Toroidal2dVertexMesh::ConstructFromMeshReader(), PottsMesh< DIM >::ConstructFromMeshReader(), and TetrahedralMesh< DIM, DIM >::SolveBoundaryElementMapping().

◆ SetIndex()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractElement< ELEMENT_DIM, SPACE_DIM >::SetIndex ( unsigned  index)

Set the index of this element in the mesh.

Parameters
indexthe new index

Definition at line 147 of file AbstractElement.cpp.

◆ SetOwnership()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractElement< ELEMENT_DIM, SPACE_DIM >::SetOwnership ( bool  ownership)

Set whether the current process owns this element.

Parameters
ownershipwhether the current process now owns this element

Definition at line 159 of file AbstractElement.cpp.

Referenced by AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::SetElementOwnerships().

◆ UpdateNode()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
virtual void AbstractElement< ELEMENT_DIM, SPACE_DIM >::UpdateNode ( const unsigned rIndex,
Node< SPACE_DIM > *  pNode 
)
pure virtual

Member Data Documentation

◆ mIndex

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned AbstractElement< ELEMENT_DIM, SPACE_DIM >::mIndex
protected

The index of this element within the mesh

Definition at line 63 of file AbstractElement.hpp.

◆ mIsDeleted

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
bool AbstractElement< ELEMENT_DIM, SPACE_DIM >::mIsDeleted
protected

Whether this element has been deleted, and hence whether its location in the mesh can be re-used.

Definition at line 69 of file AbstractElement.hpp.

◆ mNodes

◆ mOwnership

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
bool AbstractElement< ELEMENT_DIM, SPACE_DIM >::mOwnership
protected

Whether the current process owns this element.

Definition at line 72 of file AbstractElement.hpp.

◆ mpElementAttributes

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
ElementAttributes<ELEMENT_DIM, SPACE_DIM>* AbstractElement< ELEMENT_DIM, SPACE_DIM >::mpElementAttributes
protected

A pointer to an ElementAttributes object associated with this element.

Definition at line 75 of file AbstractElement.hpp.


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