#include <TetrahedralMesh.hpp>
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< unsigned > | GetContainingElementIndices (const ChastePoint< SPACE_DIM > &rTestPoint) |
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 |
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< 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 |
Classes | |
class | EdgeIterator |
Definition at line 57 of file TetrahedralMesh.hpp.
TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::TetrahedralMesh | ( | ) | [inline] |
Constructor.
Definition at line 59 of file TetrahedralMesh.cpp.
References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::Clear().
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 72 of file TetrahedralMesh.hpp.
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 823 of file TetrahedralMesh.cpp.
References AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mNodes.
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 830 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 >::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 837 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryElements.
Referenced by NonCachedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetWeightedDirectionForBoundaryElement().
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
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 965 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().
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
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 1016 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().
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::InitialiseTriangulateIo | ( | triangulateio & | mesherIo | ) | [inline, protected] |
Convenience method to tidy up a triangleio data structure before use
mesherIo | is a triangulateio class |
Definition at line 930 of file TetrahedralMesh.cpp.
Referenced by MutableMesh< ELEMENT_DIM, SPACE_DIM >::ReMesh().
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::FreeTriangulateIo | ( | triangulateio & | mesherIo | ) | [inline, protected] |
Convenience method to tidy up a triangleio data structure after use
mesherIo | is a triangulateio class |
Definition at line 947 of file TetrahedralMesh.cpp.
Referenced by MutableMesh< ELEMENT_DIM, SPACE_DIM >::ReMesh().
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 65 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 QuadraticMesh< DIM >::ConstructFromLinearMeshReader(), QuadraticMesh< DIM >::ConstructFromMeshReader(), CellBasedSimulationWithPdes< DIM >::CreateCoarsePdeMesh(), and HoneycombMeshGenerator::HoneycombMeshGenerator().
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 176 of file TetrahedralMesh.cpp.
References EXCEPTION, PetscTools::GetMyRank(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), PetscTools::GetNumProcs(), and AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mpDistributedVectorFactory.
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 223 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorBegin(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorEnd(), and AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumBoundaryElements().
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 272 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 MeshBasedCellPopulation< DIM >::WriteCellPopulationVolumeResultsToFile().
double TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetSurfaceArea | ( | ) | [inline] |
Return the surface area of the mesh.
Definition at line 287 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetBoundaryElementIteratorBegin(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetBoundaryElementIteratorEnd(), and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryElementJacobianDeterminants.
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 817 of file TetrahedralMesh.cpp.
References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshJacobianCachedData().
Referenced by MutableMesh< ELEMENT_DIM, SPACE_DIM >::RescaleMeshFromBoundaryNode().
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 316 of file TetrahedralMesh.cpp.
References RandomNumberGenerator::Instance(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mNodes, and RandomNumberGenerator::Shuffle().
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).
perm | is a vector containing the new indices |
Definition at line 327 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.
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
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 356 of file TetrahedralMesh.cpp.
References EXCEPTION, AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumElements(), and AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElements.
Referenced by FineCoarseMeshPair< DIM >::ComputeFineElementAndWeightForGivenPoint(), CellBasedSimulationWithPdes< DIM >::FindCoarseElementContainingCell(), CellBasedSimulationWithPdes< DIM >::InitialiseCoarsePdeMesh(), and AveragedSourcePde< DIM >::SetupSourceTerms().
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..
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 403 of file TetrahedralMesh.cpp.
References EXCEPTION, AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumElements(), and AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElements.
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.
rTestPoint | reference to the point |
Definition at line 448 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElements.
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.
rTestPoint | reference to the point | |
testElements | a set of elements (element indices) to look in |
Definition at line 475 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElements.
Referenced by FineCoarseMeshPair< DIM >::ComputeFineElementAndWeightForGivenPoint().
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.
rTestPoint | reference to the point |
Definition at line 508 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElements.
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 524 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 MutableMesh< ELEMENT_DIM, SPACE_DIM >::Clear(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ImportFromMesher(), and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::TetrahedralMesh().
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 547 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 613 of file TetrahedralMesh.cpp.
References EXCEPTION, and AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mNodes.
Referenced by DiscreteSystemForceCalculator::CalculateFtAndFn(), and DiscreteSystemForceCalculator::GetSamplingAngles().
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::UnflagAllElements | ( | ) | [inline] |
Unflag all elements in the mesh.
Definition at line 642 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorBegin(), and AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorEnd().
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 653 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorBegin(), and AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorEnd().
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshJacobianCachedData | ( | ) | [inline, virtual] |
Update mElementJacobians, mElementWeightedDirections and mBoundaryElementWeightedDirections.
Reimplemented in NonCachedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >.
Definition at line 844 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 >::ImportFromMesher(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshMesh(), and MutableMesh< ELEMENT_DIM, SPACE_DIM >::ReMesh().
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 893 of file TetrahedralMesh.cpp.
References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementJacobianDeterminants, and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementJacobians.
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 902 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().
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 912 of file TetrahedralMesh.cpp.
References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementJacobianDeterminants, and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementWeightedDirections.
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 921 of file TetrahedralMesh.cpp.
References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryElementJacobianDeterminants, and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryElementWeightedDirections.
TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::EdgeIterator TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::EdgesBegin | ( | ) | [inline] |
Definition at line 800 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElement(), and AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumAllElements().
Referenced by FineCoarseMeshPair< DIM >::SetUpBoxes(), MeshBasedCellPopulation< DIM >::SpringsBegin(), and MeshBasedCellPopulationWithGhostNodes< DIM >::UpdateGhostPositions().
TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::EdgeIterator TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::EdgesEnd | ( | ) | [inline] |
Definition at line 811 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(), FineCoarseMeshPair< DIM >::SetUpBoxes(), MeshBasedCellPopulation< DIM >::SpringsEnd(), and MeshBasedCellPopulationWithGhostNodes< DIM >::UpdateGhostPositions().
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 64 of file TetrahedralMesh.hpp.
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 101 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().
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 104 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_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 107 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 110 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_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 113 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<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 116 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().