ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkMesh.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
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  * http://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 <vector>
38 #include <set>
39 
40 namespace itk
41 {
42 
103 template<
104  typename TPixelType,
105  unsigned int VDimension = 3,
106  typename TMeshTraits = DefaultStaticMeshTraits< TPixelType, VDimension, VDimension >
107  >
108 class ITK_TEMPLATE_EXPORT Mesh:public PointSet< TPixelType, VDimension, TMeshTraits >
109 {
110 public:
111  ITK_DISALLOW_COPY_AND_ASSIGN(Mesh);
112 
114  using Self = Mesh;
118 
120 
122  itkNewMacro(Self);
123 
125  itkTypeMacro(Mesh, PointSet);
126 
128  using MeshTraits = TMeshTraits;
129  using PixelType = typename MeshTraits::PixelType;
130  using CellPixelType = typename MeshTraits::CellPixelType;
131 
133  static constexpr unsigned int PointDimension = TMeshTraits::PointDimension;
134  static constexpr unsigned int MaxTopologicalDimension = TMeshTraits::MaxTopologicalDimension;
135 
141  CellsAllocatedDynamicallyCellByCell } CellsAllocationMethodType;
142 
144  using CoordRepType = typename MeshTraits::CoordRepType;
145  using InterpolationWeightType = typename MeshTraits::InterpolationWeightType;
146  using PointIdentifier = typename MeshTraits::PointIdentifier;
147  using CellIdentifier = typename MeshTraits::CellIdentifier;
148  using CellFeatureIdentifier = typename MeshTraits::CellFeatureIdentifier;
149  using PointHashType = typename MeshTraits::PointHashType;
151  using PointsContainer = typename MeshTraits::PointsContainer;
152  using CellTraits = typename MeshTraits::CellTraits;
153  using CellsContainer = typename MeshTraits::CellsContainer;
154  using PointCellLinksContainer = typename MeshTraits::PointCellLinksContainer;
155  using CellLinksContainer = typename MeshTraits::CellLinksContainer;
156  using PointDataContainer = typename MeshTraits::PointDataContainer;
157  using CellDataContainer = typename MeshTraits::CellDataContainer;
158 
160  using BoundingBoxType = BoundingBox< PointIdentifier, Self::PointDimension,
162 
164  using PointsContainerPointer = typename PointsContainer::Pointer;
165  using CellsContainerPointer = typename CellsContainer::Pointer;
166  using CellsContainerConstPointer = typename CellsContainer::ConstPointer;
167  using CellLinksContainerPointer = typename CellLinksContainer::Pointer;
168  using PointDataContainerPointer = typename PointDataContainer::Pointer;
169  using CellDataContainerPointer = typename CellDataContainer::Pointer;
170  using CellDataContainerConstPointer = typename CellDataContainer::ConstPointer;
172  using CellLinksContainerConstPointer = typename CellLinksContainer::ConstPointer;
173 
175  using PointsContainerConstIterator = typename PointsContainer::ConstIterator;
176  using PointsContainerIterator = typename PointsContainer::Iterator;
177  using CellsContainerConstIterator = typename CellsContainer::ConstIterator;
178  using CellsContainerIterator = typename CellsContainer::Iterator;
179  using CellLinksContainerIterator = typename CellLinksContainer::ConstIterator;
180  using PointDataContainerIterator = typename PointDataContainer::ConstIterator;
181  using CellDataContainerIterator = typename CellDataContainer::ConstIterator;
182  using PointCellLinksContainerIterator = typename PointCellLinksContainer::const_iterator;
183 
186 
190 
193 
205  {
206 public:
209 
212  BoundaryAssignmentIdentifier() = default;
214  CellFeatureIdentifier featureId):
215  m_CellId(cellId), m_FeatureId(featureId) {}
216 
219 
222 
225  bool operator<(const Self & r) const
226  {
227  return ( ( m_CellId < r.m_CellId )
228  || ( ( m_CellId == r.m_CellId ) && ( m_FeatureId < r.m_FeatureId ) ) );
229  }
230 
233  bool operator==(const Self & r) const
234  {
235  return ( ( m_CellId == r.m_CellId ) && ( m_FeatureId == r.m_FeatureId ) );
236  }
237  }; // End Class: Mesh::BoundaryAssignmentIdentifier
239 
247  using BoundaryAssignmentsContainer =
250  using BoundaryAssignmentsContainerVector = std::vector<BoundaryAssignmentsContainerPointer>;
251 
252 protected:
253 
257 
263 
268 
279 
280 public:
282  CellIdentifier GetNumberOfCells() const;
283 
284  void PassStructure(Self *inputMesh);
285 
286  void Initialize() override;
287 
289  void CopyInformation(const DataObject *data) override;
290 
291  void Graft(const DataObject *data) override;
292 
295  const BoundingBoxType * GetBoundingBox() const;
296 
301  void SetCellLinks(CellLinksContainer *);
302 
303  CellLinksContainer * GetCellLinks();
304 
305  const CellLinksContainer * GetCellLinks() const;
306 
309  void SetCells(CellsContainer *);
310 
311  CellsContainer * GetCells();
312 
313  const CellsContainer * GetCells() const;
314 
319  void SetCellData(CellDataContainer *);
320 
321  CellDataContainer * GetCellData();
322 
323  const CellDataContainer * GetCellData() const;
324 
325 #if !defined( ITK_WRAPPING_PARSER )
326 
334  void SetBoundaryAssignments(int dimension,
336 
337 
338  BoundaryAssignmentsContainerPointer GetBoundaryAssignments(int dimension);
339 
340  const BoundaryAssignmentsContainerPointer GetBoundaryAssignments(
341  int dimension) const;
342 #endif
343 
349  void SetCell(CellIdentifier, CellAutoPointer &);
350  bool GetCell(CellIdentifier, CellAutoPointer &) const;
351 
354  void SetCellData(CellIdentifier, CellPixelType);
355  bool GetCellData(CellIdentifier, CellPixelType *) const;
357 
370  void SetBoundaryAssignment(int dimension, CellIdentifier cellId,
371  CellFeatureIdentifier featureId,
372  CellIdentifier boundaryId);
373 
382  bool GetBoundaryAssignment(int dimension, CellIdentifier cellId,
383  CellFeatureIdentifier featureId,
384  CellIdentifier *boundaryId) const;
385 
386  bool RemoveBoundaryAssignment(int dimension, CellIdentifier cellId,
387  CellFeatureIdentifier featureId);
388 
390  CellFeatureCount GetNumberOfCellBoundaryFeatures(int dimension,
391  CellIdentifier) const;
392 
395  bool GetCellBoundaryFeature(int dimension, CellIdentifier,
397 
402  CellIdentifier GetCellBoundaryFeatureNeighbors(
403  int dimension, CellIdentifier, CellFeatureIdentifier,
404  std::set< CellIdentifier > *cellSet);
405 
410  CellIdentifier GetCellNeighbors(CellIdentifier cellId,
411  std::set< CellIdentifier > *cellSet);
412 
420  bool GetAssignedCellBoundaryIfOneExists(int dimension, CellIdentifier,
422  CellAutoPointer &) const;
423 
426  void BuildCellLinks() const;
427 
431  virtual void Accept(CellMultiVisitorType *mv) const;
432 
437  itkSetMacro(CellsAllocationMethod, CellsAllocationMethodType);
438  itkGetConstReferenceMacro(CellsAllocationMethod, CellsAllocationMethodType);
440 
441 protected:
443  Mesh();
444  ~Mesh() override;
445  void PrintSelf(std::ostream & os, Indent indent) const override;
447 
451  void ReleaseCellsMemory();
452 
456 
457 private:
458  CellsAllocationMethodType m_CellsAllocationMethod;
459 }; // End Class: Mesh
460 } // end namespace itk
461 
462 #ifndef ITK_MANUAL_INSTANTIATION
463 #ifndef ITK_WRAPPING_PARSER
464 #include "itkMesh.hxx"
465 #endif
466 #endif
467 
468 #endif
CellsContainerPointer m_CellsContainer
Definition: itkMesh.h:256
typename BoundingBoxType::Pointer BoundingBoxPointer
Definition: itkMesh.h:171
typename MeshTraits::PointCellLinksContainer PointCellLinksContainer
Definition: itkMesh.h:154
std::vector< BoundaryAssignmentsContainerPointer > BoundaryAssignmentsContainerVector
Definition: itkMesh.h:250
typename MeshTraits::CellDataContainer CellDataContainer
Definition: itkMesh.h:157
bool operator<(const Self &r) const
Definition: itkMesh.h:225
CellFeatureIdentifier CellFeatureCount
Definition: itkMesh.h:185
typename PointsContainer::Pointer PointsContainerPointer
Definition: itkPointSet.h:116
typename MeshTraits::CellTraits CellTraits
Definition: itkMesh.h:152
typename MeshTraits::PixelType PixelType
Definition: itkPointSet.h:103
typename CellsContainer::ConstIterator CellsContainerConstIterator
Definition: itkMesh.h:177
typename CellDataContainer::ConstIterator CellDataContainerIterator
Definition: itkMesh.h:181
A wrapper of the STL &quot;map&quot; container.
SelfAutoPointer CellAutoPointer
typename CellLinksContainer::Pointer CellLinksContainerPointer
Definition: itkMesh.h:167
typename PointsContainer::ConstIterator PointsContainerConstIterator
Definition: itkPointSet.h:122
CellDataContainerPointer m_CellDataContainer
Definition: itkMesh.h:262
A visitor that can visit different cell types in a mesh. CellInterfaceVisitor instances can be regist...
typename CellLinksContainer::ConstPointer CellLinksContainerConstPointer
Definition: itkMesh.h:172
typename BoundaryAssignmentsContainer::Pointer BoundaryAssignmentsContainerPointer
Definition: itkMesh.h:249
Implements the N-dimensional mesh structure.
Definition: itkMesh.h:108
An abstract interface for cells.
typename MeshTraits::CellLinksContainer CellLinksContainer
Definition: itkMesh.h:155
typename MeshTraits::PointType PointType
Definition: itkPointSet.h:108
bool operator==(const Self &r) const
Definition: itkMesh.h:233
typename MeshTraits::PointsContainer PointsContainer
Definition: itkPointSet.h:109
BoundingBoxPointer m_BoundingBox
Definition: itkMesh.h:455
typename MeshTraits::CellFeatureIdentifier CellFeatureIdentifier
Definition: itkMesh.h:148
A superclass of the N-dimensional mesh structure; supports point (geometric coordinate and attribute)...
Definition: itkPointSet.h:84
BoundaryAssignmentsContainerVector m_BoundaryAssignmentsContainers
Definition: itkMesh.h:278
typename CellType::CellAutoPointer CellAutoPointer
Definition: itkMesh.h:189
typename MeshTraits::InterpolationWeightType InterpolationWeightType
Definition: itkMesh.h:145
typename PointCellLinksContainer::const_iterator PointCellLinksContainerIterator
Definition: itkMesh.h:182
typename PointsContainer::Iterator PointsContainerIterator
Definition: itkPointSet.h:123
BoundaryAssignmentIdentifier(CellIdentifier cellId, CellFeatureIdentifier featureId)
Definition: itkMesh.h:213
typename CellsContainer::Pointer CellsContainerPointer
Definition: itkMesh.h:165
CellFeatureIdentifier m_FeatureId
Definition: itkMesh.h:221
typename MeshTraits::PointDataContainer PointDataContainer
Definition: itkPointSet.h:110
CellsAllocationMethodType m_CellsAllocationMethod
Definition: itkMesh.h:458
typename CellsContainer::ConstPointer CellsContainerConstPointer
Definition: itkMesh.h:166
typename MeshTraits::PointHashType PointHashType
Definition: itkMesh.h:149
typename MeshTraits::PointIdentifier PointIdentifier
Definition: itkPointSet.h:107
Control indentation during Print() invocation.
Definition: itkIndent.h:49
typename CellLinksContainer::ConstIterator CellLinksContainerIterator
Definition: itkMesh.h:179
CellLinksContainerPointer m_CellLinksContainer
Definition: itkMesh.h:267
typename CellType::MultiVisitor CellMultiVisitorType
Definition: itkMesh.h:192
typename PointDataContainer::ConstIterator PointDataContainerIterator
Definition: itkPointSet.h:124
typename CellDataContainer::ConstPointer CellDataContainerConstPointer
Definition: itkMesh.h:170
typename MeshTraits::CellPixelType CellPixelType
Definition: itkMesh.h:130
typename MeshTraits::CellIdentifier CellIdentifier
Definition: itkMesh.h:147
Base class for most ITK classes.
Definition: itkObject.h:60
typename CellDataContainer::Pointer CellDataContainerPointer
Definition: itkMesh.h:169
typename PointDataContainer::Pointer PointDataContainerPointer
Definition: itkPointSet.h:118
typename MeshTraits::CoordRepType CoordRepType
Definition: itkPointSet.h:106
Base class for all data objects in ITK.
Represent and compute information about bounding boxes.
typename MeshTraits::CellsContainer CellsContainer
Definition: itkMesh.h:153
typename CellsContainer::Iterator CellsContainerIterator
Definition: itkMesh.h:178