ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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{
60template <typename TInputImage, typename TOutputImage>
61class ITK_TEMPLATE_EXPORT IsoContourDistanceImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
62{
63public:
64 ITK_DISALLOW_COPY_AND_MOVE(IsoContourDistanceImageFilter);
65
71
73 itkNewMacro(Self);
74
76 itkOverrideGetNameOfClassMacro(IsoContourDistanceImageFilter);
77
79 using typename Superclass::InputImageType;
80 using typename Superclass::OutputImageType;
81
84 static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
85 static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
86
89 using PixelType = typename OutputImageType::PixelType;
90 using InputPixelType = typename InputImageType::PixelType;
92
93 using OutputImageRegionType = typename OutputImageType::RegionType;
94
95 using InputSizeType = typename InputImageType::SizeType;
96 using SizeType = typename OutputImageType::SizeType;
97
98 using InputIndexType = typename InputImageType::IndexType;
99 using IndexType = typename OutputImageType::IndexType;
100
101 using InputSpacingType = typename InputImageType::SpacingType;
102
107 using RegionType = typename NarrowBandType::RegionType;
110
114 itkSetMacro(LevelSetValue, PixelRealType);
115 itkGetConstMacro(LevelSetValue, PixelRealType);
120 itkSetMacro(FarValue, PixelType);
121 itkGetConstMacro(FarValue, PixelType);
126 itkSetMacro(NarrowBanding, bool);
127 itkGetConstMacro(NarrowBanding, bool);
128 itkBooleanMacro(NarrowBanding);
131 void
133
136 {
137 return m_NarrowBand;
138 }
139
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
149protected:
151 ~IsoContourDistanceImageFilter() override = default;
152 void
153 PrintSelf(std::ostream & os, Indent indent) const override;
154
155 void
156 ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread, ThreadIdType threadId) override;
157
158 void
160 {
161 itkExceptionMacro("This class requires threadId so it must use classic multi-threading model");
162 }
163
164 void
165 GenerateData() override;
166
169
170 void
171 ThreadedGenerateDataFull(const OutputImageRegionType & outputRegionForThread, ThreadIdType threadId);
172
173 void
174 ThreadedGenerateDataBand(const OutputImageRegionType & outputRegionForThread, ThreadIdType threadId);
175
177 void
179
180 void
182
183 void
185
188
189 void
192 unsigned int center,
193 const std::vector<OffsetValueType> & stride);
194
195private:
198
200
203 std::vector<RegionType> m_NarrowBandRegion{};
204
205 std::mutex m_Mutex{};
206};
207} // namespace itk
208
209#ifndef ITK_MANUAL_INSTANTIATION
210# include "itkIsoContourDistanceImageFilter.hxx"
211#endif
212
213#endif
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
Base class for all data objects in ITK.
Control indentation during Print() invocation.
Definition itkIndent.h:50
void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId) override
void BeforeThreadedGenerateData() override
~IsoContourDistanceImageFilter() override=default
static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION ThreaderFullCallback(void *arg)
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
void SetNarrowBand(NarrowBandType *ptr)
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.
typename NodeContainerType::const_iterator ConstIterator
SmartPointer< Self > Pointer
typename NodeContainerType::iterator Iterator
Defines iteration of a local N-dimensional neighborhood of pixels across an itk::Image.
Implements transparent reference counting.
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned int ThreadIdType
itk::ITK_THREAD_RETURN_TYPE ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION