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

#include <AbstractTetrahedralMesh.hpp>

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

Classes

class  ElementIterator
 

Public Types

typedef std::vector< BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > * >::const_iterator BoundaryElementIterator
 
- Public Types inherited from AbstractMesh< ELEMENT_DIM, SPACE_DIM >
typedef std::vector< Node< SPACE_DIM > * >::const_iterator BoundaryNodeIterator
 

Public Member Functions

ElementIterator GetElementIteratorBegin (bool skipDeletedElements=true)
 
ElementIterator GetElementIteratorEnd ()
 
 AbstractTetrahedralMesh ()
 
virtual ~AbstractTetrahedralMesh ()
 
virtual unsigned GetNumElements () const
 
virtual unsigned GetNumLocalElements () const
 
virtual unsigned GetNumBoundaryElements () const
 
virtual unsigned GetNumLocalBoundaryElements () const
 
unsigned GetNumAllElements () const
 
unsigned GetNumAllBoundaryElements () const
 
virtual unsigned GetNumCableElements () const
 
virtual unsigned GetNumVertices () const
 
virtual unsigned GetMaximumNodeIndex ()
 
Element< ELEMENT_DIM, SPACE_DIM > * GetElement (unsigned index) const
 
BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > * GetBoundaryElement (unsigned index) const
 
virtual void ConstructFromMeshReader (AbstractMeshReader< ELEMENT_DIM, SPACE_DIM > &rMeshReader)=0
 
void ConstructFromMesh (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > &rOtherMesh)
 
BoundaryElementIterator GetBoundaryElementIteratorBegin () const
 
BoundaryElementIterator GetBoundaryElementIteratorEnd () const
 
virtual void GetInverseJacobianForElement (unsigned elementIndex, c_matrix< double, SPACE_DIM, ELEMENT_DIM > &rJacobian, double &rJacobianDeterminant, c_matrix< double, ELEMENT_DIM, SPACE_DIM > &rInverseJacobian) const
 
virtual void GetWeightedDirectionForBoundaryElement (unsigned elementIndex, c_vector< double, SPACE_DIM > &rWeightedDirection, double &rJacobianDeterminant) const
 
void CheckOutwardNormals ()
 
virtual void ConstructLinearMesh (unsigned width)
 
virtual void ConstructRectangularMesh (unsigned width, unsigned height, bool stagger=true)
 
virtual void ConstructCuboid (unsigned width, unsigned height, unsigned depth)
 
void ConstructRegularSlabMesh (double spaceStep, double width, double height=0, double depth=0)
 
void ConstructRegularSlabMeshWithDimensionSplit (unsigned dimension, double spaceStep, double width, double height=0, double depth=0)
 
virtual bool CalculateDesignatedOwnershipOfBoundaryElement (unsigned faceIndex)
 
virtual bool CalculateDesignatedOwnershipOfElement (unsigned elementIndex)
 
unsigned CalculateMaximumNodeConnectivityPerProcess () const
 
virtual void GetHaloNodeIndices (std::vector< unsigned > &rHaloIndices) const
 
void CalculateNodeExchange (std::vector< std::vector< unsigned > > &rNodesToSendPerProcess, std::vector< std::vector< unsigned > > &rNodesToReceivePerProcess)
 
virtual c_vector< double, 2 > CalculateMinMaxEdgeLengths ()
 
unsigned GetContainingElementIndex (const ChastePoint< SPACE_DIM > &rTestPoint, bool strict=false, std::set< unsigned > testElements=std::set< unsigned >(), bool onlyTryWithTestElements=false)
 
unsigned GetNearestElementIndexFromTestElements (const ChastePoint< SPACE_DIM > &rTestPoint, std::set< unsigned > testElements)
 
- Public Member Functions inherited from AbstractMesh< ELEMENT_DIM, SPACE_DIM >
NodeIterator GetNodeIteratorBegin (bool skipDeletedNodes=true)
 
NodeIterator GetNodeIteratorEnd ()
 
 AbstractMesh ()
 
virtual ~AbstractMesh ()
 
virtual unsigned GetNumNodes () const
 
unsigned GetNumBoundaryNodes () const
 
virtual unsigned GetNumAllNodes () const
 
unsigned GetNumNodeAttributes () const
 
Node< SPACE_DIM > * GetNode (unsigned index) const
 
virtual Node< SPACE_DIM > * GetNodeOrHaloNode (unsigned index) const
 
Node< SPACE_DIM > * GetNodeFromPrePermutationIndex (unsigned index) const
 
virtual void ReadNodesPerProcessorFile (const std::string &rNodesPerProcessorFile)
 
virtual DistributedVectorFactoryGetDistributedVectorFactory ()
 
virtual void SetDistributedVectorFactory (DistributedVectorFactory *pFactory)
 
virtual void PermuteNodes ()
 
BoundaryNodeIterator GetBoundaryNodeIteratorBegin () const
 
BoundaryNodeIterator GetBoundaryNodeIteratorEnd () const
 
std::string GetMeshFileBaseName () const
 
bool IsMeshOnDisk () const
 
const std::vector< unsigned > & rGetNodePermutation () const
 
virtual c_vector< double, SPACE_DIM > GetVectorFromAtoB (const c_vector< double, SPACE_DIM > &rLocationA, const c_vector< double, SPACE_DIM > &rLocationB)
 
double GetDistanceBetweenNodes (unsigned indexA, unsigned indexB)
 
virtual double GetWidth (const unsigned &rDimension) const
 
virtual ChasteCuboid< SPACE_DIM > CalculateBoundingBox () const
 
virtual unsigned GetNearestNodeIndex (const ChastePoint< SPACE_DIM > &rTestPoint)
 
virtual void Scale (const double xFactor=1.0, const double yFactor=1.0, const double zFactor=1.0)
 
virtual void Translate (const c_vector< double, SPACE_DIM > &rDisplacement)
 
void Translate (const double xMovement=0.0, const double yMovement=0.0, const double zMovement=0.0)
 
virtual void Rotate (c_matrix< double, SPACE_DIM, SPACE_DIM > rotationMatrix)
 
void Rotate (c_vector< double, 3 > axis, double angle)
 
void RotateX (const double theta)
 
void RotateY (const double theta)
 
void RotateZ (const double theta)
 
void Rotate (double theta)
 
virtual void RefreshMesh ()
 
bool IsMeshChanging () const
 
unsigned CalculateMaximumContainingElementsPerProcess () const
 
void SetMeshHasChangedSinceLoading ()
 

Protected Member Functions

void SetElementOwnerships ()
 
- Protected Member Functions inherited from AbstractMesh< ELEMENT_DIM, SPACE_DIM >
ChasteCuboid< SPACE_DIM > CalculateBoundingBox (const std::vector< Node< SPACE_DIM > * > &rNodes) const
 

Protected Attributes

bool mMeshIsLinear
 
std::vector< Element< ELEMENT_DIM, SPACE_DIM > * > mElements
 
std::vector< BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > * > mBoundaryElements
 
- Protected Attributes inherited from AbstractMesh< ELEMENT_DIM, SPACE_DIM >
std::vector< Node< SPACE_DIM > * > mNodes
 
std::vector< Node< SPACE_DIM > * > mBoundaryNodes
 
DistributedVectorFactorympDistributedVectorFactory
 
std::vector< unsignedmNodePermutation
 
std::string mMeshFileBaseName
 
bool mMeshChangesDuringSimulation
 

Private Member Functions

virtual unsigned SolveElementMapping (unsigned index) const =0
 
virtual unsigned SolveBoundaryElementMapping (unsigned index) const =0
 
template<class Archive >
void save (Archive &archive, const unsigned int version) const
 
template<class Archive >
void load (Archive &archive, const unsigned int version)
 

Friends

class AbstractConductivityTensors< ELEMENT_DIM, SPACE_DIM >
 
class CentroidWriter
 
class boost::serialization::access
 

Detailed Description

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

Abstract base class for all tetrahedral meshes (inherits from AbstractMesh).

Definition at line 71 of file AbstractTetrahedralMesh.hpp.

Member Typedef Documentation

◆ BoundaryElementIterator

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
typedef std::vector<BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>*>::const_iterator AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::BoundaryElementIterator

Definition of boundary element Iterator type.

Definition at line 298 of file AbstractTetrahedralMesh.hpp.

Constructor & Destructor Documentation

◆ AbstractTetrahedralMesh()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::AbstractTetrahedralMesh ( )

Constructor.

Definition at line 66 of file AbstractTetrahedralMesh.cpp.

◆ ~AbstractTetrahedralMesh()

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

Virtual destructor, since this class has virtual methods.

Definition at line 72 of file AbstractTetrahedralMesh.cpp.

Member Function Documentation

◆ CalculateDesignatedOwnershipOfBoundaryElement()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
bool AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::CalculateDesignatedOwnershipOfBoundaryElement ( unsigned  faceIndex)
virtual

Determine whether or not the current process owns node 0 of this boundary element (tie breaker to determine which process writes to file for when two or more share ownership of a face).

Parameters
faceIndexis the global index of the face
Returns
true if this process is designated owner

Reimplemented in DistributedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >, and DistributedTetrahedralMesh< DIM, DIM >.

Definition at line 684 of file AbstractTetrahedralMesh.cpp.

◆ CalculateDesignatedOwnershipOfElement()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
bool AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::CalculateDesignatedOwnershipOfElement ( unsigned  elementIndex)
virtual

Determine whether or not the current process owns node 0 of this element (tie breaker to determine which process writes to file for when two or more share ownership of an element).

Parameters
elementIndexis the global index of the element
Returns
true if this process is designated owner

Reimplemented in DistributedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >, and DistributedTetrahedralMesh< DIM, DIM >.

Definition at line 701 of file AbstractTetrahedralMesh.cpp.

Referenced by AbstractFunctionalCalculator< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::Calculate().

◆ CalculateMaximumNodeConnectivityPerProcess()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::CalculateMaximumNodeConnectivityPerProcess ( ) const
Returns
Iterates through local nodes and finds two representative nodes with a maximum number of containing elements for all locally owned nodes. The two representative nodes are the most connected nodes at the boundary and the most connected node in the interior. This is because a well-element-connected boundary node is likely to be connected to more nodes since fewer of them will overlap.) At these representative node the node connectivity (number of nodes in forward star) is determined.

Useful for determining FEM matrix fill.

Definition at line 718 of file AbstractTetrahedralMesh.cpp.

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

◆ CalculateMinMaxEdgeLengths()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
c_vector< double, 2 > AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::CalculateMinMaxEdgeLengths ( )
virtual

Computes the minimum and maximum lengths of the edges in the mesh. Overridden in Distributed case

Todo:
Should be const
Returns
The minimum and maximum edge lengths in the mesh

Reimplemented in DistributedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >, and DistributedTetrahedralMesh< DIM, DIM >.

Definition at line 968 of file AbstractTetrahedralMesh.cpp.

Referenced by DistributedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::CalculateMinMaxEdgeLengths(), and FineCoarseMeshPair< DIM >::SetUpBoxes().

◆ CalculateNodeExchange()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::CalculateNodeExchange ( std::vector< std::vector< unsigned > > &  rNodesToSendPerProcess,
std::vector< std::vector< unsigned > > &  rNodesToReceivePerProcess 
)

Get the nodes which will need to be exchanged between remote processes. If we have an element which node indices outside the local [mLo, mHi) region then we know that those nodes will need to be recieved from a remote process, while those inside the range [mLo, mHi) will need to be sent

Parameters
rNodesToSendPerProcess(output) a vector which will be of size GetNumProcs() where each internal vector except i=GetMyRank() contains an ordered list of indices of nodes to send to process i
rNodesToReceivePerProcess(output) a vector which will be of size GetNumProcs() for information to receive for process i
Todo:
#2426 This is where we should check that we DO HAVE a Tetrahedral (not Distributed mesh)

Definition at line 885 of file AbstractTetrahedralMesh.cpp.

References PetscTools::GetNumProcs().

◆ CheckOutwardNormals()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::CheckOutwardNormals ( )

Check whether mesh has outward-facing normals.

This throws a suitable exception if an inward facing normal is found.

(Returns nothing)

Definition at line 186 of file AbstractTetrahedralMesh.cpp.

References Node< SPACE_DIM >::ContainingElementsBegin(), Node< SPACE_DIM >::ContainingElementsEnd(), EXCEPTION, AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex(), and Node< SPACE_DIM >::rGetLocation().

◆ ConstructCuboid()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ConstructCuboid ( unsigned  width,
unsigned  height,
unsigned  depth 
)
virtual

Construct a 3D cuboid grid on [0,width]x[0,height]x[0,depth].

Parameters
widthwidth of the mesh (in the x-direction)
heightheight of the mesh (in the y-direction)
depthdepth of the mesh (in the z-direction).

In this method the width is also THE NUMBER OF ELEMENTS IN THE x-direction, and similarly with the y and z directions.

Overridden in DistributedTetrahedralMesh

Reimplemented in QuadraticMesh< DIM >, DistributedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >, and DistributedTetrahedralMesh< DIM, DIM >.

Definition at line 384 of file AbstractTetrahedralMesh.cpp.

Referenced by QuadraticMesh< DIM >::ConstructCuboid(), and DistributedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ConstructCuboid().

◆ ConstructFromMesh()

◆ ConstructFromMeshReader()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
virtual void AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ConstructFromMeshReader ( AbstractMeshReader< ELEMENT_DIM, SPACE_DIM > &  rMeshReader)
pure virtual

◆ ConstructLinearMesh()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ConstructLinearMesh ( unsigned  width)
virtual

Construct a 1D linear grid on [0,width]

ELEMENT_DIM must be equal to 1. If SPACE_DIM > 1 then the y & z default to 0.0 for every node.

Parameters
widthwidth of the mesh (in the x-direction)

In this method the width is also THE NUMBER OF ELEMENTS IN THE x-direction.

Overridden in DistributedTetrahedralMesh

Reimplemented in QuadraticMesh< DIM >, DistributedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >, and DistributedTetrahedralMesh< DIM, DIM >.

Definition at line 245 of file AbstractTetrahedralMesh.cpp.

Referenced by QuadraticMesh< DIM >::ConstructLinearMesh(), and DistributedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ConstructLinearMesh().

◆ ConstructRectangularMesh()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ConstructRectangularMesh ( unsigned  width,
unsigned  height,
bool  stagger = true 
)
virtual

Construct a 2D rectangular grid on [0,width]x[0,height].

Diagonals can be staggered so that there is no preferred diffusion propagation direction.

Parameters
widthwidth of the mesh (in the x-direction)
heightheight of the mesh (in the y-direction)
staggerwhether the mesh should 'jumble' up the elements (defaults to true)

In this method the width is also THE NUMBER OF ELEMENTS IN THE x-direction, and similarly with the y direction.

Overridden in DistributedTetrahedralMesh

Reimplemented in QuadraticMesh< DIM >, DistributedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >, and DistributedTetrahedralMesh< DIM, DIM >.

Definition at line 277 of file AbstractTetrahedralMesh.cpp.

Referenced by QuadraticMesh< DIM >::ConstructRectangularMesh(), and DistributedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ConstructRectangularMesh().

◆ ConstructRegularSlabMesh()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ConstructRegularSlabMesh ( double  spaceStep,
double  width,
double  height = 0,
double  depth = 0 
)

Create a 1D mesh on [0, width], 2D mesh on [0, width]x[0 height] with staggering or 3D mesh on [0, width]x[0 height]x[0 depth with a given axis-aligned space step. If SPACE_DIM > ELEMENT_DIM then the y & z default to 0.0 for every node.

Parameters
spaceStepThe axis-aligned space step
widthThe width (x-dimension)
heightThe height (y-dimension - ignored if ELEMENT_DIM is 1D)
depthThe depth (z-dimension -ignored in 1D and 2D)

Definition at line 569 of file AbstractTetrahedralMesh.cpp.

Referenced by CardiacElectroMechProbRegularGeom< DIM >::CardiacElectroMechProbRegularGeom(), and CuboidMeshConstructor< ELEMENT_DIM, SPACE_DIM >::Construct().

◆ ConstructRegularSlabMeshWithDimensionSplit()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ConstructRegularSlabMeshWithDimensionSplit ( unsigned  dimension,
double  spaceStep,
double  width,
double  height = 0,
double  depth = 0 
)

This is a wrapper method to ConstructRegularSlabMesh() which is useful for parallel distributed meshes. By default slabs are split across processes in the top dimension (y in 2d, z in 3d) but it may be more useful to split in x or y.

Create a 1D mesh on [0, width], 2D mesh on [0, width]x[0 height] with staggering or 3D mesh on [0, width]x[0 height]x[0 depth with a given axis-aligned space step. If SPACE_DIM > ELEMENT_DIM then the y & z default to 0.0 for every node.

Parameters
dimensionThe dimension/axis to be split over the processes. When dimension=0 the split is on x, dimension=1 indicates split on y etc. If dimension >= SPACE_DIM then an exception is thrown.
spaceStepThe axis-aligned space step
widthThe width (x-dimension)
heightThe height (y-dimension - ignored if ELEMENT_DIM is 1D)
depthThe depth (z-dimension -ignored in 1D and 2D)

Definition at line 617 of file AbstractTetrahedralMesh.cpp.

◆ GetBoundaryElement()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > * AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetBoundaryElement ( unsigned  index) const

Get the boundary element with a given index in the mesh.

Parameters
indexthe global index of the boundary element
Returns
a pointer to the boundary element.

Definition at line 148 of file AbstractTetrahedralMesh.cpp.

Referenced by AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ConstructFromMesh(), BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::MergeFromArchive(), and Cylindrical2dMesh::ReconstructCylindricalMesh().

◆ GetBoundaryElementIteratorBegin()

◆ GetBoundaryElementIteratorEnd()

◆ GetContainingElementIndex()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetContainingElementIndex ( const ChastePoint< SPACE_DIM > &  rTestPoint,
bool  strict = false,
std::set< unsigned testElements = std::set<unsigned>(),
bool  onlyTryWithTestElements = false 
)

Return the element index for the first element that contains a test point

Parameters
rTestPointreference to the point
strictShould the element returned contain the point in the interior and not on an edge/face/vertex (default = not strict)
testElementsa set of guesses for the element (a set of element indices), to be checked first for potential efficiency improvements. (default = empty set)
onlyTryWithTestElementsDo not continue with other elements after trying the with testElements (for cases where you know the testPoint must be in the set of test elements or maybe outside the mesh).
Returns
element index

Definition at line 992 of file AbstractTetrahedralMesh.cpp.

References EXCEPTION.

Referenced by AveragedSourceEllipticPde< DIM >::SetupSourceTerms(), AveragedSourceParabolicPde< DIM >::SetupSourceTerms(), and VolumeDependentAveragedSourceEllipticPde< DIM >::SetupSourceTerms().

◆ GetElement()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
Element< ELEMENT_DIM, SPACE_DIM > * AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElement ( unsigned  index) const

◆ GetElementIteratorBegin()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ElementIterator AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorBegin ( bool  skipDeletedElements = true)
inline
Returns
an iterator to the first element in the mesh.
Parameters
skipDeletedElementswhether to include deleted element

Definition at line 725 of file AbstractTetrahedralMesh.hpp.

Referenced by AbstractCorrectionTermAssembler< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::AbstractCorrectionTermAssembler(), QuadraturePointsGroup< DIM >::QuadraturePointsGroup(), QuadraticMeshHelper< DIM >::AddInternalNodesToElements(), AbstractFunctionalCalculator< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::Calculate(), AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::DoAssemble(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ExportToMesher(), Toroidal2dMesh::GenerateVectorsOfElementsStraddlingCylindricalPeriodicBoundaries(), Cylindrical2dMesh::GenerateVectorsOfElementsStraddlingPeriodicBoundaries(), Toroidal2dMesh::GenerateVectorsOfElementsStraddlingToroidalPeriodicBoundaries(), PapillaryFibreCalculator::GetRadiusVectors(), Cylindrical2dMesh::ReconstructCylindricalMesh(), Toroidal2dMesh::ReconstructCylindricalMesh(), Toroidal2dMesh::ReconstructToroidalMesh(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshJacobianCachedData(), Toroidal2dMesh::ReMesh(), PapillaryFibreCalculator::SmoothStructureTensors(), Cylindrical2dMesh::UseTheseElementsToDecideMeshing(), CellPopulationElementWriter< ELEMENT_DIM, SPACE_DIM >::Visit(), VtkMeshWriter< ELEMENT_DIM, SPACE_DIM >::WriteFilesUsingMesh(), and XdmfMeshWriter< ELEMENT_DIM, SPACE_DIM >::WriteFilesUsingMesh().

◆ GetElementIteratorEnd()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ElementIterator AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorEnd ( )
inline

◆ GetHaloNodeIndices()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetHaloNodeIndices ( std::vector< unsigned > &  rHaloIndices) const
virtual

Utility method to give the functionality of iterating through the halo nodes of a process. Will return an empty std::vector (i.e. no halo nodes) unless overridden by distributed derived classes.

Parameters
rHaloIndicesA vector to fill with the global indices of the nodes which are locally halos

Reimplemented in DistributedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >, and DistributedTetrahedralMesh< DIM, DIM >.

Definition at line 829 of file AbstractTetrahedralMesh.cpp.

◆ GetInverseJacobianForElement()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetInverseJacobianForElement ( unsigned  elementIndex,
c_matrix< double, SPACE_DIM, ELEMENT_DIM > &  rJacobian,
double rJacobianDeterminant,
c_matrix< double, ELEMENT_DIM, SPACE_DIM > &  rInverseJacobian 
) const
virtual

Compute the inverse Jacobian for a given element in the mesh.

Parameters
elementIndexindex of an element
rJacobianthe Jacobian matrix
rJacobianDeterminantthe determinant of the Jacobian matrix
rInverseJacobianthe inverse Jacobian matrix

Reimplemented in NonCachedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >, TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >, TetrahedralMesh< 3, 3 >, and TetrahedralMesh< DIM, DIM >.

Definition at line 167 of file AbstractTetrahedralMesh.cpp.

Referenced by AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::AssembleOnElement().

◆ GetMaximumNodeIndex()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetMaximumNodeIndex ( )
virtual
Returns
the largest index of nodes on this process. Overwritten in subclasses and used for setting up NodeMaps
the largest node index on this process.

Reimplemented in NodesOnlyMesh< SPACE_DIM >, NodesOnlyMesh< 2 >, and NodesOnlyMesh< DIM >.

Definition at line 135 of file AbstractTetrahedralMesh.cpp.

◆ GetNearestElementIndexFromTestElements()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNearestElementIndexFromTestElements ( const ChastePoint< SPACE_DIM > &  rTestPoint,
std::set< unsigned testElements 
)

As with GetNearestElementIndex() except only searches in the given set of elements.

Parameters
rTestPointreference to the point
testElementsa set of elements (element indices) to look in
Returns
element index

Definition at line 1039 of file AbstractTetrahedralMesh.cpp.

References EXCEPT_IF_NOT.

◆ GetNumAllBoundaryElements()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumAllBoundaryElements ( ) const
Returns
the total number of boundary elements (including those marked as deleted).

Definition at line 117 of file AbstractTetrahedralMesh.cpp.

Referenced by Cylindrical2dMesh::ReconstructCylindricalMesh(), and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshJacobianCachedData().

◆ GetNumAllElements()

◆ GetNumBoundaryElements()

◆ GetNumCableElements()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumCableElements ( ) const
virtual
Returns
the number of cable elements that are actually in use.

This will always return zero until overridden in the MixedDimensionMesh class

Reimplemented in MixedDimensionMesh< ELEMENT_DIM, SPACE_DIM >, and MixedDimensionMesh< ELEMENT_DIM, ELEMENT_DIM >.

Definition at line 123 of file AbstractTetrahedralMesh.cpp.

◆ GetNumElements()

◆ GetNumLocalBoundaryElements()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumLocalBoundaryElements ( ) const
virtual
Returns
the number of local boundary elements that are in use on this process (only over-ridden when the mesh is distributed).

Reimplemented in DistributedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >, and DistributedTetrahedralMesh< DIM, DIM >.

Definition at line 111 of file AbstractTetrahedralMesh.cpp.

Referenced by QuadraticMeshHelper< DIM >::AddInternalNodesToBoundaryElements().

◆ GetNumLocalElements()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumLocalElements ( ) const
virtual
Returns
the number of local elements that are in use on this process (only over-ridden when the mesh is distributed).

Reimplemented in DistributedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >, and DistributedTetrahedralMesh< DIM, DIM >.

Definition at line 93 of file AbstractTetrahedralMesh.cpp.

Referenced by QuadraticMeshHelper< DIM >::AddInternalNodesToElements().

◆ GetNumVertices()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumVertices ( ) const
virtual
Returns
the number of vertices (nodes which are also corners of elements). For a linear mesh all nodes are vertices, so this method is a synonym for GetNumNodes. However, it is over-ridden in quadratic meshes.

Reimplemented in QuadraticMesh< DIM >.

Definition at line 129 of file AbstractTetrahedralMesh.cpp.

Referenced by CmguiDeformedSolutionsWriter< DIM >::CmguiDeformedSolutionsWriter().

◆ GetWeightedDirectionForBoundaryElement()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetWeightedDirectionForBoundaryElement ( unsigned  elementIndex,
c_vector< double, SPACE_DIM > &  rWeightedDirection,
double rJacobianDeterminant 
) const
virtual

Compute the weighted direction for a given boundary element.

Parameters
elementIndexindex of an element
rWeightedDirectionthe weighted direction vector
rJacobianDeterminantthe determinant of the Jacobian matrix

Reimplemented in NonCachedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >, TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >, TetrahedralMesh< 3, 3 >, and TetrahedralMesh< DIM, DIM >.

Definition at line 177 of file AbstractTetrahedralMesh.cpp.

◆ load()

◆ save()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
template<class Archive >
void AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::save ( Archive &  archive,
const unsigned int  version 
) const
inlineprivate

◆ SetElementOwnerships()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::SetElementOwnerships ( )
protectedvirtual

Sets the ownership of each element according to which nodes are owned by the process.

Information on node ownership comes from the distributed vector factory and an element is "owned" if one or more of its nodes are owned

Reimplemented from AbstractMesh< ELEMENT_DIM, SPACE_DIM >.

Reimplemented in DistributedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >, and DistributedTetrahedralMesh< DIM, DIM >.

Definition at line 45 of file AbstractTetrahedralMesh.cpp.

References AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), and AbstractElement< ELEMENT_DIM, SPACE_DIM >::SetOwnership().

◆ SolveBoundaryElementMapping()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
virtual unsigned AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::SolveBoundaryElementMapping ( unsigned  index) const
privatepure virtual

Pure virtual solve boundary element mapping method. For a boundary element with a given global index, get the local index used by this process. Overridden in TetrahedralMesh and DistributedTetrahedralMesh classes.

Parameters
indexthe global index of the boundary element
Returns
local index

Implemented in DistributedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >, DistributedTetrahedralMesh< DIM, DIM >, TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >, TetrahedralMesh< 3, 3 >, and TetrahedralMesh< DIM, DIM >.

◆ SolveElementMapping()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
virtual unsigned AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::SolveElementMapping ( unsigned  index) const
privatepure virtual

Pure virtual solve element mapping method. For an element with a given global index, get the local index used by this process. Overridden in TetrahedralMesh and DistributedTetrahedralMesh classes.

Parameters
indexthe global index of the element
Returns
local index

Implemented in DistributedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >, DistributedTetrahedralMesh< DIM, DIM >, TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >, TetrahedralMesh< 3, 3 >, and TetrahedralMesh< DIM, DIM >.

Friends And Related Symbol Documentation

◆ AbstractConductivityTensors< ELEMENT_DIM, SPACE_DIM >

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
friend class AbstractConductivityTensors< ELEMENT_DIM, SPACE_DIM >
friend

Definition at line 800 of file AbstractTetrahedralMesh.hpp.

◆ boost::serialization::access

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
friend class boost::serialization::access
friend

Needed for serialization.

Definition at line 104 of file AbstractTetrahedralMesh.hpp.

◆ CentroidWriter

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
friend class CentroidWriter
friend

Definition at line 74 of file AbstractTetrahedralMesh.hpp.

Member Data Documentation

◆ mBoundaryElements

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
std::vector<BoundaryElement<ELEMENT_DIM-1, SPACE_DIM> *> AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryElements
protected

◆ mElements

◆ mMeshIsLinear

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
bool AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mMeshIsLinear
protected

Most tet meshes are linear (set to true). Set to false in quadratics.

Definition at line 79 of file AbstractTetrahedralMesh.hpp.

Referenced by QuadraticMesh< DIM >::QuadraticMesh(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::load(), and AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::save().


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