ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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
24#include <stack>
25
26namespace itk
27{
75
76template <typename TInputImage, typename TOutputImage, typename TFunction1, typename TFunction2>
77class ITK_TEMPLATE_EXPORT ValuedRegionalExtremaImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
78{
79public:
80 ITK_DISALLOW_COPY_AND_MOVE(ValuedRegionalExtremaImageFilter);
81
87
89 using InputImageType = TInputImage;
90 using OutputImageType = TOutputImage;
91 using InputImagePointer = typename InputImageType::Pointer;
92 using InputImageConstPointer = typename InputImageType::ConstPointer;
93 using InputImageRegionType = typename InputImageType::RegionType;
94 using InputImagePixelType = typename InputImageType::PixelType;
95 using ISizeType = typename InputImageType::SizeType;
96 using OutputImagePointer = typename OutputImageType::Pointer;
97 using OutputImageConstPointer = typename OutputImageType::ConstPointer;
98 using OutputImageRegionType = typename OutputImageType::RegionType;
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);
121
125 itkSetMacro(MarkerValue, typename TInputImage::PixelType);
126 itkGetConstReferenceMacro(MarkerValue, typename TInputImage::PixelType);
128
132 itkGetConstMacro(Flat, bool);
133
136
137protected:
140 void
141 PrintSelf(std::ostream & os, Indent indent) const override;
142
146 void
148
150 void
151 EnlargeOutputRequestedRegion(DataObject * itkNotUsed(output)) override;
152
153 void
154 GenerateData() override;
155
156private:
157 typename TInputImage::PixelType m_MarkerValue{};
158
159 bool m_FullyConnected{ false };
160 bool m_Flat{ false };
161
162 using OutIndexType = typename OutputImageType::IndexType;
163 using InIndexType = typename InputImageType::IndexType;
166 using IndexStack = std::stack<OutIndexType>;
167}; // end of class
168} // end namespace itk
169
170#ifndef ITK_MANUAL_INSTANTIATION
171# include "itkValuedRegionalExtremaImageFilter.hxx"
172#endif
173
174#endif
Const version of ShapedNeighborhoodIterator, defining iteration of a local N-dimensional neighborhood...
Base class for all data objects in ITK.
Control indentation during Print() invocation.
Definition itkIndent.h:50
A neighborhood iterator which can take on an arbitrary shape.
Implements transparent reference counting.
typename OutputImageType::ConstPointer OutputImageConstPointer
typename InputImageType::RegionType InputImageRegionType
ShapedNeighborhoodIterator< OutputImageType > NOutputIterator
typename OutputImageType::PixelType OutputImagePixelType
ConstShapedNeighborhoodIterator< InputImageType > ConstInputIterator
typename InputImageType::ConstPointer InputImageConstPointer
ImageToImageFilter< TInputImage, TOutputImage > Superclass
typename OutputImageType::RegionType OutputImageRegionType
~ValuedRegionalExtremaImageFilter() override=default
void PrintSelf(std::ostream &os, Indent indent) const override
void EnlargeOutputRequestedRegion(DataObject *output) override
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....