ITK  6.0.0
Insight Toolkit
Public Types | Public Member Functions | Static Public Member Functions | Protected Attributes | List of all members
itk::SimplexMesh< TPixelType, VDimension, TMeshTraits > Class Template Reference

#include <itkSimplexMesh.h>

Detailed Description

template<typename TPixelType, unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
class itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >

The class represents a 2-simplex mesh.

A simplex mesh can be used for deformable model segmentation of 3D image data. To create a simplex mesh one needs a triangle mesh, which can be converted to using the class itkTriangleMeshToSimplexMeshFilter. The back filtering (from simplex to triangle mesh)is done through a itkSimplexMeshToTriangleMeshFilter.

Author
Thomas Boettger. Division Medical and Biological Informatics, German Cancer Research Center, Heidelberg.

Definition at line 47 of file itkSimplexMesh.h.

+ Inheritance diagram for itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >:
+ Collaboration diagram for itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >:

Public Types

using CellAutoPointer = typename CellType::CellAutoPointer
 
using ConstPointer = SmartPointer< const Self >
 
using CovariantVectorType = CovariantVector< typename VectorType::ValueType, 3 >
 
using GeometryMapConstIterator = typename GeometryMapType::ConstIterator
 
using GeometryMapIterator = typename GeometryMapType::Iterator
 
using GeometryMapPointer = typename GeometryMapType::Pointer
 
using GeometryMapType = itk::MapContainer< SizeValueType, SimplexMeshGeometry * >
 
using IndexArray = typename SimplexMeshGeometry::IndexArray
 
using LineType = itk::LineCell< CellType >
 
using MeshClassCellsAllocationMethodEnum = itk::MeshEnums::MeshClassCellsAllocationMethod
 
using MeshTraits = TMeshTraits
 
using NeighborListType = std::vector< SizeValueType >
 
using NeighborSetIterator = typename NeighborSetType::iterator
 
using NeighborSetType = std::set< SizeValueType >
 
using PixelType = typename MeshTraits::PixelType
 
using Pointer = SmartPointer< Self >
 
using PointIdentifier = typename TMeshTraits::PointIdentifier
 
using PointsContainer = typename MeshTraits::PointsContainer
 
using PointsContainerIterator = typename Superclass::PointsContainer::Iterator
 
using PointType = typename TMeshTraits::PointType
 
using Self = SimplexMesh
 
using Superclass = Mesh< TPixelType, VDimension, TMeshTraits >
 
using VectorType = typename PointType::VectorType
 
- Public Types inherited from itk::Mesh< TPixelType, 3, DefaultStaticMeshTraits< TPixelType, 3, 3, TPixelType, TPixelType, TPixelType > >
using BoundaryAssignmentsContainer = MapContainer< BoundaryAssignmentIdentifier, CellIdentifier >
 
using BoundaryAssignmentsContainerPointer = typename BoundaryAssignmentsContainer::Pointer
 
using BoundaryAssignmentsContainerVector = std::vector< BoundaryAssignmentsContainerPointer >
 
using BoundingBoxPointer = typename BoundingBoxType::Pointer
 
using BoundingBoxType = BoundingBox< PointIdentifier, Self::PointDimension, CoordRepType, PointsContainer >
 
using CellAutoPointer = typename CellType::CellAutoPointer
 
using CellDataContainer = typename MeshTraits::CellDataContainer
 
using CellDataContainerConstPointer = typename CellDataContainer::ConstPointer
 
using CellDataContainerIterator = typename CellDataContainer::ConstIterator
 
using CellDataContainerPointer = typename CellDataContainer::Pointer
 
using CellFeatureCount = CellFeatureIdentifier
 
using CellFeatureIdentifier = typename MeshTraits::CellFeatureIdentifier
 
using CellIdentifier = typename MeshTraits::CellIdentifier
 
using CellLinksContainer = typename MeshTraits::CellLinksContainer
 
using CellLinksContainerConstPointer = typename CellLinksContainer::ConstPointer
 
using CellLinksContainerIterator = typename CellLinksContainer::ConstIterator
 
using CellLinksContainerPointer = typename CellLinksContainer::Pointer
 
using CellMultiVisitorType = typename CellType::MultiVisitor
 
using CellPixelType = typename MeshTraits::CellPixelType
 
using CellsContainer = typename MeshTraits::CellsContainer
 
using CellsContainerConstIterator = typename CellsContainer::ConstIterator
 
using CellsContainerConstPointer = typename CellsContainer::ConstPointer
 
using CellsContainerIterator = typename CellsContainer::Iterator
 
using CellsContainerPointer = typename CellsContainer::Pointer
 
using CellsVectorContainer = typename itk::VectorContainer< IdentifierType >
 
using CellsVectorContainerPointer = typename CellsVectorContainer::Pointer
 
using CellTraits = typename MeshTraits::CellTraits
 
using CellType = CellInterface< CellPixelType, CellTraits >
 
using ConstPointer = SmartPointer< const Self >
 
using CoordRepType = typename MeshTraits::CoordRepType
 
using InterpolationWeightType = typename MeshTraits::InterpolationWeightType
 
using MeshClassCellsAllocationMethodEnum = MeshEnums::MeshClassCellsAllocationMethod
 
using MeshTraits = DefaultStaticMeshTraits< TPixelType, 3, 3, TPixelType, TPixelType, TPixelType >
 
using OutputHexahedronCellType = itk::HexahedronCell< CellType >
 
using OutputLineCellType = itk::LineCell< CellType >
 
using OutputPolygonCellType = itk::PolygonCell< CellType >
 
using OutputPolyLineCellType = itk::PolyLineCell< CellType >
 
using OutputQuadraticEdgeCellType = itk::QuadraticEdgeCell< CellType >
 
using OutputQuadraticTriangleCellType = itk::QuadraticTriangleCell< CellType >
 
using OutputQuadrilateralCellType = itk::QuadrilateralCell< CellType >
 
using OutputTetrahedronCellType = itk::TetrahedronCell< CellType >
 
using OutputTriangleCellType = itk::TriangleCell< CellType >
 
using OutputVertexCellType = itk::VertexCell< CellType >
 
using PixelType = typename MeshTraits::PixelType
 
using PointCellLinksContainer = typename MeshTraits::PointCellLinksContainer
 
using PointCellLinksContainerIterator = typename PointCellLinksContainer::const_iterator
 
using PointDataContainer = typename MeshTraits::PointDataContainer
 
using PointDataContainerIterator = typename PointDataContainer::ConstIterator
 
using PointDataContainerPointer = typename PointDataContainer::Pointer
 
using Pointer = SmartPointer< Self >
 
using PointHashType = typename MeshTraits::PointHashType
 
using PointIdentifier = typename MeshTraits::PointIdentifier
 
using PointsContainer = typename MeshTraits::PointsContainer
 
using PointsContainerConstIterator = typename PointsContainer::ConstIterator
 
using PointsContainerIterator = typename PointsContainer::Iterator
 
using PointsContainerPointer = typename PointsContainer::Pointer
 
using PointType = typename MeshTraits::PointType
 
using Self = Mesh
 
using Superclass = PointSet< TPixelType, VDimension, DefaultStaticMeshTraits< TPixelType, 3, 3, TPixelType, TPixelType, TPixelType > >
 
- Public Types inherited from itk::PointSet< TPixelType, VDimension, DefaultStaticMeshTraits< TPixelType, 3, 3, TPixelType, TPixelType, TPixelType > >
using ConstPointer = SmartPointer< const Self >
 
using MeshTraits = DefaultStaticMeshTraits< TPixelType, 3, 3, TPixelType, TPixelType, TPixelType >
 
using PixelType = typename MeshTraits::PixelType
 
using PointDataContainer = typename MeshTraits::PointDataContainer
 
using PointDataContainerConstPointer = typename PointDataContainer::ConstPointer
 
using PointDataContainerIterator = typename PointDataContainer::ConstIterator
 
using PointDataContainerPointer = typename PointDataContainer::Pointer
 
using Pointer = SmartPointer< Self >
 
using PointIdentifier = typename MeshTraits::PointIdentifier
 
using PointType = typename MeshTraits::PointType
 
using Self = PointSet
 
using Superclass = PointSetBase< typename TMeshTraits::PointsContainer >
 
- Public Types inherited from itk::PointSetBase< TMeshTraits::PointsContainer >
using ConstPointer = SmartPointer< const Self >
 
using CoordRepType = typename PointType::CoordRepType
 
using Pointer = SmartPointer< Self >
 
using PointIdentifier = typename TPointsContainer::ElementIdentifier
 
using PointsContainer = TMeshTraits::PointsContainer
 
using PointsContainerConstIterator = typename PointsContainer::ConstIterator
 
using PointsContainerConstPointer = typename PointsContainer::ConstPointer
 
using PointsContainerIterator = typename PointsContainer::Iterator
 
using PointsContainerPointer = typename PointsContainer::Pointer
 
using PointsVectorContainer = typename itk::VectorContainer< PointIdentifier, CoordRepType >
 
using PointsVectorContainerPointer = typename PointsVectorContainer::Pointer
 
using PointType = typename TPointsContainer::Element
 
using RegionType = long
 
using Self = PointSetBase
 
using Superclass = DataObject
 
- Public Types inherited from itk::DataObject
using ConstPointer = SmartPointer< const Self >
 
using DataObjectIdentifierType = std::string
 
using DataObjectPointerArraySizeType = std::vector< Pointer >::size_type
 
using Pointer = SmartPointer< Self >
 
using Self = DataObject
 
using Superclass = Object
 
- Public Types inherited from itk::Object
using ConstPointer = SmartPointer< const Self >
 
using Pointer = SmartPointer< Self >
 
using Self = Object
 
using Superclass = LightObject
 
- Public Types inherited from itk::LightObject
using ConstPointer = SmartPointer< const Self >
 
using Pointer = SmartPointer< Self >
 
using Self = LightObject
 

Public Member Functions

CellIdentifier AddEdge (PointIdentifier startPointId, PointIdentifier endPointId)
 
CellIdentifier AddFace (CellAutoPointer &cellPointer)
 
void AddNeighbor (PointIdentifier pointIdx, PointIdentifier neighborIdx)
 
CovariantVectorType ComputeNormal (PointIdentifier idx) const
 
void CopyInformation (const DataObject *data) override
 
PointType GetBarycentricCoordinates (PointIdentifier idx) const
 
double GetDistance (PointIdentifier idx) const
 
virtual const GeometryMapPointerGetGeometryData () const
 
virtual CellIdentifier GetLastCellId () const
 
double GetMeanCurvature (PointIdentifier idx) const
 
const char * GetNameOfClass () const override
 
IndexArray GetNeighbors (PointIdentifier idx) const
 
NeighborListTypeGetNeighbors (PointIdentifier idx, unsigned int radius, NeighborListType *list=nullptr) const
 
double GetPhi (PointIdentifier idx) const
 
double GetRadius (PointIdentifier idx) const
 
PointType GetReferenceMetrics (PointIdentifier idx) const
 
CellIdentifier ReplaceFace (CellIdentifier replaceIndex, CellAutoPointer &cellPointer)
 
void ReplaceNeighbor (PointIdentifier pointIdx, PointIdentifier oldIdx, PointIdentifier newIdx)
 
void SetBarycentricCoordinates (PointIdentifier idx, PointType value)
 
void SetDistance (PointIdentifier idx, double value)
 
virtual void SetGeometryData (GeometryMapPointer _arg)
 
void SetGeometryData (PointIdentifier pointId, SimplexMeshGeometry *)
 
virtual void SetLastCellId (CellIdentifier _arg)
 
void SetMeanCurvature (PointIdentifier idx, double value)
 
void SetPhi (PointIdentifier idx, double value)
 
void SetRadius (PointIdentifier idx, double value)
 
void SetReferenceMetrics (PointIdentifier idx, PointType value)
 
void SwapNeighbors (PointIdentifier pointIdx, PointIdentifier firstIdx, PointIdentifier secondIdx)
 
- Public Member Functions inherited from itk::Mesh< TPixelType, 3, DefaultStaticMeshTraits< TPixelType, 3, 3, TPixelType, TPixelType, TPixelType > >
virtual void Accept (CellMultiVisitorType *mv) const
 
void BuildCellLinks () const
 
void CopyInformation (const DataObject *data) override
 
void DeleteUnusedCellData ()
 
bool GetAssignedCellBoundaryIfOneExists (int dimension, CellIdentifier, CellFeatureIdentifier, CellAutoPointer &) const
 
bool GetBoundaryAssignment (int dimension, CellIdentifier cellId, CellFeatureIdentifier featureId, CellIdentifier *boundaryId) const
 
BoundaryAssignmentsContainerPointer GetBoundaryAssignments (int dimension)
 
const BoundaryAssignmentsContainerPointer GetBoundaryAssignments (int dimension) const
 
const BoundingBoxTypeGetBoundingBox () const
 
bool GetCell (CellIdentifier, CellAutoPointer &) const
 
bool GetCellBoundaryFeature (int dimension, CellIdentifier, CellFeatureIdentifier, CellAutoPointer &) const
 
CellIdentifier GetCellBoundaryFeatureNeighbors (int dimension, CellIdentifier, CellFeatureIdentifier, std::set< CellIdentifier > *cellSet)
 
CellDataContainerGetCellData ()
 
const CellDataContainerGetCellData () const
 
bool GetCellData (CellIdentifier, CellPixelType *) const
 
CellLinksContainerGetCellLinks ()
 
const CellLinksContainerGetCellLinks () const
 
CellIdentifier GetCellNeighbors (CellIdentifier cellId, std::set< CellIdentifier > *cellSet)
 
CellsContainerGetCells ()
 
const CellsContainerGetCells () const
 
virtual CellsVectorContainerGetCellsArray ()
 
const char * GetNameOfClass () const override
 
CellFeatureCount GetNumberOfCellBoundaryFeatures (int dimension, CellIdentifier) const
 
CellIdentifier GetNumberOfCells () const
 
void Graft (const DataObject *data) override
 
void Initialize () override
 
void PassStructure (Self *inputMesh)
 
bool RemoveBoundaryAssignment (int dimension, CellIdentifier cellId, CellFeatureIdentifier featureId)
 
void SetBoundaryAssignment (int dimension, CellIdentifier cellId, CellFeatureIdentifier featureId, CellIdentifier boundaryId)
 
void SetBoundaryAssignments (int dimension, BoundaryAssignmentsContainer *)
 
void SetCell (CellIdentifier, CellAutoPointer &)
 
void SetCellData (CellDataContainer *)
 
void SetCellData (CellIdentifier, CellPixelType)
 
void SetCellLinks (CellLinksContainer *)
 
void SetCells (CellsContainer *)
 
virtual void SetCellsArray (CellsVectorContainer *)
 
virtual void SetCellsArray (CellsVectorContainer *, int cellType)
 
virtual void SetCellsAllocationMethod (MeshClassCellsAllocationMethodEnum _arg)
 
virtual const MeshClassCellsAllocationMethodEnumGetCellsAllocationMethod () const
 
- Public Member Functions inherited from itk::PointSet< TPixelType, VDimension, DefaultStaticMeshTraits< TPixelType, 3, 3, TPixelType, TPixelType, TPixelType > >
const char * GetNameOfClass () const override
 
PointDataContainerGetPointData ()
 
const PointDataContainerGetPointData () const
 
bool GetPointData (PointIdentifier, PixelType *) const
 
void Graft (const DataObject *data) override
 
void Initialize () override
 
void SetPointData (PointDataContainer *)
 
void SetPointData (PointIdentifier, PixelType)
 
- Public Member Functions inherited from itk::PointSetBase< TMeshTraits::PointsContainer >
void CopyInformation (const DataObject *data) override
 
virtual RegionType GetBufferedRegion () const
 
virtual RegionType GetMaximumNumberOfRegions () const
 
const char * GetNameOfClass () const override
 
PointIdentifier GetNumberOfPoints () const
 
PointType GetPoint (PointIdentifier) const
 
bool GetPoint (PointIdentifier, PointType *) const
 
PointsContainerGetPoints ()
 
const PointsContainerGetPoints () const
 
virtual RegionType GetRequestedRegion () const
 
void Initialize () override
 
void PassStructure (Self *inputPointSet)
 
bool RequestedRegionIsOutsideOfTheBufferedRegion () override
 
virtual void SetBufferedRegion (const RegionType &region)
 
void SetPoint (PointIdentifier, PointType)
 
void SetPoints (PointsContainer *)
 
void SetPoints (PointsVectorContainer *)
 
void SetPointsByCoordinates (const std::vector< CoordRepType > &coordinates)
 
void SetRequestedRegion (const DataObject *data) override
 
virtual void SetRequestedRegion (const RegionType &region)
 
void SetRequestedRegionToLargestPossibleRegion () override
 
void UpdateOutputInformation () override
 
bool VerifyRequestedRegion () override
 
- Public Member Functions inherited from itk::DataObject
virtual void CopyInformation (const DataObject *)
 
virtual void DataHasBeenGenerated ()
 
void DisconnectPipeline ()
 
bool GetDataReleased () const
 
const char * GetNameOfClass () const override
 
virtual const bool & GetReleaseDataFlag () const
 
SmartPointer< ProcessObjectGetSource () const
 
DataObjectPointerArraySizeType GetSourceOutputIndex () const
 
const DataObjectIdentifierTypeGetSourceOutputName () const
 
virtual ModifiedTimeType GetUpdateMTime () const
 
virtual void Graft (const DataObject *)
 
virtual void Initialize ()
 
virtual void PrepareForNewData ()
 
virtual void PropagateRequestedRegion ()
 
void ReleaseData ()
 
virtual void ReleaseDataFlagOn ()
 
virtual bool RequestedRegionIsOutsideOfTheBufferedRegion ()
 
virtual void ResetPipeline ()
 
void SetReleaseDataFlag (bool flag)
 
virtual void SetRequestedRegion (const DataObject *)
 
virtual void SetRequestedRegionToLargestPossibleRegion ()
 
bool ShouldIReleaseData () const
 
virtual void Update ()
 
virtual void UpdateOutputData ()
 
virtual void UpdateOutputInformation ()
 
void UpdateSource () const
 
virtual bool VerifyRequestedRegion ()
 
void SetPipelineMTime (ModifiedTimeType time)
 
virtual const ModifiedTimeTypeGetPipelineMTime () const
 
virtual void SetRealTimeStamp (RealTimeStamp _arg)
 
virtual const RealTimeStampGetRealTimeStamp () const
 
- Public Member Functions inherited from itk::Object
unsigned long AddObserver (const EventObject &event, Command *cmd) const
 
unsigned long AddObserver (const EventObject &event, std::function< void(const EventObject &)> function) const
 
LightObject::Pointer CreateAnother () const override
 
virtual void DebugOff () const
 
virtual void DebugOn () const
 
CommandGetCommand (unsigned long tag)
 
bool GetDebug () const
 
MetaDataDictionaryGetMetaDataDictionary ()
 
const MetaDataDictionaryGetMetaDataDictionary () const
 
virtual ModifiedTimeType GetMTime () const
 
const char * GetNameOfClass () const override
 
virtual const TimeStampGetTimeStamp () const
 
bool HasObserver (const EventObject &event) const
 
void InvokeEvent (const EventObject &)
 
void InvokeEvent (const EventObject &) const
 
virtual void Modified () const
 
void Register () const override
 
void RemoveAllObservers ()
 
void RemoveObserver (unsigned long tag) const
 
void SetDebug (bool debugFlag) const
 
void SetReferenceCount (int) override
 
void UnRegister () const noexcept override
 
void SetMetaDataDictionary (const MetaDataDictionary &rhs)
 
void SetMetaDataDictionary (MetaDataDictionary &&rrhs)
 
virtual void SetObjectName (std::string _arg)
 
virtual const std::string & GetObjectName () const
 
- Public Member Functions inherited from itk::LightObject
Pointer Clone () const
 
virtual Pointer CreateAnother () const
 
virtual void Delete ()
 
virtual const char * GetNameOfClass () const
 
virtual int GetReferenceCount () const
 
void Print (std::ostream &os, Indent indent=0) const
 
virtual void Register () const
 
virtual void SetReferenceCount (int)
 
virtual void UnRegister () const noexcept
 

Static Public Member Functions

static Pointer New ()
 
- Static Public Member Functions inherited from itk::Mesh< TPixelType, 3, DefaultStaticMeshTraits< TPixelType, 3, 3, TPixelType, TPixelType, TPixelType > >
static Pointer New ()
 
- Static Public Member Functions inherited from itk::PointSet< TPixelType, VDimension, DefaultStaticMeshTraits< TPixelType, 3, 3, TPixelType, TPixelType, TPixelType > >
static Pointer New ()
 
- Static Public Member Functions inherited from itk::DataObject
static bool GetGlobalReleaseDataFlag ()
 
static void GlobalReleaseDataFlagOff ()
 
static void GlobalReleaseDataFlagOn ()
 
static Pointer New ()
 
static void SetGlobalReleaseDataFlag (bool val)
 
- Static Public Member Functions inherited from itk::Object
static bool GetGlobalWarningDisplay ()
 
static void GlobalWarningDisplayOff ()
 
static void GlobalWarningDisplayOn ()
 
static Pointer New ()
 
static void SetGlobalWarningDisplay (bool val)
 
- Static Public Member Functions inherited from itk::LightObject
static void BreakOnError ()
 
static Pointer New ()
 

Protected Member Functions

 SimplexMesh ()
 
 ~SimplexMesh () override
 
void PrintSelf (std::ostream &os, Indent indent) const override
 
- Protected Member Functions inherited from itk::Mesh< TPixelType, 3, DefaultStaticMeshTraits< TPixelType, 3, 3, TPixelType, TPixelType, TPixelType > >
void ReleaseCellsMemory ()
 
 Mesh ()
 
 ~Mesh () override
 
void PrintSelf (std::ostream &os, Indent indent) const override
 
- Protected Member Functions inherited from itk::PointSet< TPixelType, VDimension, DefaultStaticMeshTraits< TPixelType, 3, 3, TPixelType, TPixelType, TPixelType > >
 PointSet ()=default
 
 ~PointSet () override=default
 
- Protected Member Functions inherited from itk::PointSetBase< TMeshTraits::PointsContainer >
 PointSetBase ()=default
 
void PrintSelf (std::ostream &os, Indent indent) const override
 
 ~PointSetBase () override=0
 
- Protected Member Functions inherited from itk::DataObject
 DataObject ()
 
void PrintSelf (std::ostream &os, Indent indent) const override
 
virtual void PropagateResetPipeline ()
 
 ~DataObject () override
 
- Protected Member Functions inherited from itk::Object
 Object ()
 
bool PrintObservers (std::ostream &os, Indent indent) const
 
void PrintSelf (std::ostream &os, Indent indent) const override
 
virtual void SetTimeStamp (const TimeStamp &timeStamp)
 
 ~Object () override
 
- Protected Member Functions inherited from itk::LightObject
virtual LightObject::Pointer InternalClone () const
 
 LightObject ()
 
virtual void PrintHeader (std::ostream &os, Indent indent) const
 
virtual void PrintSelf (std::ostream &os, Indent indent) const
 
virtual void PrintTrailer (std::ostream &os, Indent indent) const
 
virtual ~LightObject ()
 

Protected Attributes

GeometryMapPointer m_GeometryData {}
 
CellIdentifier m_LastCellId {}
 
- Protected Attributes inherited from itk::Mesh< TPixelType, 3, DefaultStaticMeshTraits< TPixelType, 3, 3, TPixelType, TPixelType, TPixelType > >
CellsVectorContainerPointer cellOutputVectorContainer
 
BoundaryAssignmentsContainerVector m_BoundaryAssignmentsContainers
 
BoundingBoxPointer m_BoundingBox
 
CellDataContainerPointer m_CellDataContainer
 
CellLinksContainerPointer m_CellLinksContainer
 
CellsContainerPointer m_CellsContainer
 
- Protected Attributes inherited from itk::PointSet< TPixelType, VDimension, DefaultStaticMeshTraits< TPixelType, 3, 3, TPixelType, TPixelType, TPixelType > >
PointDataContainerPointer m_PointDataContainer
 
- Protected Attributes inherited from itk::PointSetBase< TMeshTraits::PointsContainer >
RegionType m_BufferedRegion
 
RegionType m_MaximumNumberOfRegions
 
RegionType m_NumberOfRegions
 
PointsContainerPointer m_PointsContainer
 
RegionType m_RequestedNumberOfRegions
 
RegionType m_RequestedRegion
 
- Protected Attributes inherited from itk::LightObject
std::atomic< int > m_ReferenceCount {}
 

Additional Inherited Members

- Static Public Attributes inherited from itk::Mesh< TPixelType, 3, DefaultStaticMeshTraits< TPixelType, 3, 3, TPixelType, TPixelType, TPixelType > >
static constexpr unsigned int MaxTopologicalDimension
 
static constexpr unsigned int PointDimension
 
- Static Public Attributes inherited from itk::PointSetBase< TMeshTraits::PointsContainer >
static constexpr unsigned int PointDimension
 

Member Typedef Documentation

◆ CellAutoPointer

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
using itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::CellAutoPointer = typename CellType::CellAutoPointer

Definition at line 92 of file itkSimplexMesh.h.

◆ ConstPointer

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
using itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::ConstPointer = SmartPointer<const Self>

Standard type alias.

Definition at line 62 of file itkSimplexMesh.h.

◆ CovariantVectorType

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
using itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::CovariantVectorType = CovariantVector<typename VectorType::ValueType, 3>

Definition at line 86 of file itkSimplexMesh.h.

◆ GeometryMapConstIterator

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
using itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::GeometryMapConstIterator = typename GeometryMapType::ConstIterator

Definition at line 106 of file itkSimplexMesh.h.

◆ GeometryMapIterator

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
using itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::GeometryMapIterator = typename GeometryMapType::Iterator

iterator definition for iterating over a geometry map

Definition at line 105 of file itkSimplexMesh.h.

◆ GeometryMapPointer

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
using itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::GeometryMapPointer = typename GeometryMapType::Pointer

smartpointer def for the geometry map

Definition at line 102 of file itkSimplexMesh.h.

◆ GeometryMapType

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
using itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::GeometryMapType = itk::MapContainer<SizeValueType, SimplexMeshGeometry *>

map containing a SimplexMeshGeometry data object for each mesh point

Definition at line 99 of file itkSimplexMesh.h.

◆ IndexArray

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
using itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::IndexArray = typename SimplexMeshGeometry::IndexArray

definition for array of indices.

Definition at line 65 of file itkSimplexMesh.h.

◆ LineType

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
using itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::LineType = itk::LineCell<CellType>

Definition at line 95 of file itkSimplexMesh.h.

◆ MeshClassCellsAllocationMethodEnum

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
using itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::MeshClassCellsAllocationMethodEnum = itk::MeshEnums::MeshClassCellsAllocationMethod

Definition at line 109 of file itkSimplexMesh.h.

◆ MeshTraits

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
using itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::MeshTraits = TMeshTraits

Hold on to the type information specified by the template parameters.

Definition at line 118 of file itkSimplexMesh.h.

◆ NeighborListType

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
using itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::NeighborListType = std::vector<SizeValueType>

Definition at line 74 of file itkSimplexMesh.h.

◆ NeighborSetIterator

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
using itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::NeighborSetIterator = typename NeighborSetType::iterator

Definition at line 71 of file itkSimplexMesh.h.

◆ NeighborSetType

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
using itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::NeighborSetType = std::set<SizeValueType>

definition for a set of neighbor indices

Definition at line 68 of file itkSimplexMesh.h.

◆ PixelType

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
using itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::PixelType = typename MeshTraits::PixelType

Definition at line 119 of file itkSimplexMesh.h.

◆ Pointer

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
using itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::Pointer = SmartPointer<Self>

Standard type alias.

Definition at line 59 of file itkSimplexMesh.h.

◆ PointIdentifier

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
using itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::PointIdentifier = typename TMeshTraits::PointIdentifier

Definition at line 80 of file itkSimplexMesh.h.

◆ PointsContainer

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
using itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::PointsContainer = typename MeshTraits::PointsContainer

Definition at line 120 of file itkSimplexMesh.h.

◆ PointsContainerIterator

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
using itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::PointsContainerIterator = typename Superclass::PointsContainer::Iterator

Definition at line 122 of file itkSimplexMesh.h.

◆ PointType

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
using itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::PointType = typename TMeshTraits::PointType

Definition at line 77 of file itkSimplexMesh.h.

◆ Self

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
using itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::Self = SimplexMesh

Standard type alias.

Definition at line 53 of file itkSimplexMesh.h.

◆ Superclass

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
using itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::Superclass = Mesh<TPixelType, VDimension, TMeshTraits>

Standard type alias.

Definition at line 56 of file itkSimplexMesh.h.

◆ VectorType

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
using itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::VectorType = typename PointType::VectorType

Definition at line 83 of file itkSimplexMesh.h.

Constructor & Destructor Documentation

◆ SimplexMesh()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::SimplexMesh ( )
protected

Constructor for use by New() method.

◆ ~SimplexMesh()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::~SimplexMesh ( )
overrideprotected

Constructor for use by New() method.

Member Function Documentation

◆ AddEdge()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
CellIdentifier itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::AddEdge ( PointIdentifier  startPointId,
PointIdentifier  endPointId 
)

Add a new edge to the simplex mesh by specifying the ids of the start and end point of the edge Note: This can destroy the simplex mesh structure! Better use the simplex mesh modification or creation filters

◆ AddFace()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
CellIdentifier itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::AddFace ( CellAutoPointer cellPointer)

Add a new simplex mesh cell to the mesh by passing an AutoPointer of a previously created simplex mesh cell

Note: This can destroy the simplex mesh structure! You should use the simplex mesh modification or creation filters.

◆ AddNeighbor()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
void itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::AddNeighbor ( PointIdentifier  pointIdx,
PointIdentifier  neighborIdx 
)

Add a neighbor to a point. Note: This can destroy the simplex mesh topology! Better use the simplex mesh creation filters.

◆ ComputeNormal()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
CovariantVectorType itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::ComputeNormal ( PointIdentifier  idx) const

compute the normal vector in the specified mesh point

◆ CopyInformation()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
void itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::CopyInformation ( const DataObject data)
overridevirtual

copy all necessary information from passed object to the mesh

Reimplemented from itk::DataObject.

◆ GetBarycentricCoordinates()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
PointType itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::GetBarycentricCoordinates ( PointIdentifier  idx) const

Set the barycentric coordinates for a specified point

◆ GetDistance()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
double itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::GetDistance ( PointIdentifier  idx) const

Get the distance to the foot point for the specified point

◆ GetGeometryData()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
virtual const GeometryMapPointer & itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::GetGeometryData ( ) const
virtual

returns the current map of geometrydata

◆ GetLastCellId()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
virtual CellIdentifier itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::GetLastCellId ( ) const
virtual

Set the id value valid for new cells

◆ GetMeanCurvature()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
double itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::GetMeanCurvature ( PointIdentifier  idx) const

Get the mean curvature for the specified point

◆ GetNameOfClass()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
const char * itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::GetNameOfClass ( ) const
overridevirtual
See also
LightObject::GetNameOfClass()

Reimplemented from itk::DataObject.

◆ GetNeighbors() [1/2]

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
IndexArray itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::GetNeighbors ( PointIdentifier  idx) const

Get the three direct neighbors of a point

◆ GetNeighbors() [2/2]

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
NeighborListType * itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::GetNeighbors ( PointIdentifier  idx,
unsigned int  radius,
NeighborListType list = nullptr 
) const

Get all neighbor points with a specified radius

◆ GetPhi()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
double itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::GetPhi ( PointIdentifier  idx) const

Get the simplex angle for the specified point

◆ GetRadius()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
double itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::GetRadius ( PointIdentifier  idx) const

Get the circum circles radius for the specified point

◆ GetReferenceMetrics()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
PointType itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::GetReferenceMetrics ( PointIdentifier  idx) const

Return the reference metrics for the specified point

◆ New()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
static Pointer itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::New ( )
static

Method for creation through the object factory.

◆ PrintSelf()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
void itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::PrintSelf ( std::ostream &  os,
Indent  indent 
) const
overrideprotectedvirtual

Constructor for use by New() method.

Reimplemented from itk::DataObject.

◆ ReplaceFace()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
CellIdentifier itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::ReplaceFace ( CellIdentifier  replaceIndex,
CellAutoPointer cellPointer 
)

Replaces the cell specified by replaceIndex with the new cell passed by its AutoPointer

◆ ReplaceNeighbor()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
void itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::ReplaceNeighbor ( PointIdentifier  pointIdx,
PointIdentifier  oldIdx,
PointIdentifier  newIdx 
)

Replace a neighbor of a specific point by a new one

◆ SetBarycentricCoordinates()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
void itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::SetBarycentricCoordinates ( PointIdentifier  idx,
PointType  value 
)

Set the geometry data for a specified point

◆ SetDistance()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
void itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::SetDistance ( PointIdentifier  idx,
double  value 
)

Set the distance to the foot point for the specified point

◆ SetGeometryData() [1/2]

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
virtual void itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::SetGeometryData ( GeometryMapPointer  _arg)
virtual

set the map of geometrydata to the new pointer

◆ SetGeometryData() [2/2]

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
void itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::SetGeometryData ( PointIdentifier  pointId,
SimplexMeshGeometry  
)

Set the geometry data for a specified point

◆ SetLastCellId()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
virtual void itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::SetLastCellId ( CellIdentifier  _arg)
virtual

Get the first free id for new cells

◆ SetMeanCurvature()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
void itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::SetMeanCurvature ( PointIdentifier  idx,
double  value 
)

Set the mean curvature for the specified point

◆ SetPhi()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
void itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::SetPhi ( PointIdentifier  idx,
double  value 
)

Set the simplex angle for the specified point

◆ SetRadius()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
void itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::SetRadius ( PointIdentifier  idx,
double  value 
)

Set the circum circles radius for the specified point

◆ SetReferenceMetrics()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
void itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::SetReferenceMetrics ( PointIdentifier  idx,
PointType  value 
)

Set the reference metrics for a specified point

◆ SwapNeighbors()

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
void itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::SwapNeighbors ( PointIdentifier  pointIdx,
PointIdentifier  firstIdx,
PointIdentifier  secondIdx 
)

Swap the order of two neighbors

Member Data Documentation

◆ m_GeometryData

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
GeometryMapPointer itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::m_GeometryData {}
protected

The map stores a SimplexMeshGeometry object for each mesh point

Definition at line 300 of file itkSimplexMesh.h.

◆ m_LastCellId

template<typename TPixelType , unsigned int VDimension = 3, typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension, TPixelType, TPixelType, TPixelType>>
CellIdentifier itk::SimplexMesh< TPixelType, VDimension, TMeshTraits >::m_LastCellId {}
protected

The last cell id is the index which is used for insertion of new cells. It increases during mesh creation. This is done because one cannot rely on the size of the map or the highest index when cells are removed.

Definition at line 308 of file itkSimplexMesh.h.


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