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

#include <AbstractMesh.hpp>

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

Classes

class  NodeIterator
 

Public Types

typedef std::vector< Node< SPACE_DIM > * >::const_iterator BoundaryNodeIterator
 

Public Member Functions

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

virtual void SetElementOwnerships ()
 
ChasteCuboid< SPACE_DIM > CalculateBoundingBox (const std::vector< Node< SPACE_DIM > * > &rNodes) const
 

Protected Attributes

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 SolveNodeMapping (unsigned index) const =0
 
template<class Archive >
void serialize (Archive &archive, const unsigned int version)
 

Friends

class TestDistributedTetrahedralMesh
 
template<unsigned A_DIMENSION>
class NodesOnlyMesh
 
template<unsigned A_DIMENSION>
class QuadraticMeshHelper
 
class boost::serialization::access
 

Detailed Description

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

Abstract base class for all meshes.

Definition at line 60 of file AbstractMesh.hpp.

Member Typedef Documentation

◆ BoundaryNodeIterator

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
typedef std::vector<Node<SPACE_DIM>*>::const_iterator AbstractMesh< ELEMENT_DIM, SPACE_DIM >::BoundaryNodeIterator

Definition of boundary node Iterator type.

Definition at line 148 of file AbstractMesh.hpp.

Constructor & Destructor Documentation

◆ AbstractMesh()

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

Constructor.

Definition at line 44 of file AbstractMesh.cpp.

◆ ~AbstractMesh()

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

Virtual destructor, since this class has virtual methods.

Definition at line 52 of file AbstractMesh.cpp.

Member Function Documentation

◆ CalculateBoundingBox() [1/2]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
ChasteCuboid< SPACE_DIM > AbstractMesh< ELEMENT_DIM, SPACE_DIM >::CalculateBoundingBox ( ) const
virtual

Calculate the bounding box (width extremes for all dimensions of the mesh). Overridden in distributed case

Returns
The minimum and maximum co-ordinates of any node in each dimension

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

Definition at line 277 of file AbstractMesh.cpp.

Referenced by DistributedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::CalculateBoundingBox(), and Cylindrical2dMesh::UpdateTopAndBottom().

◆ CalculateBoundingBox() [2/2]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
ChasteCuboid< SPACE_DIM > AbstractMesh< ELEMENT_DIM, SPACE_DIM >::CalculateBoundingBox ( const std::vector< Node< SPACE_DIM > * > &  rNodes) const
protected

Calculate a bounding box from a set of nodes. A generalised version of the public CalculateBoundingBox method.

Parameters
rNodesthe list of nodes to calculate the bounding box for.
Returns
the bounding box in the form of a ChasteCuboid.
Todo:
#1322 use a const version of NodeIterator here

Definition at line 227 of file AbstractMesh.cpp.

Referenced by Electrodes< DIM >::Electrodes(), and FineCoarseMeshPair< DIM >::SetUpBoxes().

◆ CalculateMaximumContainingElementsPerProcess()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned AbstractMesh< ELEMENT_DIM, SPACE_DIM >::CalculateMaximumContainingElementsPerProcess ( ) const
Returns
Iterates through local nodes and finds the maximum number of containing elements for all locally owned nodes Useful for determining FEM matrix fill.

Definition at line 491 of file AbstractMesh.cpp.

Referenced by LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::InitialiseForSolve(), and AbstractTetrahedralMeshWriter< ELEMENT_DIM, SPACE_DIM >::WriteNclFile().

◆ GetBoundaryNodeIteratorBegin()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
AbstractMesh< ELEMENT_DIM, SPACE_DIM >::BoundaryNodeIterator AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetBoundaryNodeIteratorBegin ( ) const

◆ GetBoundaryNodeIteratorEnd()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
AbstractMesh< ELEMENT_DIM, SPACE_DIM >::BoundaryNodeIterator AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetBoundaryNodeIteratorEnd ( ) const
Returns
a pointer to *one past* the last boundary node in the mesh (for consistency with STL iterators).

Definition at line 175 of file AbstractMesh.cpp.

Referenced by PapillaryFibreCalculator::GetRadiusVectorForOneElement().

◆ GetDistanceBetweenNodes()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
double AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetDistanceBetweenNodes ( unsigned  indexA,
unsigned  indexB 
)

Return the distance between two nodes.

This method calls GetVectorFromAtoB(), which is overridden in some daughter classes (e.g. Cylindrical2dMesh).

Parameters
indexAa node index
indexBa node index
Returns
distance between two nodes.

Definition at line 212 of file AbstractMesh.cpp.

Referenced by HeterotypicBoundaryLengthWriter< ELEMENT_DIM, SPACE_DIM >::Visit().

◆ GetDistributedVectorFactory()

◆ GetMeshFileBaseName()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
std::string AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetMeshFileBaseName ( ) const
Returns
mMeshFileBaseName.

Definition at line 181 of file AbstractMesh.cpp.

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

◆ GetNearestNodeIndex()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNearestNodeIndex ( const ChastePoint< SPACE_DIM > &  rTestPoint)
virtual

GetNearestNodeIndex iterates through all nodes in the mesh and returns the global index with the smallest distance to the provided point.

This method is overridden in the distributed case to return the global node index.

This method uses GetVectorFromAtoB distance and hence may return a correct solution in non-Euclidean space, but only if this method is overridden in a subclass (see e.g. Cylindrical2dMesh for an example of this).

Parameters
rTestPointreference to the point
Returns
node index

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

Definition at line 285 of file AbstractMesh.cpp.

Referenced by DistributedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNearestNodeIndex().

◆ GetNode()

◆ GetNodeFromPrePermutationIndex()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
Node< SPACE_DIM > * AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNodeFromPrePermutationIndex ( unsigned  index) const

Get the node with a given index in the mesh, prior to any node permutation being applied. For non-permuted meshes, this will have the same effect as GetNode.

This method is intended for use by the archiving code, to enable checkpoint migration, so that we can load the correct cells and boundary conditions after the mesh has been re-partitioned.

If unsure, use GetNode in preference to this method!

Parameters
indexthe global index of the node prior to a permutation being applied
Returns
a pointer to the node

Definition at line 106 of file AbstractMesh.cpp.

Referenced by BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::MergeFromArchive().

◆ GetNodeIteratorBegin()

◆ GetNodeIteratorEnd()

◆ GetNodeOrHaloNode()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
Node< SPACE_DIM > * AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNodeOrHaloNode ( unsigned  index) const
virtual

Get the node with a given index in the mesh (synonym of GetNode() unless overridden in a distributed mesh).

Parameters
indexthe global index of the node
Returns
a pointer to the node.

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

Definition at line 100 of file AbstractMesh.cpp.

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

◆ GetNumAllNodes()

◆ GetNumBoundaryNodes()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumBoundaryNodes ( ) const
Returns
the number of boundary nodes in the mesh.

Definition at line 72 of file AbstractMesh.cpp.

Referenced by BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, 3 >::DefineConstantDirichletOnMeshBoundary().

◆ GetNumNodeAttributes()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumNodeAttributes ( ) const
Returns
the number of attributes on each node. Note, this implementation assumes that all nodes have the same number of attributes so that the first node in the container is indicative of the others.

Definition at line 84 of file AbstractMesh.cpp.

◆ GetNumNodes()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumNodes ( ) const
virtual
Returns
the number of nodes that are actually in use.

Overridden in MutableMesh and DistributedTetrahedralMesh.

Reimplemented in PottsMesh< DIM >, PottsMesh< SPACE_DIM >, DistributedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >, DistributedTetrahedralMesh< DIM, DIM >, MutableMesh< ELEMENT_DIM, SPACE_DIM >, MutableMesh< 2, 2 >, MutableMesh< DIM, DIM >, MutableMesh< SPACE_DIM, SPACE_DIM >, NodesOnlyMesh< SPACE_DIM >, NodesOnlyMesh< 2 >, NodesOnlyMesh< DIM >, ImmersedBoundaryMesh< ELEMENT_DIM, SPACE_DIM >, ImmersedBoundaryMesh< 2, 2 >, ImmersedBoundaryMesh< DIM, DIM >, MutableVertexMesh< ELEMENT_DIM, SPACE_DIM >, MutableVertexMesh< 2, 2 >, MutableVertexMesh< DIM, DIM >, VertexMesh< ELEMENT_DIM, SPACE_DIM >, VertexMesh< 2, 2 >, and VertexMesh< DIM, DIM >.

Definition at line 66 of file AbstractMesh.cpp.

Referenced by AbstractContinuumMechanicsSolver< DIM >::AbstractContinuumMechanicsSolver(), CardiacElectroMechProbRegularGeom< DIM >::CardiacElectroMechProbRegularGeom(), CmguiDeformedSolutionsWriter< DIM >::CmguiDeformedSolutionsWriter(), Hdf5ToTxtConverter< ELEMENT_DIM, SPACE_DIM >::Hdf5ToTxtConverter(), Hdf5ToVtkConverter< ELEMENT_DIM, SPACE_DIM >::Hdf5ToVtkConverter(), StreeterFibreGenerator< SPACE_DIM >::StreeterFibreGenerator(), VoltageInterpolaterOntoMechanicsMesh< DIM >::VoltageInterpolaterOntoMechanicsMesh(), AbstractFunctionalCalculator< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::Calculate(), CuboidMeshConstructor< ELEMENT_DIM, SPACE_DIM >::Construct(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ConstructFromMesh(), NodePartitioner< ELEMENT_DIM, SPACE_DIM >::DumbPartitioning(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ExportToMesher(), CryptCellsGenerator< CELL_CYCLE_MODEL >::Generate(), NonlinearElasticityTools< DIM >::GetNodesByComponentValue(), LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::PrepareForSetupLinearSystem(), AbstractTetrahedralMeshWriter< ELEMENT_DIM, SPACE_DIM >::WriteFilesUsingMesh(), XdmfMeshWriter< ELEMENT_DIM, SPACE_DIM >::WriteFilesUsingMesh(), AbstractTetrahedralMeshWriter< ELEMENT_DIM, SPACE_DIM >::WriteNclFile(), DiscreteSystemForceCalculator::WriteResultsToFile(), and LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::WriteVtkResultsToFile().

◆ GetVectorFromAtoB()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
c_vector< double, SPACE_DIM > AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetVectorFromAtoB ( const c_vector< double, SPACE_DIM > &  rLocationA,
const c_vector< double, SPACE_DIM > &  rLocationB 
)
virtual
Returns
a vector between two points in space.

This method is overridden in some daughter classes (e.g. Cylindrical2dMesh).

Parameters
rLocationAa c_vector of coordinates
rLocationBa c_vector of coordinates
Returns
c_vector from location A to location B.

Reimplemented in ImmersedBoundaryMesh< ELEMENT_DIM, SPACE_DIM >, ImmersedBoundaryMesh< 2, 2 >, ImmersedBoundaryMesh< DIM, DIM >, PeriodicNodesOnlyMesh< SPACE_DIM >, VertexMesh< ELEMENT_DIM, SPACE_DIM >, VertexMesh< 2, 2 >, and VertexMesh< DIM, DIM >.

Definition at line 204 of file AbstractMesh.cpp.

Referenced by VertexMesh< ELEMENT_DIM, SPACE_DIM >::GetVectorFromAtoB().

◆ GetWidth()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
double AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetWidth ( const unsigned rDimension) const
virtual

Calculate the 'width' of any dimension of the mesh.

This method is overridden in some daughter classes (e.g. Cylindrical2dMesh).

Parameters
rDimensiona dimension (0,1 or 2)
Returns
The maximum distance between any nodes in this dimension.

Reimplemented in NodesOnlyMesh< SPACE_DIM >, NodesOnlyMesh< 2 >, NodesOnlyMesh< DIM >, Cylindrical2dMesh, Cylindrical2dNodesOnlyMesh, PeriodicNodesOnlyMesh< SPACE_DIM >, Toroidal2dMesh, Cylindrical2dVertexMesh, and Toroidal2dVertexMesh.

Definition at line 220 of file AbstractMesh.cpp.

Referenced by NodesOnlyMesh< SPACE_DIM >::GetWidth(), Cylindrical2dMesh::GetWidth(), Cylindrical2dNodesOnlyMesh::GetWidth(), PeriodicNodesOnlyMesh< SPACE_DIM >::GetWidth(), and Cylindrical2dVertexMesh::GetWidth().

◆ IsMeshChanging()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
bool AbstractMesh< ELEMENT_DIM, SPACE_DIM >::IsMeshChanging ( ) const
Returns
Whether the mesh changes (used in archiving).

Definition at line 485 of file AbstractMesh.cpp.

◆ IsMeshOnDisk()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
bool AbstractMesh< ELEMENT_DIM, SPACE_DIM >::IsMeshOnDisk ( ) const

Get whether this mesh was read from file.

Returns
whether this mesh was read from file

Definition at line 192 of file AbstractMesh.cpp.

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

◆ PermuteNodes()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractMesh< ELEMENT_DIM, SPACE_DIM >::PermuteNodes ( )
virtual

Permute the nodes so that they appear in a different order in mNodes (and their mIndex's are altered accordingly).

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

Definition at line 162 of file AbstractMesh.cpp.

References NEVER_REACHED.

◆ ReadNodesPerProcessorFile()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractMesh< ELEMENT_DIM, SPACE_DIM >::ReadNodesPerProcessorFile ( const std::string &  rNodesPerProcessorFile)
virtual

Read in the number of nodes per processor from file.

Parameters
rNodesPerProcessorFilethe name of the file

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

Definition at line 120 of file AbstractMesh.cpp.

References NEVER_REACHED.

◆ RefreshMesh()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractMesh< ELEMENT_DIM, SPACE_DIM >::RefreshMesh ( )
virtual

This method allows the mesh properties to be re-calculated after one or more nodes have been moved.

Reimplemented in TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >, TetrahedralMesh< 3, 3 >, TetrahedralMesh< DIM, DIM >, Cylindrical2dNodesOnlyMesh, PeriodicNodesOnlyMesh< SPACE_DIM >, and Toroidal2dMesh.

Definition at line 480 of file AbstractMesh.cpp.

Referenced by Cylindrical2dNodesOnlyMesh::RefreshMesh(), and PeriodicNodesOnlyMesh< SPACE_DIM >::RefreshMesh().

◆ rGetNodePermutation()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
const std::vector< unsigned > & AbstractMesh< ELEMENT_DIM, SPACE_DIM >::rGetNodePermutation ( ) const
Returns
mNodePermutation.

When empty (most meshes) there is no node permutation When non-empty (parallel distributed meshes) then for a given original_index mNodePermutation[original_index] holds the new assigned index of that node in memory

Definition at line 198 of file AbstractMesh.cpp.

Referenced by QuadraticMeshHelper< DIM >::AddInternalNodesToElements(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::save(), and AbstractTetrahedralMeshWriter< ELEMENT_DIM, SPACE_DIM >::WriteNclFile().

◆ Rotate() [1/3]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractMesh< ELEMENT_DIM, SPACE_DIM >::Rotate ( c_matrix< double, SPACE_DIM, SPACE_DIM >  rotationMatrix)
virtual

Do a general mesh rotation with a positive determinant orthonormal rotation matrix. This is the rotation method that actually does the work. Should be overridden when the child class has halo nodes.

Parameters
rotationMatrixis a Ublas rotation matrix of the correct form

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

Definition at line 385 of file AbstractMesh.cpp.

◆ Rotate() [2/3]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractMesh< ELEMENT_DIM, SPACE_DIM >::Rotate ( c_vector< double, 3 >  axis,
double  angle 
)

Do an angle axis rotation.

Parameters
axisis the axis of rotation (does not need to be normalised)
angleis the angle of rotation in radians

Definition at line 399 of file AbstractMesh.cpp.

◆ Rotate() [3/3]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractMesh< ELEMENT_DIM, SPACE_DIM >::Rotate ( double  theta)

Rotating a 2D mesh equates that rotation around the z-axis.

Parameters
thetais the angle of rotation in radians

Definition at line 474 of file AbstractMesh.cpp.

◆ RotateX()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractMesh< ELEMENT_DIM, SPACE_DIM >::RotateX ( const double  theta)

Rotate the mesh about the x-axis.

Parameters
thetais the angle of rotation in radians

Definition at line 424 of file AbstractMesh.cpp.

◆ RotateY()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractMesh< ELEMENT_DIM, SPACE_DIM >::RotateY ( const double  theta)

Rotate the mesh about the y-axis.

Parameters
thetais the angle of rotation in radians

Definition at line 440 of file AbstractMesh.cpp.

◆ RotateZ()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractMesh< ELEMENT_DIM, SPACE_DIM >::RotateZ ( const double  theta)

Rotate the mesh about the z-axis.

Parameters
thetais the angle of rotation in radians

Definition at line 457 of file AbstractMesh.cpp.

References EXCEPTION.

◆ Scale()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractMesh< ELEMENT_DIM, SPACE_DIM >::Scale ( const double  xFactor = 1.0,
const double  yFactor = 1.0,
const double  zFactor = 1.0 
)
virtual

Scale the mesh.

Parameters
xFactoris the scale in the x-direction (defaults to 1.0)
yFactoris the scale in the y-direction (defaults to 1.0)
zFactoris the scale in the z-direction (defaults to 1.0)

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

Definition at line 322 of file AbstractMesh.cpp.

Referenced by DistributedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::Scale(), and Cylindrical2dVertexMesh::Scale().

◆ serialize()

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

Serialize the mesh.

Parameters
archivethe archive
versionthe current version of this class

Definition at line 87 of file AbstractMesh.hpp.

References ProcessSpecificArchive< Archive >::Get(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mMeshChangesDuringSimulation, and AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mpDistributedVectorFactory.

◆ SetDistributedVectorFactory()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractMesh< ELEMENT_DIM, SPACE_DIM >::SetDistributedVectorFactory ( DistributedVectorFactory pFactory)
virtual

Set method for mpDistributedVectorFactory. Must be called before the mesh is used for anything. This only actually impacts the DistributedTetrahedralMesh subclass, in which the supplied factory is then used to specify the node distribution among the processes.

Parameters
pFactorya factory to use for this mesh

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

Definition at line 147 of file AbstractMesh.cpp.

References EXCEPTION, PetscTools::GetNumProcs(), and DistributedVectorFactory::GetNumProcs().

Referenced by AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::load(), and DistributedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::SetDistributedVectorFactory().

◆ SetElementOwnerships()

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

◆ SetMeshHasChangedSinceLoading()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractMesh< ELEMENT_DIM, SPACE_DIM >::SetMeshHasChangedSinceLoading ( )

Set whether the mesh has been modified since it was read from file. This prevents the archiving code just blithely storing the original, unmodified, mesh.

Definition at line 506 of file AbstractMesh.cpp.

Referenced by ImmersedBoundaryCellPopulation< DIM >::GetTetrahedralMeshForPdeModifier(), and VertexBasedCellPopulation< DIM >::GetTetrahedralMeshForPdeModifier().

◆ SolveNodeMapping()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
virtual unsigned AbstractMesh< ELEMENT_DIM, SPACE_DIM >::SolveNodeMapping ( unsigned  index) const
privatepure virtual

◆ Translate() [1/2]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractMesh< ELEMENT_DIM, SPACE_DIM >::Translate ( const c_vector< double, SPACE_DIM > &  rDisplacement)
virtual

Translate the mesh given the displacement vector. This is the translation method that actually does the work. Should be overridden when the child class has halo nodes.

Parameters
rDisplacementis a translation vector of the correct size

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

Definition at line 371 of file AbstractMesh.cpp.

◆ Translate() [2/2]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractMesh< ELEMENT_DIM, SPACE_DIM >::Translate ( const double  xMovement = 0.0,
const double  yMovement = 0.0,
const double  zMovement = 0.0 
)

Translate the mesh given the coordinate displacements separately.

Parameters
xMovementis the x-displacement (defaults to 0.0)
yMovementis the y-displacement (defaults to 0.0)
zMovementis the z-displacement (defaults to 0.0)

Definition at line 344 of file AbstractMesh.cpp.

Friends And Related Symbol Documentation

◆ boost::serialization::access

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

Needed for serialization.

Definition at line 79 of file AbstractMesh.hpp.

◆ NodesOnlyMesh

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
template<unsigned A_DIMENSION>
friend class NodesOnlyMesh
friend

Definition at line 63 of file AbstractMesh.hpp.

◆ QuadraticMeshHelper

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
template<unsigned A_DIMENSION>
friend class QuadraticMeshHelper
friend

Definition at line 64 of file AbstractMesh.hpp.

◆ TestDistributedTetrahedralMesh

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

Definition at line 62 of file AbstractMesh.hpp.

Member Data Documentation

◆ mBoundaryNodes

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
std::vector<Node<SPACE_DIM> *> AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryNodes
protected

Vector of pointers to boundary nodes in the mesh.

Definition at line 99 of file AbstractMesh.hpp.

Referenced by QuadraticMeshHelper< DIM >::AddNodeToBoundaryElement(), and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ImportFromMesher().

◆ mMeshChangesDuringSimulation

◆ mMeshFileBaseName

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
std::string AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mMeshFileBaseName
protected

If the mesh is constructed from file using a MeshReader, this member variable stores the base name of these files.

Definition at line 120 of file AbstractMesh.hpp.

◆ mNodePermutation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
std::vector<unsigned> AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mNodePermutation
protected

Vector containing node permutation information.

When empty (most meshes) there is no node permutation When non-empty (parallel distributed meshes) then for a given original_index mNodePermutation[original_index] holds the new assigned index of that node in memory

Definition at line 114 of file AbstractMesh.hpp.

◆ mNodes

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
std::vector<Node<SPACE_DIM> *> AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mNodes
protected

◆ mpDistributedVectorFactory

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
DistributedVectorFactory* AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mpDistributedVectorFactory
protected

DistributedVectorFactory capable of reproducing the given number of nodes owned by each processor.

Definition at line 105 of file AbstractMesh.hpp.

Referenced by AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::load(), and AbstractMesh< ELEMENT_DIM, SPACE_DIM >::serialize().


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