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"
47#include "itkPolygonCell.h"
48#include "itkTetrahedronCell.h"
49#include "itkHexahedronCell.h"
52
53
54namespace itk
55{
56
123template <typename TPixelType,
124 unsigned int VDimension = 3,
125 typename TMeshTraits = DefaultStaticMeshTraits<TPixelType, VDimension, VDimension>>
126class ITK_TEMPLATE_EXPORT Mesh : public PointSet<TPixelType, VDimension, TMeshTraits>
127{
128public:
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;
147 using PixelType = typename MeshTraits::PixelType;
148 using CellPixelType = typename MeshTraits::CellPixelType;
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
171 using CoordRepType = typename MeshTraits::CoordRepType;
172 using InterpolationWeightType = typename MeshTraits::InterpolationWeightType;
173 using PointIdentifier = typename MeshTraits::PointIdentifier;
174 using CellIdentifier = typename MeshTraits::CellIdentifier;
175 using CellFeatureIdentifier = typename MeshTraits::CellFeatureIdentifier;
176 using PointHashType = typename MeshTraits::PointHashType;
178 using PointsContainer = typename MeshTraits::PointsContainer;
179 using CellTraits = typename MeshTraits::CellTraits;
180 using CellsContainer = typename MeshTraits::CellsContainer;
181 using PointCellLinksContainer = typename MeshTraits::PointCellLinksContainer;
182 using CellLinksContainer = typename MeshTraits::CellLinksContainer;
183 using PointDataContainer = typename MeshTraits::PointDataContainer;
184 using CellDataContainer = typename MeshTraits::CellDataContainer;
185
189
192
203
205 using PointsContainerConstIterator = typename PointsContainer::ConstIterator;
206 using PointsContainerIterator = typename PointsContainer::Iterator;
207 using CellsContainerConstIterator = typename CellsContainer::ConstIterator;
208 using CellsContainerIterator = typename CellsContainer::Iterator;
209 using CellLinksContainerIterator = typename CellLinksContainer::ConstIterator;
210 using PointDataContainerIterator = typename PointDataContainer::ConstIterator;
211 using CellDataContainerIterator = typename CellDataContainer::ConstIterator;
212 using PointCellLinksContainerIterator = typename PointCellLinksContainer::const_iterator;
213
216
230
233
245 {
246 public:
249
254 : m_CellId(cellId)
255 , m_FeatureId(featureId)
256 {}
257
260
263
266 bool
267 operator<(const Self & r) const
268 {
269 return ((m_CellId < r.m_CellId) || ((m_CellId == r.m_CellId) && (m_FeatureId < r.m_FeatureId)));
270 }
271
274 bool
275 operator==(const Self & r) const
276 {
277 return ((m_CellId == r.m_CellId) && (m_FeatureId == r.m_FeatureId));
278 }
279 }; // End Class: Mesh::BoundaryAssignmentIdentifier
291 using BoundaryAssignmentsContainerVector = std::vector<BoundaryAssignmentsContainerPointer>;
292
293protected:
296 CellsContainerPointer m_CellsContainer{};
297
303 CellDataContainerPointer m_CellDataContainer{};
304
308 mutable CellLinksContainerPointer m_CellLinksContainer{};
309
319 BoundaryAssignmentsContainerVector m_BoundaryAssignmentsContainers{};
320
321public:
323 CellIdentifier
325
329 void
330 PassStructure(Self * inputMesh);
331
335 void
336 Initialize() override;
337
339 void
340 CopyInformation(const DataObject * data) override;
341
342 void
343 Graft(const DataObject * data) override;
344
347 const BoundingBoxType *
349
354 void
356
360
362 const CellLinksContainer *
364
367 void
369
375 virtual void
377
382 virtual void
384
388 virtual CellsVectorContainer *
390
394
396 const CellsContainer *
397 GetCells() const;
398
403 void
405
409
411 const CellDataContainer *
412 GetCellData() const;
413
417 void
419
420#if !defined(ITK_WRAPPING_PARSER)
429 void
431
435
438 GetBoundaryAssignments(int dimension) const;
439#endif
440
446 void
448
455 bool
457
463
470 bool
472
485 void
487 CellIdentifier cellId,
488 CellFeatureIdentifier featureId,
489 CellIdentifier boundaryId);
490
499 bool
501 CellIdentifier cellId,
502 CellFeatureIdentifier featureId,
503 CellIdentifier * boundaryId) const;
504
508 bool
510
516
525 bool
527
536 std::set<CellIdentifier> * cellSet);
537
543 GetCellNeighbors(CellIdentifier cellId, std::set<CellIdentifier> * cellSet);
544
552 bool
554
557 void
559
565 virtual void
567
572 itkSetMacro(CellsAllocationMethod, MeshClassCellsAllocationMethodEnum);
573 itkGetConstReferenceMacro(CellsAllocationMethod, MeshClassCellsAllocationMethodEnum);
576protected:
579 ~Mesh() override;
580 void
581 PrintSelf(std::ostream & os, Indent indent) const override;
588 void
590
593 BoundingBoxPointer m_BoundingBox{};
594
595private:
596 MeshClassCellsAllocationMethodEnum m_CellsAllocationMethod{};
597
599 void
600 CreateCell(int cellType, CellAutoPointer &);
601}; // End Class: Mesh
602
604extern ITKMesh_EXPORT std::ostream &
605 operator<<(std::ostream & out, const MeshEnums::MeshClassCellsAllocationMethod value);
606} // end namespace itk
607
608#ifndef ITK_MANUAL_INSTANTIATION
609# ifndef ITK_WRAPPING_PARSER
610# include "itkMesh.hxx"
611# endif
612#endif
613
614#endif
Represent and compute information about bounding boxes.
A visitor that can visit different cell types in a mesh. CellInterfaceVisitor instances can be regist...
An abstract interface for cells.
SelfAutoPointer CellAutoPointer
Base class for all data objects in ITK.
Represents a hexahedron (cuboid) for a Mesh.
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Represents a line segment for a Mesh.
Definition: itkLineCell.h:41
A wrapper of the STL "map" container.
BoundaryAssignmentIdentifier(CellIdentifier cellId, CellFeatureIdentifier featureId)
Definition: itkMesh.h:253
bool operator==(const Self &r) const
Definition: itkMesh.h:275
CellFeatureIdentifier m_FeatureId
Definition: itkMesh.h:262
Implements the N-dimensional mesh structure.
Definition: itkMesh.h:127
typename CellsContainer::ConstPointer CellsContainerConstPointer
Definition: itkMesh.h:196
CellDataContainer * GetCellData()
typename MeshTraits::CellDataContainer CellDataContainer
Definition: itkMesh.h:184
bool GetBoundaryAssignment(int dimension, CellIdentifier cellId, CellFeatureIdentifier featureId, CellIdentifier *boundaryId) const
typename PointsContainer::Pointer PointsContainerPointer
Definition: itkMesh.h:194
typename CellLinksContainer::Pointer CellLinksContainerPointer
Definition: itkMesh.h:197
typename MeshTraits::CellFeatureIdentifier CellFeatureIdentifier
Definition: itkMesh.h:175
typename CellDataContainer::ConstIterator CellDataContainerIterator
Definition: itkMesh.h:211
void PrintSelf(std::ostream &os, Indent indent) const override
typename MeshTraits::PointCellLinksContainer PointCellLinksContainer
Definition: itkMesh.h:181
CellsContainer * GetCells()
void ReleaseCellsMemory()
CellLinksContainer * GetCellLinks()
void SetBoundaryAssignment(int dimension, CellIdentifier cellId, CellFeatureIdentifier featureId, CellIdentifier boundaryId)
virtual void SetCellsArray(CellsVectorContainer *)
BoundaryAssignmentsContainerPointer GetBoundaryAssignments(int dimension)
CellFeatureIdentifier CellFeatureCount
Definition: itkMesh.h:215
typename BoundingBoxType::Pointer BoundingBoxPointer
Definition: itkMesh.h:201
bool GetCellData(CellIdentifier, CellPixelType *) const
virtual void SetCellsArray(CellsVectorContainer *, int cellType)
void SetCellData(CellIdentifier, CellPixelType)
void SetCellLinks(CellLinksContainer *)
typename PointsContainer::Iterator PointsContainerIterator
Definition: itkMesh.h:206
typename MeshTraits::PointHashType PointHashType
Definition: itkMesh.h:176
typename itk::VectorContainer< IdentifierType > CellsVectorContainer
Definition: itkMesh.h:187
typename CellDataContainer::ConstPointer CellDataContainerConstPointer
Definition: itkMesh.h:200
typename MeshTraits::PixelType PixelType
Definition: itkMesh.h:147
typename MeshTraits::PointIdentifier PointIdentifier
Definition: itkMesh.h:173
void BuildCellLinks() const
const BoundingBoxType * GetBoundingBox() const
typename CellsContainer::Pointer CellsContainerPointer
Definition: itkMesh.h:195
CellIdentifier GetNumberOfCells() const
typename CellDataContainer::Pointer CellDataContainerPointer
Definition: itkMesh.h:199
typename CellLinksContainer::ConstIterator CellLinksContainerIterator
Definition: itkMesh.h:209
typename CellsVectorContainer::Pointer CellsVectorContainerPointer
Definition: itkMesh.h:188
typename PointCellLinksContainer::const_iterator PointCellLinksContainerIterator
Definition: itkMesh.h:212
typename MeshTraits::CellIdentifier CellIdentifier
Definition: itkMesh.h:174
CellIdentifier GetCellNeighbors(CellIdentifier cellId, std::set< CellIdentifier > *cellSet)
void PassStructure(Self *inputMesh)
typename PointDataContainer::Pointer PointDataContainerPointer
Definition: itkMesh.h:198
typename MeshTraits::CellPixelType CellPixelType
Definition: itkMesh.h:148
CellIdentifier GetCellBoundaryFeatureNeighbors(int dimension, CellIdentifier, CellFeatureIdentifier, std::set< CellIdentifier > *cellSet)
TMeshTraits MeshTraits
Definition: itkMesh.h:146
bool GetCell(CellIdentifier, CellAutoPointer &) const
const BoundaryAssignmentsContainerPointer GetBoundaryAssignments(int dimension) const
typename CellsContainer::ConstIterator CellsContainerConstIterator
Definition: itkMesh.h:207
const CellDataContainer * GetCellData() const
const CellLinksContainer * GetCellLinks() const
typename MeshTraits::CellsContainer CellsContainer
Definition: itkMesh.h:180
void Initialize() override
void SetCells(CellsContainer *)
void CreateCell(int cellType, CellAutoPointer &)
virtual CellsVectorContainer * GetCellsArray()
bool RemoveBoundaryAssignment(int dimension, CellIdentifier cellId, CellFeatureIdentifier featureId)
void SetBoundaryAssignments(int dimension, BoundaryAssignmentsContainer *)
typename CellType::CellAutoPointer CellAutoPointer
Definition: itkMesh.h:229
typename CellType::MultiVisitor CellMultiVisitorType
Definition: itkMesh.h:232
void Graft(const DataObject *data) override
void SetCellData(CellDataContainer *)
virtual void Accept(CellMultiVisitorType *mv) const
bool GetAssignedCellBoundaryIfOneExists(int dimension, CellIdentifier, CellFeatureIdentifier, CellAutoPointer &) const
CellsVectorContainerPointer cellOutputVectorContainer
Definition: itkMesh.h:298
typename MeshTraits::PointsContainer PointsContainer
Definition: itkMesh.h:178
typename MeshTraits::CellLinksContainer CellLinksContainer
Definition: itkMesh.h:182
typename MeshTraits::InterpolationWeightType InterpolationWeightType
Definition: itkMesh.h:172
typename MeshTraits::CoordRepType CoordRepType
Definition: itkMesh.h:171
void DeleteUnusedCellData()
typename PointDataContainer::ConstIterator PointDataContainerIterator
Definition: itkMesh.h:210
typename CellsContainer::Iterator CellsContainerIterator
Definition: itkMesh.h:208
std::vector< BoundaryAssignmentsContainerPointer > BoundaryAssignmentsContainerVector
Definition: itkMesh.h:291
typename CellLinksContainer::ConstPointer CellLinksContainerConstPointer
Definition: itkMesh.h:202
typename BoundaryAssignmentsContainer::Pointer BoundaryAssignmentsContainerPointer
Definition: itkMesh.h:290
void CopyInformation(const DataObject *data) override
typename MeshTraits::CellTraits CellTraits
Definition: itkMesh.h:179
const CellsContainer * GetCells() const
bool GetCellBoundaryFeature(int dimension, CellIdentifier, CellFeatureIdentifier, CellAutoPointer &) const
~Mesh() override
void SetCell(CellIdentifier, CellAutoPointer &)
typename MeshTraits::PointDataContainer PointDataContainer
Definition: itkMesh.h:183
typename MeshTraits::PointType PointType
Definition: itkMesh.h:177
typename PointsContainer::ConstIterator PointsContainerConstIterator
Definition: itkMesh.h:205
CellFeatureCount GetNumberOfCellBoundaryFeatures(int dimension, CellIdentifier) const
Base class for most ITK classes.
Definition: itkObject.h:62
A superclass of the N-dimensional mesh structure; supports point (geometric coordinate and attribute)...
Definition: itkPointSet.h:82
Represents a series of connected line segments for a Mesh.
Represents a polygon in a Mesh.
Represents a second order line segment for a Mesh.
Represents a second order triangular patch for a Mesh.
Represents a quadrilateral for a Mesh.
TetrahedronCell represents a tetrahedron for a Mesh.
Represents a single vertex for a Mesh.
Definition: itkVertexCell.h:38
Define a front-end to the STL "vector" container that conforms to the IndexedContainerInterface.
SmartPointer< const Self > ConstPointer
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
ITKCommon_EXPORT std::ostream & operator<<(std::ostream &out, typename AnatomicalOrientation::CoordinateEnum value)
bool operator<(const Index< VDimension > &one, const Index< VDimension > &two)
Definition: itkIndex.h:557