ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkAttributeMorphologyBaseImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
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  * http://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 
24 #define PAMI
25 
26 namespace itk
27 {
62 template< typename TInputImage, typename TOutputImage, typename TAttribute, typename TFunction >
63 class ITK_TEMPLATE_EXPORT AttributeMorphologyBaseImageFilter:
64  public ImageToImageFilter< TInputImage, TOutputImage >
65 {
66 public:
72 
76  using InputImagePointer = typename Superclass::InputImagePointer;
77 
82  using OutputPixelType = typename TOutputImage::PixelType;
83  using OutputInternalPixelType = typename TOutputImage::InternalPixelType;
84  using InputPixelType = typename TInputImage::PixelType;
85  using InputInternalPixelType = typename TInputImage::InternalPixelType;
87  using OffsetType = typename TInputImage::OffsetType;
88  using SizeType = typename TInputImage::SizeType;
89 
90  static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
91 
95  using InputImageType = TInputImage;
96  using OutputImageType = TOutputImage;
97 // using IndexType = typename TInputImage::IndexType;
98 // using SizeType = typename TInputImage::SizeType;
100  using ListType = std::list< IndexType >;
101  using AttributeType = TAttribute;
102 
108 
113 
117  itkNewMacro(Self);
118 
125  itkSetMacro(FullyConnected, bool);
126  itkGetConstReferenceMacro(FullyConnected, bool);
127  itkBooleanMacro(FullyConnected);
129 
135  itkSetMacro(Lambda, AttributeType);
136  itkGetConstMacro(Lambda, AttributeType);
138 
139 protected:
141  {
142  m_FullyConnected = false;
143  m_AttributeValuePerPixel = 1;
144  m_Lambda = 0;
145  }
146 
149  void PrintSelf(std::ostream & os, Indent indent) const override;
150 
154  void GenerateData() override;
155 
159  void GenerateInputRequestedRegion() override;
160 
165  void EnlargeOutputRequestedRegion( DataObject * itkNotUsed(output) ) override;
166 
168 
169 private:
170 
173 
174  // some constants used several times in the code
175  static constexpr OffsetValueType INACTIVE = -1;
176  static constexpr OffsetValueType ACTIVE = -2;
177  static constexpr OffsetValueType ROOT = -3;
178 
179  // Just used for area/volume openings at the moment
181 
182  using OffsetVecType = std::vector< OffsetType >;
183  // offset in the linear array.
184  using OffsetDirectVecType = std::vector< OffsetValueType >;
185 
186  void SetupOffsetVec(OffsetDirectVecType & PosOffsets, OffsetVecType & Offsets);
187 
189  {
190 public:
193  };
194 
195  GreyAndPos * m_SortPixels;
197 #ifndef PAMI
198  bool *m_Processed;
199 #endif
200  // This is a bit ugly, but I can't see an easy way around
202 
204  {
205 public:
206  TFunction m_TFunction;
207  bool operator()(GreyAndPos const & l, GreyAndPos const & r) const
208  {
209  if ( m_TFunction(l.Val, r.Val) )
210  {
211  return true;
212  }
213  if ( l.Val == r.Val )
214  {
215  return ( l.Pos < r.Pos );
216  }
217  return false;
218  }
219  };
220 
221 #ifdef PAMI
222  // version from PAMI. Note - using the AuxData array rather than the
223  // parent array to store area
225  {
226  m_Parent[x] = ACTIVE;
227  m_AuxData[x] = m_AttributeValuePerPixel;
228  }
229 
231  {
232  if ( m_Parent[x] >= 0 )
233  {
234  m_Parent[x] = FindRoot(m_Parent[x]);
235  return ( m_Parent[x] );
236  }
237  else
238  {
239  return ( x );
240  }
241  }
242 
244  {
245  return ( ( m_Raw[x] == m_Raw[y] ) || ( m_AuxData[x] < m_Lambda ) );
246  }
247 
249  {
250  OffsetValueType r = FindRoot(n);
251 
252  if ( r != p )
253  {
254  if ( Criterion(r, p) )
255  {
256  m_AuxData[p] += m_AuxData[r];
257  m_Parent[r] = p;
258  }
259  else
260  {
261  m_AuxData[p] = m_Lambda;
262  }
263  }
264  }
265 
266 #else
267  // version from ISMM paper
268  void MakeSet(OffsetValueType x)
269  {
270  m_Parent[x] = ACTIVE;
271  m_AuxData[x] = m_AttributeValuePerPixel;
272  }
273 
274  void Link(OffsetValueType x, OffsetValueType y)
275  {
276  if ( ( m_Parent[y] == ACTIVE ) && ( m_Parent[x] == ACTIVE ) )
277  {
278  // should be a call to MergeAuxData
279  m_AuxData[y] = m_AuxData[x] + m_AuxData[y];
280  m_AuxData[x] = -m_AttributeValuePerPixel;
281  }
282  else if ( m_Parent[x] == ACTIVE )
283  {
284  m_AuxData[x] = -m_AttributeValuePerPixel;
285  }
286  else
287  {
288  m_AuxData[y] = -m_AttributeValuePerPixel;
289  m_Parent[y] = INACTIVE;
290  }
291  m_Parent[x] = y;
292  }
293 
294  OffsetValueType FindRoot(OffsetValueType x)
295  {
296  if ( m_Parent[x] >= 0 )
297  {
298  m_Parent[x] = FindRoot(m_Parent[x]);
299  return ( m_Parent[x] );
300  }
301  else
302  {
303  return ( x );
304  }
305  }
306 
307  bool Equiv(OffsetValueType x, OffsetValueType y)
308  {
309  return ( ( m_Raw[x] == m_Raw[y] ) || ( m_Parent[x] == ACTIVE ) );
310  }
311 
312  void Union(OffsetValueType n, OffsetValueType p)
313  {
314  OffsetValueType r = FindRoot(n);
315 
316  if ( r != p )
317  {
318  if ( Equiv(r, p) )
319  {
320  Link(r, p);
321  }
322  else if ( m_Parent[p] == ACTIVE )
323  {
324  m_Parent[p] = INACTIVE;
325  m_AuxData[p] = -m_AttributeValuePerPixel;
326  }
327  }
328  }
329 
330 #endif
331 };
332 } // end namespace itk
333 
334 #ifndef ITK_MANUAL_INSTANTIATION
335 #include "itkAttributeMorphologyBaseImageFilter.hxx"
336 #endif
337 
338 #endif
bool Criterion(OffsetValueType x, OffsetValueType y)
bool operator()(GreyAndPos const &l, GreyAndPos const &r) const
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Base class for all process objects that output image data.
void Union(OffsetValueType n, OffsetValueType p)
typename InputImageType::Pointer InputImagePointer
TOutputImage OutputImageType
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
signed long OffsetValueType
Definition: itkIntTypes.h:94
Base class for all data objects in ITK.