Chaste Commit::ca8ccdedf819b6e02855bc0e8e6f50bdecbc5208
|
#include <VertexMesh.hpp>
Classes | |
class | VertexElementIterator |
Public Member Functions | |
VertexElementIterator | GetElementIteratorBegin (bool skipDeletedElements=true) |
VertexElementIterator | GetElementIteratorEnd () |
VertexMesh (std::vector< Node< SPACE_DIM > * > nodes, std::vector< VertexElement< ELEMENT_DIM, SPACE_DIM > * > vertexElements) | |
VertexMesh (std::vector< Node< SPACE_DIM > * > nodes, std::vector< VertexElement< ELEMENT_DIM - 1, SPACE_DIM > * > faces, std::vector< VertexElement< ELEMENT_DIM, SPACE_DIM > * > vertexElements) | |
VertexMesh (TetrahedralMesh< 2, 2 > &rMesh, bool isPeriodic=false, bool isBounded=false, bool scaleBoundByEdgeLength=true, double maxDelaunayEdgeLength=DBL_MAX, bool offsetNewBoundaryNodes=false) | |
Alternative 2D 'Voronoi' constructor. | |
VertexMesh (TetrahedralMesh< 3, 3 > &rMesh) | |
VertexMesh () | |
virtual | ~VertexMesh () |
unsigned | GetNumEdges () const |
Edge< SPACE_DIM > * | GetEdge (unsigned index) const |
const EdgeHelper< SPACE_DIM > & | rGetEdgeHelper () const |
virtual unsigned | GetNumNodes () const |
virtual unsigned | GetNumElements () const |
unsigned | GetNumAllElements () const |
virtual unsigned | GetNumFaces () const |
VertexElement< ELEMENT_DIM, SPACE_DIM > * | GetElement (unsigned index) const |
VertexElement< ELEMENT_DIM - 1, SPACE_DIM > * | GetFace (unsigned index) const |
virtual c_vector< double, SPACE_DIM > | GetCentroidOfElement (unsigned index) |
void | ConstructFromMeshReader (AbstractMeshReader< ELEMENT_DIM, SPACE_DIM > &rMeshReader) |
virtual void | Clear () |
unsigned | GetDelaunayNodeIndexCorrespondingToVoronoiElementIndex (unsigned elementIndex) |
unsigned | GetVoronoiElementIndexCorrespondingToDelaunayNodeIndex (unsigned nodeIndex) |
unsigned | GetRosetteRankOfElement (unsigned index) |
virtual c_vector< double, SPACE_DIM > | GetVectorFromAtoB (const c_vector< double, SPACE_DIM > &rLocationA, const c_vector< double, SPACE_DIM > &rLocationB) |
virtual double | GetVolumeOfElement (unsigned index) |
virtual double | GetSurfaceAreaOfElement (unsigned index) |
c_vector< double, SPACE_DIM > | GetAreaGradientOfElementAtNode (VertexElement< ELEMENT_DIM, SPACE_DIM > *pElement, unsigned localIndex) |
c_vector< double, SPACE_DIM > | GetPreviousEdgeGradientOfElementAtNode (VertexElement< ELEMENT_DIM, SPACE_DIM > *pElement, unsigned localIndex) |
c_vector< double, SPACE_DIM > | GetNextEdgeGradientOfElementAtNode (VertexElement< ELEMENT_DIM, SPACE_DIM > *pElement, unsigned localIndex) |
c_vector< double, SPACE_DIM > | GetPerimeterGradientOfElementAtNode (VertexElement< ELEMENT_DIM, SPACE_DIM > *pElement, unsigned localIndex) |
virtual c_vector< double, 3 > | CalculateMomentsOfElement (unsigned index) |
double | GetEdgeLength (unsigned elementIndex1, unsigned elementIndex2) |
double | GetElongationShapeFactorOfElement (unsigned elementIndex) |
double | CalculateUnitNormalToFaceWithArea (VertexElement< ELEMENT_DIM - 1, SPACE_DIM > *pFace, c_vector< double, SPACE_DIM > &rNormal) |
virtual double | CalculateAreaOfFace (VertexElement< ELEMENT_DIM - 1, SPACE_DIM > *pFace) |
c_vector< double, SPACE_DIM > | GetShortAxisOfElement (unsigned index) |
std::set< unsigned > | GetNeighbouringNodeIndices (unsigned nodeIndex) |
std::set< unsigned > | GetNeighbouringNodeNotAlsoInElement (unsigned nodeIndex, unsigned elemIndex) |
std::set< unsigned > | GetNeighbouringElementIndices (unsigned elementIndex) |
virtual VertexMesh< ELEMENT_DIM, SPACE_DIM > * | GetMeshForVtk () |
bool | IsNearExistingNodes (c_vector< double, SPACE_DIM > newNodeLocation, std::vector< Node< SPACE_DIM > * > nodesToCheck, double minClearance) |
Public Member Functions inherited from AbstractMesh< ELEMENT_DIM, SPACE_DIM > | |
NodeIterator | GetNodeIteratorBegin (bool skipDeletedNodes=true) |
NodeIterator | GetNodeIteratorEnd () |
AbstractMesh () | |
virtual | ~AbstractMesh () |
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 DistributedVectorFactory * | GetDistributedVectorFactory () |
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 |
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 | |
unsigned | SolveNodeMapping (unsigned index) const |
unsigned | SolveElementMapping (unsigned index) const |
unsigned | SolveBoundaryElementMapping (unsigned index) const |
void | GenerateEdgesFromElements (std::vector< VertexElement< ELEMENT_DIM, SPACE_DIM > * > &rElements) |
void | GenerateVerticesFromElementCircumcentres (TetrahedralMesh< ELEMENT_DIM, SPACE_DIM > &rMesh) |
bool | ElementIncludesPoint (const c_vector< double, SPACE_DIM > &rTestPoint, unsigned elementIndex) |
unsigned | GetLocalIndexForElementEdgeClosestToPoint (const c_vector< double, SPACE_DIM > &rTestPoint, unsigned elementIndex) |
template<class Archive > | |
void | save (Archive &archive, const unsigned int version) const |
template<class Archive > | |
void | load (Archive &archive, const unsigned int version) |
Protected Member Functions inherited from AbstractMesh< ELEMENT_DIM, SPACE_DIM > | |
virtual void | SetElementOwnerships () |
ChasteCuboid< SPACE_DIM > | CalculateBoundingBox (const std::vector< Node< SPACE_DIM > * > &rNodes) const |
Protected Attributes | |
std::vector< VertexElement< ELEMENT_DIM, SPACE_DIM > * > | mElements |
std::vector< VertexElement< ELEMENT_DIM - 1, SPACE_DIM > * > | mFaces |
EdgeHelper< SPACE_DIM > | mEdgeHelper |
std::map< unsigned, unsigned > | mVoronoiElementIndexMap |
TetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * | mpDelaunayMesh |
Protected Attributes inherited from AbstractMesh< ELEMENT_DIM, SPACE_DIM > | |
std::vector< Node< SPACE_DIM > * > | mNodes |
std::vector< Node< SPACE_DIM > * > | mBoundaryNodes |
DistributedVectorFactory * | mpDistributedVectorFactory |
std::vector< unsigned > | mNodePermutation |
std::string | mMeshFileBaseName |
bool | mMeshChangesDuringSimulation |
Friends | |
class | TestVertexMesh |
class | boost::serialization::access |
Additional Inherited Members | |
Public Types inherited from AbstractMesh< ELEMENT_DIM, SPACE_DIM > | |
typedef std::vector< Node< SPACE_DIM > * >::const_iterator | BoundaryNodeIterator |
A vertex-based mesh class, in which elements may contain different numbers of nodes. This is facilitated by the VertexElement class.
This class has two applications in the cell_based code.
First, VertexMesh is used as a member of the MeshBasedCellPopulation class to represent a Voronoi tessellation, the dual to a Delaunay mesh, which allows the shapes of cells to be visualised in simulations of a class of off-lattice cell centre-based models.
Second, VertexMesh serves as a parent class for MutableVertexMesh, which is used as a member of the VertexBasedCellPopulation class to represent the junctional network of cells that forms the basis of simulations of off-lattice vertex-based models.
Definition at line 77 of file VertexMesh.hpp.
VertexMesh< ELEMENT_DIM, SPACE_DIM >::VertexMesh | ( | std::vector< Node< SPACE_DIM > * > | nodes, |
std::vector< VertexElement< ELEMENT_DIM, SPACE_DIM > * > | vertexElements | ||
) |
Default constructor.
nodes | vector of pointers to nodes |
vertexElements | vector of pointers to VertexElements |
Definition at line 45 of file VertexMesh.cpp.
References VertexMesh< ELEMENT_DIM, SPACE_DIM >::Clear(), VertexMesh< ELEMENT_DIM, SPACE_DIM >::GenerateEdgesFromElements(), VertexElement< ELEMENT_DIM, SPACE_DIM >::GetFace(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetIndex(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNode(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), VertexMesh< ELEMENT_DIM, SPACE_DIM >::mElements, VertexMesh< ELEMENT_DIM, SPACE_DIM >::mFaces, AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mMeshChangesDuringSimulation, and AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mNodes.
VertexMesh< ELEMENT_DIM, SPACE_DIM >::VertexMesh | ( | std::vector< Node< SPACE_DIM > * > | nodes, |
std::vector< VertexElement< ELEMENT_DIM - 1, SPACE_DIM > * > | faces, | ||
std::vector< VertexElement< ELEMENT_DIM, SPACE_DIM > * > | vertexElements | ||
) |
Constructor.
nodes | vector of pointers to nodes |
faces | vector of pointer to VertexElements |
vertexElements | vector of pointers to VertexElement<3,3>s |
Definition at line 110 of file VertexMesh.cpp.
References VertexMesh< ELEMENT_DIM, SPACE_DIM >::Clear().
VertexMesh< ELEMENT_DIM, SPACE_DIM >::VertexMesh | ( | TetrahedralMesh< 2, 2 > & | rMesh, |
bool | isPeriodic = false , |
||
bool | isBounded = false , |
||
bool | scaleBoundByEdgeLength = true , |
||
double | maxDelaunayEdgeLength = DBL_MAX , |
||
bool | offsetNewBoundaryNodes = false |
||
) |
Alternative 2D 'Voronoi' constructor.
This VertexMesh constructor is currently only defined for 2D meshes.
Creates a Voronoi tessellation of a given tetrahedral mesh, which must be Delaunay (see TetrahedralMesh::CheckIsVoronoi).
rMesh | a tetrahedral mesh |
isPeriodic | a boolean that indicates whether the mesh is periodic or not. Defaults to false. |
isBounded | a boolean to indicate whether to bound the voronoi tesselation. Defaults to false. |
scaleBoundByEdgeLength | whether to scale the distance bounding nodes are placed from the mesh. Defaults to true. |
maxDelaunayEdgeLength | the maximum edge length in the mesh. Edges longer than this are ignored in boundary calculation. Defaults to DBL_MAX so there ia no max length. |
offsetNewBoundaryNodes | whether to add new node towards the centre of the boundary edges or not. Defaults to false. |
VertexMesh< ELEMENT_DIM, SPACE_DIM >::VertexMesh | ( | TetrahedralMesh< 3, 3 > & | rMesh | ) |
Alternative 3D 'Voronoi' constructor. Creates a Voronoi tessellation of a given tetrahedral mesh, which must be Delaunay (see TetrahedralMesh::CheckIsVoronoi).
rMesh | a tetrahedral mesh |
VertexMesh< ELEMENT_DIM, SPACE_DIM >::VertexMesh | ( | ) |
Default constructor for use by serializer.
Definition at line 676 of file VertexMesh.cpp.
|
virtual |
Destructor.
Definition at line 684 of file VertexMesh.cpp.
|
virtual |
Get the area of a given face in 3D. Uses CalculateUnitNormalToFaceWithArea
This needs to be overridden in daughter classes for non-Euclidean metrics.
pFace | a face in the mesh |
Definition at line 1814 of file VertexMesh.cpp.
References NEVER_REACHED.
|
virtual |
Compute the second moments and product moment of area for a given 2D element about its centroid. These are:
I_xx, the second moment of area about an axis through the centroid of the element parallel to the x-axis;
I_yy, the second moment of area about an axis through the centroid of the element parallel to the y-axis;
and I_xy, product moment of area through the centroid of the element.
Formulae for these quantities may be found e.g. in the following reference:
Mechanics of Materials James M. Gere (Author), Barry J. Goodno. Cengage Learning; 8th edition (January 1, 2012)
This method is used within GetShortAxisOfElement() to compute the direction of the shortest principal axis passing through the centroid, or 'short axis', of the element.
Note that by definition, the second moments of area must be non-negative, while the product moment of area may not be.
index | the global index of a specified vertex element |
Definition at line 1554 of file VertexMesh.cpp.
References AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeLocation(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), and NEVER_REACHED.
double VertexMesh< ELEMENT_DIM, SPACE_DIM >::CalculateUnitNormalToFaceWithArea | ( | VertexElement< ELEMENT_DIM - 1, SPACE_DIM > * | pFace, |
c_vector< double, SPACE_DIM > & | rNormal | ||
) |
Compute the unit normal vector to a given face in 3D. This is achieved by calculating scaled normal, which is the effective sum of signed areas of triangle forming the face. Note: this may return the outward or inward normal, depending on the face chirality.
pFace | a face in the mesh |
rNormal | vector in which to return the unit normal |
Definition at line 1776 of file VertexMesh.cpp.
References NEVER_REACHED, and VectorProduct().
|
virtual |
Delete mNodes, mFaces and mElements.
Reimplemented in MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >, MutableVertexMesh< 2, 2 >, and MutableVertexMesh< DIM, DIM >.
Definition at line 788 of file VertexMesh.cpp.
Referenced by VertexMesh< ELEMENT_DIM, SPACE_DIM >::VertexMesh(), VertexMesh< ELEMENT_DIM, SPACE_DIM >::VertexMesh(), and VertexMesh< 2, 2 >::GenerateVerticesFromElementCircumcentres().
void VertexMesh< ELEMENT_DIM, SPACE_DIM >::ConstructFromMeshReader | ( | AbstractMeshReader< ELEMENT_DIM, SPACE_DIM > & | rMeshReader | ) |
Construct the mesh using a MeshReader.
rMeshReader | the mesh reader |
Referenced by VertexMesh< ELEMENT_DIM, SPACE_DIM >::load().
|
protected |
Test whether a given point lies inside a given element.
We use a winding number test, which counts the number of times the polygon associated with the element winds around the given point. The point is outside only when this "winding number" vanishes; otherwise, the point is inside.
One must decide whether a point on the polygon's boundary is inside or outside: we adopt the standard convention that a point on a left or bottom edge is inside, and a point on a right or top edge is outside. This way, if two distinct polygons share a common boundary segment, then a point on that segment will be in one polygon or the other, but not both at the same time.
rTestPoint | the point to test |
elementIndex | global index of the element in the mesh |
Definition at line 1366 of file VertexMesh.cpp.
References AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeLocation(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), and NEVER_REACHED.
|
protected |
Build edges from elements. Populates edges in EdgeHelper class
rElements | from which edges are built |
Get Doxygen to ignore, since it's confused by explicit instantiation of templated methods
Get Doxygen to ignore, since it's confused by explicit instantiation of templated methods
Definition at line 569 of file VertexMesh.cpp.
Referenced by MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >::MutableVertexMesh(), and VertexMesh< ELEMENT_DIM, SPACE_DIM >::VertexMesh().
|
protected |
Populate mNodes with locations corresponding to the element circumcentres of a given TetrahedralMesh. Used by 'Voronoi' constructors.
rMesh | a tetrahedral mesh |
Definition at line 599 of file VertexMesh.cpp.
Referenced by Cylindrical2dVertexMesh::Cylindrical2dVertexMesh(), and Toroidal2dVertexMesh::Toroidal2dVertexMesh().
c_vector< double, SPACE_DIM > VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetAreaGradientOfElementAtNode | ( | VertexElement< ELEMENT_DIM, SPACE_DIM > * | pElement, |
unsigned | localIndex | ||
) |
Compute the area gradient of a 2D element at one of its nodes.
N.B. This calls GetVectorFromAtoB(), which can be overridden in daughter classes for non-Euclidean metrics.
pElement | pointer to a specified vertex element |
localIndex | local index of a node in this element |
Definition at line 1672 of file VertexMesh.cpp.
References NEVER_REACHED.
Referenced by FarhadifarForce< DIM >::AddForceContribution(), and NagaiHondaForce< DIM >::AddForceContribution().
|
virtual |
Compute the centroid of an element.
A formula for the centroid of a plane polygon may be found e.g. in the following reference:
Mechanics of Materials James M. Gere (Author), Barry J. Goodno. Cengage Learning; 8th edition (January 1, 2012)
This needs to be overridden in daughter classes for non-Euclidean metrics.
index | the global index of a specified vertex element |
Definition at line 852 of file VertexMesh.cpp.
References AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeLocation(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), and NEVER_REACHED.
unsigned VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetDelaunayNodeIndexCorrespondingToVoronoiElementIndex | ( | unsigned | elementIndex | ) |
elementIndex | global index of an element in the Voronoi mesh |
Definition at line 712 of file VertexMesh.cpp.
References UNSIGNED_UNSET.
Referenced by Cylindrical2dVertexMesh::GetMeshForVtk(), Toroidal2dVertexMesh::GetMeshForVtk(), CellPopulationAreaWriter< ELEMENT_DIM, SPACE_DIM >::Visit(), and VoronoiDataWriter< ELEMENT_DIM, SPACE_DIM >::Visit().
Edge< SPACE_DIM > * VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetEdge | ( | unsigned | index | ) | const |
Fetches an edge.
index | global index of the edge |
Definition at line 587 of file VertexMesh.cpp.
double VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetEdgeLength | ( | unsigned | elementIndex1, |
unsigned | elementIndex2 | ||
) |
elementIndex1 | index of an element in the mesh |
elementIndex2 | index of an element in the mesh |
Definition at line 624 of file VertexMesh.cpp.
Referenced by HeterotypicBoundaryLengthWriter< ELEMENT_DIM, SPACE_DIM >::Visit().
VertexElement< ELEMENT_DIM, SPACE_DIM > * VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetElement | ( | unsigned | index | ) | const |
index | the global index of a specified vertex element. |
Definition at line 838 of file VertexMesh.cpp.
Referenced by CellPopulationElementWriter< ELEMENT_DIM, SPACE_DIM >::Visit().
|
inline |
skipDeletedElements | whether to include deleted element |
Definition at line 732 of file VertexMesh.hpp.
Referenced by FarhadifarForce< DIM >::AddForceContribution(), NagaiHondaForce< DIM >::AddForceContribution(), WelikyOsterForce< DIM >::AddForceContribution(), Cylindrical2dVertexMesh::GetMeshForVtk(), Toroidal2dVertexMesh::GetMeshForVtk(), VertexMeshWriter< ELEMENT_DIM, SPACE_DIM >::MakeVtkMesh(), CellPopulationAreaWriter< ELEMENT_DIM, SPACE_DIM >::Visit(), and VoronoiDataWriter< ELEMENT_DIM, SPACE_DIM >::Visit().
|
inline |
Definition at line 739 of file VertexMesh.hpp.
Referenced by FarhadifarForce< DIM >::AddForceContribution(), NagaiHondaForce< DIM >::AddForceContribution(), WelikyOsterForce< DIM >::AddForceContribution(), Cylindrical2dVertexMesh::GetMeshForVtk(), Toroidal2dVertexMesh::GetMeshForVtk(), VertexMeshWriter< ELEMENT_DIM, SPACE_DIM >::MakeVtkMesh(), CellPopulationAreaWriter< ELEMENT_DIM, SPACE_DIM >::Visit(), and VoronoiDataWriter< ELEMENT_DIM, SPACE_DIM >::Visit().
double VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetElongationShapeFactorOfElement | ( | unsigned | elementIndex | ) |
Get the elongation shape factor of a given element. This is defined as the square root of the ratio of the two second moments of the element around its principal axes.
elementIndex | index of an element in the mesh |
Definition at line 659 of file VertexMesh.cpp.
VertexElement< ELEMENT_DIM - 1, SPACE_DIM > * VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetFace | ( | unsigned | index | ) | const |
index | the global index of a specified face. |
Definition at line 845 of file VertexMesh.cpp.
|
protected |
Get the local index of a given element which is the start vertex of the edge of the element that the overlapping point rTestPoint is closest to.
rTestPoint | the point to test |
elementIndex | global index of the element in the mesh |
Definition at line 1488 of file VertexMesh.cpp.
References AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeLocation(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), and NEVER_REACHED.
|
virtual |
Return a pointer to the vertex mesh
This method may be overridden in daughter classes for non-Euclidean metrics. This can then be used when writing to VTK.
Reimplemented in Cylindrical2dVertexMesh, and Toroidal2dVertexMesh.
Definition at line 1027 of file VertexMesh.cpp.
Referenced by VertexMeshWriter< ELEMENT_DIM, SPACE_DIM >::WriteVtkUsingMesh().
std::set< unsigned > VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetNeighbouringElementIndices | ( | unsigned | elementIndex | ) |
Given an element, find a set containing the indices of its neighbouring elements.
elementIndex | global index of the element |
Definition at line 993 of file VertexMesh.cpp.
References AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNode(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), and Node< SPACE_DIM >::rGetContainingElementIndices().
Referenced by IsolatedLabelledCellKiller< DIM >::CheckAndLabelCellsForApoptosisOrDeath(), and HeterotypicBoundaryLengthWriter< ELEMENT_DIM, SPACE_DIM >::Visit().
std::set< unsigned > VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetNeighbouringNodeIndices | ( | unsigned | nodeIndex | ) |
Given a node, find a set containing the indices of its neighbouring nodes.
nodeIndex | global index of the node |
Definition at line 923 of file VertexMesh.cpp.
std::set< unsigned > VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetNeighbouringNodeNotAlsoInElement | ( | unsigned | nodeIndex, |
unsigned | elemIndex | ||
) |
Given a node and one of its containing elements, find a set containing the indices of those neighbouring node(s) that are NOT also in the element.
Note that we allow for more than one such index, since there is no reason a priori to assume that each node is contained by exactly three elements.
nodeIndex | global index of the node |
elemIndex | global index of the element |
Definition at line 953 of file VertexMesh.cpp.
References EXCEPTION, AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex(), and AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNumNodes().
c_vector< double, SPACE_DIM > VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetNextEdgeGradientOfElementAtNode | ( | VertexElement< ELEMENT_DIM, SPACE_DIM > * | pElement, |
unsigned | localIndex | ||
) |
Compute the gradient of the edge of a 2D element starting at its nodes.
N.B. This calls GetVectorFromAtoB(), which can be overridden in daughter classes for non-Euclidean metrics.
pElement | pointer to a specified vertex element |
localIndex | local index of a node in this element |
Definition at line 1730 of file VertexMesh.cpp.
References NEVER_REACHED.
Referenced by FarhadifarForce< DIM >::AddForceContribution(), and NagaiHondaForce< DIM >::AddForceContribution().
unsigned VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetNumAllElements | ( | ) | const |
Definition at line 826 of file VertexMesh.cpp.
unsigned VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetNumEdges | ( | ) | const |
Gets the number of edges in the mesh
Definition at line 581 of file VertexMesh.cpp.
|
virtual |
Reimplemented in MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >, MutableVertexMesh< 2, 2 >, and MutableVertexMesh< DIM, DIM >.
Definition at line 820 of file VertexMesh.cpp.
|
virtual |
Definition at line 832 of file VertexMesh.cpp.
|
virtual |
Reimplemented from AbstractMesh< ELEMENT_DIM, SPACE_DIM >.
Reimplemented in MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >, MutableVertexMesh< 2, 2 >, and MutableVertexMesh< DIM, DIM >.
Definition at line 814 of file VertexMesh.cpp.
Referenced by VertexMeshWriter< ELEMENT_DIM, SPACE_DIM >::MakeVtkMesh(), and VertexMeshWriter< ELEMENT_DIM, SPACE_DIM >::WriteFilesUsingMesh().
c_vector< double, SPACE_DIM > VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetPerimeterGradientOfElementAtNode | ( | VertexElement< ELEMENT_DIM, SPACE_DIM > * | pElement, |
unsigned | localIndex | ||
) |
Compute the gradient of the perimeter of a 2D element at its nodes. This returns the sum of GetPreviousEdgeGradientAtNode() and GetNextEdgeGradientAtNode().
pElement | pointer to a specified vertex element |
localIndex | local index of a node in this element |
Definition at line 1754 of file VertexMesh.cpp.
References NEVER_REACHED.
c_vector< double, SPACE_DIM > VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetPreviousEdgeGradientOfElementAtNode | ( | VertexElement< ELEMENT_DIM, SPACE_DIM > * | pElement, |
unsigned | localIndex | ||
) |
Compute the gradient of the edge of a 2D element ending at its nodes.
N.B. This calls GetVectorFromAtoB(), which can be overridden in daughter classes for non-Euclidean metrics.
pElement | pointer to a specified vertex element |
localIndex | local index of a node in this element |
Definition at line 1701 of file VertexMesh.cpp.
References NEVER_REACHED.
unsigned VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetRosetteRankOfElement | ( | unsigned | index | ) |
Get the "rosette rank" of an element.
This is defined as the maximum number of elements shared by any node in the specified element.
index | the global index of a specified vertex element |
Definition at line 764 of file VertexMesh.cpp.
References AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNode(), and AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNumNodes().
c_vector< double, SPACE_DIM > VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetShortAxisOfElement | ( | unsigned | index | ) |
Compute the direction of the shortest principal axis passing through the centroid, or 'short axis', of a given element. This is the eigenvector associated with the eigenvalue of largest magnitude of the inertia matrix
J = ( I_xx -I_xy ) ( -I_xy I_yy )
whose entries are computed by calling the method CalculateMomentsOfElement().
Note that if the nodes owned by the element are supplied in clockwise rather than anticlockwise manner, or if this arises when any periodicity is enforced, then the sign of each moment may be incorrect change. This means that we need to consider the eigenvalue of largest magnitude rather than largest value when computing the short axis of the element.
If the element is a regular polygon then the eigenvalues of the inertia tensor are equal: in this case we return a random unit vector.
This method is only implemented in 2D at present.
index | the global index of a specified vertex element |
Definition at line 1613 of file VertexMesh.cpp.
References RandomNumberGenerator::Instance(), NEVER_REACHED, and RandomNumberGenerator::ranf().
Referenced by ShortAxisVertexBasedDivisionRule< SPACE_DIM >::CalculateCellDivisionVector().
|
virtual |
Compute the surface area (or perimeter in 2D) of an element.
This needs to be overridden in daughter classes for non-Euclidean metrics.
index | the global index of a specified vertex element |
Definition at line 1331 of file VertexMesh.cpp.
References VertexElement< ELEMENT_DIM, SPACE_DIM >::GetFace(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex(), VertexElement< ELEMENT_DIM, SPACE_DIM >::GetNumFaces(), and AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNumNodes().
Referenced by FarhadifarForce< DIM >::AddForceContribution(), NagaiHondaForce< DIM >::AddForceContribution(), WelikyOsterForce< DIM >::AddForceContribution(), and VoronoiDataWriter< ELEMENT_DIM, SPACE_DIM >::Visit().
|
virtual |
Overridden GetVectorFromAtoB() method. Returns a vector between two points in space.
If the mesh is being used to represent a Voronoi tessellation, and mpDelaunayMesh is not NULL, then use that to compute GetVectorFromAtoB.
rLocationA | a c_vector of coordinates |
rLocationB | a c_vector of coordinates |
Reimplemented from AbstractMesh< ELEMENT_DIM, SPACE_DIM >.
Definition at line 1258 of file VertexMesh.cpp.
References AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetVectorFromAtoB().
Referenced by WelikyOsterForce< DIM >::AddForceContribution().
|
virtual |
Get the volume (or area in 2D, or length in 1D) of an element.
This needs to be overridden in daughter classes for non-Euclidean metrics.
index | the global index of a specified vertex element |
Definition at line 1274 of file VertexMesh.cpp.
References VertexElement< ELEMENT_DIM, SPACE_DIM >::GetFace(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeLocation(), VertexElement< ELEMENT_DIM, SPACE_DIM >::GetNumFaces(), and AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNumNodes().
Referenced by FarhadifarForce< DIM >::AddForceContribution(), NagaiHondaForce< DIM >::AddForceContribution(), WelikyOsterForce< DIM >::AddForceContribution(), Toroidal2dVertexMesh::ConstructFromMeshReader(), CellPopulationAreaWriter< ELEMENT_DIM, SPACE_DIM >::Visit(), and VoronoiDataWriter< ELEMENT_DIM, SPACE_DIM >::Visit().
unsigned VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetVoronoiElementIndexCorrespondingToDelaunayNodeIndex | ( | unsigned | nodeIndex | ) |
nodeIndex | global index of a node in the Delaunay mesh |
Definition at line 738 of file VertexMesh.cpp.
References EXCEPTION, and UNSIGNED_UNSET.
bool VertexMesh< ELEMENT_DIM, SPACE_DIM >::IsNearExistingNodes | ( | c_vector< double, SPACE_DIM > | newNodeLocation, |
std::vector< Node< SPACE_DIM > * > | nodesToCheck, | ||
double | minClearance | ||
) |
Helper method to determine if a point in space is near to any existing nodes. Used when making bounded Voronoi Tesselations
newNodeLocation | the location of the proposed node |
nodesToCheck | the nodes to check for proximity to the proposed node |
minClearance | the minimum clearance |
Definition at line 1033 of file VertexMesh.cpp.
Referenced by Cylindrical2dVertexMesh::Cylindrical2dVertexMesh(), and Toroidal2dVertexMesh::Toroidal2dVertexMesh().
|
inlineprotected |
Load a mesh by using VertexMeshReader and the location in ArchiveLocationInfo.
archive | the archive |
version | the current version of this class |
Definition at line 213 of file VertexMesh.hpp.
References VertexMesh< ELEMENT_DIM, SPACE_DIM >::ConstructFromMeshReader(), ArchiveLocationInfo::GetArchiveDirectory(), and ArchiveLocationInfo::GetMeshFilename().
const EdgeHelper< SPACE_DIM > & VertexMesh< ELEMENT_DIM, SPACE_DIM >::rGetEdgeHelper | ( | ) | const |
Fetches EdgeHelper.
Definition at line 593 of file VertexMesh.cpp.
|
inlineprotected |
Archive the VertexMesh and its member variables. Note that this will write out a VertexMeshWriter file to wherever ArchiveLocationInfo has specified.
archive | the archive |
version | the current version of this class |
Definition at line 196 of file VertexMesh.hpp.
References ArchiveLocationInfo::GetArchiveRelativePath(), ArchiveLocationInfo::GetMeshFilename(), and VertexMeshWriter< ELEMENT_DIM, SPACE_DIM >::WriteFilesUsingMesh().
|
protected |
Solve boundary element mapping method. This overridden method is required as it is pure virtual in the base class.
index | the global index of the boundary element |
Definition at line 704 of file VertexMesh.cpp.
|
protected |
Solve element mapping method. This overridden method is required as it is pure virtual in the base class.
index | the global index of the element |
Definition at line 697 of file VertexMesh.cpp.
|
protectedvirtual |
Solve node mapping method. This overridden method is required as it is pure virtual in the base class.
index | the global index of the node |
Implements AbstractMesh< ELEMENT_DIM, SPACE_DIM >.
Definition at line 690 of file VertexMesh.cpp.
Needed for serialization.
Definition at line 186 of file VertexMesh.hpp.
Definition at line 79 of file VertexMesh.hpp.
|
protected |
Object that owns and manages edges
Definition at line 90 of file VertexMesh.hpp.
|
protected |
Vector of pointers to VertexElements.
Definition at line 84 of file VertexMesh.hpp.
Referenced by Cylindrical2dVertexMesh::Cylindrical2dVertexMesh(), MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >::MutableVertexMesh(), Toroidal2dVertexMesh::Toroidal2dVertexMesh(), VertexMesh< ELEMENT_DIM, SPACE_DIM >::VertexElementIterator::VertexElementIterator(), VertexMesh< ELEMENT_DIM, SPACE_DIM >::VertexMesh(), Toroidal2dVertexMesh::ConstructFromMeshReader(), VertexMesh< 2, 2 >::ElementIncludesPoint(), VertexMesh< 2, 2 >::GetFace(), and VertexMesh< 2, 2 >::GetNumFaces().
|
protected |
Vector of pointers to VertexElements.
Definition at line 87 of file VertexMesh.hpp.
Referenced by MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >::MutableVertexMesh(), VertexMesh< ELEMENT_DIM, SPACE_DIM >::VertexMesh(), and VertexMesh< 2, 2 >::SolveElementMapping().
|
protected |
Delaunay tetrahedral mesh that is used only when the vertex mesh is used to represent a Voronoi tessellation. A pointer to the Delaunay mesh is required in this case because the Delaunay mesh may be a subclass of TetrahedralMesh, which overrides methods such as GetVectorFromAtoB().
Definition at line 108 of file VertexMesh.hpp.
Referenced by Cylindrical2dVertexMesh::Cylindrical2dVertexMesh(), Toroidal2dVertexMesh::Toroidal2dVertexMesh(), VertexMesh< 2, 2 >::GenerateVerticesFromElementCircumcentres(), Cylindrical2dVertexMesh::GetMeshForVtk(), and Toroidal2dVertexMesh::GetMeshForVtk().
|
protected |
Map that is used only when the vertex mesh is used to represent a Voronoi tessellation, the dual to a Delaunay tetrahedral mesh. The map consists of pairs (index1, index2), where index1 denotes the global index of a node in the Delaunay mesh and index2 denotes the global index of the corresponding element in the Voronoi mesh.
Definition at line 99 of file VertexMesh.hpp.