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
22#include "itkNarrowBand.h"
24#include "itkNumericTraits.h"
25#include <mutex>
26
27namespace itk
28{
58template <typename TInputImage, typename TOutputImage>
59class ITK_TEMPLATE_EXPORT IsoContourDistanceImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
60{
61public:
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
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>));
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
151protected:
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
171
172 void
173 ThreadedGenerateDataFull(const OutputImageRegionType & outputRegionForThread, ThreadIdType threadId);
174
175 void
176 ThreadedGenerateDataBand(const OutputImageRegionType & outputRegionForThread, ThreadIdType threadId);
177
179 void
181
182 void
184
185 void
187
190
191 void
194 unsigned int center,
195 const std::vector<OffsetValueType> & stride);
196
197private:
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
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
Base class for all data objects in ITK.
Base class for all process objects that output image data.
typename OutputImageType::RegionType OutputImageRegionType
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Compute an approximate distance from an interpolated isocontour to the close grid points.
void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId) override
typename NarrowBandType::Iterator BandIterator
void BeforeThreadedGenerateData() override
~IsoContourDistanceImageFilter() override=default
typename InputImageType::SizeType InputSizeType
typename NarrowBandType::RegionType RegionType
typename OutputImageType::SizeType SizeType
static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION ThreaderFullCallback(void *arg)
typename OutputImageType::IndexType IndexType
typename InputImageType::SpacingType InputSpacingType
typename NarrowBandType::ConstIterator ConstBandIterator
void ThreadedGenerateDataBand(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId)
void DynamicThreadedGenerateData(const OutputImageRegionType &) override
void ThreadedGenerateDataFull(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId)
void PrintSelf(std::ostream &os, Indent indent) const override
typename NumericTraits< InputPixelType >::RealType PixelRealType
typename NarrowBandType::Pointer NarrowBandPointer
void SetNarrowBand(NarrowBandType *ptr)
typename InputImageType::PixelType InputPixelType
typename InputImageType::IndexType InputIndexType
typename OutputImageType::PixelType PixelType
void EnlargeOutputRequestedRegion(DataObject *) override
void GenerateInputRequestedRegion() override
void ComputeValue(const InputNeighbordIteratorType &inNeigIt, OutputNeighborhoodIteratorType &outNeigIt, unsigned int center, const std::vector< OffsetValueType > &stride)
Narrow Band class.
Definition: itkNarrowBand.h:52
typename NodeContainerType::const_iterator ConstIterator
Definition: itkNarrowBand.h:70
typename NodeContainerType::iterator Iterator
Definition: itkNarrowBand.h:71
Defines iteration of a local N-dimensional neighborhood of pixels across an itk::Image.
Define additional traits for native types such as int or float.
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
SmartPointer< Self > Pointer
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned int ThreadIdType
Definition: itkIntTypes.h:102
itk::ITK_THREAD_RETURN_TYPE ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION