ITK  6.0.0
Insight Toolkit
itkMesh.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright NumFOCUS
4  *
5  * Licensed under the Apache License, Version 2.0 (the "License");
6  * you may not use this file except in compliance with the License.
7  * You may obtain a copy of the License at
8  *
9  * https://www.apache.org/licenses/LICENSE-2.0.txt
10  *
11  * Unless required by applicable law or agreed to in writing, software
12  * distributed under the License is distributed on an "AS IS" BASIS,
13  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14  * See the License for the specific language governing permissions and
15  * limitations under the License.
16  *
17  *=========================================================================*/
18 /*=========================================================================
19  *
20  * Portions of this file are subject to the VTK Toolkit Version 3 copyright.
21  *
22  * Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
23  *
24  * For complete copyright, license and disclaimer of warranty information
25  * please refer to the NOTICE file at the top of the ITK source tree.
26  *
27  *=========================================================================*/
28 #ifndef itkMesh_h
29 #define itkMesh_h
30 
31 
32 #include "itkPointSet.h"
33 
34 #include "itkBoundingBox.h"
35 #include "itkCellInterface.h"
36 #include "itkMapContainer.h"
37 #include "itkCommonEnums.h"
38 #include "ITKMeshExport.h"
39 #include <vector>
40 #include <set>
41 #include "itkVectorContainer.h"
42 #include "itkVertexCell.h"
43 #include "itkLineCell.h"
44 #include "itkPolyLineCell.h"
45 #include "itkTriangleCell.h"
46 #include "itkQuadrilateralCell.h"
47 #include "itkPolygonCell.h"
48 #include "itkTetrahedronCell.h"
49 #include "itkHexahedronCell.h"
50 #include "itkQuadraticEdgeCell.h"
52 
53 
54 namespace itk
55 {
56 
123 template <typename TPixelType,
124  unsigned int VDimension = 3,
125  typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension>>
126 class ITK_TEMPLATE_EXPORT Mesh : public PointSet<TPixelType, VDimension, TMeshTraits>
127 {
128 public:
129  ITK_DISALLOW_COPY_AND_MOVE(Mesh);
130 
132  using Self = Mesh;
136 
137  using typename Superclass::RegionType;
138 
140  itkNewMacro(Self);
141 
143  itkOverrideGetNameOfClassMacro(Mesh);
144 
146  using MeshTraits = TMeshTraits;
150 
152  static constexpr unsigned int PointDimension = TMeshTraits::PointDimension;
153  static constexpr unsigned int MaxTopologicalDimension = TMeshTraits::MaxTopologicalDimension;
154 
155 #if !defined(ITK_LEGACY_REMOVE)
156  using CellsAllocationMethodType = MeshClassCellsAllocationMethodEnum;
158  // We need to expose the enum values at the class level
159  // for backwards compatibility
160  static constexpr CellsAllocationMethodType CellsAllocationMethodUndefined =
161  MeshClassCellsAllocationMethodEnum::CellsAllocationMethodUndefined;
162  static constexpr CellsAllocationMethodType CellsAllocatedAsStaticArray =
164  static constexpr CellsAllocationMethodType CellsAllocatedAsADynamicArray =
166  static constexpr CellsAllocationMethodType CellsAllocatedDynamicallyCellByCell =
168 #endif
169 
172 #ifndef ITK_FUTURE_LEGACY_REMOVE
173  using CoordRepType ITK_FUTURE_DEPRECATED(
174  "ITK 6 discourages using `CoordRepType`. Please use `CoordinateType` instead!") = CoordinateType;
175 #endif
189 
193 
196 
207 
209  using PointsContainerConstIterator = typename PointsContainer::ConstIterator;
210  using PointsContainerIterator = typename PointsContainer::Iterator;
211  using CellsContainerConstIterator = typename CellsContainer::ConstIterator;
212  using CellsContainerIterator = typename CellsContainer::Iterator;
213  using CellLinksContainerIterator = typename CellLinksContainer::ConstIterator;
214  using PointDataContainerIterator = typename PointDataContainer::ConstIterator;
215  using CellDataContainerIterator = typename CellDataContainer::ConstIterator;
216  using PointCellLinksContainerIterator = typename PointCellLinksContainer::const_iterator;
217 
220 
234 
237 
249  {
250  public:
253 
256  BoundaryAssignmentIdentifier() = default;
258  : m_CellId(cellId)
259  , m_FeatureId(featureId)
260  {}
261 
264 
267 
270  bool
271  operator<(const Self & r) const
272  {
273  return ((m_CellId < r.m_CellId) || ((m_CellId == r.m_CellId) && (m_FeatureId < r.m_FeatureId)));
274  }
275 
278  bool
279  operator==(const Self & r) const
280  {
281  return ((m_CellId == r.m_CellId) && (m_FeatureId == r.m_FeatureId));
282  }
283  }; // End Class: Mesh::BoundaryAssignmentIdentifier
295  using BoundaryAssignmentsContainerVector = std::vector<BoundaryAssignmentsContainerPointer>;
296 
297 protected:
300  CellsContainerPointer m_CellsContainer{};
301 
307  CellDataContainerPointer m_CellDataContainer{};
308 
312  mutable CellLinksContainerPointer m_CellLinksContainer{};
313 
323  BoundaryAssignmentsContainerVector m_BoundaryAssignmentsContainers{};
324 
325 public:
327  CellIdentifier
328  GetNumberOfCells() const;
329 
333  void
334  PassStructure(Self * inputMesh);
335 
339  void
340  Initialize() override;
341 
343  void
344  CopyInformation(const DataObject * data) override;
345 
346  void
347  Graft(const DataObject * data) override;
348 
351  const BoundingBoxType *
352  GetBoundingBox() const;
353 
358  void
359  SetCellLinks(CellLinksContainer *);
360 
362  CellLinksContainer *
363  GetCellLinks();
364 
366  const CellLinksContainer *
367  GetCellLinks() const;
368 
371  void
372  SetCells(CellsContainer *);
373 
379  virtual void
380  SetCellsArray(CellsVectorContainer *);
381 
386  virtual void
387  SetCellsArray(CellsVectorContainer *, int cellType);
388 
392  virtual CellsVectorContainer *
393  GetCellsArray();
394 
396  CellsContainer *
397  GetCells();
398 
400  const CellsContainer *
401  GetCells() const;
402 
407  void
408  SetCellData(CellDataContainer *);
409 
411  CellDataContainer *
412  GetCellData();
413 
415  const CellDataContainer *
416  GetCellData() const;
417 
421  void
422  DeleteUnusedCellData();
423 
424 #if !defined(ITK_WRAPPING_PARSER)
425 
433  void
434  SetBoundaryAssignments(int dimension, BoundaryAssignmentsContainer *);
435 
437  BoundaryAssignmentsContainerPointer
438  GetBoundaryAssignments(int dimension);
439 
441  const BoundaryAssignmentsContainerPointer
442  GetBoundaryAssignments(int dimension) const;
443 #endif
444 
450  void
451  SetCell(CellIdentifier, CellAutoPointer &);
452 
459  bool
460  GetCell(CellIdentifier, CellAutoPointer &) const;
461 
466  void SetCellData(CellIdentifier, CellPixelType);
467 
474  bool
475  GetCellData(CellIdentifier, CellPixelType *) const;
476 
489  void
490  SetBoundaryAssignment(int dimension,
491  CellIdentifier cellId,
492  CellFeatureIdentifier featureId,
493  CellIdentifier boundaryId);
494 
503  bool
504  GetBoundaryAssignment(int dimension,
505  CellIdentifier cellId,
506  CellFeatureIdentifier featureId,
507  CellIdentifier * boundaryId) const;
508 
512  bool
513  RemoveBoundaryAssignment(int dimension, CellIdentifier cellId, CellFeatureIdentifier featureId);
514 
518  CellFeatureCount
519  GetNumberOfCellBoundaryFeatures(int dimension, CellIdentifier) const;
520 
529  bool
530  GetCellBoundaryFeature(int dimension, CellIdentifier, CellFeatureIdentifier, CellAutoPointer &) const;
531 
536  CellIdentifier
537  GetCellBoundaryFeatureNeighbors(int dimension,
538  CellIdentifier,
539  CellFeatureIdentifier,
540  std::set<CellIdentifier> * cellSet);
541 
546  CellIdentifier
547  GetCellNeighbors(CellIdentifier cellId, std::set<CellIdentifier> * cellSet);
548 
556  bool
557  GetAssignedCellBoundaryIfOneExists(int dimension, CellIdentifier, CellFeatureIdentifier, CellAutoPointer &) const;
558 
561  void
562  BuildCellLinks() const;
563 
569  virtual void
570  Accept(CellMultiVisitorType * mv) const;
571 
576  itkSetMacro(CellsAllocationMethod, MeshClassCellsAllocationMethodEnum);
577  itkGetConstReferenceMacro(CellsAllocationMethod, MeshClassCellsAllocationMethodEnum);
580 protected:
582  Mesh();
583  ~Mesh() override;
584  void
585  PrintSelf(std::ostream & os, Indent indent) const override;
592  void
593  ReleaseCellsMemory();
594 
597  BoundingBoxPointer m_BoundingBox{};
598 
599 private:
600  MeshClassCellsAllocationMethodEnum m_CellsAllocationMethod{};
601 
603  void
604  CreateCell(int cellType, CellAutoPointer &);
605 }; // End Class: Mesh
606 
608 extern ITKMesh_EXPORT std::ostream &
609  operator<<(std::ostream & out, const MeshEnums::MeshClassCellsAllocationMethod value);
610 } // end namespace itk
611 
612 #ifndef ITK_MANUAL_INSTANTIATION
613 # ifndef ITK_WRAPPING_PARSER
614 # include "itkMesh.hxx"
615 # endif
616 #endif
617 
618 #endif
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::PointCellLinksContainer
typename MeshTraits::PointCellLinksContainer PointCellLinksContainer
Definition: itkMesh.h:185
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::CellTraits
typename MeshTraits::CellTraits CellTraits
Definition: itkMesh.h:183
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::CellMultiVisitorType
typename CellType::MultiVisitor CellMultiVisitorType
Definition: itkMesh.h:236
itkHexahedronCell.h
itk::DefaultDynamicMeshTraits::CellTraits
itkMakeCellTraitsMacro CellTraits
Definition: itkDefaultDynamicMeshTraits.h:113
itk::VertexCell
Represents a single vertex for a Mesh.
Definition: itkVertexCell.h:37
ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAddImageFilter.h:94
itk::QuadrilateralCell
Represents a quadrilateral for a Mesh.
Definition: itkQuadrilateralCell.h:37
itk::Mesh::BoundaryAssignmentIdentifier::operator==
bool operator==(const Self &r) const
Definition: itkMesh.h:279
itk::Mesh::BoundaryAssignmentIdentifier
Definition: itkMesh.h:248
itk::PolyLineCell
Represents a series of connected line segments for a Mesh.
Definition: itkPolyLineCell.h:36
itkQuadrilateralCell.h
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::CellDataContainerIterator
typename CellDataContainer::ConstIterator CellDataContainerIterator
Definition: itkMesh.h:215
itk::PointSet< TCoordinate, VDimension, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::PointDataContainerPointer
typename PointDataContainer::Pointer PointDataContainerPointer
Definition: itkPointSet.h:108
itk::operator<
bool operator<(const Index< VDimension > &one, const Index< VDimension > &two)
Definition: itkIndex.h:566
itk::detail::VectorContainer
Define a front-end to the STL "vector" container that conforms to the IndexedContainerInterface.
Definition: itkVectorContainer.h:51
itk::PointSet
A superclass of the N-dimensional mesh structure; supports point (geometric coordinate and attribute)...
Definition: itkPointSet.h:81
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::CellLinksContainer
typename MeshTraits::CellLinksContainer CellLinksContainer
Definition: itkMesh.h:186
itk::DefaultDynamicMeshTraits::CoordinateType
TCoordinate CoordinateType
Definition: itkDefaultDynamicMeshTraits.h:72
itk::PointSetBase< DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > ::PointsContainer >::CoordinateType
typename PointType::CoordinateType CoordinateType
Definition: itkPointSetBase.h:74
itk::DefaultDynamicMeshTraits::PointIdentifier
IdentifierType PointIdentifier
Definition: itkDefaultDynamicMeshTraits.h:85
itk::MeshEnums::MeshClassCellsAllocationMethod::CellsAllocatedDynamicallyCellByCell
itk::Mesh::cellOutputVectorContainer
CellsVectorContainerPointer cellOutputVectorContainer
Definition: itkMesh.h:302
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::BoundingBoxPointer
typename BoundingBoxType::Pointer BoundingBoxPointer
Definition: itkMesh.h:205
itk::DefaultDynamicMeshTraits::InterpolationWeightType
TInterpolationWeight InterpolationWeightType
Definition: itkDefaultDynamicMeshTraits.h:77
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::CellLinksContainerPointer
typename CellLinksContainer::Pointer CellLinksContainerPointer
Definition: itkMesh.h:201
itkLineCell.h
itk::DefaultDynamicMeshTraits::PointCellLinksContainer
std::set< CellIdentifier > PointCellLinksContainer
Definition: itkDefaultDynamicMeshTraits.h:126
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::CellsContainerPointer
typename CellsContainer::Pointer CellsContainerPointer
Definition: itkMesh.h:199
itk::operator<<
ITKCommon_EXPORT std::ostream & operator<<(std::ostream &out, typename AnatomicalOrientation::CoordinateEnum value)
itk::PolygonCell
Represents a polygon in a Mesh.
Definition: itkPolygonCell.h:49
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::CellsVectorContainerPointer
typename CellsVectorContainer::Pointer CellsVectorContainerPointer
Definition: itkMesh.h:192
itk::PointSetBase< DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > ::PointsContainer >::PointsContainerIterator
typename PointsContainer::Iterator PointsContainerIterator
Definition: itkPointSetBase.h:95
itk::PointSetBase< DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > ::PointsContainer >::PointsContainer
DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > ::PointsContainer PointsContainer
Definition: itkPointSetBase.h:80
itk::CellInterface::MultiVisitor
A visitor that can visit different cell types in a mesh. CellInterfaceVisitor instances can be regist...
Definition: itkCellInterface.h:160
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::CellDataContainerConstPointer
typename CellDataContainer::ConstPointer CellDataContainerConstPointer
Definition: itkMesh.h:204
itk::MapContainer
A wrapper of the STL "map" container.
Definition: itkMapContainer.h:45
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::CellAutoPointer
typename CellType::CellAutoPointer CellAutoPointer
Definition: itkMesh.h:233
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::CellLinksContainerConstPointer
typename CellLinksContainer::ConstPointer CellLinksContainerConstPointer
Definition: itkMesh.h:206
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::CellIdentifier
typename MeshTraits::CellIdentifier CellIdentifier
Definition: itkMesh.h:178
itk::QuadraticEdgeCell
Represents a second order line segment for a Mesh.
Definition: itkQuadraticEdgeCell.h:35
itk::TetrahedronCell
TetrahedronCell represents a tetrahedron for a Mesh.
Definition: itkTetrahedronCell.h:36
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::CellsContainer
typename MeshTraits::CellsContainer CellsContainer
Definition: itkMesh.h:184
itk::PointSet< TCoordinate, VDimension, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::PointIdentifier
typename MeshTraits::PointIdentifier PointIdentifier
Definition: itkPointSet.h:103
itk::DefaultDynamicMeshTraits::CellPixelType
TCellPixelType CellPixelType
Definition: itkDefaultDynamicMeshTraits.h:71
itk::DefaultDynamicMeshTraits::CellFeatureIdentifier
IdentifierType CellFeatureIdentifier
Definition: itkDefaultDynamicMeshTraits.h:94
itk::DefaultDynamicMeshTraits::CellIdentifier
IdentifierType CellIdentifier
Definition: itkDefaultDynamicMeshTraits.h:89
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itkPolyLineCell.h
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::CellsContainerConstPointer
typename CellsContainer::ConstPointer CellsContainerConstPointer
Definition: itkMesh.h:200
itk::CellInterface::CellAutoPointer
SelfAutoPointer CellAutoPointer
Definition: itkCellInterface.h:138
itk::PointSet< TCoordinate, VDimension, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::PixelType
typename MeshTraits::PixelType PixelType
Definition: itkPointSet.h:100
itk::HexahedronCell
Represents a hexahedron (cuboid) for a Mesh.
Definition: itkHexahedronCell.h:45
itk::TriangleCell
Definition: itkTriangleCell.h:46
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::CellsVectorContainer
typename itk::VectorContainer< IdentifierType > CellsVectorContainer
Definition: itkMesh.h:191
itk::Mesh::BoundaryAssignmentIdentifier::m_FeatureId
CellFeatureIdentifier m_FeatureId
Definition: itkMesh.h:266
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::CellsContainerConstIterator
typename CellsContainer::ConstIterator CellsContainerConstIterator
Definition: itkMesh.h:211
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::InterpolationWeightType
typename MeshTraits::InterpolationWeightType InterpolationWeightType
Definition: itkMesh.h:176
itk::LineCell
Represents a line segment for a Mesh.
Definition: itkLineCell.h:40
itk::PointSet< TCoordinate, VDimension, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::PointDataContainerIterator
typename PointDataContainer::ConstIterator PointDataContainerIterator
Definition: itkPointSet.h:112
itkVertexCell.h
itk::DefaultDynamicMeshTraits::PixelType
TPixelType PixelType
Definition: itkDefaultDynamicMeshTraits.h:70
itk::MeshEnums::MeshClassCellsAllocationMethod
MeshClassCellsAllocationMethod
Definition: itkCommonEnums.h:192
itk::CellInterface
An abstract interface for cells.
Definition: itkCellInterface.h:97
itkVectorContainer.h
itk::MeshEnums::MeshClassCellsAllocationMethod::CellsAllocatedAsADynamicArray
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::CellFeatureIdentifier
typename MeshTraits::CellFeatureIdentifier CellFeatureIdentifier
Definition: itkMesh.h:179
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::CellLinksContainerIterator
typename CellLinksContainer::ConstIterator CellLinksContainerIterator
Definition: itkMesh.h:213
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itkCommonEnums.h
itk::PointSetBase< DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > ::PointsContainer >::PointsContainerConstIterator
typename PointsContainer::ConstIterator PointsContainerConstIterator
Definition: itkPointSetBase.h:94
itk::Mesh
Implements the N-dimensional mesh structure.
Definition: itkMesh.h:126
itkPolygonCell.h
itkQuadraticEdgeCell.h
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::CellPixelType
typename MeshTraits::CellPixelType CellPixelType
Definition: itkMesh.h:148
itk::PointSet< TCoordinate, VDimension, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::PointDataContainer
typename MeshTraits::PointDataContainer PointDataContainer
Definition: itkPointSet.h:105
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::PointCellLinksContainerIterator
typename PointCellLinksContainer::const_iterator PointCellLinksContainerIterator
Definition: itkMesh.h:216
itk::Object
Base class for most ITK classes.
Definition: itkObject.h:61
itkCellInterface.h
itk::Point
A templated class holding a geometric point in n-Dimensional space.
Definition: itkPoint.h:53
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::BoundaryAssignmentsContainerVector
std::vector< BoundaryAssignmentsContainerPointer > BoundaryAssignmentsContainerVector
Definition: itkMesh.h:295
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::CellsContainerIterator
typename CellsContainer::Iterator CellsContainerIterator
Definition: itkMesh.h:212
itk::PointSet< TCoordinate, VDimension, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::PointType
typename MeshTraits::PointType PointType
Definition: itkPointSet.h:104
itkPointSet.h
itk::BoundingBox
Represent and compute information about bounding boxes.
Definition: itkBoundingBox.h:70
AddImageFilter
Definition: itkAddImageFilter.h:81
itkMapContainer.h
itkQuadraticTriangleCell.h
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::CellFeatureCount
CellFeatureIdentifier CellFeatureCount
Definition: itkMesh.h:219
itk::PointSetBase< DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > ::PointsContainer >::PointsContainerPointer
typename PointsContainer::Pointer PointsContainerPointer
Definition: itkPointSetBase.h:90
itk::Mesh::BoundaryAssignmentIdentifier::m_CellId
CellIdentifier m_CellId
Definition: itkMesh.h:263
itk::QuadraticTriangleCell
Represents a second order triangular patch for a Mesh.
Definition: itkQuadraticTriangleCell.h:36
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::CellDataContainerPointer
typename CellDataContainer::Pointer CellDataContainerPointer
Definition: itkMesh.h:203
itkTriangleCell.h
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::BoundaryAssignmentsContainerPointer
typename BoundaryAssignmentsContainer::Pointer BoundaryAssignmentsContainerPointer
Definition: itkMesh.h:294
itk::DefaultDynamicMeshTraits
A simple structure that holds type information for a mesh and its cells.
Definition: itkDefaultDynamicMeshTraits.h:63
itk::MeshEnums::MeshClassCellsAllocationMethod::CellsAllocatedAsStaticArray
itk::Mesh::BoundaryAssignmentIdentifier::BoundaryAssignmentIdentifier
BoundaryAssignmentIdentifier(CellIdentifier cellId, CellFeatureIdentifier featureId)
Definition: itkMesh.h:257
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::PointHashType
typename MeshTraits::PointHashType PointHashType
Definition: itkMesh.h:180
itkTetrahedronCell.h
itk::Mesh< TCoordinate, 2, DefaultDynamicMeshTraits< TCoordinate, 2, 2, TCoordinate > >::CellDataContainer
typename MeshTraits::CellDataContainer CellDataContainer
Definition: itkMesh.h:188
itkBoundingBox.h
itk::DataObject
Base class for all data objects in ITK.
Definition: itkDataObject.h:293