Chaste
Release::2018.1
|
#include <TetrahedralMesh.hpp>
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 | GetContainingElementIndexWithInitialGuess (const ChastePoint< SPACE_DIM > &rTestPoint, unsigned startingElementGuess, bool strict=false) |
unsigned | GetNearestElementIndex (const ChastePoint< SPACE_DIM > &rTestPoint) |
std::vector< unsigned > | GetContainingElementIndices (const ChastePoint< SPACE_DIM > &rTestPoint) |
virtual void | Clear () |
double | GetAngleBetweenNodes (unsigned indexA, unsigned indexB) |
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 () |
Public Member Functions inherited from AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > | |
ElementIterator | GetElementIteratorBegin (bool skipDeletedElements=true) |
ElementIterator | GetElementIteratorEnd () |
AbstractTetrahedralMesh () | |
virtual | ~AbstractTetrahedralMesh () |
virtual unsigned | GetNumElements () const |
virtual unsigned | GetNumLocalElements () const |
virtual unsigned | GetNumBoundaryElements () const |
virtual unsigned | GetNumLocalBoundaryElements () const |
unsigned | GetNumAllElements () const |
unsigned | GetNumAllBoundaryElements () const |
virtual unsigned | GetNumCableElements () const |
virtual unsigned | GetNumVertices () const |
virtual unsigned | GetMaximumNodeIndex () |
Element< ELEMENT_DIM, SPACE_DIM > * | GetElement (unsigned index) const |
BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > * | GetBoundaryElement (unsigned index) const |
void | ConstructFromMesh (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > &rOtherMesh) |
BoundaryElementIterator | GetBoundaryElementIteratorBegin () const |
BoundaryElementIterator | GetBoundaryElementIteratorEnd () const |
void | CheckOutwardNormals () |
virtual void | ConstructLinearMesh (unsigned width) |
virtual void | ConstructRectangularMesh (unsigned width, unsigned height, bool stagger=true) |
virtual void | ConstructCuboid (unsigned width, unsigned height, unsigned depth) |
void | ConstructRegularSlabMesh (double spaceStep, double width, double height=0, double depth=0) |
void | ConstructRegularSlabMeshWithDimensionSplit (unsigned dimension, double spaceStep, double width, double height=0, double depth=0) |
virtual bool | CalculateDesignatedOwnershipOfBoundaryElement (unsigned faceIndex) |
virtual bool | CalculateDesignatedOwnershipOfElement (unsigned elementIndex) |
unsigned | CalculateMaximumNodeConnectivityPerProcess () const |
virtual void | GetHaloNodeIndices (std::vector< unsigned > &rHaloIndices) const |
void | CalculateNodeExchange (std::vector< std::vector< unsigned > > &rNodesToSendPerProcess, std::vector< std::vector< unsigned > > &rNodesToReceivePerProcess) |
virtual c_vector< double, 2 > | CalculateMinMaxEdgeLengths () |
unsigned | GetContainingElementIndex (const ChastePoint< SPACE_DIM > &rTestPoint, bool strict=false, std::set< unsigned > testElements=std::set< unsigned >(), bool onlyTryWithTestElements=false) |
unsigned | GetNearestElementIndexFromTestElements (const ChastePoint< SPACE_DIM > &rTestPoint, std::set< unsigned > testElements) |
Public Member Functions inherited from AbstractMesh< ELEMENT_DIM, SPACE_DIM > | |
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 DistributedVectorFactory * | GetDistributedVectorFactory () |
virtual void | SetDistributedVectorFactory (DistributedVectorFactory *pFactory) |
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) |
bool | IsMeshChanging () const |
unsigned | CalculateMaximumContainingElementsPerProcess () const |
void | SetMeshHasChangedSinceLoading () |
Protected Member Functions | |
unsigned | SolveNodeMapping (unsigned index) const |
unsigned | SolveElementMapping (unsigned index) const |
unsigned | SolveBoundaryElementMapping (unsigned index) const |
template<class MESHER_IO > | |
void | ExportToMesher (NodeMap &map, MESHER_IO &mesherInput, int *elementList=nullptr) |
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 Member Functions inherited from AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > | |
void | SetElementOwnerships () |
Protected Member Functions inherited from AbstractMesh< ELEMENT_DIM, SPACE_DIM > | |
ChasteCuboid< SPACE_DIM > | CalculateBoundingBox (const std::vector< Node< SPACE_DIM > * > &rNodes) 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 |
Protected Attributes inherited from AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > | |
bool | mMeshIsLinear |
std::vector< Element < ELEMENT_DIM, SPACE_DIM > * > | mElements |
std::vector< BoundaryElement < ELEMENT_DIM-1, SPACE_DIM > * > | mBoundaryElements |
Protected Attributes inherited from AbstractMesh< ELEMENT_DIM, SPACE_DIM > | |
std::vector< Node< SPACE_DIM > * > | mNodes |
std::vector< Node< SPACE_DIM > * > | mBoundaryNodes |
DistributedVectorFactory * | mpDistributedVectorFactory |
std::vector< unsigned > | mNodePermutation |
std::string | mMeshFileBaseName |
bool | mMeshChangesDuringSimulation |
Private Member Functions | |
template<class Archive > | |
void | serialize (Archive &archive, const unsigned int version) |
Friends | |
class | TestTetrahedralMesh |
class | TestCryptSimulation2d |
class | boost::serialization::access |
Additional Inherited Members | |
Public Types inherited from AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > | |
typedef std::vector < BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > * >::const_iterator | BoundaryElementIterator |
Public Types inherited from AbstractMesh< ELEMENT_DIM, SPACE_DIM > | |
typedef std::vector< Node < SPACE_DIM > * >::const_iterator | BoundaryNodeIterator |
Forward declaration for triangle helper methods (used in MutableMesh QuadraticMesh) A concrete tetrahedral mesh class.
Definition at line 64 of file TetrahedralMesh.hpp.
TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::TetrahedralMesh | ( | ) |
Constructor.
Definition at line 62 of file TetrahedralMesh.cpp.
bool TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::CheckIsConforming | ( | ) |
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 236 of file TetrahedralMesh.cpp.
|
virtual |
Clear all the data in the mesh.
Reimplemented in NodesOnlyMesh< SPACE_DIM >, NodesOnlyMesh< DIM >, NodesOnlyMesh< 2 >, MutableMesh< ELEMENT_DIM, SPACE_DIM >, MutableMesh< ELEMENT_DIM, ELEMENT_DIM >, MutableMesh< 2, 2 >, MutableMesh< DIM, DIM >, and MutableMesh< SPACE_DIM, SPACE_DIM >.
Definition at line 459 of file TetrahedralMesh.cpp.
Referenced by MutableMesh< ELEMENT_DIM, SPACE_DIM >::Clear(), and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ImportFromMesher().
|
virtual |
Construct the mesh using a MeshReader.
rMeshReader | the mesh reader |
Implements AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >.
Definition at line 68 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(), AbstractMeshReader< ELEMENT_DIM, SPACE_DIM >::HasNodePermutation(), ElementData::NodeIndices, AbstractMeshReader< ELEMENT_DIM, SPACE_DIM >::Reset(), and AbstractElement< ELEMENT_DIM, SPACE_DIM >::SetAttribute().
Referenced by QuadraticMesh< DIM >::ConstructFromLinearMeshReader(), QuadraticMesh< DIM >::ConstructFromMeshReader(), NodesOnlyMesh< SPACE_DIM >::ConstructFromMeshReader(), VertexBasedCellPopulation< SPACE_DIM >::GetTetrahedralMeshForPdeModifier(), and HoneycombMeshGenerator::HoneycombMeshGenerator().
TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::EdgeIterator TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::EdgesBegin | ( | ) |
Definition at line 625 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElement(), and AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumAllElements().
TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::EdgeIterator TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::EdgesEnd | ( | ) |
Definition at line 636 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumAllElements().
Referenced by MeshBasedCellPopulationWithGhostNodes< DIM >::ApplyGhostForces(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::EdgeIterator::GetNodeA(), and MeshBasedCellPopulation< ELEMENT_DIM, SPACE_DIM >::SpringIterator::SpringIterator().
|
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 791 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().
|
protected |
Convenience method to tidy up a triangleio data structure after use
mesherIo | is a triangulateio class |
Definition at line 772 of file TetrahedralMesh.cpp.
double TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetAngleBetweenNodes | ( | unsigned | indexA, |
unsigned | indexB | ||
) |
indexA | a node index |
indexB | a node index |
Definition at line 482 of file TetrahedralMesh.cpp.
References EXCEPTION.
Referenced by DiscreteSystemForceCalculator::CalculateFtAndFn(), and DiscreteSystemForceCalculator::GetSamplingAngles().
unsigned TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetContainingElementIndexWithInitialGuess | ( | const ChastePoint< SPACE_DIM > & | rTestPoint, |
unsigned | startingElementGuess, | ||
bool | strict = false |
||
) |
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 372 of file TetrahedralMesh.cpp.
References EXCEPTION.
std::vector< unsigned > TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetContainingElementIndices | ( | const ChastePoint< SPACE_DIM > & | rTestPoint | ) |
rTestPoint | reference to the point |
Definition at line 444 of file TetrahedralMesh.cpp.
|
virtual |
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 727 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().
|
virtual |
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 718 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 | ( | const ChastePoint< SPACE_DIM > & | rTestPoint | ) |
"Closest" means that the minimum interpolation weights for the testPoint are maximised for this element.
rTestPoint | reference to the point |
Definition at line 417 of file TetrahedralMesh.cpp.
References EXCEPT_IF_NOT.
double TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetSurfaceArea | ( | ) |
Definition at line 302 of file TetrahedralMesh.cpp.
double TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetVolume | ( | ) |
Definition at line 287 of file TetrahedralMesh.cpp.
Referenced by CellPopulationAreaWriter< ELEMENT_DIM, SPACE_DIM >::Visit().
|
virtual |
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 746 of file TetrahedralMesh.cpp.
References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryElementJacobianDeterminants, and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryElementWeightedDirections.
|
virtual |
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 737 of file TetrahedralMesh.cpp.
References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementJacobianDeterminants, and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElementWeightedDirections.
|
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 842 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().
|
protected |
Convenience method to tidy up a triangleio data structure before use
mesherIo | is a triangulateio class |
Definition at line 754 of file TetrahedralMesh.cpp.
|
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 331 of file TetrahedralMesh.cpp.
References RandomNumberGenerator::Instance(), and RandomNumberGenerator::Shuffle().
void TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::PermuteNodes | ( | const std::vector< unsigned > & | perm | ) |
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 343 of file TetrahedralMesh.cpp.
|
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 190 of file TetrahedralMesh.cpp.
References EXCEPTION, PetscTools::GetMyRank(), and PetscTools::GetNumProcs().
|
virtual |
Update mElementJacobians, mElementWeightedDirections and mBoundaryElementWeightedDirections.
Reimplemented in NonCachedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >.
Definition at line 669 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 >::ImportFromMesher(), and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshMesh().
|
virtual |
Overridden RefreshMesh method. This method calls RefreshJacobianCachedData.
Reimplemented from AbstractMesh< ELEMENT_DIM, SPACE_DIM >.
Reimplemented in Cylindrical2dNodesOnlyMesh.
Definition at line 642 of file TetrahedralMesh.cpp.
References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshJacobianCachedData().
Referenced by Cylindrical2dNodesOnlyMesh::RefreshMesh().
|
inlineprivate |
Serialize the mesh.
archive | the archive |
version | the current version of this class |
Definition at line 78 of file TetrahedralMesh.hpp.
|
protectedvirtual |
Overridden solve boundary element mapping method.
index | the global index of the boundary element |
Implements AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >.
Definition at line 662 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryElements.
|
protectedvirtual |
Overridden solve element mapping method.
index | the global index of the element |
Implements AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >.
Definition at line 655 of file TetrahedralMesh.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mElements.
|
protectedvirtual |
Overridden solve node mapping method.
index | the global index of the node |
Implements AbstractMesh< ELEMENT_DIM, SPACE_DIM >.
Reimplemented in NodesOnlyMesh< DIM >, and NodesOnlyMesh< 2 >.
Definition at line 648 of file TetrahedralMesh.cpp.
References AbstractMesh< ELEMENT_DIM, SPACE_DIM >::mNodes.
|
friend |
Needed for serialization.
Definition at line 70 of file TetrahedralMesh.hpp.
|
protected |
Vector storing the determinant of the Jacobian matrix for each boundary element in the mesh.
Definition at line 125 of file TetrahedralMesh.hpp.
Referenced by TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetWeightedDirectionForBoundaryElement(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshJacobianCachedData(), and Cylindrical2dMesh::ReMesh().
|
protected |
Vector storing the weighted direction for each boundary element in the mesh.
Definition at line 122 of file TetrahedralMesh.hpp.
Referenced by TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetWeightedDirectionForBoundaryElement(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshJacobianCachedData(), and Cylindrical2dMesh::ReMesh().
|
protected |
Vector storing the inverse Jacobian matrix for each element in the mesh.
Definition at line 116 of file TetrahedralMesh.hpp.
Referenced by TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetInverseJacobianForElement(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshJacobianCachedData(), and Cylindrical2dMesh::UseTheseElementsToDecideMeshing().
|
protected |
Vector storing the Jacobian determinant for each element in the mesh.
Definition at line 119 of file TetrahedralMesh.hpp.
Referenced by TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetInverseJacobianForElement(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetJacobianForElement(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetWeightedDirectionForElement(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshJacobianCachedData(), and Cylindrical2dMesh::UseTheseElementsToDecideMeshing().
|
protected |
Vector storing the Jacobian matrix for each element in the mesh.
Definition at line 113 of file TetrahedralMesh.hpp.
Referenced by TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetInverseJacobianForElement(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetJacobianForElement(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshJacobianCachedData(), and Cylindrical2dMesh::UseTheseElementsToDecideMeshing().
|
protected |
Vector storing the weighted direction for each element in the mesh.
Definition at line 110 of file TetrahedralMesh.hpp.
Referenced by TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetWeightedDirectionForElement(), and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshJacobianCachedData().