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 
21 #include "itkImageToImageFilter.h"
22 
23 namespace itk
24 {
54 template <typename TInputImage, typename TOutputImage>
55 class ITK_TEMPLATE_EXPORT HMinimaImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
56 {
57 public:
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
109  itkConceptMacro(InputEqualityComparableCheck, (Concept::EqualityComparable<InputImagePixelType>));
112  // End concept checking
113 #endif
114 
115 protected:
117  ~HMinimaImageFilter() override = default;
118  void
119  PrintSelf(std::ostream & os, Indent indent) const override;
120 
124  void
125  GenerateInputRequestedRegion() override;
126 
128  void
129  EnlargeOutputRequestedRegion(DataObject * itkNotUsed(output)) override;
130 
133  void
134  GenerateData() override;
135 
136 private:
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
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAddImageFilter.h:94
itk::Concept::OStreamWritable
Definition: itkConceptChecking.h:638
itk::ImageSource::OutputImagePointer
typename OutputImageType::Pointer OutputImagePointer
Definition: itkImageSource.h:91
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::ImageToImageFilter::InputImagePixelType
typename InputImageType::PixelType InputImagePixelType
Definition: itkImageToImageFilter.h:133
itk::ImageToImageFilter
Base class for filters that take an image as input and produce an image as output.
Definition: itkImageToImageFilter.h:108
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::ImageToImageFilter::InputImagePointer
typename InputImageType::Pointer InputImagePointer
Definition: itkImageToImageFilter.h:130
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::ImageToImageFilter::InputImageType
TInputImage InputImageType
Definition: itkImageToImageFilter.h:129
itkImageToImageFilter.h
itk::HMinimaImageFilter
Suppress local minima whose depth below the baseline is less than h.
Definition: itkHMinimaImageFilter.h:55
itk::ImageSource::OutputImageRegionType
typename OutputImageType::RegionType OutputImageRegionType
Definition: itkImageSource.h:92
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::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:139
itk::ImageSource::OutputImagePixelType
typename OutputImageType::PixelType OutputImagePixelType
Definition: itkImageSource.h:93
itk::Concept::Convertible
Definition: itkConceptChecking.h:216
itk::ImageToImageFilter::InputImageRegionType
typename InputImageType::RegionType InputImageRegionType
Definition: itkImageToImageFilter.h:132
itk::HMinimaImageFilter::OutputImageConstPointer
typename OutputImageType::ConstPointer OutputImageConstPointer
Definition: itkHMinimaImageFilter.h:74
itk::Concept::EqualityComparable
Definition: itkConceptChecking.h:306
itk::ImageToImageFilter::InputImageConstPointer
typename InputImageType::ConstPointer InputImageConstPointer
Definition: itkImageToImageFilter.h:131
itk::ImageSource::OutputImageType
TOutputImage OutputImageType
Definition: itkImageSource.h:90
itk::DataObject
Base class for all data objects in ITK.
Definition: itkDataObject.h:293