ITK  6.0.0
Insight Toolkit
itkValuedRegionalExtremaImageFilter.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 itkValuedRegionalExtremaImageFilter_h
19 #define itkValuedRegionalExtremaImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
24 #include <stack>
25 
26 namespace itk
27 {
76 template <typename TInputImage, typename TOutputImage, typename TFunction1, typename TFunction2>
77 class ITK_TEMPLATE_EXPORT ValuedRegionalExtremaImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
78 {
79 public:
80  ITK_DISALLOW_COPY_AND_MOVE(ValuedRegionalExtremaImageFilter);
81 
87 
89  using InputImageType = TInputImage;
90  using OutputImageType = TOutputImage;
94  using InputImagePixelType = typename InputImageType::PixelType;
99  using OutputImagePixelType = typename OutputImageType::PixelType;
100 
102  static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
103  static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
104 
106  itkNewMacro(Self);
107 
109  itkOverrideGetNameOfClassMacro(ValuedRegionalExtremaImageFilter);
110 
117  itkSetMacro(FullyConnected, bool);
118  itkGetConstReferenceMacro(FullyConnected, bool);
119  itkBooleanMacro(FullyConnected);
125  itkSetMacro(MarkerValue, typename TInputImage::PixelType);
126  itkGetConstReferenceMacro(MarkerValue, typename TInputImage::PixelType);
132  itkGetConstMacro(Flat, bool);
133 
134 #ifdef ITK_USE_CONCEPT_CHECKING
135  // Begin concept checking
138  // End concept checking
139 #endif
140 
141 protected:
143  ~ValuedRegionalExtremaImageFilter() override = default;
144  void
145  PrintSelf(std::ostream & os, Indent indent) const override;
146 
150  void
151  GenerateInputRequestedRegion() override;
152 
154  void
155  EnlargeOutputRequestedRegion(DataObject * itkNotUsed(output)) override;
156 
157  void
158  GenerateData() override;
159 
160 private:
161  typename TInputImage::PixelType m_MarkerValue{};
162 
163  bool m_FullyConnected{ false };
164  bool m_Flat{ false };
165 
170  using IndexStack = std::stack<OutIndexType>;
171 }; // end of class
172 } // end namespace itk
173 
174 #ifndef ITK_MANUAL_INSTANTIATION
175 # include "itkValuedRegionalExtremaImageFilter.hxx"
176 #endif
177 
178 #endif
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::ValuedRegionalExtremaImageFilter< TInputImage, TOutputImage, std::less< TInputImage::PixelType >, std::less< TOutputImage::PixelType > >::ISizeType
typename InputImageType::SizeType ISizeType
Definition: itkValuedRegionalExtremaImageFilter.h:95
ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAddImageFilter.h:94
itk::Concept::HasNumericTraits
Definition: itkConceptChecking.h:716
itk::ImageSource::OutputImagePointer
typename OutputImageType::Pointer OutputImagePointer
Definition: itkImageSource.h:91
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::ShapedNeighborhoodIterator
A neighborhood iterator which can take on an arbitrary shape.
Definition: itkShapedNeighborhoodIterator.h:150
itk::ValuedRegionalExtremaImageFilter< TInputImage, TOutputImage, std::less< TInputImage::PixelType >, std::less< TOutputImage::PixelType > >::IndexStack
std::stack< OutIndexType > IndexStack
Definition: itkValuedRegionalExtremaImageFilter.h:170
itkConstantBoundaryCondition.h
itk::ImageToImageFilter::InputImagePixelType
typename InputImageType::PixelType InputImagePixelType
Definition: itkImageToImageFilter.h:133
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
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::ValuedRegionalExtremaImageFilter
Uses a flooding algorithm to set all voxels that are not a regional extrema to the max or min of the ...
Definition: itkValuedRegionalExtremaImageFilter.h:77
itk::ConstShapedNeighborhoodIterator
Const version of ShapedNeighborhoodIterator, defining iteration of a local N-dimensional neighborhood...
Definition: itkConstShapedNeighborhoodIterator.h:72
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
itkShapedNeighborhoodIterator.h
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::ValuedRegionalExtremaImageFilter< TInputImage, TOutputImage, std::less< TInputImage::PixelType >, std::less< TOutputImage::PixelType > >::InIndexType
typename InputImageType::IndexType InIndexType
Definition: itkValuedRegionalExtremaImageFilter.h:167
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::HasPixelTraits
Definition: itkConceptChecking.h:753
itk::ValuedRegionalExtremaImageFilter< TInputImage, TOutputImage, std::less< TInputImage::PixelType >, std::less< TOutputImage::PixelType > >::OutIndexType
typename OutputImageType::IndexType OutIndexType
Definition: itkValuedRegionalExtremaImageFilter.h:166
itk::ImageToImageFilter::InputImageRegionType
typename InputImageType::RegionType InputImageRegionType
Definition: itkImageToImageFilter.h:132
itk::ImageToImageFilter::InputImageConstPointer
typename InputImageType::ConstPointer InputImageConstPointer
Definition: itkImageToImageFilter.h:131
itk::ValuedRegionalExtremaImageFilter< TInputImage, TOutputImage, std::less< TInputImage::PixelType >, std::less< TOutputImage::PixelType > >::OutputImageConstPointer
typename OutputImageType::ConstPointer OutputImageConstPointer
Definition: itkValuedRegionalExtremaImageFilter.h:97
itk::ImageSource::OutputImageType
TOutputImage OutputImageType
Definition: itkImageSource.h:90
itk::DataObject
Base class for all data objects in ITK.
Definition: itkDataObject.h:293