ITK  5.4.0
Insight Toolkit
itkAttributeMorphologyBaseImageFilter.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 itkAttributeMorphologyBaseImageFilter_h
19 #define itkAttributeMorphologyBaseImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 #include <vector>
23 #include <memory> // For unique_ptr.
24 
25 namespace itk
26 {
53 template <typename TInputImage, typename TOutputImage, typename TAttribute, typename TFunction>
54 class ITK_TEMPLATE_EXPORT AttributeMorphologyBaseImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
55 {
56 public:
62 
66  using typename Superclass::InputImagePointer;
67 
72  using OutputPixelType = typename TOutputImage::PixelType;
73  using OutputInternalPixelType = typename TOutputImage::InternalPixelType;
74  using InputPixelType = typename TInputImage::PixelType;
75  using InputInternalPixelType = typename TInputImage::InternalPixelType;
77  using OffsetType = typename TInputImage::OffsetType;
78  using SizeType = typename TInputImage::SizeType;
79 
80  static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
81 
85  using InputImageType = TInputImage;
86  using OutputImageType = TOutputImage;
87  // using IndexType = typename TInputImage::IndexType;
88  // using SizeType = typename TInputImage::SizeType;
90  using ListType = std::list<IndexType>;
91  using AttributeType = TAttribute;
92 
98 
102  itkOverrideGetNameOfClassMacro(AttributeMorphologyBaseImageFilter);
103 
107  itkNewMacro(Self);
108 
115  itkSetMacro(FullyConnected, bool);
116  itkGetConstReferenceMacro(FullyConnected, bool);
117  itkBooleanMacro(FullyConnected);
125  itkSetMacro(Lambda, AttributeType);
126  itkGetConstMacro(Lambda, AttributeType);
129 protected:
131  {
132  m_FullyConnected = false;
133  m_AttributeValuePerPixel = 1;
134  m_Lambda = 0;
135  }
136 
137  ~AttributeMorphologyBaseImageFilter() override = default;
139  void
140  PrintSelf(std::ostream & os, Indent indent) const override;
141 
145  void
146  GenerateData() override;
147 
151  void
152  GenerateInputRequestedRegion() override;
153 
158  void
159  EnlargeOutputRequestedRegion(DataObject * itkNotUsed(output)) override;
160 
161  AttributeType m_AttributeValuePerPixel{};
162 
163 private:
164  bool m_FullyConnected{};
165  AttributeType m_Lambda{};
166 
167  // some constants used several times in the code
168  static constexpr OffsetValueType INACTIVE = -1;
169  static constexpr OffsetValueType ACTIVE = -2;
170 
171  // Just used for area/volume openings at the moment
172  std::unique_ptr<AttributeType[]> m_AuxData;
173 
174  using OffsetVecType = std::vector<OffsetType>;
175  // offset in the linear array.
176  using OffsetDirectVecType = std::vector<OffsetValueType>;
177 
178  void
179  SetupOffsetVec(OffsetDirectVecType & PosOffsets, OffsetVecType & Offsets);
180 
181  // m_SortPixels contains offsets into the raw image
182  // it is sorted with a stable sort by grey level as the
183  // first step in the algorithm. The sorting step avoids
184  // the need to explicitly locate regional extrema.
185  std::unique_ptr<OffsetValueType[]> m_SortPixels;
186  std::unique_ptr<OffsetValueType[]> m_Parent;
187 
188  // This is a bit ugly, but I can't see an easy way around
189  std::unique_ptr<InputPixelType[]> m_Raw;
190 
192  {
193  public:
194  TFunction m_TFunction;
195  // buf contains the raw data, which is what
196  // we want to sort by. i.e. the first value in
197  // the sorted buffer will be the location of the
198  // largest or smallest pixel.
200  bool
201  operator()(OffsetValueType const & l, OffsetValueType const & r) const
202  {
203  return (m_TFunction(buf[l], buf[r]));
204  }
205  };
206 
207  CompareOffsetType m_CompareOffset{};
208  // version from PAMI. Note - using the AuxData array rather than the
209  // parent array to store area
210  void
212  {
213  m_Parent[x] = ACTIVE;
214  m_AuxData[x] = m_AttributeValuePerPixel;
215  }
216 
219  {
220  if (m_Parent[x] >= 0)
221  {
222  m_Parent[x] = FindRoot(m_Parent[x]);
223  return (m_Parent[x]);
224  }
225  else
226  {
227  return (x);
228  }
229  }
230 
231  bool
233  {
234  return ((m_Raw[x] == m_Raw[y]) || (m_AuxData[x] < m_Lambda));
235  }
236 
237  void
239  {
240  OffsetValueType r = FindRoot(n);
241 
242  if (r != p)
243  {
244  if (Criterion(r, p))
245  {
246  m_AuxData[p] += m_AuxData[r];
247  m_Parent[r] = p;
248  }
249  else
250  {
251  m_AuxData[p] = m_Lambda;
252  }
253  }
254  }
255 };
256 } // end namespace itk
257 
258 #ifndef ITK_MANUAL_INSTANTIATION
259 # include "itkAttributeMorphologyBaseImageFilter.hxx"
260 #endif
261 
262 #endif
itk::AttributeMorphologyBaseImageFilter::m_AuxData
std::unique_ptr< AttributeType[]> m_AuxData
Definition: itkAttributeMorphologyBaseImageFilter.h:172
itk::AttributeMorphologyBaseImageFilter< TInputImage, TOutputImage, TAttribute, std::greater< TInputImage::PixelType > >::OutputPixelType
typename TOutputImage::PixelType OutputPixelType
Definition: itkAttributeMorphologyBaseImageFilter.h:72
itk::AttributeMorphologyBaseImageFilter< TInputImage, TOutputImage, TAttribute, std::greater< TInputImage::PixelType > >::OffsetDirectVecType
std::vector< OffsetValueType > OffsetDirectVecType
Definition: itkAttributeMorphologyBaseImageFilter.h:176
itk::AttributeMorphologyBaseImageFilter< TInputImage, TOutputImage, TAttribute, std::greater< TInputImage::PixelType > >::OffsetType
typename TInputImage::OffsetType OffsetType
Definition: itkAttributeMorphologyBaseImageFilter.h:77
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::AttributeMorphologyBaseImageFilter< TInputImage, TOutputImage, TAttribute, std::greater< TInputImage::PixelType > >::OffsetVecType
std::vector< OffsetType > OffsetVecType
Definition: itkAttributeMorphologyBaseImageFilter.h:174
itk::AttributeMorphologyBaseImageFilter::CompareOffsetType::buf
InputPixelType * buf
Definition: itkAttributeMorphologyBaseImageFilter.h:199
itk::AttributeMorphologyBaseImageFilter::Criterion
bool Criterion(OffsetValueType x, OffsetValueType y)
Definition: itkAttributeMorphologyBaseImageFilter.h:232
itk::AttributeMorphologyBaseImageFilter::m_SortPixels
std::unique_ptr< OffsetValueType[]> m_SortPixels
Definition: itkAttributeMorphologyBaseImageFilter.h:185
itk::AttributeMorphologyBaseImageFilter< TInputImage, TOutputImage, TAttribute, std::greater< TInputImage::PixelType > >::InputPixelType
typename TInputImage::PixelType InputPixelType
Definition: itkAttributeMorphologyBaseImageFilter.h:74
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::AttributeMorphologyBaseImageFilter::AttributeMorphologyBaseImageFilter
AttributeMorphologyBaseImageFilter(const Self &)
Definition: itkAttributeMorphologyBaseImageFilter.h:138
itk::AttributeMorphologyBaseImageFilter< TInputImage, TOutputImage, TAttribute, std::greater< TInputImage::PixelType > >::InputInternalPixelType
typename TInputImage::InternalPixelType InputInternalPixelType
Definition: itkAttributeMorphologyBaseImageFilter.h:75
itk::AttributeMorphologyBaseImageFilter< TInputImage, TOutputImage, TAttribute, std::greater< TInputImage::PixelType > >::SizeType
typename TInputImage::SizeType SizeType
Definition: itkAttributeMorphologyBaseImageFilter.h:78
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::AttributeMorphologyBaseImageFilter< TInputImage, TOutputImage, TAttribute, std::greater< TInputImage::PixelType > >::IndexType
typename TInputImage::IndexType IndexType
Definition: itkAttributeMorphologyBaseImageFilter.h:76
itk::AttributeMorphologyBaseImageFilter::CompareOffsetType::m_TFunction
TFunction m_TFunction
Definition: itkAttributeMorphologyBaseImageFilter.h:194
itk::AttributeMorphologyBaseImageFilter
Morphological opening by attributes.
Definition: itkAttributeMorphologyBaseImageFilter.h:54
itk::ImageToImageFilter::InputImageType
TInputImage InputImageType
Definition: itkImageToImageFilter.h:129
itk::AttributeMorphologyBaseImageFilter::Union
void Union(OffsetValueType n, OffsetValueType p)
Definition: itkAttributeMorphologyBaseImageFilter.h:238
itk::AttributeMorphologyBaseImageFilter::MakeSet
void MakeSet(OffsetValueType x)
Definition: itkAttributeMorphologyBaseImageFilter.h:211
itkImageToImageFilter.h
itk::OffsetValueType
long OffsetValueType
Definition: itkIntTypes.h:94
itk::AttributeMorphologyBaseImageFilter< TInputImage, TOutputImage, TAttribute, std::greater< TInputImage::PixelType > >::RegionType
typename TOutputImage::RegionType RegionType
Definition: itkAttributeMorphologyBaseImageFilter.h:89
itk::AttributeMorphologyBaseImageFilter::FindRoot
OffsetValueType FindRoot(OffsetValueType x)
Definition: itkAttributeMorphologyBaseImageFilter.h:218
itk::AttributeMorphologyBaseImageFilter< TInputImage, TOutputImage, TAttribute, std::greater< TInputImage::PixelType > >::AttributeType
TAttribute AttributeType
Definition: itkAttributeMorphologyBaseImageFilter.h:91
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:139
itk::AttributeMorphologyBaseImageFilter< TInputImage, TOutputImage, TAttribute, std::greater< TInputImage::PixelType > >::OutputInternalPixelType
typename TOutputImage::InternalPixelType OutputInternalPixelType
Definition: itkAttributeMorphologyBaseImageFilter.h:73
itk::AttributeMorphologyBaseImageFilter::m_Raw
std::unique_ptr< InputPixelType[]> m_Raw
Definition: itkAttributeMorphologyBaseImageFilter.h:189
itk::AttributeMorphologyBaseImageFilter::CompareOffsetType
Definition: itkAttributeMorphologyBaseImageFilter.h:191
itk::AttributeMorphologyBaseImageFilter< TInputImage, TOutputImage, TAttribute, std::greater< TInputImage::PixelType > >::ListType
std::list< IndexType > ListType
Definition: itkAttributeMorphologyBaseImageFilter.h:90
itk::AttributeMorphologyBaseImageFilter::CompareOffsetType::operator()
bool operator()(OffsetValueType const &l, OffsetValueType const &r) const
Definition: itkAttributeMorphologyBaseImageFilter.h:201
itk::ImageSource::OutputImageType
TOutputImage OutputImageType
Definition: itkImageSource.h:90
itk::DataObject
Base class for all data objects in ITK.
Definition: itkDataObject.h:293
itk::AttributeMorphologyBaseImageFilter::m_Parent
std::unique_ptr< OffsetValueType[]> m_Parent
Definition: itkAttributeMorphologyBaseImageFilter.h:186
itk::AttributeMorphologyBaseImageFilter::AttributeMorphologyBaseImageFilter
AttributeMorphologyBaseImageFilter()
Definition: itkAttributeMorphologyBaseImageFilter.h:130