ITK  6.0.0
Insight Toolkit
itkIsolatedConnectedImageFilter.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 itkIsolatedConnectedImageFilter_h
19#define itkIsolatedConnectedImageFilter_h
20
22
23namespace itk
24{
71template <typename TInputImage, typename TOutputImage>
72class ITK_TEMPLATE_EXPORT IsolatedConnectedImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
73{
74public:
75 ITK_DISALLOW_COPY_AND_MOVE(IsolatedConnectedImageFilter);
76
82
84 itkNewMacro(Self);
85
87 itkOverrideGetNameOfClassMacro(IsolatedConnectedImageFilter);
88 using InputImageType = TInputImage;
92 using InputImagePixelType = typename InputImageType::PixelType;
97 using OutputImageType = TOutputImage;
100 using OutputImagePixelType = typename OutputImageType::PixelType;
101
102 using SeedsContainerType = std::vector<IndexType>;
103
105
106 void
107 PrintSelf(std::ostream & os, Indent indent) const override;
108
111 void
112 AddSeed1(const IndexType & seed);
113
121 void
122 SetSeed1(const IndexType & seed);
123
125 void
127
129 void
130 AddSeed2(const IndexType & seed);
131
138 void
139 SetSeed2(const IndexType & seed);
140
142 void
144
146 virtual const SeedsContainerType &
147 GetSeeds1() const;
148 virtual const SeedsContainerType &
149 GetSeeds2() const;
154 itkSetMacro(Lower, InputImagePixelType);
155 itkGetConstReferenceMacro(Lower, InputImagePixelType);
160 itkSetMacro(Upper, InputImagePixelType);
161 itkGetConstReferenceMacro(Upper, InputImagePixelType);
164#if !defined(ITK_LEGACY_REMOVE)
169 void
170 SetUpperValueLimit(InputImagePixelType upperValue)
171 {
172 this->SetUpper(upperValue);
173 }
174
175 InputImagePixelType
176 GetUpperValueLimit()
177 {
178 return this->GetUpper();
179 }
180#endif
181
184 itkSetMacro(IsolatedValueTolerance, InputImagePixelType);
185 itkGetConstReferenceMacro(IsolatedValueTolerance, InputImagePixelType);
191 itkSetMacro(ReplaceValue, OutputImagePixelType);
192 itkGetConstReferenceMacro(ReplaceValue, OutputImagePixelType);
196 itkGetConstReferenceMacro(IsolatedValue, InputImagePixelType);
197
200 itkSetMacro(FindUpperThreshold, bool);
201 itkBooleanMacro(FindUpperThreshold);
202 itkGetConstReferenceMacro(FindUpperThreshold, bool);
207 itkGetConstReferenceMacro(ThresholdingFailed, bool);
208
209#ifdef ITK_USE_CONCEPT_CHECKING
210 // Begin concept checking
212 // End concept checking
213#endif
214
215protected:
217 ~IsolatedConnectedImageFilter() override = default;
220
223
224 OutputImagePixelType m_ReplaceValue{};
225
226 InputImagePixelType m_IsolatedValue{};
227 InputImagePixelType m_IsolatedValueTolerance{};
228
229 bool m_FindUpperThreshold{};
230 bool m_ThresholdingFailed{};
231
232 // Override since the filter needs all the data for the algorithm
233 void
235
236 // Override since the filter produces the entire dataset
237 void
239
240 void
241 GenerateData() override;
242};
243} // end namespace itk
244
245#ifndef ITK_MANUAL_INSTANTIATION
246# include "itkIsolatedConnectedImageFilter.hxx"
247#endif
248
249#endif
Base class for all data objects in ITK.
Base class for all process objects that output image data.
typename OutputImageType::PixelType OutputImagePixelType
typename OutputImageType::RegionType OutputImageRegionType
TOutputImage OutputImageType
typename OutputImageType::Pointer OutputImagePointer
Base class for filters that take an image as input and produce an image as output.
typename InputImageType::ConstPointer InputImageConstPointer
typename InputImageType::Pointer InputImagePointer
typename InputImageType::PixelType InputImagePixelType
typename InputImageType::RegionType InputImageRegionType
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Label pixels that are connected to one set of seeds but not another.
typename NumericTraits< InputImagePixelType >::RealType InputRealType
void AddSeed1(const IndexType &seed)
typename InputImageType::IndexType IndexType
void EnlargeOutputRequestedRegion(DataObject *output) override
virtual const SeedsContainerType & GetSeeds1() const
void SetSeed1(const IndexType &seed)
void SetSeed2(const IndexType &seed)
virtual const SeedsContainerType & GetSeeds2() const
void GenerateInputRequestedRegion() override
typename InputImageType::SizeType SizeType
void PrintSelf(std::ostream &os, Indent indent) const override
void AddSeed2(const IndexType &seed)
~IsolatedConnectedImageFilter() override=default
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< const Self > ConstPointer
SmartPointer< Self > Pointer
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....