ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkPointSet.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 itkPointSet_h
29 #define itkPointSet_h
30 
31 #include "itkDataObject.h"
33 #include <vector>
34 #include <set>
35 
36 namespace itk
37 {
38 
79 template<
80  typename TPixelType,
81  unsigned int VDimension = 3,
82  typename TMeshTraits = DefaultStaticMeshTraits< TPixelType, VDimension, VDimension >
83  >
84 class ITK_TEMPLATE_EXPORT PointSet:public DataObject
85 {
86 public:
87  ITK_DISALLOW_COPY_AND_ASSIGN(PointSet);
88 
90  using Self = PointSet;
94 
96  itkNewMacro(Self);
97 
99  itkTypeMacro(PointSet, Object);
100 
102  using MeshTraits = TMeshTraits;
104 
111 
113  static constexpr unsigned int PointDimension = TMeshTraits::PointDimension;
114 
116  using PointsContainerPointer = typename PointsContainer::Pointer;
117  using PointsContainerConstPointer = typename PointsContainer::ConstPointer;
118  using PointDataContainerPointer = typename PointDataContainer::Pointer;
119  using PointDataContainerConstPointer = typename PointDataContainer::ConstPointer;
120 
122  using PointsContainerConstIterator = typename PointsContainer::ConstIterator;
123  using PointsContainerIterator = typename PointsContainer::Iterator;
124  using PointDataContainerIterator = typename PointDataContainer::ConstIterator;
125 
127  using RegionType = long;
128 
131  itkGetConstMacro(MaximumNumberOfRegions, RegionType);
132 
133 protected:
137 
143 
144 public:
146  void PassStructure(Self *inputPointSet);
147 
148  void Initialize() override;
149 
150  PointIdentifier GetNumberOfPoints() const;
151 
155  void SetPoints(PointsContainer *);
156 
157  PointsContainer * GetPoints();
158 
159  const PointsContainer * GetPoints() const;
160 
161  void SetPointData(PointDataContainer *);
162 
163  PointDataContainer * GetPointData();
164 
165  const PointDataContainer * GetPointData() const;
166 
169  void SetPoint(PointIdentifier, PointType);
170  bool GetPoint(PointIdentifier, PointType *) const;
171  PointType GetPoint(PointIdentifier) const;
173 
176  void SetPointData(PointIdentifier, PixelType);
177  bool GetPointData(PointIdentifier, PixelType *) const;
179 
181  void UpdateOutputInformation() override;
182 
183  void SetRequestedRegionToLargestPossibleRegion() override;
184 
185  void CopyInformation(const DataObject *data) override;
186 
187  void Graft(const DataObject *data) override;
188 
189  bool RequestedRegionIsOutsideOfTheBufferedRegion() override;
190 
191  bool VerifyRequestedRegion() override;
192 
197  void SetRequestedRegion(const DataObject *data) override;
198 
200  virtual void SetRequestedRegion(const RegionType & region);
201 
202  itkGetConstMacro(RequestedRegion, RegionType);
203 
205  virtual void SetBufferedRegion(const RegionType & region);
206 
207  itkGetConstMacro(BufferedRegion, RegionType);
208 
209 protected:
211  PointSet();
212  ~PointSet() override = default;
213  void PrintSelf(std::ostream & os, Indent indent) const override;
215 
216  // If the RegionType is ITK_UNSTRUCTURED_REGION, then the following
217  // variables represent the maximum number of region that the data
218  // object can be broken into, which region out of how many is
219  // currently in the buffered region, and the number of regions and
220  // the specific region requested for the update. Data objects that
221  // do not support any division of the data can simply leave the
222  // MaximumNumberOfRegions as 1. The RequestedNumberOfRegions and
223  // RequestedRegion are used to define the currently requested
224  // region. The LargestPossibleRegion is always requested region = 0
225  // and number of regions = 1;
231 }; // End Class: PointSet
232 } // end namespace itk
233 
234 #ifndef ITK_MANUAL_INSTANTIATION
235 #include "itkPointSet.hxx"
236 #endif
237 
238 /*
239 #ifndef ITK_MANUAL_INSTANTIATION
240 #include "itkPointSet.hxx"
241 #endif
242 */
243 
244 #endif
RegionType m_BufferedRegion
Definition: itkPointSet.h:229
A wrapper of the STL &quot;map&quot; container.
RegionType m_RequestedRegion
Definition: itkPointSet.h:230
PointsContainerPointer m_PointsContainer
Definition: itkPointSet.h:131
class ITK_FORWARD_EXPORT DataObject
Definition: itkDataObject.h:41
RegionType m_MaximumNumberOfRegions
Definition: itkPointSet.h:226
PointDataContainerPointer m_PointDataContainer
Definition: itkPointSet.h:142
RegionType m_NumberOfRegions
Definition: itkPointSet.h:227
A superclass of the N-dimensional mesh structure; supports point (geometric coordinate and attribute)...
Definition: itkPointSet.h:84
RegionType m_RequestedNumberOfRegions
Definition: itkPointSet.h:228
Control indentation during Print() invocation.
Definition: itkIndent.h:49
A simple structure that holds type information for a mesh and its cells.
Base class for most ITK classes.
Definition: itkObject.h:60
Base class for all data objects in ITK.