#include <TetrahedralMesh.hpp>
Inherits AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >.
Inherited by MutableMesh< ELEMENT_DIM, SPACE_DIM >, MutableMesh< 2, 2 >, MutableMesh< DIM, DIM >, and NonCachedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >.
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 | 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) |
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) |
void | RefreshMesh () |
void | PermuteNodes () |
void | PermuteNodesWithMetisBinaries (unsigned numProcs) |
void | PermuteNodes (std::vector< unsigned > &perm) |
unsigned | GetContainingElementIndex (ChastePoint< SPACE_DIM > testPoint, bool strict=false, std::set< unsigned > testElements=std::set< unsigned >(), bool onlyTryWithTestElements=false) |
unsigned | GetContainingElementIndexWithInitialGuess (ChastePoint< SPACE_DIM > testPoint, unsigned startingElementGuess, bool strict=false) |
unsigned | GetNearestElementIndex (ChastePoint< SPACE_DIM > testPoint) |
unsigned | GetNearestElementIndexFromTestElements (ChastePoint< SPACE_DIM > testPoint, std::set< unsigned > testElements) |
todo #1299 refactor common code with GetNearestElementIndex | |
std::vector< unsigned > | GetContainingElementIndices (ChastePoint< SPACE_DIM > testPoint) |
virtual void | Clear () |
std::set< unsigned > | CalculateBoundaryOfFlaggedRegion () |
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 |
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< double > | mElementJacobianDeterminants |
std::vector< c_vector< double, SPACE_DIM > > | mBoundaryElementWeightedDirections |
std::vector< double > | mBoundaryElementJacobianDeterminants |
Private Member Functions | |
template<class Archive > | |
void | serialize (Archive &archive, const unsigned int version) |
Friends | |
class | TestTetrahedralMesh |
class | TestCryptSimulation2d |
class | boost::serialization::access |
A concrete tetrahedral mesh class.
Definition at line 54 of file TetrahedralMesh.hpp.
TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::TetrahedralMesh | ( | ) | [inline] |
Constructor.
Definition at line 50 of file TetrahedralMesh.cpp.
References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::Clear().
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 818 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorBegin(), and AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorEnd().
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.
Definition at line 214 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorBegin(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorEnd(), and AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumBoundaryElements().
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::Clear | ( | ) | [inline, virtual] |
Clear all the data in the mesh.
Reimplemented in MutableMesh< ELEMENT_DIM, SPACE_DIM >, MutableMesh< 2, 2 >, and MutableMesh< DIM, DIM >.
Definition at line 795 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryElements, AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryNodes, AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElements, and AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mNodes.
Referenced by TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::TetrahedralMesh().
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ConstructFromMeshReader | ( | AbstractMeshReader< ELEMENT_DIM, SPACE_DIM > & | rMeshReader | ) | [inline, virtual] |
Construct the mesh using a MeshReader.
rMeshReader | the mesh reader |
Implements AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >.
Definition at line 56 of file TetrahedralMesh.cpp.
References 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 >::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 HoneycombMeshGenerator::HoneycombMeshGenerator(), and MutableMesh< ELEMENT_DIM, SPACE_DIM >::ReMesh().
TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::EdgeIterator TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::EdgesBegin | ( | ) | [inline] |
Definition at line 1051 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElement(), and AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumAllElements().
Referenced by VoronoiTessellation< DIM >::Initialise(), FineCoarseMeshPair< DIM >::SetUpBoxesOnFineMesh(), MeshBasedTissue< DIM >::SpringsBegin(), MeshBasedTissueWithGhostNodes< DIM >::UpdateGhostPositions(), and VertexMesh< ELEMENT_DIM, SPACE_DIM >::VertexMesh().
TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::EdgeIterator TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::EdgesEnd | ( | ) | [inline] |
Definition at line 1062 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumAllElements().
Referenced by TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::EdgeIterator::GetNodeA(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::EdgeIterator::GetNodeB(), VoronoiTessellation< DIM >::Initialise(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::EdgeIterator::operator++(), FineCoarseMeshPair< DIM >::SetUpBoxesOnFineMesh(), MeshBasedTissue< DIM >::SpringsEnd(), MeshBasedTissueWithGhostNodes< DIM >::UpdateGhostPositions(), and VertexMesh< ELEMENT_DIM, SPACE_DIM >::VertexMesh().
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::FlagElementsNotContainingNodes | ( | const std::set< unsigned > | rNodes | ) | [inline] |
Flag all elements not containing ANY of the given nodes
rNodes | set of nodes to check for |
Definition at line 924 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorBegin(), and AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorEnd().
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].
indexA | a node index | |
indexB | a node index |
Definition at line 884 of file TetrahedralMesh.cpp.
References AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mNodes.
unsigned TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetContainingElementIndex | ( | ChastePoint< SPACE_DIM > | testPoint, | |
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
testPoint | ||
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 621 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumElements(), and AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElements.
Referenced by FineCoarseMeshPair< DIM >::ComputeFineElementsAndWeightsForCoarseQuadPoints().
unsigned TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetContainingElementIndexWithInitialGuess | ( | ChastePoint< SPACE_DIM > | testPoint, | |
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..
testPoint | ||
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 669 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumElements(), and AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElements.
std::vector< unsigned > TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetContainingElementIndices | ( | ChastePoint< SPACE_DIM > | testPoint | ) | [inline] |
Return all element indices for elements that are known to contain a test point.
testPoint | the point |
Definition at line 779 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElements.
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.
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 1150 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 VoronoiTessellation< DIM >::GenerateVerticesFromElementCircumcentres(), and VertexMesh< ELEMENT_DIM, SPACE_DIM >::GenerateVerticesFromElementCircumcentres().
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.
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 1141 of file TetrahedralMesh.cpp.
References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementJacobianDeterminants, and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementJacobians.
unsigned TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNearestElementIndex | ( | ChastePoint< SPACE_DIM > | testPoint | ) | [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.
testPoint | the point |
Definition at line 714 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElements.
unsigned TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNearestElementIndexFromTestElements | ( | ChastePoint< SPACE_DIM > | testPoint, | |
std::set< unsigned > | testElements | |||
) | [inline] |
todo #1299 refactor common code with GetNearestElementIndex
As with GetNearestElementIndex() except only searches in the given set of elements.
testPoint | the point | |
testElements | a set of elements (element indices) to look in |
Definition at line 747 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElements.
Referenced by FineCoarseMeshPair< DIM >::ComputeFineElementsAndWeightsForCoarseQuadPoints().
double TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetSurfaceArea | ( | ) | [inline] |
Return the surface area of the mesh.
Definition at line 278 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetBoundaryElementIteratorBegin(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetBoundaryElementIteratorEnd(), and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryElementJacobianDeterminants.
double TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetVolume | ( | ) | [inline] |
Return the volume of the mesh, calculated by adding the determinant of each element and dividing by n!, where n is the element dimension.
Definition at line 263 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorBegin(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorEnd(), and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementJacobianDeterminants.
Referenced by MeshBasedTissue< DIM >::WriteTissueAreaResultsToFile().
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.
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 1169 of file TetrahedralMesh.cpp.
References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryElementJacobianDeterminants, and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryElementWeightedDirections.
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.
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 1160 of file TetrahedralMesh.cpp.
References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementJacobianDeterminants, and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementWeightedDirections.
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::PermuteNodes | ( | 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).
perm | is a vector containing the new indices |
Definition at line 462 of file TetrahedralMesh.cpp.
References AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumAllNodes(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), and AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mNodes.
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::PermuteNodes | ( | ) | [inline, virtual] |
Permute the nodes 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 440 of file TetrahedralMesh.cpp.
References RandomNumberGenerator::Instance(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mNodes, and RandomNumberGenerator::randMod().
Referenced by TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::PermuteNodesWithMetisBinaries().
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::PermuteNodesWithMetisBinaries | ( | unsigned | numProcs | ) | [inline] |
Permute the nodes so that they appear in a different order in mNodes (and their mIndex's are altered accordingly) using Metis binaries.
numProcs | Number of processors (e.g. number of partitions) |
Definition at line 487 of file TetrahedralMesh.cpp.
References PetscTools::AmMaster(), PetscTools::Barrier(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumAllElements(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumAllNodes(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumElements(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), OutputFileHandler::GetOutputDirectoryFullPath(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElements, OutputFileHandler::OpenOutputFile(), and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::PermuteNodes().
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ReadNodesPerProcessorFile | ( | const std::string & | rNodesPerProcessorFile | ) | [inline, virtual] |
Read in the number of nodes per processor from file.
rNodesPerProcessorFile | the name of the file |
Reimplemented from AbstractMesh< ELEMENT_DIM, SPACE_DIM >.
Definition at line 167 of file TetrahedralMesh.cpp.
References PetscTools::GetMyRank(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), PetscTools::GetNumProcs(), and AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mpDistributedVectorFactory.
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshJacobianCachedData | ( | ) | [inline, virtual] |
Update mElementJacobians, mElementWeightedDirections and mBoundaryElementWeightedDirections.
Reimplemented in NonCachedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >.
Definition at line 1095 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetBoundaryElementIteratorBegin(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetBoundaryElementIteratorEnd(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorBegin(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorEnd(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumAllBoundaryElements(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumAllElements(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryElementJacobianDeterminants, TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryElementWeightedDirections, TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementInverseJacobians, TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementJacobianDeterminants, TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementJacobians, and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementWeightedDirections.
Referenced by TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ConstructFromMeshReader(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshMesh(), and MutableMesh< ELEMENT_DIM, SPACE_DIM >::ReMesh().
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 1068 of file TetrahedralMesh.cpp.
References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshJacobianCachedData().
Referenced by MutableMesh< ELEMENT_DIM, SPACE_DIM >::RescaleMeshFromBoundaryNode(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::Rotate(), and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::Translate().
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::Rotate | ( | double | theta | ) | [inline] |
Rotating a 2D mesh equates that rotation around the z-axis.
theta | is the angle of rotation in radians |
Definition at line 434 of file TetrahedralMesh.cpp.
References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RotateZ().
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::Rotate | ( | c_vector< double, 3 > | axis, | |
double | angle | |||
) | [inline] |
Do an angle axis rotation.
axis | is the axis of rotation (does not need to be normalised) | |
angle | is the angle of rotation in radians |
Definition at line 357 of file TetrahedralMesh.cpp.
References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::Rotate().
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::Rotate | ( | c_matrix< double, SPACE_DIM, SPACE_DIM > | rotationMatrix | ) | [inline] |
Do a general mesh rotation with a positive determinant orthonormal rotation matrix. This is the rotation method that actually does the work.
rotationMatrix | is a Ublas rotation matrix of the correct form |
Definition at line 342 of file TetrahedralMesh.cpp.
References AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumAllNodes(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mNodes, and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshMesh().
Referenced by TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::Rotate(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RotateX(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RotateY(), and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RotateZ().
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RotateX | ( | const double | theta | ) | [inline] |
Rotate the mesh about the x-axis.
theta | is the angle of rotation in radians |
Definition at line 382 of file TetrahedralMesh.cpp.
References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::Rotate().
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RotateY | ( | const double | theta | ) | [inline] |
Rotate the mesh about the y-axis.
theta | is the angle of rotation in radians |
Definition at line 398 of file TetrahedralMesh.cpp.
References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::Rotate().
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RotateZ | ( | const double | theta | ) | [inline] |
Rotate the mesh about the z-axis.
theta | is the angle of rotation in radians |
Definition at line 416 of file TetrahedralMesh.cpp.
References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::Rotate().
Referenced by TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::Rotate().
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::serialize | ( | Archive & | archive, | |
const unsigned int | version | |||
) | [inline, private] |
Serialize the mesh.
archive | the archive | |
version | the current version of this class |
Reimplemented from AbstractMesh< ELEMENT_DIM, SPACE_DIM >.
Reimplemented in MutableMesh< ELEMENT_DIM, SPACE_DIM >, NonCachedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >, QuadraticMesh< DIM >, Cylindrical2dMesh, MutableMesh< 2, 2 >, and MutableMesh< DIM, DIM >.
Definition at line 69 of file TetrahedralMesh.hpp.
unsigned TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::SolveBoundaryElementMapping | ( | unsigned | index | ) | const [inline, protected, virtual] |
Overridden solve boundary element mapping method.
index | the global index of the boundary element |
Implements AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >.
Definition at line 1088 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryElements.
Referenced by NonCachedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetWeightedDirectionForBoundaryElement().
unsigned TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::SolveElementMapping | ( | unsigned | index | ) | const [inline, protected, virtual] |
Overridden solve element mapping method.
index | the global index of the element |
Implements AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >.
Definition at line 1081 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElements.
Referenced by NonCachedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetInverseJacobianForElement().
unsigned TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::SolveNodeMapping | ( | unsigned | index | ) | const [inline, protected, virtual] |
Overridden solve node mapping method.
index | the global index of the node |
Implements AbstractMesh< ELEMENT_DIM, SPACE_DIM >.
Definition at line 1074 of file TetrahedralMesh.cpp.
References AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mNodes.
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::Translate | ( | const double | xMovement = 0.0 , |
|
const double | yMovement = 0.0 , |
|||
const double | zMovement = 0.0 | |||
) | [inline] |
Translate the mesh given the coordinate displacements separately.
xMovement | is the x-displacement (defaults to 0.0) | |
yMovement | is the y-displacement (defaults to 0.0) | |
zMovement | is the z-displacement (defaults to 0.0) |
Definition at line 307 of file TetrahedralMesh.cpp.
References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::Translate().
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::Translate | ( | const c_vector< double, SPACE_DIM > & | rDisplacement | ) | [inline] |
Translate the mesh given the displacement vector. This is the translation method that actually does the work.
rDisplacement | is a translation vector of the correct size |
Definition at line 328 of file TetrahedralMesh.cpp.
References AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumAllNodes(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mNodes, and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshMesh().
Referenced by HoneycombMeshGenerator::GetCircularMesh(), and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::Translate().
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::UnflagAllElements | ( | ) | [inline] |
Unflag all elements in the mesh.
Definition at line 913 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorBegin(), and AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorEnd().
friend class boost::serialization::access [friend] |
Needed for serialization.
Reimplemented from AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >.
Reimplemented in MutableMesh< ELEMENT_DIM, SPACE_DIM >, NonCachedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >, QuadraticMesh< DIM >, Cylindrical2dMesh, MutableMesh< 2, 2 >, and MutableMesh< DIM, DIM >.
Definition at line 61 of file TetrahedralMesh.hpp.
std::vector<double> TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryElementJacobianDeterminants [protected] |
Vector storing the determinant of the Jacobian matrix for each boundary element in the mesh.
Definition at line 113 of file TetrahedralMesh.hpp.
Referenced by TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetSurfaceArea(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetWeightedDirectionForBoundaryElement(), MutableMesh< ELEMENT_DIM, SPACE_DIM >::MoveMergeNode(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshJacobianCachedData(), MutableMesh< ELEMENT_DIM, SPACE_DIM >::ReIndex(), and Cylindrical2dMesh::ReMesh().
std::vector< c_vector<double, SPACE_DIM> > TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryElementWeightedDirections [protected] |
Vector storing the weighted direction for each boundary element in the mesh.
Definition at line 110 of file TetrahedralMesh.hpp.
Referenced by TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetWeightedDirectionForBoundaryElement(), MutableMesh< ELEMENT_DIM, SPACE_DIM >::MoveMergeNode(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshJacobianCachedData(), MutableMesh< ELEMENT_DIM, SPACE_DIM >::ReIndex(), and Cylindrical2dMesh::ReMesh().
std::vector< c_matrix<double, ELEMENT_DIM, SPACE_DIM> > TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementInverseJacobians [protected] |
Vector storing the inverse Jacobian matrix for each element in the mesh.
Definition at line 104 of file TetrahedralMesh.hpp.
Referenced by MutableMesh< ELEMENT_DIM, SPACE_DIM >::CheckIsVoronoi(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetInverseJacobianForElement(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshJacobianCachedData(), MutableMesh< ELEMENT_DIM, SPACE_DIM >::ReIndex(), and Cylindrical2dMesh::UseTheseElementsToDecideMeshing().
std::vector<double> TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementJacobianDeterminants [protected] |
Vector storing the Jacobian determinant for each element in the mesh.
Definition at line 107 of file TetrahedralMesh.hpp.
Referenced by TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetInverseJacobianForElement(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetJacobianForElement(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetVolume(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetWeightedDirectionForElement(), MutableMesh< ELEMENT_DIM, SPACE_DIM >::MoveMergeNode(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshJacobianCachedData(), MutableMesh< ELEMENT_DIM, SPACE_DIM >::ReIndex(), and Cylindrical2dMesh::UseTheseElementsToDecideMeshing().
std::vector< c_matrix<double, SPACE_DIM, ELEMENT_DIM> > TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementJacobians [protected] |
Vector storing the Jacobian matrix for each element in the mesh.
Definition at line 101 of file TetrahedralMesh.hpp.
Referenced by MutableMesh< ELEMENT_DIM, SPACE_DIM >::CheckIsVoronoi(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetInverseJacobianForElement(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetJacobianForElement(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshJacobianCachedData(), MutableMesh< ELEMENT_DIM, SPACE_DIM >::ReIndex(), and Cylindrical2dMesh::UseTheseElementsToDecideMeshing().
std::vector< c_vector<double, SPACE_DIM> > TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementWeightedDirections [protected] |
Vector storing the weighted direction for each element in the mesh.
Definition at line 98 of file TetrahedralMesh.hpp.
Referenced by TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetWeightedDirectionForElement(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshJacobianCachedData(), MutableMesh< ELEMENT_DIM, SPACE_DIM >::ReIndex(), and MutableMesh< ELEMENT_DIM, SPACE_DIM >::SetNode().