ITK  4.9.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 Mesh:public PointSet< TPixelType, VDimension, TMeshTraits >
109 {
110 public:
112  typedef Mesh Self;
116 
118 
120  itkNewMacro(Self);
121 
123  itkTypeMacro(Mesh, PointSet);
124 
126  typedef TMeshTraits MeshTraits;
127  typedef typename MeshTraits::PixelType PixelType;
128  typedef typename MeshTraits::CellPixelType CellPixelType;
129 
131  itkStaticConstMacro(PointDimension, unsigned int,
132  TMeshTraits::PointDimension);
133  itkStaticConstMacro(MaxTopologicalDimension, unsigned int,
134  TMeshTraits::MaxTopologicalDimension);
136 
143 
145  typedef typename MeshTraits::CoordRepType CoordRepType;
146  typedef typename MeshTraits::InterpolationWeightType InterpolationWeightType;
147  typedef typename MeshTraits::PointIdentifier PointIdentifier;
148  typedef typename MeshTraits::CellIdentifier CellIdentifier;
149  typedef typename MeshTraits::CellFeatureIdentifier CellFeatureIdentifier;
150  typedef typename MeshTraits::PointHashType PointHashType;
151  typedef typename MeshTraits::PointType PointType;
152  typedef typename MeshTraits::PointsContainer PointsContainer;
153  typedef typename MeshTraits::CellTraits CellTraits;
154  typedef typename MeshTraits::CellsContainer CellsContainer;
155  typedef typename MeshTraits::PointCellLinksContainer PointCellLinksContainer;
156  typedef typename MeshTraits::CellLinksContainer CellLinksContainer;
157  typedef typename MeshTraits::PointDataContainer PointDataContainer;
158  typedef typename MeshTraits::CellDataContainer CellDataContainer;
159 
161  typedef BoundingBox< PointIdentifier, itkGetStaticConstMacro(PointDimension),
163 
165  typedef typename PointsContainer::Pointer PointsContainerPointer;
166  typedef typename CellsContainer::Pointer CellsContainerPointer;
167  typedef typename CellsContainer::ConstPointer CellsContainerConstPointer;
168  typedef typename CellLinksContainer::Pointer CellLinksContainerPointer;
169  typedef typename PointDataContainer::Pointer PointDataContainerPointer;
170  typedef typename CellDataContainer::Pointer CellDataContainerPointer;
171  typedef typename CellDataContainer::ConstPointer CellDataContainerConstPointer;
173  typedef typename CellLinksContainer::ConstPointer CellLinksContainerConstPointer;
174 
176  typedef typename PointsContainer::ConstIterator PointsContainerConstIterator;
177  typedef typename PointsContainer::Iterator PointsContainerIterator;
178  typedef typename CellsContainer::ConstIterator CellsContainerConstIterator;
179  typedef typename CellsContainer::Iterator CellsContainerIterator;
180  typedef typename CellLinksContainer::ConstIterator CellLinksContainerIterator;
181  typedef typename PointDataContainer::ConstIterator PointDataContainerIterator;
182  typedef typename CellDataContainer::ConstIterator CellDataContainerIterator;
183  typedef typename PointCellLinksContainer::const_iterator PointCellLinksContainerIterator;
184 
187 
191 
194 
206  {
207 public:
210 
215  CellFeatureIdentifier featureId):
216  m_CellId(cellId), m_FeatureId(featureId) {}
218 
221 
224 
227  bool operator<(const Self & r) const
228  {
229  return ( ( m_CellId < r.m_CellId )
230  || ( ( m_CellId == r.m_CellId ) && ( m_FeatureId < r.m_FeatureId ) ) );
231  }
232 
235  bool operator==(const Self & r) const
236  {
237  return ( ( m_CellId == r.m_CellId ) && ( m_FeatureId == r.m_FeatureId ) );
238  }
239  }; // End Class: Mesh::BoundaryAssignmentIdentifier
241 
253  typedef std::vector< BoundaryAssignmentsContainerPointer >
255 
256 protected:
257 
261 
267 
272 
283 
284 public:
287 
288  void PassStructure(Self *inputMesh);
289 
290  virtual void Initialize() ITK_OVERRIDE;
291 
293  virtual void CopyInformation(const DataObject *data) ITK_OVERRIDE;
294 
295  virtual void Graft(const DataObject *data) ITK_OVERRIDE;
296 
299  const BoundingBoxType * GetBoundingBox() const;
300 
306 
307  CellLinksContainer * GetCellLinks();
308 
309  const CellLinksContainer * GetCellLinks() const;
310 
313  void SetCells(CellsContainer *);
314 
315  CellsContainer * GetCells();
316 
317  const CellsContainer * GetCells() const;
318 
324 
325  CellDataContainer * GetCellData();
326 
327  const CellDataContainer * GetCellData() const;
328 
329 #if !defined( ITK_WRAPPING_PARSER )
330 
338  void SetBoundaryAssignments(int dimension,
340 
341 
343 
345  int dimension) const;
346 #endif
347 
354  bool GetCell(CellIdentifier, CellAutoPointer &) const;
355 
361 
374  void SetBoundaryAssignment(int dimension, CellIdentifier cellId,
375  CellFeatureIdentifier featureId,
376  CellIdentifier boundaryId);
377 
386  bool GetBoundaryAssignment(int dimension, CellIdentifier cellId,
387  CellFeatureIdentifier featureId,
388  CellIdentifier *boundaryId) const;
389 
390  bool RemoveBoundaryAssignment(int dimension, CellIdentifier cellId,
391  CellFeatureIdentifier featureId);
392 
395  CellIdentifier) const;
396 
399  bool GetCellBoundaryFeature(int dimension, CellIdentifier,
401 
407  int dimension, CellIdentifier, CellFeatureIdentifier,
408  std::set< CellIdentifier > *cellSet);
409 
415  std::set< CellIdentifier > *cellSet);
416 
426  CellAutoPointer &) const;
427 
430  void BuildCellLinks() const;
431 
435  virtual void Accept(CellMultiVisitorType *mv) const;
436 
441  itkSetMacro(CellsAllocationMethod, CellsAllocationMethodType);
442  itkGetConstReferenceMacro(CellsAllocationMethod, CellsAllocationMethodType);
444 
445 protected:
447  Mesh();
448  ~Mesh();
449  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
451 
455  void ReleaseCellsMemory();
456 
460 
461 private:
462  Mesh(const Self &) ITK_DELETE_FUNCTION;
463  void operator=(const Self &) ITK_DELETE_FUNCTION;
464 
466 }; // End Class: Mesh
467 } // end namespace itk
468 
469 #ifndef ITK_MANUAL_INSTANTIATION
470 #ifndef ITK_WRAPPING_PARSER
471 #include "itkMesh.hxx"
472 #endif
473 #endif
474 
475 #endif
static const unsigned int PointDimension
Definition: itkMesh.h:132
CellsContainerPointer m_CellsContainer
Definition: itkMesh.h:260
TMeshTraits MeshTraits
Definition: itkMesh.h:123
PointsContainer::ConstIterator PointsContainerConstIterator
Definition: itkMesh.h:176
BoundingBox< PointIdentifier, itkGetStaticConstMacro(PointDimension), CoordRepType, PointsContainer > BoundingBoxType
Definition: itkMesh.h:162
const BoundingBoxType * GetBoundingBox() const
bool operator<(const Self &r) const
Definition: itkMesh.h:227
CellIdentifier GetCellNeighbors(CellIdentifier cellId, std::set< CellIdentifier > *cellSet)
MeshTraits::PointDataContainer PointDataContainer
Definition: itkMesh.h:157
SelfAutoPointer CellAutoPointer
MeshTraits::PointCellLinksContainer PointCellLinksContainer
Definition: itkMesh.h:155
virtual void CopyInformation(const DataObject *data) override
CellFeatureCount GetNumberOfCellBoundaryFeatures(int dimension, CellIdentifier) const
CellsContainer::Iterator CellsContainerIterator
Definition: itkMesh.h:179
PointCellLinksContainer::const_iterator PointCellLinksContainerIterator
Definition: itkMesh.h:183
A wrapper of the STL &quot;map&quot; container.
BoundaryAssignmentIdentifier Self
Definition: itkMesh.h:209
MeshTraits::PointsContainer PointsContainer
Definition: itkMesh.h:152
MeshTraits::PixelType PixelType
Definition: itkMesh.h:127
MeshTraits::CellIdentifier CellIdentifier
Definition: itkMesh.h:148
MapContainer< BoundaryAssignmentIdentifier, CellIdentifier > BoundaryAssignmentsContainer
Definition: itkMesh.h:250
void SetCells(CellsContainer *)
MeshTraits::PointType PointType
Definition: itkMesh.h:151
CellDataContainerPointer m_CellDataContainer
Definition: itkMesh.h:266
A visitor that can visit different cell types in a mesh. CellInterfaceVisitor instances can be regist...
PointDataContainer::Pointer PointDataContainerPointer
Definition: itkMesh.h:169
void SetBoundaryAssignment(int dimension, CellIdentifier cellId, CellFeatureIdentifier featureId, CellIdentifier boundaryId)
bool GetCellBoundaryFeature(int dimension, CellIdentifier, CellFeatureIdentifier, CellAutoPointer &) const
Implements the N-dimensional mesh structure.
Definition: itkMesh.h:108
CellDataContainer::Pointer CellDataContainerPointer
Definition: itkMesh.h:170
bool GetCell(CellIdentifier, CellAutoPointer &) const
Superclass::RegionType RegionType
Definition: itkMesh.h:117
An abstract interface for cells.
CellsContainer::ConstIterator CellsContainerConstIterator
Definition: itkMesh.h:178
std::vector< BoundaryAssignmentsContainerPointer > BoundaryAssignmentsContainerVector
Definition: itkMesh.h:254
MeshTraits::CellDataContainer CellDataContainer
Definition: itkMesh.h:158
void SetCellData(CellDataContainer *)
BoundaryAssignmentsContainerPointer GetBoundaryAssignments(int dimension)
CellLinksContainer::ConstIterator CellLinksContainerIterator
Definition: itkMesh.h:180
BoundingBoxType::Pointer BoundingBoxPointer
Definition: itkMesh.h:172
Mesh Self
Definition: itkMesh.h:112
CellFeatureIdentifier CellFeatureCount
Definition: itkMesh.h:186
CellLinksContainer::ConstPointer CellLinksContainerConstPointer
Definition: itkMesh.h:173
void SetCellLinks(CellLinksContainer *)
bool operator==(const Self &r) const
Definition: itkMesh.h:235
CellLinksContainer::Pointer CellLinksContainerPointer
Definition: itkMesh.h:168
PointsContainer::Pointer PointsContainerPointer
Definition: itkMesh.h:165
MeshTraits::InterpolationWeightType InterpolationWeightType
Definition: itkMesh.h:146
BoundingBoxPointer m_BoundingBox
Definition: itkMesh.h:459
PointDataContainer::ConstIterator PointDataContainerIterator
Definition: itkMesh.h:181
CellDataContainer * GetCellData()
MeshTraits::CoordRepType CoordRepType
Definition: itkMesh.h:145
PointSet< TPixelType, VDimension, TMeshTraits > Superclass
Definition: itkMesh.h:113
bool GetAssignedCellBoundaryIfOneExists(int dimension, CellIdentifier, CellFeatureIdentifier, CellAutoPointer &) const
CellType::MultiVisitor CellMultiVisitorType
Definition: itkMesh.h:193
A superclass of the N-dimensional mesh structure; supports point (geometric coordinate and attribute)...
Definition: itkPointSet.h:84
BoundaryAssignmentsContainerVector m_BoundaryAssignmentsContainers
Definition: itkMesh.h:282
CellType::CellAutoPointer CellAutoPointer
Definition: itkMesh.h:190
MeshTraits::PointsContainer PointsContainer
Definition: itkPointSet.h:107
MeshTraits::PointIdentifier PointIdentifier
Definition: itkMesh.h:147
void PassStructure(Self *inputMesh)
virtual void Graft(const DataObject *data) override
void ReleaseCellsMemory()
virtual void Accept(CellMultiVisitorType *mv) const
MeshTraits::CellLinksContainer CellLinksContainer
Definition: itkMesh.h:156
CellsContainer::Pointer CellsContainerPointer
Definition: itkMesh.h:166
MeshTraits::CellFeatureIdentifier CellFeatureIdentifier
Definition: itkMesh.h:149
BoundaryAssignmentIdentifier(CellIdentifier cellId, CellFeatureIdentifier featureId)
Definition: itkMesh.h:214
bool RemoveBoundaryAssignment(int dimension, CellIdentifier cellId, CellFeatureIdentifier featureId)
SmartPointer< const Self > ConstPointer
Definition: itkMesh.h:115
CellsContainer * GetCells()
CellLinksContainer * GetCellLinks()
CellsAllocationMethodType
Definition: itkMesh.h:139
void BuildCellLinks() const
CellsContainer::ConstPointer CellsContainerConstPointer
Definition: itkMesh.h:167
void PrintSelf(std::ostream &os, Indent indent) const override
CellFeatureIdentifier m_FeatureId
Definition: itkMesh.h:223
BoundaryAssignmentsContainer::Pointer BoundaryAssignmentsContainerPointer
Definition: itkMesh.h:252
MeshTraits::CellsContainer CellsContainer
Definition: itkMesh.h:154
CellIdentifier GetCellBoundaryFeatureNeighbors(int dimension, CellIdentifier, CellFeatureIdentifier, std::set< CellIdentifier > *cellSet)
CellsAllocationMethodType m_CellsAllocationMethod
Definition: itkMesh.h:465
virtual void Initialize() override
CellInterface< CellPixelType, CellTraits > CellType
Definition: itkMesh.h:189
MeshTraits::CellPixelType CellPixelType
Definition: itkMesh.h:128
Control indentation during Print() invocation.
Definition: itkIndent.h:49
SmartPointer< Self > Pointer
Definition: itkMesh.h:114
CellLinksContainerPointer m_CellLinksContainer
Definition: itkMesh.h:271
void SetBoundaryAssignments(int dimension, BoundaryAssignmentsContainer *)
void SetCell(CellIdentifier, CellAutoPointer &)
bool GetBoundaryAssignment(int dimension, CellIdentifier cellId, CellFeatureIdentifier featureId, CellIdentifier *boundaryId) const
CellDataContainer::ConstPointer CellDataContainerConstPointer
Definition: itkMesh.h:171
MeshTraits::CellTraits CellTraits
Definition: itkMesh.h:153
MeshTraits::PointHashType PointHashType
Definition: itkMesh.h:150
CellIdentifier GetNumberOfCells() const
CellDataContainer::ConstIterator CellDataContainerIterator
Definition: itkMesh.h:182
static const unsigned int MaxTopologicalDimension
Definition: itkMesh.h:134
Base class for all data objects in ITK.
Represent and compute information about bounding boxes.
PointsContainer::Iterator PointsContainerIterator
Definition: itkMesh.h:177