ITK  6.0.0
Insight Toolkit
itkHMinimaImageFilter.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 itkHMinimaImageFilter_h
19#define itkHMinimaImageFilter_h
20
22
23namespace itk
24{
54template <typename TInputImage, typename TOutputImage>
55class ITK_TEMPLATE_EXPORT HMinimaImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
56{
57public:
58 ITK_DISALLOW_COPY_AND_MOVE(HMinimaImageFilter);
59
65
67 using InputImageType = TInputImage;
71 using InputImagePixelType = typename InputImageType::PixelType;
72 using OutputImageType = TOutputImage;
76 using OutputImagePixelType = typename OutputImageType::PixelType;
77
79 static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
80 static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
81
83 itkNewMacro(Self);
84
86 itkOverrideGetNameOfClassMacro(HMinimaImageFilter);
87
92 itkSetMacro(Height, InputImagePixelType);
93 itkGetConstMacro(Height, InputImagePixelType);
102 itkSetMacro(FullyConnected, bool);
103 itkGetConstReferenceMacro(FullyConnected, bool);
104 itkBooleanMacro(FullyConnected);
107#ifdef ITK_USE_CONCEPT_CHECKING
108 // Begin concept checking
112 // End concept checking
113#endif
114
115protected:
117 ~HMinimaImageFilter() override = default;
118 void
119 PrintSelf(std::ostream & os, Indent indent) const override;
120
124 void
126
128 void
129 EnlargeOutputRequestedRegion(DataObject * itkNotUsed(output)) override;
130
133 void
134 GenerateData() override;
135
136private:
138 unsigned long m_NumberOfIterationsUsed{ 1 };
139 bool m_FullyConnected{ false };
140}; // end of class
141} // end namespace itk
142
143#ifndef ITK_MANUAL_INSTANTIATION
144# include "itkHMinimaImageFilter.hxx"
145#endif
146
147#endif
Base class for all data objects in ITK.
Suppress local minima whose depth below the baseline is less than h.
void EnlargeOutputRequestedRegion(DataObject *) override
void PrintSelf(std::ostream &os, Indent indent) const override
void GenerateInputRequestedRegion() override
typename OutputImageType::ConstPointer OutputImageConstPointer
void GenerateData() override
~HMinimaImageFilter() override=default
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
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....