ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkContourExtractor2DImageFilter.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 #ifndef itkContourExtractor2DImageFilter_h
19 #define itkContourExtractor2DImageFilter_h
20 
21 #include "itkImageToPathFilter.h"
23 #include "itkConceptChecking.h"
24 #include <unordered_map>
25 #include <deque>
26 #include <list>
27 
28 namespace itk
29 {
95 template< typename TInputImage >
96 class ITK_TEMPLATE_EXPORT ContourExtractor2DImageFilter:
97  public ImageToPathFilter< TInputImage, PolyLineParametricPath< 2 > >
98 {
99 public:
100  ITK_DISALLOW_COPY_AND_ASSIGN(ContourExtractor2DImageFilter);
101 
103  static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
104 
106  using InputImageType = TInputImage;
108 
114 
116  itkNewMacro(Self);
117 
120 
122  using InputImagePointer = typename InputImageType::Pointer;
123  using InputPixelType = typename InputImageType::PixelType;
125  using InputOffsetType = typename InputImageType::OffsetType;
130 
133 
134  using VertexListConstPointer = typename VertexListType::ConstPointer;
135 
138  itkSetMacro(ReverseContourOrientation, bool);
139  itkGetConstReferenceMacro(ReverseContourOrientation, bool);
140  itkBooleanMacro(ReverseContourOrientation);
142 
146  itkSetMacro(VertexConnectHighPixels, bool);
147  itkGetConstReferenceMacro(VertexConnectHighPixels, bool);
148  itkBooleanMacro(VertexConnectHighPixels);
150 
153  void SetRequestedRegion(const InputRegionType region);
154 
155  itkGetConstReferenceMacro(RequestedRegion, InputRegionType);
156  void ClearRequestedRegion();
157 
160  itkSetMacro(ContourValue, InputRealType);
161  itkGetConstReferenceMacro(ContourValue, InputRealType);
163 
164 #ifdef ITK_USE_CONCEPT_CHECKING
165  // Begin concept checking
166  itkConceptMacro( DimensionShouldBe2,
168  itkConceptMacro( InputPixelTypeComparable,
170  itkConceptMacro( InputHasPixelTraitsCheck,
172  itkConceptMacro( InputHasNumericTraitsCheck,
174  // End concept checking
175 #endif
176 
177 protected:
178 
180  ~ContourExtractor2DImageFilter() override = default;
181  void PrintSelf(std::ostream & os, Indent indent) const override;
182 
183  void GenerateData() override;
184 
188  void GenerateInputRequestedRegion() override;
189 
190 private:
191  VertexType InterpolateContourPosition(InputPixelType fromValue,
192  InputPixelType toValue,
193  InputIndexType fromIndex,
194  InputOffsetType toOffset);
195 
196  void AddSegment(const VertexType from, const VertexType to);
197 
198  void FillOutputs();
199 
206 
207  // Represent each contour as deque of vertices to facilitate addition of
208  // nodes at beginning or end. At the end of the processing, we will copy
209  // the contour into a PolyLineParametricPath.
210  // We subclass the deque to store an additional bit of information: an
211  // identification number for each growing contour. We use this number so
212  // that when it becomes necessary to merge two growing contours, we can
213  // merge the newer one into the older one. This helps because then we can
214  // guarantee that the output contour list is ordered from left to right,
215  // top to bottom (in terms of the first pixel of the contour encountered
216  // by the marching squares). Currently we make no guarantees that this
217  // pixel is the first pixel in the contour list, just that the contours
218  // are so ordered in the output. Ensuring this latter condition (first
219  // pixel traversed = first pixel in contour) would be possible by either
220  // changing the merging rules, which would make the contouring operation
221  //slower, or by storing additional data as to which pixel was first.
222  class ContourType:public std::deque< VertexType >
223  {
224 public:
225  unsigned int m_ContourNumber;
226  };
227 
228  // Store all the growing contours in a list. We may need to delete contours
229  // from anywhere in the sequence (when we merge them together), so we need to
230  // use a list instead of a vector or similar.
231  using ContourContainer = std::list< ContourType >;
232  using ContourRef = typename ContourContainer::iterator;
233 
234  // declare the hash function we are using for the hash_map.
235  struct VertexHash {
236  using CoordinateType = typename VertexType::CoordRepType;
237  inline SizeValueType operator()(const VertexType & k) const
238  {
239  // Xor the hashes of the vertices together, after multiplying the
240  // first by some number, so that identical (x,y) vertex indices
241  // don't all hash to the same bucket. This is a decent if not
242  // optimal hash.
243  const SizeValueType hashVertex1 = this->float_hash(k[0] * 0xbeef);
244  const SizeValueType hashVertex2 = this->float_hash(k[1]);
245  const SizeValueType hashValue = hashVertex1 ^ hashVertex2;
246 
247  return hashValue;
248  }
249 
250  // Define hash function for floats. Based on method from
251  // http://www.brpreiss.com/books/opus4/html/page217.html
252  inline SizeValueType float_hash(const CoordinateType & k) const
253  {
254  if ( k == 0 )
255  {
256  return 0;
257  }
258  int exponent;
259  CoordinateType mantissa = std::frexp(k, &exponent);
260  auto value = static_cast< SizeValueType >( std::fabs(mantissa) );
261  value = ( 2 * value - 1 ) * ~0U;
262  return value;
263  }
264  };
265 
266  // We use a hash to associate the endpoints of each contour with the
267  // contour itself. This makes it easy to look up which contour we should add
268  // a new arc to.
269  // We can't store the contours themselves in the hashtable because we
270  // need to have two tables (one to hash from beginpoint -> contour and one
271  // for endpoint -> contour), and sometimes will remove a contour from the
272  // tables (if it has been closed or merged with another contour). So in the
273  // hash table we store a reference to the contour. Because sometimes we will
274  // need to merge contours, we need to be able to quickly remove contours
275  // from our list when they have been merged into another. Thus, we store
276  // an iterator pointing to the contour in the list.
277 
278  using VertexToContourMap = std::unordered_map< VertexType, ContourRef, VertexHash >;
279  using VertexMapIterator = typename VertexToContourMap::iterator;
280  using VertexContourRefPair = typename VertexToContourMap::value_type;
281 
282  // The contours we find in the image are stored here
284 
285  // And indexed by their beginning and ending points here
288 };
289 } // end namespace itk
290 
291 #ifndef ITK_MANUAL_INSTANTIATION
292 #include "itkContourExtractor2DImageFilter.hxx"
293 #endif
294 
295 #endif
typename OutputPathType::VertexType VertexType
typename VertexToContourMap::value_type VertexContourRefPair
Represent a path of line segments through ND Space.
Light weight base class for most itk classes.
Define numeric traits for std::vector.
typename NumericTraits< InputPixelType >::RealType InputRealType
SizeValueType float_hash(const CoordinateType &k) const
unsigned long SizeValueType
Definition: itkIntTypes.h:83
Base class for filters that take an image as input and produce an path as output. ...
typename OutputPathType::VertexListType VertexListType
typename InputImageType::Pointer InputImagePointer
SizeValueType operator()(const VertexType &k) const
typename InputImageType::RegionType InputRegionType
typename VertexToContourMap::iterator VertexMapIterator
typename ContourContainer::iterator ContourRef
typename InputImageType::OffsetType InputOffsetType
A templated class holding a point in n-Dimensional image space.
std::unordered_map< VertexType, ContourRef, VertexHash > VertexToContourMap
typename OutputPathType::Pointer OutputPathPointer
Define a front-end to the STL &quot;vector&quot; container that conforms to the IndexedContainerInterface.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
typename InputImageType::IndexType InputIndexType
Computes a list of PolyLineParametricPath objects from the contours in a 2D image.
#define itkConceptMacro(name, concept)
typename InputImageType::PixelType InputPixelType
typename VertexListType::ConstPointer VertexListConstPointer