TetrahedralMesh< ELEMENT_DIM, SPACE_DIM > Class Template Reference

#include <TetrahedralMesh.hpp>

Inherits AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >.

Inherited by MutableMesh< ELEMENT_DIM, SPACE_DIM >, MutableMesh< 2, 2 >, MutableMesh< DIM, DIM >, MutableMesh< SPACE_DIM, SPACE_DIM >, and NonCachedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >.

Collaboration diagram for TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >:
Collaboration graph
[legend]

List of all members.

Classes

class  EdgeIterator

Public Member Functions

 TetrahedralMesh ()
void ConstructFromMeshReader (AbstractMeshReader< ELEMENT_DIM, SPACE_DIM > &rMeshReader)
void ReadNodesPerProcessorFile (const std::string &rNodesPerProcessorFile)
bool CheckIsConforming ()
double GetVolume ()
double GetSurfaceArea ()
void RefreshMesh ()
void PermuteNodes ()
void PermuteNodes (const std::vector< unsigned > &perm)
unsigned GetContainingElementIndex (const ChastePoint< SPACE_DIM > &rTestPoint, bool strict=false, std::set< unsigned > testElements=std::set< unsigned >(), bool onlyTryWithTestElements=false)
unsigned GetContainingElementIndexWithInitialGuess (const ChastePoint< SPACE_DIM > &rTestPoint, unsigned startingElementGuess, bool strict=false)
unsigned GetNearestElementIndex (const ChastePoint< SPACE_DIM > &rTestPoint)
unsigned GetNearestElementIndexFromTestElements (const ChastePoint< SPACE_DIM > &rTestPoint, std::set< unsigned > testElements)
std::vector< unsignedGetContainingElementIndices (const ChastePoint< SPACE_DIM > &rTestPoint)
virtual void Clear ()
std::set< unsignedCalculateBoundaryOfFlaggedRegion ()
double GetAngleBetweenNodes (unsigned indexA, unsigned indexB)
void UnflagAllElements ()
void FlagElementsNotContainingNodes (const std::set< unsigned > rNodes)
virtual void RefreshJacobianCachedData ()
virtual void GetJacobianForElement (unsigned elementIndex, c_matrix< double, SPACE_DIM, SPACE_DIM > &rJacobian, double &rJacobianDeterminant) 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 GetWeightedDirectionForElement (unsigned elementIndex, c_vector< double, SPACE_DIM > &rWeightedDirection, double &rJacobianDeterminant) const
virtual void GetWeightedDirectionForBoundaryElement (unsigned elementIndex, c_vector< double, SPACE_DIM > &rWeightedDirection, double &rJacobianDeterminant) const
EdgeIterator EdgesBegin ()
EdgeIterator EdgesEnd ()

Protected Member Functions

unsigned SolveNodeMapping (unsigned index) const
unsigned SolveElementMapping (unsigned index) const
unsigned SolveBoundaryElementMapping (unsigned index) const
template<class MESHER_IO >
void ExportToMesher (NodeMap &map, MESHER_IO &mesherInput, int *elementList=NULL)
template<class MESHER_IO >
void ImportFromMesher (MESHER_IO &mesherOutput, unsigned numberOfElements, int *elementList, unsigned numberOfFaces, int *faceList, int *edgeMarkerList)
void InitialiseTriangulateIo (triangulateio &mesherIo)
void FreeTriangulateIo (triangulateio &mesherIo)

Protected Attributes

std::vector< c_vector< double,
SPACE_DIM > > 
mElementWeightedDirections
std::vector< c_matrix< double,
SPACE_DIM, ELEMENT_DIM > > 
mElementJacobians
std::vector< c_matrix< double,
ELEMENT_DIM, SPACE_DIM > > 
mElementInverseJacobians
std::vector< doublemElementJacobianDeterminants
std::vector< c_vector< double,
SPACE_DIM > > 
mBoundaryElementWeightedDirections
std::vector< doublemBoundaryElementJacobianDeterminants

Private Member Functions

template<class Archive >
void serialize (Archive &archive, const unsigned int version)

Friends

class TestTetrahedralMesh
class TestCryptSimulation2d
class boost::serialization::access

Detailed Description

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

Forward declaration for triangle helper methods (used in MutableMesh QuadraticMesh) A concrete tetrahedral mesh class.

Definition at line 57 of file TetrahedralMesh.hpp.


Constructor & Destructor Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::TetrahedralMesh (  )  [inline]

Constructor.

Definition at line 53 of file TetrahedralMesh.cpp.

References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::Clear().


Member Function Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
std::set< unsigned > TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::CalculateBoundaryOfFlaggedRegion (  )  [inline]

Return the set of nodes which are on the boundary of the flagged region(s).

Definition at line 550 of file TetrahedralMesh.cpp.

References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorBegin(), and AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorEnd().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
bool TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::CheckIsConforming (  )  [inline]

Check whether mesh is conforming Conforming (defn.): the intersection of two elements should be either the empty set, a vertex, an edge or a face.

It may be possible to construct non-conforming meshes which contain internal faces owned by only one element: two coplanar triangular faces of two elements form a square, but the same square on the adjacent pair of elements is formed by splitting the diagonal the other way.

Returns:
false if there are any orphaned internal faces

Definition at line 225 of file TetrahedralMesh.cpp.

References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorBegin(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorEnd(), and AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumBoundaryElements().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::Clear (  )  [inline, virtual]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ConstructFromMeshReader ( AbstractMeshReader< ELEMENT_DIM, SPACE_DIM > &  rMeshReader  )  [inline, virtual]

Construct the mesh using a MeshReader.

Parameters:
rMeshReader the mesh reader

Implements AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >.

Definition at line 59 of file TetrahedralMesh.cpp.

References Node< SPACE_DIM >::AddNodeAttribute(), ElementData::AttributeValue, AbstractMeshReader< ELEMENT_DIM, SPACE_DIM >::GetMeshFileBaseName(), AbstractMeshReader< ELEMENT_DIM, SPACE_DIM >::GetNextElementData(), AbstractMeshReader< ELEMENT_DIM, SPACE_DIM >::GetNextFaceData(), AbstractMeshReader< ELEMENT_DIM, SPACE_DIM >::GetNextNode(), AbstractMeshReader< ELEMENT_DIM, SPACE_DIM >::GetNodeAttributes(), AbstractMeshReader< ELEMENT_DIM, SPACE_DIM >::GetNumElementAttributes(), AbstractMeshReader< ELEMENT_DIM, SPACE_DIM >::GetNumElements(), AbstractMeshReader< ELEMENT_DIM, SPACE_DIM >::GetNumFaceAttributes(), AbstractMeshReader< ELEMENT_DIM, SPACE_DIM >::GetNumFaces(), AbstractMeshReader< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryElements, AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryNodes, AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElements, AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mMeshFileBaseName, AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mNodes, ElementData::NodeIndices, TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshJacobianCachedData(), AbstractMeshReader< ELEMENT_DIM, SPACE_DIM >::Reset(), and AbstractElement< ELEMENT_DIM, SPACE_DIM >::SetRegion().

Referenced by CylindricalHoneycombMeshGenerator::CylindricalHoneycombMeshGenerator(), and HoneycombMeshGenerator::HoneycombMeshGenerator().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::EdgeIterator TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::EdgesBegin (  )  [inline]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::EdgeIterator TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::EdgesEnd (  )  [inline]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
template<class MESHER_IO >
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ExportToMesher ( NodeMap map,
MESHER_IO &  mesherInput,
int *  elementList = NULL 
) [inline, protected]

Export the mesh (currently only the nodes) to an external mesher This is determined at compile time when the MESHER_IO template is instantiated to either

  • triangulateio (for triangle remesher in 2D)
  • tetgenio (for tetgen remesher in 3D) Since conditional compilation is not allowed, care must be taken to only use common data members in this method
    Parameters:
    map is a NodeMap which associates the indices of nodes in the old mesh with indices of nodes in the new mesh. This should be created with the correct size (NumAllNodes)
    mesherInput is a triangulateio or tetgenio class (decided at compile time) Note that only nodes are exported and thus any late mesh is based on the convex hull
    elementList is a pointer to either mesherInput.trianglelist or mesherInput.tetrahedronlist (used when we are not remeshing, but converting the form of an existing mesh) Note that this should have been re-malloced and mesherInput.numberoftriangles or mesherInput.numberoftetrahedra should be set prior to the call when elementList is non-NULL

Definition at line 968 of file TetrahedralMesh.cpp.

References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorBegin(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorEnd(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumAllNodes(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mNodes, NodeMap::SetDeleted(), and NodeMap::SetNewIndex().

Referenced by MutableMesh< ELEMENT_DIM, SPACE_DIM >::ReMesh().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::FlagElementsNotContainingNodes ( const std::set< unsigned rNodes  )  [inline]

Flag all elements not containing ANY of the given nodes

Parameters:
rNodes set of nodes to check for

Definition at line 658 of file TetrahedralMesh.cpp.

References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorBegin(), and AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorEnd().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::FreeTriangulateIo ( triangulateio &  mesherIo  )  [inline, protected]

Convenience method to tidy up a triangleio data structure after use

Parameters:
mesherIo is a triangulateio class

Definition at line 949 of file TetrahedralMesh.cpp.

Referenced by MutableMesh< ELEMENT_DIM, SPACE_DIM >::ReMesh().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
double TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetAngleBetweenNodes ( unsigned  indexA,
unsigned  indexB 
) [inline]

Calcuate the angle between the node at indexB and the x axis about the node at indexA. The angle returned is in the range (-pi,pi].

Parameters:
indexA a node index
indexB a node index

Definition at line 618 of file TetrahedralMesh.cpp.

References EXCEPTION, and AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mNodes.

Referenced by DiscreteSystemForceCalculator::CalculateFtAndFn(), and DiscreteSystemForceCalculator::GetSamplingAngles().

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

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

Parameters:
rTestPoint reference to the point
strict Should the element returned contain the point in the interior and not on an edge/face/vertex (default = not strict)
testElements a set of guesses for the element (a set of element indices), to be checked first for potential efficiency improvements. (default = empty set)
onlyTryWithTestElements Do 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).

Definition at line 362 of file TetrahedralMesh.cpp.

References EXCEPTION, AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumElements(), and AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElements.

Referenced by FineCoarseMeshPair< DIM >::ComputeCoarseElementForGivenPoint(), FineCoarseMeshPair< DIM >::ComputeFineElementAndWeightForGivenPoint(), CellBasedPdeHandler< DIM >::FindCoarseElementContainingCell(), CellBasedPdeHandler< DIM >::InitialiseCellPdeElementMap(), VolumeDependentAveragedSourcePde< DIM >::SetupSourceTerms(), and AveragedSourcePde< DIM >::SetupSourceTerms().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetContainingElementIndexWithInitialGuess ( const ChastePoint< SPACE_DIM > &  rTestPoint,
unsigned  startingElementGuess,
bool  strict = false 
) [inline]

Return the element index for the first element that contains a test point. Like GetContainingElementIndex but uses the user given element (M say) as the first element checked, and then checks M+1,M+2,..,Ne,0,1..

Parameters:
rTestPoint reference to the point
startingElementGuess Which element to try first.
strict Should the element returned contain the point in the interior and not on an edge/face/vertex (default = not strict)

Definition at line 409 of file TetrahedralMesh.cpp.

References EXCEPTION, AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumElements(), and AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElements.

Referenced by CellBasedPdeHandler< DIM >::UpdateCellPdeElementMap().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
std::vector< unsigned > TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetContainingElementIndices ( const ChastePoint< SPACE_DIM > &  rTestPoint  )  [inline]

Return all element indices for elements that are known to contain a test point.

Parameters:
rTestPoint reference to the point

Definition at line 512 of file TetrahedralMesh.cpp.

References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElements.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void TetrahedralMesh< 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 [inline, virtual]

Get the Jacobian matrix, its inverse and its determinant for a given element.

Parameters:
elementIndex index of the element in the mesh
rJacobian the Jacobian matrix
rJacobianDeterminant the determinant of the Jacobian matrix
rInverseJacobian the inverse Jacobian matrix

Reimplemented from AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >.

Reimplemented in NonCachedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >.

Definition at line 904 of file TetrahedralMesh.cpp.

References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementInverseJacobians, TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementJacobianDeterminants, and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementJacobians.

Referenced by VertexMesh< ELEMENT_DIM, SPACE_DIM >::GenerateVerticesFromElementCircumcentres(), and CellwiseDataGradient< DIM >::SetupGradients().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetJacobianForElement ( unsigned  elementIndex,
c_matrix< double, SPACE_DIM, SPACE_DIM > &  rJacobian,
double rJacobianDeterminant 
) const [inline, virtual]

Get the Jacobian matrix and its determinant for a given element.

Parameters:
elementIndex index of the element in the mesh
rJacobian the Jacobian matrix
rJacobianDeterminant the determinant of the Jacobian matrix

Reimplemented in NonCachedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >.

Definition at line 895 of file TetrahedralMesh.cpp.

References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementJacobianDeterminants, and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementJacobians.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNearestElementIndex ( const ChastePoint< SPACE_DIM > &  rTestPoint  )  [inline]

Return the element index for an element is closest to the testPoint.

"Closest" means that the minimum interpolation weights for the testPoint are maximised for this element.

Parameters:
rTestPoint reference to the point

Definition at line 455 of file TetrahedralMesh.cpp.

References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElements.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNearestElementIndexFromTestElements ( const ChastePoint< SPACE_DIM > &  rTestPoint,
std::set< unsigned testElements 
) [inline]

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

Parameters:
rTestPoint reference to the point
testElements a set of elements (element indices) to look in

Definition at line 481 of file TetrahedralMesh.cpp.

References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElements.

Referenced by FineCoarseMeshPair< DIM >::ComputeCoarseElementForGivenPoint(), and FineCoarseMeshPair< DIM >::ComputeFineElementAndWeightForGivenPoint().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
double TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetSurfaceArea (  )  [inline]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
double TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetVolume (  )  [inline]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetWeightedDirectionForBoundaryElement ( unsigned  elementIndex,
c_vector< double, SPACE_DIM > &  rWeightedDirection,
double rJacobianDeterminant 
) const [inline, virtual]

Get the weighted direction and the determinant of the Jacobian for a given boundary element.

Parameters:
elementIndex index of the element in the mesh
rWeightedDirection the weighted direction
rJacobianDeterminant the determinant of the Jacobian matrix

Reimplemented from AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >.

Reimplemented in NonCachedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >.

Definition at line 923 of file TetrahedralMesh.cpp.

References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryElementJacobianDeterminants, and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryElementWeightedDirections.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetWeightedDirectionForElement ( unsigned  elementIndex,
c_vector< double, SPACE_DIM > &  rWeightedDirection,
double rJacobianDeterminant 
) const [inline, virtual]

Get the weighted direction and the determinant of the Jacobian for a given element.

Parameters:
elementIndex index of the element in the mesh
rWeightedDirection the weighted direction
rJacobianDeterminant the determinant of the Jacobian matrix

Reimplemented in NonCachedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >.

Definition at line 914 of file TetrahedralMesh.cpp.

References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementJacobianDeterminants, and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementWeightedDirections.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
template<class MESHER_IO >
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ImportFromMesher ( MESHER_IO &  mesherOutput,
unsigned  numberOfElements,
int *  elementList,
unsigned  numberOfFaces,
int *  faceList,
int *  edgeMarkerList 
) [inline, protected]

Import the mesh from an external mesher This is determined at compile time when the MESHER_IO template is instantiated to either

  • triangulateio (for triangle remesher in 2D)
  • tetgenio (for tetgen remesher in 3D) Since conditional compilation is not allowed, care must be taken to only use common data members in this method (or add arguments
    Todo:
    #1545 ...)
    Parameters:
    mesherOutput is a triangulateio or tetgenio class (decided at compile time)
    numberOfElements is a copy of either mesherOutput.numberoftriangles or mesherOutput.numberoftetrahedra
    elementList is a pointer to either mesherOutput.trianglelist or mesherOutput.tetrahedronlist
    numberOfFaces is a copy of either mesherOutput.edges or mesherOutput.numberoftrifaces
    faceList is a pointer to either mesherOutput.edgelist or mesherOutput.trifacelist
    edgeMarkerList is a pointer to either mesherOutput.edgemarkerlist or NULL

Definition at line 1019 of file TetrahedralMesh.cpp.

References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::Clear(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryElements, AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryNodes, AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElements, AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mNodes, and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshJacobianCachedData().

Referenced by MutableMesh< ELEMENT_DIM, SPACE_DIM >::ReMesh().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::InitialiseTriangulateIo ( triangulateio &  mesherIo  )  [inline, protected]

Convenience method to tidy up a triangleio data structure before use

Parameters:
mesherIo is a triangulateio class

Definition at line 931 of file TetrahedralMesh.cpp.

Referenced by MutableMesh< ELEMENT_DIM, SPACE_DIM >::ReMesh().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::PermuteNodes ( const std::vector< unsigned > &  perm  )  [inline]

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

Parameters:
perm is a vector containing the new indices

Definition at line 333 of file TetrahedralMesh.cpp.

References AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumAllNodes(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mNodes, and AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mNodesPermutation.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::PermuteNodes (  )  [inline, virtual]

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

Reimplemented from AbstractMesh< ELEMENT_DIM, SPACE_DIM >.

Definition at line 321 of file TetrahedralMesh.cpp.

References RandomNumberGenerator::Instance(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mNodes, and RandomNumberGenerator::Shuffle().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ReadNodesPerProcessorFile ( const std::string &  rNodesPerProcessorFile  )  [inline, virtual]

Read in the number of nodes per processor from file.

Parameters:
rNodesPerProcessorFile the name of the file

Reimplemented from AbstractMesh< ELEMENT_DIM, SPACE_DIM >.

Definition at line 180 of file TetrahedralMesh.cpp.

References EXCEPTION, PetscTools::GetMyRank(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), PetscTools::GetNumProcs(), and AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mpDistributedVectorFactory.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshJacobianCachedData (  )  [inline, virtual]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshMesh (  )  [inline, virtual]

Overridden RefreshMesh method. This method calls RefreshJacobianCachedData.

Reimplemented from AbstractMesh< ELEMENT_DIM, SPACE_DIM >.

Definition at line 819 of file TetrahedralMesh.cpp.

References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshJacobianCachedData().

Referenced by MutableMesh< ELEMENT_DIM, SPACE_DIM >::RescaleMeshFromBoundaryNode().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
template<class Archive >
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::serialize ( Archive &  archive,
const unsigned int  version 
) [inline, private]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::SolveBoundaryElementMapping ( unsigned  index  )  const [inline, protected, virtual]

Overridden solve boundary element mapping method.

Parameters:
index the global index of the boundary element

Implements AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >.

Definition at line 839 of file TetrahedralMesh.cpp.

References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryElements.

Referenced by NonCachedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetWeightedDirectionForBoundaryElement().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::SolveElementMapping ( unsigned  index  )  const [inline, protected, virtual]

Overridden solve element mapping method.

Parameters:
index the global index of the element

Implements AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >.

Definition at line 832 of file TetrahedralMesh.cpp.

References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElements.

Referenced by NonCachedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetInverseJacobianForElement().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::SolveNodeMapping ( unsigned  index  )  const [inline, protected, virtual]

Overridden solve node mapping method.

Parameters:
index the global index of the node

Implements AbstractMesh< ELEMENT_DIM, SPACE_DIM >.

Definition at line 825 of file TetrahedralMesh.cpp.

References AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mNodes.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::UnflagAllElements (  )  [inline]

Friends And Related Function Documentation

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

Member Data Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
std::vector<double> TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryElementJacobianDeterminants [protected]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
std::vector< c_vector<double, SPACE_DIM> > TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryElementWeightedDirections [protected]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
std::vector< c_matrix<double, ELEMENT_DIM, SPACE_DIM> > TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementInverseJacobians [protected]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
std::vector<double> TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementJacobianDeterminants [protected]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
std::vector< c_matrix<double, SPACE_DIM, ELEMENT_DIM> > TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementJacobians [protected]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
std::vector< c_vector<double, SPACE_DIM> > TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementWeightedDirections [protected]

The documentation for this class was generated from the following files:
Generated on Thu Dec 22 13:07:57 2011 for Chaste by  doxygen 1.6.3