Chaste Release::3.1
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > Class Template Reference

#include <AbstractTetrahedralMesh.hpp>

Inheritance diagram for AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >:
Collaboration diagram for AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >:

List of all members.

Classes

class  ElementIterator

Public Types

typedef std::vector
< BoundaryElement< ELEMENT_DIM-1,
SPACE_DIM >
* >::const_iterator 
BoundaryElementIterator

Public Member Functions

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
Element< ELEMENT_DIM, SPACE_DIM > * GetElement (unsigned index) const
BoundaryElement< ELEMENT_DIM-1,
SPACE_DIM > * 
GetBoundaryElement (unsigned index) const
virtual void ConstructFromMeshReader (AbstractMeshReader< ELEMENT_DIM, SPACE_DIM > &rMeshReader)=0
void ConstructFromMesh (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > &rOtherMesh)
BoundaryElementIterator GetBoundaryElementIteratorBegin () const
BoundaryElementIterator GetBoundaryElementIteratorEnd () 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 GetWeightedDirectionForBoundaryElement (unsigned elementIndex, c_vector< double, SPACE_DIM > &rWeightedDirection, double &rJacobianDeterminant) 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)
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 ()

Protected Member Functions

void SetElementOwnerships ()

Protected Attributes

bool mMeshIsLinear
std::vector< Element
< ELEMENT_DIM, SPACE_DIM > * > 
mElements
std::vector< BoundaryElement
< ELEMENT_DIM-1, SPACE_DIM > * > 
mBoundaryElements

Private Member Functions

virtual unsigned SolveElementMapping (unsigned index) const =0
virtual unsigned SolveBoundaryElementMapping (unsigned index) const =0
template<class Archive >
void save (Archive &archive, const unsigned int version) const
template<class Archive >
void load (Archive &archive, const unsigned int version)

Friends

class AbstractConductivityTensors< ELEMENT_DIM, SPACE_DIM >
class boost::serialization::access

Detailed Description

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

Abstract base class for all tetrahedral meshes (inherits from AbstractMesh).

Definition at line 71 of file AbstractTetrahedralMesh.hpp.


Member Typedef Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
typedef std::vector<BoundaryElement<ELEMENT_DIM-1, SPACE_DIM> *>::const_iterator AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::BoundaryElementIterator

Definition of boundary element Iterator type.

Definition at line 295 of file AbstractTetrahedralMesh.hpp.


Constructor & Destructor Documentation

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

Constructor.

Definition at line 65 of file AbstractTetrahedralMesh.cpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::~AbstractTetrahedralMesh ( ) [virtual]

Virtual destructor, since this class has virtual methods.

Definition at line 71 of file AbstractTetrahedralMesh.cpp.


Member Function Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
bool AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::CalculateDesignatedOwnershipOfBoundaryElement ( unsigned  faceIndex) [virtual]

Determine whether or not the current process owns node 0 of this boundary element (tie breaker to determine which process writes to file for when two or more share ownership of a face).

Parameters:
faceIndexis the global index of the face

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

Definition at line 614 of file AbstractTetrahedralMesh.cpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
bool AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::CalculateDesignatedOwnershipOfElement ( unsigned  elementIndex) [virtual]

Determine whether or not the current process owns node 0 of this element (tie breaker to determine which process writes to file for when two or more share ownership of an element).

Parameters:
elementIndexis the global index of the element

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

Definition at line 631 of file AbstractTetrahedralMesh.cpp.

Referenced by AbstractFunctionalCalculator< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::Calculate().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::CalculateMaximumNodeConnectivityPerProcess ( ) const
Returns:
Iterates through local nodes and finds the node with the/a maximum number of containing elements for all locally owned nodes. At that representative node the node connectivity (number of nodes in forward star) is determined.

Useful for determining FEM matrix fill.

Definition at line 648 of file AbstractTetrahedralMesh.cpp.

References AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
c_vector< double, 2 > AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::CalculateMinMaxEdgeLengths ( ) [virtual]

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

Todo:
Should be const
Returns:
The minimum and maximum edge lengths in the mesh

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

Definition at line 837 of file AbstractTetrahedralMesh.cpp.

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

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::CalculateNodeExchange ( std::vector< std::vector< unsigned > > &  rNodesToSendPerProcess,
std::vector< std::vector< unsigned > > &  rNodesToReceivePerProcess 
)

Get the nodes which will need to be exchanged between remote processes. If we have an element which node indices outside the local [mLo, mHi) region then we know that those nodes will need to be recieved from a remote process, while those inside the range [mLo, mHi) will need to be sent

Parameters:
rNodesToSendPerProcess(output) a vector which will be of size GetNumProcs() where each internal vector except i=GetMyRank() contains an ordered list of indices of nodes to send to process i
rNodesToReceivePerProcess(output) a vector which will be of size GetNumProcs() for information to receive for process i

Definition at line 754 of file AbstractTetrahedralMesh.cpp.

References PetscTools::GetNumProcs().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::CheckOutwardNormals ( )

Check whether mesh has outward-facing normals.

This throws a suitable exception if an inward facing normal is found.

Returns:
nothing

Definition at line 180 of file AbstractTetrahedralMesh.cpp.

References Node< SPACE_DIM >::ContainingElementsBegin(), Node< SPACE_DIM >::ContainingElementsEnd(), EXCEPTION, AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex(), and Node< SPACE_DIM >::rGetLocation().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ConstructCuboid ( unsigned  width,
unsigned  height,
unsigned  depth 
) [virtual]

Construct a 3D cuboid grid on [0,width]x[0,height]x[0,depth].

Parameters:
widthwidth of the mesh (in the x-direction)
heightheight of the mesh (in the y-direction)
depthdepth of the mesh (in the z-direction).

In this method the width is also THE NUMBER OF ELEMENTS IN THE x-direction, and similarly with the y and z directions.

Overridden in DistributedTetrahedralMesh

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

Definition at line 381 of file AbstractTetrahedralMesh.cpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
virtual void AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ConstructFromMeshReader ( AbstractMeshReader< ELEMENT_DIM, SPACE_DIM > &  rMeshReader) [pure virtual]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ConstructLinearMesh ( unsigned  width) [virtual]

Construct a 1D linear grid on [0,width]

ELEMENT_DIM must be equal to 1. If SPACE_DIM > 1 then the y & z default to 0.0 for every node.

Parameters:
widthwidth of the mesh (in the x-direction)

In this method the width is also THE NUMBER OF ELEMENTS IN THE x-direction.

Overridden in DistributedTetrahedralMesh

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

Definition at line 239 of file AbstractTetrahedralMesh.cpp.

Referenced by QuadraticMesh< DIM >::ConstructLinearMesh().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ConstructRectangularMesh ( unsigned  width,
unsigned  height,
bool  stagger = true 
) [virtual]

Construct a 2D rectangular grid on [0,width]x[0,height].

Diagonals can be staggered so that there is no preferred diffusion propagation direction.

Parameters:
widthwidth of the mesh (in the x-direction)
heightheight of the mesh (in the y-direction)
staggerwhether the mesh should 'jumble' up the elements (defaults to true)

In this method the width is also THE NUMBER OF ELEMENTS IN THE x-direction, and similarly with the y direction.

Overridden in DistributedTetrahedralMesh

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

Definition at line 271 of file AbstractTetrahedralMesh.cpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ConstructRegularSlabMesh ( double  spaceStep,
double  width,
double  height = 0,
double  depth = 0 
)

Create a 1D mesh on [0, width], 2D mesh on [0, width]x[0 height] with staggering or 3D mesh on [0, width]x[0 height]x[0 depth with a given axis-aligned space step. If SPACE_DIM > ELEMENT_DIM then the y & z default to 0.0 for every node.

Parameters:
spaceStepThe axis-aligned space step
widthThe width (x-dimension)
heightThe height (y-dimension - ignored if ELEMENT_DIM is 1D)
depthThe depth (z-dimension -ignored in 1D and 2D)

Definition at line 566 of file AbstractTetrahedralMesh.cpp.

References EXCEPTION.

Referenced by CardiacElectroMechProbRegularGeom< DIM >::CardiacElectroMechProbRegularGeom(), CuboidMeshConstructor< ELEMENT_DIM, SPACE_DIM >::Construct(), and CellBasedPdeHandler< DIM >::UseCoarsePdeMesh().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > * AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetBoundaryElement ( unsigned  index) const

Get the boundary element with a given index in the mesh.

Parameters:
indexthe global index of the boundary element
Returns:
a pointer to the boundary element.

Definition at line 141 of file AbstractTetrahedralMesh.cpp.

Referenced by AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ConstructFromMesh(), BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::MergeFromArchive(), and Cylindrical2dMesh::ReconstructCylindricalMesh().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
Element< ELEMENT_DIM, SPACE_DIM > * AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElement ( unsigned  index) const
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ElementIterator AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorBegin ( bool  skipDeletedElements = true) [inline]

Get an iterator to the first element in the mesh.

Parameters:
skipDeletedElementswhether to include deleted element

Definition at line 671 of file AbstractTetrahedralMesh.hpp.

Referenced by AbstractCorrectionTermAssembler< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::AbstractCorrectionTermAssembler(), QuadraticMeshHelper< DIM >::AddInternalNodesToElements(), BidomainProblem< DIM >::AnalyseMeshForBath(), IncompressibleNonlinearElasticitySolver< DIM >::AssembleSystem(), CompressibleNonlinearElasticitySolver< DIM >::AssembleSystem(), AbstractFunctionalCalculator< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::Calculate(), AdaptiveTetrahedralMesh::ConstructFromDistributedMesh(), ExtendedBidomainTissue< SPACE_DIM >::CreateIntracellularConductivityTensorSecondCell(), AbstractContinuumMechanicsAssembler< DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >::DoAssemble(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ExportToMesher(), Cylindrical2dMesh::GenerateVectorsOfElementsStraddlingPeriodicBoundaries(), PapillaryFibreCalculator::GetRadiusVectors(), AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::Initialise(), Cylindrical2dMesh::ReconstructCylindricalMesh(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshJacobianCachedData(), AbstractContinuumMechanicsSolver< DIM >::RemovePressureDummyValuesThroughLinearInterpolation(), PapillaryFibreCalculator::SmoothStructureTensors(), Cylindrical2dMesh::UseTheseElementsToDecideMeshing(), AbstractNonlinearElasticitySolver< DIM >::WriteCurrentDeformationGradients(), and VtkMeshWriter< ELEMENT_DIM, SPACE_DIM >::WriteFilesUsingMesh().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetHaloNodeIndices ( std::vector< unsigned > &  rHaloIndices) const [virtual]

Utility method to give the functionality of iterating through the halo nodes of a process. Will return an empty std::vector (i.e. no halo nodes) unless overridden by distributed derived classes.

Parameters:
rHaloIndicesA vector to fill with the global indices of the nodes which are locally halos

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

Definition at line 698 of file AbstractTetrahedralMesh.cpp.

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

Compute the inverse Jacobian for a given element in the mesh.

Parameters:
elementIndexindex of an element
rJacobianthe Jacobian matrix
rJacobianDeterminantthe determinant of the Jacobian matrix
rInverseJacobianthe inverse Jacobian matrix

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

Definition at line 160 of file AbstractTetrahedralMesh.cpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumAllBoundaryElements ( ) const

Get the total number of boundary elements (including those marked as deleted).

Definition at line 116 of file AbstractTetrahedralMesh.cpp.

Referenced by Cylindrical2dMesh::ReconstructCylindricalMesh(), and TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::RefreshJacobianCachedData().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumCableElements ( ) const [virtual]

Get the number of cable elements that are actually in use.

This will always return zero until overridden in the MixedDimensionMesh class

Reimplemented in MixedDimensionMesh< ELEMENT_DIM, SPACE_DIM >.

Definition at line 123 of file AbstractTetrahedralMesh.cpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumLocalBoundaryElements ( ) const [virtual]

Get the number of local boundary elements that are in use on this process (only over-ridden when the mesh is distributed).

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

Definition at line 110 of file AbstractTetrahedralMesh.cpp.

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

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumLocalElements ( ) const [virtual]

Get the number of local elements that are in use on this process (only over-ridden when the mesh is distributed).

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

Definition at line 92 of file AbstractTetrahedralMesh.cpp.

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

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
unsigned AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumVertices ( ) const [virtual]

Get the number of vertices (nodes which are also corners of elements). For a linear mesh all nodes are vertices, so this method is a synonym for GetNumNodes. However, it is over-ridden in quadratic meshes.

Reimplemented in QuadraticMesh< DIM >.

Definition at line 128 of file AbstractTetrahedralMesh.cpp.

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

Compute the weighted direction for a given boundary element.

Parameters:
elementIndexindex of an element
rWeightedDirectionthe weighted direction vector
rJacobianDeterminantthe determinant of the Jacobian matrix

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

Definition at line 170 of file AbstractTetrahedralMesh.cpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
template<class Archive >
void AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::load ( Archive &  archive,
const unsigned int  version 
) [inline, private]

Loads a mesh by using TrianglesMeshReader and the location in ArchiveLocationInfo.

Parameters:
archive
version

Todo:
#1199 make this work for everything else...

Definition at line 187 of file AbstractTetrahedralMesh.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
template<class Archive >
void AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::save ( Archive &  archive,
const unsigned int  version 
) const [inline, private]

Archive the AbstractTetrahedralMesh. Note that this will write out a TrianglesMeshWriter file to wherever ArchiveLocationInfo has specified.

If the mesh is MutableMesh (or a subclass) the file is written by examining the current mesh.

If the mesh is not mutable then the file is a copy of the original file the mesh was read from.

Parameters:
archivethe archive
versionthe current version of this class

Definition at line 115 of file AbstractTetrahedralMesh.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::SetElementOwnerships ( ) [protected, virtual]

Sets the ownership of each element according to which nodes are owned by the process.

Information on node ownership comes from the distributed vector factory and an element is "owned" if one or more of its nodes are owned

Reimplemented from AbstractMesh< ELEMENT_DIM, SPACE_DIM >.

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

Definition at line 44 of file AbstractTetrahedralMesh.cpp.

References AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), and AbstractElement< ELEMENT_DIM, SPACE_DIM >::SetOwnership().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
virtual unsigned AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::SolveBoundaryElementMapping ( unsigned  index) const [private, pure virtual]

Pure virtual solve boundary element mapping method. For a boundary element with a given global index, get the local index used by this process. Overridden in TetrahedralMesh and DistributedTetrahedralMesh classes.

Parameters:
indexthe global index of the boundary element

Implemented in DistributedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >, TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >, DistributedTetrahedralMesh< DIM, DIM >, TetrahedralMesh< DIM, DIM >, and TetrahedralMesh< 3, 3 >.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
virtual unsigned AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::SolveElementMapping ( unsigned  index) const [private, pure virtual]

Pure virtual solve element mapping method. For an element with a given global index, get the local index used by this process. Overridden in TetrahedralMesh and DistributedTetrahedralMesh classes.

Parameters:
indexthe global index of the element

Implemented in DistributedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >, TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >, DistributedTetrahedralMesh< DIM, DIM >, TetrahedralMesh< DIM, DIM >, and TetrahedralMesh< 3, 3 >.


Friends And Related Function Documentation


Member Data Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
std::vector<BoundaryElement<ELEMENT_DIM-1, SPACE_DIM> *> AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mBoundaryElements [protected]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
bool AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::mMeshIsLinear [protected]

Most tet meshes are linear (set to true). Set to false in quadratics.

Definition at line 79 of file AbstractTetrahedralMesh.hpp.

Referenced by AbstractTetrahedralMesh< SPACE_DIM, SPACE_DIM >::load(), QuadraticMesh< DIM >::QuadraticMesh(), and AbstractTetrahedralMesh< SPACE_DIM, SPACE_DIM >::save().


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