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;
140 using VertexListConstPointer = typename VertexListType::ConstPointer;
142 using VertexValueType = typename VertexType::ValueType;
143
146
153
156 itkSetMacro(ReverseContourOrientation, bool);
157 itkGetConstReferenceMacro(ReverseContourOrientation, bool);
158 itkBooleanMacro(ReverseContourOrientation);
160
164 itkSetMacro(VertexConnectHighPixels, bool);
165 itkGetConstReferenceMacro(VertexConnectHighPixels, bool);
166 itkBooleanMacro(VertexConnectHighPixels);
168
170 itkSetMacro(LabelContours, bool);
171 itkGetConstReferenceMacro(LabelContours, bool);
172 itkBooleanMacro(LabelContours);
174
177 void
179
180 itkGetConstReferenceMacro(RequestedRegion, InputRegionType);
181 void
183
186 itkSetMacro(ContourValue, InputRealType);
187 itkGetConstReferenceMacro(ContourValue, InputRealType);
189
191 itkConceptMacro(InputPixelTypeComparable, (Concept::Comparable<InputPixelType>));
192 itkConceptMacro(InputHasPixelTraitsCheck, (Concept::HasPixelTraits<InputPixelType>));
193 itkConceptMacro(InputHasNumericTraitsCheck, (Concept::HasNumericTraits<InputPixelType>));
194
195protected:
197 ~ContourExtractor2DImageFilter() override = default;
198 void
199 PrintSelf(std::ostream & os, Indent indent) const override;
200
201 void
202 GenerateData() override;
203
207 void
209
210private:
211 using LabelsContainer = std::vector<InputPixelType>;
212 using LabelsConstIterator = typename LabelsContainer::const_iterator;
213 using LabelsIterator = typename LabelsContainer::iterator;
214
221
222 // Represent each contour as deque of vertices to facilitate addition of
223 // nodes at beginning or end. At the end of the processing, we will copy
224 // the contour into a PolyLineParametricPath.
225 // We subclass the deque to store an additional bit of information: an
226 // identification number for each growing contour. We use this number so
227 // that when it becomes necessary to merge two growing contours, we can
228 // merge the newer one into the older one. This helps because then we can
229 // guarantee that the output contour list is ordered from left to right,
230 // top to bottom (in terms of the first pixel of the contour encountered
231 // by the marching squares). Currently we make no guarantees that this
232 // pixel is the first pixel in the contour list, just that the contours
233 // are so ordered in the output. Ensuring this latter condition (first
234 // pixel traversed = first pixel in contour) would be possible by either
235 // changing the merging rules, which would make the contouring operation
236 // slower, or by storing additional data as to which pixel was first.
237 class ContourType : public std::deque<VertexType>
238 {
239 public:
240 unsigned int m_ContourNumber;
241 };
242
243 // Store all the growing contours in a list. We may need to delete contours
244 // from anywhere in the sequence (when we merge them together), so we need to
245 // use a list instead of a vector or similar.
246 using ContourConstIterator = typename ContourType::const_iterator;
247 using ContourContainerType = std::list<ContourType>;
248 using ContourContainerIterator = typename ContourContainerType::iterator;
249 using ContourContainerConstIterator = typename ContourContainerType::const_iterator;
250
251 // We use a hash to associate the endpoints of each contour with the
252 // contour itself. This makes it easy to look up which contour we should add
253 // a new arc to.
254
255 // We can't store the contours themselves in the hashtable because we
256 // need to have two tables (one to hash from beginpoint -> contour and one
257 // for endpoint -> contour), and sometimes will remove a contour from the
258 // tables (if it has been closed or merged with another contour). Because
259 // sometimes we will need to merge contours, we need to be able to quickly
260 // remove contours from our list when they have been merged into
261 // another. Thus, we store an iterator pointing to the contour in the list.
262
264 {
265 using CoordinateType = typename VertexType::CoordinateType;
266 inline size_t
267 operator()(const VertexType & v) const noexcept
268 {
269 return std::hash<CoordinateType>{}(v[0]) ^ (std::hash<CoordinateType>{}(v[1]) << 1);
270 }
271 };
272 using VertexToContourContainerIteratorMap = std::unordered_map<VertexType, ContourContainerIterator, VertexHash>;
273 using VertexToContourContainerIteratorMapIterator = typename VertexToContourContainerIteratorMap::iterator;
274 using VertexToContourContainerIteratorMapConstIterator = typename VertexToContourContainerIteratorMap::const_iterator;
275 using VertexToContourContainerIteratorMapKeyValuePair = typename VertexToContourContainerIteratorMap::value_type;
276
277 // Subroutine to create contours for a single label.
278 void
280 const InputImageType * input,
281 const InputRegionType usableRegion,
282 SizeValueType totalNumberOfPixels,
283 ContourContainerType & contoursOutput);
284
285 // Subroutine to handle the case that the supplied values are
286 // labels, which are *not* compared to a contour value.
287 void
289
292 InputPixelType toValue,
293 InputIndexType fromIndex,
294 InputOffsetType toOffset);
295
303
304 void
305 AddSegment(const VertexType from, const VertexType to, ContourData & contourData);
306
307 void
308 FillOutputs(const std::vector<InputPixelType> & allLabels,
309 std::unordered_map<InputPixelType, ContourContainerType> & labelsContoursOutput);
310
312};
313} // end namespace itk
314
315#ifndef ITK_MANUAL_INSTANTIATION
316# include "itkContourExtractor2DImageFilter.hxx"
317#endif
318
319#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
ImageRegionIterator< InputImageType > RegionIterator
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
typename OutputPathType::VertexListType VertexListType
void AddSegment(const VertexType from, const VertexType to, ContourData &contourData)
void PrintSelf(std::ostream &os, Indent indent) const override
ImageRegionConstIterator< InputImageType > RegionConstIterator
ImageRegionIndexRange< InputImageType::ImageDimension > RegionIndexRange
typename OutputPathType::Pointer OutputPathPointer
typename InputImageType::Pointer InputImagePointer
typename LabelsContainer::const_iterator LabelsConstIterator
typename InputImageType::OffsetType InputOffsetType
ImageRegionRange< InputImageType > RegionRange
typename VertexListType::ConstPointer VertexListConstPointer
typename ContourContainerType::iterator ContourContainerIterator
typename VertexToContourContainerIteratorMap::const_iterator VertexToContourContainerIteratorMapConstIterator
typename OutputPathType::VertexType VertexType
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