ITK  6.0.0
Insight Toolkit
itkHConcaveImageFilter.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 itkHConcaveImageFilter_h
19 #define itkHConcaveImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 
23 namespace itk
24 {
44 template <typename TInputImage, typename TOutputImage>
45 class ITK_TEMPLATE_EXPORT HConcaveImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
46 {
47 public:
48  ITK_DISALLOW_COPY_AND_MOVE(HConcaveImageFilter);
49 
55 
57  using InputImageType = TInputImage;
61  using InputImagePixelType = typename InputImageType::PixelType;
62  using OutputImageType = TOutputImage;
66  using OutputImagePixelType = typename OutputImageType::PixelType;
67 
69  static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
70  static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
71 
73  itkNewMacro(Self);
74 
76  itkOverrideGetNameOfClassMacro(HConcaveImageFilter);
77 
82  itkSetMacro(Height, InputImagePixelType);
83  itkGetConstMacro(Height, InputImagePixelType);
92  itkSetMacro(FullyConnected, bool);
93  itkGetConstReferenceMacro(FullyConnected, bool);
94  itkBooleanMacro(FullyConnected);
97 #ifdef ITK_USE_CONCEPT_CHECKING
98  // Begin concept checking
102  // End concept checking
103 #endif
104 
105 protected:
107  ~HConcaveImageFilter() override = default;
108  void
109  PrintSelf(std::ostream & os, Indent indent) const override;
110 
114  void
115  GenerateInputRequestedRegion() override;
116 
118  void
119  EnlargeOutputRequestedRegion(DataObject * itkNotUsed(output)) override;
120 
123  void
124  GenerateData() override;
125 
126 private:
128  unsigned long m_NumberOfIterationsUsed{ 1 };
129  bool m_FullyConnected{ false };
130 }; // end of class
131 } // end namespace itk
132 
133 #ifndef ITK_MANUAL_INSTANTIATION
134 # include "itkHConcaveImageFilter.hxx"
135 #endif
136 
137 #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::HConcaveImageFilter::OutputImageConstPointer
typename OutputImageType::ConstPointer OutputImageConstPointer
Definition: itkHConcaveImageFilter.h:64
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::HConcaveImageFilter
Identify local minima whose depth below the baseline is greater than h.
Definition: itkHConcaveImageFilter.h:45
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::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::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