ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkContourExtractor2DImageFilter.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#ifndef itkContourExtractor2DImageFilter_h
19#define itkContourExtractor2DImageFilter_h
20
21#include <deque>
22#include <list>
23#include <unordered_map>
24#include "itkConceptChecking.h"
26#include "itkImageRegionRange.h"
28#include "itkIndexRange.h"
30
31namespace itk
32{
104
105template <typename TInputImage>
106class ITK_TEMPLATE_EXPORT ContourExtractor2DImageFilter
107 : public ImageToPathFilter<TInputImage, PolyLineParametricPath<2>>
108{
109public:
110 ITK_DISALLOW_COPY_AND_MOVE(ContourExtractor2DImageFilter);
111
113 static constexpr unsigned int InputImageDimension{ TInputImage::ImageDimension };
114
116 using InputImageType = TInputImage;
118
124
126 itkNewMacro(Self);
127
129 itkOverrideGetNameOfClassMacro(ContourExtractor2DImageFilter);
130
132 using InputImagePointer = typename InputImageType::Pointer;
133 using InputIndexType = typename InputImageType::IndexType;
134 using InputSizeType = typename InputImageType::SizeType;
135 using InputOffsetType = typename InputImageType::OffsetType;
136 using InputPixelType = typename InputImageType::PixelType;
137 using InputRegionType = typename InputImageType::RegionType;
143
146
148#ifndef ITK_FUTURE_LEGACY_REMOVE
149 using RegionIndexRange ITK_FUTURE_DEPRECATED(
150 "Please use `itk::ImageRegionIndexRange` or `itk::MakeIndexRange` directly!") =
152#endif
155
156#ifndef ITK_FUTURE_LEGACY_REMOVE
157 using RegionIterator ITK_FUTURE_DEPRECATED("Please use `itk::ImageRegionIterator<TImage>` directly!") =
159 using RegionConstIterator ITK_FUTURE_DEPRECATED("Please use `itk::ImageRegionConstIterator<TImage>` directly!") =
161#endif
162
166 itkSetMacro(ReverseContourOrientation, bool);
167 itkGetConstReferenceMacro(ReverseContourOrientation, bool);
168 itkBooleanMacro(ReverseContourOrientation);
174 itkSetMacro(VertexConnectHighPixels, bool);
175 itkGetConstReferenceMacro(VertexConnectHighPixels, bool);
176 itkBooleanMacro(VertexConnectHighPixels);
180 itkSetMacro(LabelContours, bool);
181 itkGetConstReferenceMacro(LabelContours, bool);
182 itkBooleanMacro(LabelContours);
186 void
188
189 itkGetConstReferenceMacro(RequestedRegion, InputRegionType);
190 void
192
196 itkSetMacro(ContourValue, InputRealType);
197 itkGetConstReferenceMacro(ContourValue, InputRealType);
200 itkConceptMacro(InputPixelTypeComparable, (Concept::Comparable<InputPixelType>));
201 itkConceptMacro(InputHasPixelTraitsCheck, (Concept::HasPixelTraits<InputPixelType>));
202 itkConceptMacro(InputHasNumericTraitsCheck, (Concept::HasNumericTraits<InputPixelType>));
203
204protected:
206 ~ContourExtractor2DImageFilter() override = default;
207 void
208 PrintSelf(std::ostream & os, Indent indent) const override;
209
210 void
211 GenerateData() override;
212
216 void
218
219private:
220 using LabelsContainer = std::vector<InputPixelType>;
221 using LabelsConstIterator = typename LabelsContainer::const_iterator;
222 using LabelsIterator = typename LabelsContainer::iterator;
223
230
231 // Represent each contour as deque of vertices to facilitate addition of
232 // nodes at beginning or end. At the end of the processing, we will copy
233 // the contour into a PolyLineParametricPath.
234 // We subclass the deque to store an additional bit of information: an
235 // identification number for each growing contour. We use this number so
236 // that when it becomes necessary to merge two growing contours, we can
237 // merge the newer one into the older one. This helps because then we can
238 // guarantee that the output contour list is ordered from left to right,
239 // top to bottom (in terms of the first pixel of the contour encountered
240 // by the marching squares). Currently we make no guarantees that this
241 // pixel is the first pixel in the contour list, just that the contours
242 // are so ordered in the output. Ensuring this latter condition (first
243 // pixel traversed = first pixel in contour) would be possible by either
244 // changing the merging rules, which would make the contouring operation
245 // slower, or by storing additional data as to which pixel was first.
246 class ContourType : public std::deque<VertexType>
247 {
248 public:
249 unsigned int m_ContourNumber;
250 };
251
252 // Store all the growing contours in a list. We may need to delete contours
253 // from anywhere in the sequence (when we merge them together), so we need to
254 // use a list instead of a vector or similar.
255 using ContourConstIterator = typename ContourType::const_iterator;
256 using ContourContainerType = std::list<ContourType>;
257 using ContourContainerIterator = typename ContourContainerType::iterator;
258 using ContourContainerConstIterator = typename ContourContainerType::const_iterator;
259
260 // We use a hash to associate the endpoints of each contour with the
261 // contour itself. This makes it easy to look up which contour we should add
262 // a new arc to.
263
264 // We can't store the contours themselves in the hashtable because we
265 // need to have two tables (one to hash from beginpoint -> contour and one
266 // for endpoint -> contour), and sometimes will remove a contour from the
267 // tables (if it has been closed or merged with another contour). Because
268 // sometimes we will need to merge contours, we need to be able to quickly
269 // remove contours from our list when they have been merged into
270 // another. Thus, we store an iterator pointing to the contour in the list.
271
273 {
275 inline size_t
276 operator()(const VertexType & v) const noexcept
277 {
278 return std::hash<CoordinateType>{}(v[0]) ^ (std::hash<CoordinateType>{}(v[1]) << 1);
279 }
280 };
281 using VertexToContourContainerIteratorMap = std::unordered_map<VertexType, ContourContainerIterator, VertexHash>;
282 using VertexToContourContainerIteratorMapIterator = typename VertexToContourContainerIteratorMap::iterator;
283 using VertexToContourContainerIteratorMapConstIterator = typename VertexToContourContainerIteratorMap::const_iterator;
284 using VertexToContourContainerIteratorMapKeyValuePair = typename VertexToContourContainerIteratorMap::value_type;
285
286 // Subroutine to create contours for a single label.
287 void
289 const InputImageType * input,
290 const InputRegionType usableRegion,
291 SizeValueType totalNumberOfPixels,
292 ContourContainerType & contoursOutput);
293
294 // Subroutine to handle the case that the supplied values are
295 // labels, which are *not* compared to a contour value.
296 void
298
301 InputPixelType toValue,
302 InputIndexType fromIndex,
303 InputOffsetType toOffset);
304
312
313 void
314 AddSegment(const VertexType from, const VertexType to, ContourData & contourData);
315
316 void
317 FillOutputs(const std::vector<InputPixelType> & allLabels,
318 std::unordered_map<InputPixelType, ContourContainerType> & labelsContoursOutput);
319
321};
322} // end namespace itk
323
324#ifndef ITK_MANUAL_INSTANTIATION
325# include "itkContourExtractor2DImageFilter.hxx"
326#endif
327
328#endif
typename InputImageType::IndexType InputIndexType
typename InputImageType::PixelType InputPixelType
typename ContourContainerType::const_iterator ContourContainerConstIterator
void SetRequestedRegion(const InputRegionType region)
~ContourExtractor2DImageFilter() override=default
typename VertexToContourContainerIteratorMap::value_type VertexToContourContainerIteratorMapKeyValuePair
typename VertexToContourContainerIteratorMap::iterator VertexToContourContainerIteratorMapIterator
void GenerateInputRequestedRegion() override
ImageRegionRange< const InputImageType > RegionConstRange
typename LabelsContainer::iterator LabelsIterator
VertexType InterpolateContourPosition(InputPixelType fromValue, InputPixelType toValue, InputIndexType fromIndex, InputOffsetType toOffset)
ImageToPathFilter< InputImageType, OutputPathType > Superclass
typename NumericTraits< InputPixelType >::RealType InputRealType
typename InputImageType::SizeType InputSizeType
void AddSegment(const VertexType from, const VertexType to, ContourData &contourData)
void PrintSelf(std::ostream &os, Indent indent) const override
typename InputImageType::Pointer InputImagePointer
typename LabelsContainer::const_iterator LabelsConstIterator
typename InputImageType::OffsetType InputOffsetType
ImageRegionRange< InputImageType > RegionRange
typename ContourContainerType::iterator ContourContainerIterator
typename VertexToContourContainerIteratorMap::const_iterator VertexToContourContainerIteratorMapConstIterator
typename InputImageType::RegionType InputRegionType
void FillOutputs(const std::vector< InputPixelType > &allLabels, std::unordered_map< InputPixelType, ContourContainerType > &labelsContoursOutput)
std::unordered_map< VertexType, ContourContainerIterator, VertexHash > VertexToContourContainerIteratorMap
typename ContourType::const_iterator ContourConstIterator
void CreateSingleContour(InputPixelType label, const InputImageType *input, const InputRegionType usableRegion, SizeValueType totalNumberOfPixels, ContourContainerType &contoursOutput)
A multi-dimensional iterator templated over image type that walks a region of pixels.
A multi-dimensional iterator templated over image type that walks a region of pixels.
Control indentation during Print() invocation.
Definition itkIndent.h:50
Represent a path of line segments through ND Space.
VectorContainer< unsigned int, VertexType > VertexListType
Implements transparent reference counting.
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned long SizeValueType
Definition itkIntTypes.h:86
IndexRange< VDimension, false > ImageRegionIndexRange
vcl_size_t operator()(const VertexType &v) const noexcept