ITK  6.0.0
Insight Toolkit
itkIsoContourDistanceImageFilter.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 itkIsoContourDistanceImageFilter_h
19 #define itkIsoContourDistanceImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 #include "itkNarrowBand.h"
24 #include "itkNumericTraits.h"
25 #include <mutex>
26 
27 namespace itk
28 {
58 template <typename TInputImage, typename TOutputImage>
59 class ITK_TEMPLATE_EXPORT IsoContourDistanceImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
60 {
61 public:
62  ITK_DISALLOW_COPY_AND_MOVE(IsoContourDistanceImageFilter);
63 
69 
71  itkNewMacro(Self);
72 
74  itkOverrideGetNameOfClassMacro(IsoContourDistanceImageFilter);
75 
77  using typename Superclass::InputImageType;
78  using typename Superclass::OutputImageType;
79 
82  static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
83  static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
84 
87  using PixelType = typename OutputImageType::PixelType;
88  using InputPixelType = typename InputImageType::PixelType;
90 
92 
95 
98 
99  using InputSpacingType = typename InputImageType::SpacingType;
100 
108 
111  itkSetMacro(LevelSetValue, PixelRealType);
112  itkGetConstMacro(LevelSetValue, PixelRealType);
117  itkSetMacro(FarValue, PixelType);
118  itkGetConstMacro(FarValue, PixelType);
123  itkSetMacro(NarrowBanding, bool);
124  itkGetConstMacro(NarrowBanding, bool);
125  itkBooleanMacro(NarrowBanding);
129  void
130  SetNarrowBand(NarrowBandType * ptr);
131 
134  {
135  return m_NarrowBand;
136  }
137 
138 #ifdef ITK_USE_CONCEPT_CHECKING
139  // Begin concept checking
140  itkConceptMacro(InputEqualityComparableCheck, (Concept::EqualityComparable<InputPixelType>));
141  itkConceptMacro(OutputEqualityComparableCheck, (Concept::EqualityComparable<PixelType>));
143  itkConceptMacro(DoubleConvertibleToOutputCheck, (Concept::Convertible<double, PixelType>));
144  itkConceptMacro(InputConvertibleToOutputCheck, (Concept::Convertible<InputPixelType, PixelType>));
145  itkConceptMacro(OutputAdditiveOperatorsCheck, (Concept::AdditiveOperators<PixelType>));
146  itkConceptMacro(InputOStreamWritableCheck, (Concept::OStreamWritable<InputPixelType>));
147  itkConceptMacro(OutputOStreamWritableCheck, (Concept::OStreamWritable<PixelType>));
148  // End concept checking
149 #endif
150 
151 protected:
153  ~IsoContourDistanceImageFilter() override = default;
154  void
155  PrintSelf(std::ostream & os, Indent indent) const override;
156 
157  void
158  ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread, ThreadIdType threadId) override;
159 
160  void
162  {
163  itkExceptionMacro("This class requires threadId so it must use classic multi-threading model");
164  }
165 
166  void
167  GenerateData() override;
168 
170  ThreaderFullCallback(void * arg);
171 
172  void
173  ThreadedGenerateDataFull(const OutputImageRegionType & outputRegionForThread, ThreadIdType threadId);
174 
175  void
176  ThreadedGenerateDataBand(const OutputImageRegionType & outputRegionForThread, ThreadIdType threadId);
177 
179  void
180  BeforeThreadedGenerateData() override;
181 
182  void
183  GenerateInputRequestedRegion() override;
184 
185  void
186  EnlargeOutputRequestedRegion(DataObject *) override;
187 
190 
191  void
192  ComputeValue(const InputNeighbordIteratorType & inNeigIt,
193  OutputNeighborhoodIteratorType & outNeigIt,
194  unsigned int center,
195  const std::vector<OffsetValueType> & stride);
196 
197 private:
198  PixelRealType m_LevelSetValue{};
199  PixelType m_FarValue{};
200 
201  InputSpacingType m_Spacing{};
202 
203  bool m_NarrowBanding{};
204  NarrowBandPointer m_NarrowBand{};
205  std::vector<RegionType> m_NarrowBandRegion{};
206 
207  std::mutex m_Mutex{};
208 };
209 } // namespace itk
210 
211 #ifndef ITK_MANUAL_INSTANTIATION
212 # include "itkIsoContourDistanceImageFilter.hxx"
213 #endif
214 
215 #endif
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::Concept::OStreamWritable
Definition: itkConceptChecking.h:638
itk::IsoContourDistanceImageFilter::InputSpacingType
typename InputImageType::SpacingType InputSpacingType
Definition: itkIsoContourDistanceImageFilter.h:99
itkNeighborhoodIterator.h
itk::IsoContourDistanceImageFilter::BandIterator
typename NarrowBandType::Iterator BandIterator
Definition: itkIsoContourDistanceImageFilter.h:107
itk::IsoContourDistanceImageFilter::InputSizeType
typename InputImageType::SizeType InputSizeType
Definition: itkIsoContourDistanceImageFilter.h:93
itk::NarrowBand
Narrow Band class.
Definition: itkNarrowBand.h:51
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::IsoContourDistanceImageFilter::NarrowBandPointer
typename NarrowBandType::Pointer NarrowBandPointer
Definition: itkIsoContourDistanceImageFilter.h:104
itk::SmartPointer< Self >
itk::IsoContourDistanceImageFilter::InputIndexType
typename InputImageType::IndexType InputIndexType
Definition: itkIsoContourDistanceImageFilter.h:96
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::NarrowBand::ConstIterator
typename NodeContainerType::const_iterator ConstIterator
Definition: itkNarrowBand.h:70
itk::IsoContourDistanceImageFilter
Compute an approximate distance from an interpolated isocontour to the close grid points.
Definition: itkIsoContourDistanceImageFilter.h:59
itk::IsoContourDistanceImageFilter::GetNarrowBand
NarrowBandPointer GetNarrowBand() const
Definition: itkIsoContourDistanceImageFilter.h:133
itk::Concept::SameDimension
Definition: itkConceptChecking.h:696
itk::ThreadIdType
unsigned int ThreadIdType
Definition: itkIntTypes.h:102
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::ImageToImageFilter
Base class for filters that take an image as input and produce an image as output.
Definition: itkImageToImageFilter.h:108
itk::IsoContourDistanceImageFilter::ConstBandIterator
typename NarrowBandType::ConstIterator ConstBandIterator
Definition: itkIsoContourDistanceImageFilter.h:106
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION
itk::ITK_THREAD_RETURN_TYPE ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION
Definition: itkThreadSupport.h:89
itk::IsoContourDistanceImageFilter::SizeType
typename OutputImageType::SizeType SizeType
Definition: itkIsoContourDistanceImageFilter.h:94
itk::IsoContourDistanceImageFilter::InputPixelType
typename InputImageType::PixelType InputPixelType
Definition: itkIsoContourDistanceImageFilter.h:88
itk::IsoContourDistanceImageFilter::RegionType
typename NarrowBandType::RegionType RegionType
Definition: itkIsoContourDistanceImageFilter.h:105
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itkImageToImageFilter.h
itk::NumericTraits
Define additional traits for native types such as int or float.
Definition: itkNumericTraits.h:60
itk::NeighborhoodIterator
Defines iteration of a local N-dimensional neighborhood of pixels across an itk::Image.
Definition: itkNeighborhoodIterator.h:212
itk::ImageSource::OutputImageRegionType
typename OutputImageType::RegionType OutputImageRegionType
Definition: itkImageSource.h:92
itk::Concept::AdditiveOperators
Definition: itkConceptChecking.h:360
itkConceptMacro
#define itkConceptMacro(name, concept)
Definition: itkConceptChecking.h:65
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::ConstNeighborhoodIterator
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
Definition: itkConstNeighborhoodIterator.h:51
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:139
itk::Concept::Convertible
Definition: itkConceptChecking.h:216
itk::IsoContourDistanceImageFilter::PixelType
typename OutputImageType::PixelType PixelType
Definition: itkIsoContourDistanceImageFilter.h:87
itk::IsoContourDistanceImageFilter::IndexType
typename OutputImageType::IndexType IndexType
Definition: itkIsoContourDistanceImageFilter.h:97
itkNumericTraits.h
itk::BandNode
Definition: itkNarrowBand.h:35
itk::IsoContourDistanceImageFilter::PixelRealType
typename NumericTraits< InputPixelType >::RealType PixelRealType
Definition: itkIsoContourDistanceImageFilter.h:89
itk::NarrowBand::Iterator
typename NodeContainerType::iterator Iterator
Definition: itkNarrowBand.h:71
itk::IsoContourDistanceImageFilter::DynamicThreadedGenerateData
void DynamicThreadedGenerateData(const OutputImageRegionType &) override
Definition: itkIsoContourDistanceImageFilter.h:161
itk::Concept::EqualityComparable
Definition: itkConceptChecking.h:306
itkNarrowBand.h
itk::DataObject
Base class for all data objects in ITK.
Definition: itkDataObject.h:293