ITK  4.9.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 >
64  public ImageToImageFilter< TInputImage, TOutputImage >
65 {
66 public:
72 
77 
82  typedef typename TOutputImage::PixelType OutputPixelType;
83  typedef typename TOutputImage::InternalPixelType OutputInternalPixelType;
84  typedef typename TInputImage::PixelType InputPixelType;
85  typedef typename TInputImage::InternalPixelType InputInternalPixelType;
86  typedef typename TInputImage::IndexType IndexType;
87  typedef typename TInputImage::OffsetType OffsetType;
88  typedef typename TInputImage::SizeType SizeType;
89 
90  itkStaticConstMacro(ImageDimension, unsigned int,
91  TOutputImage::ImageDimension);
92 
96  typedef TInputImage InputImageType;
97  typedef TOutputImage OutputImageType;
98 // typedef typename TInputImage::IndexType IndexType;
99 // typedef typename TInputImage::SizeType SizeType;
100  typedef typename TOutputImage::RegionType RegionType;
101  typedef std::list< IndexType > ListType;
102  typedef TAttribute AttributeType;
103 
109 
114 
118  itkNewMacro(Self);
119 
126  itkSetMacro(FullyConnected, bool);
127  itkGetConstReferenceMacro(FullyConnected, bool);
128  itkBooleanMacro(FullyConnected);
130 
136  itkSetMacro(Lambda, AttributeType);
137  itkGetConstMacro(Lambda, AttributeType);
139 
140 protected:
142  {
143  m_FullyConnected = false;
145  m_Lambda = 0;
146  }
147 
150  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
151 
155  void GenerateData() ITK_OVERRIDE;
156 
160  void GenerateInputRequestedRegion() ITK_OVERRIDE;
161 
166  void EnlargeOutputRequestedRegion( DataObject * itkNotUsed(output) ) ITK_OVERRIDE;
167 
169 
170 private:
171 
174 
175  // some constants used several times in the code
176  itkStaticConstMacro(INACTIVE, OffsetValueType, -1);
177  itkStaticConstMacro(ACTIVE, OffsetValueType, -2);
178  itkStaticConstMacro(ROOT, OffsetValueType, -3);
179 
180  // Just used for area/volume openings at the moment
182 
183  typedef std::vector< OffsetType > OffsetVecType;
184  // offset in the linear array.
185  typedef std::vector< OffsetValueType > OffsetDirectVecType;
186 
187  void SetupOffsetVec(OffsetDirectVecType & PosOffsets, OffsetVecType & Offsets);
188 
190  {
191 public:
193  OffsetValueType Pos;
194  };
195 
196  GreyAndPos * m_SortPixels;
198 #ifndef PAMI
199  bool *m_Processed;
200 #endif
201  // This is a bit ugly, but I can't see an easy way around
203 
205  {
206 public:
207  TFunction m_TFunction;
208  bool operator()(GreyAndPos const & l, GreyAndPos const & r) const
209  {
210  if ( m_TFunction(l.Val, r.Val) )
211  {
212  return true;
213  }
214  if ( l.Val == r.Val )
215  {
216  return ( l.Pos < r.Pos );
217  }
218  return false;
219  }
220  };
221 
222 #ifdef PAMI
223  // version from PAMI. Note - using the AuxData array rather than the
224  // parent array to store area
226  {
227  m_Parent[x] = ACTIVE;
229  }
230 
232  {
233  if ( m_Parent[x] >= 0 )
234  {
235  m_Parent[x] = FindRoot(m_Parent[x]);
236  return ( m_Parent[x] );
237  }
238  else
239  {
240  return ( x );
241  }
242  }
243 
245  {
246  return ( ( m_Raw[x] == m_Raw[y] ) || ( m_AuxData[x] < m_Lambda ) );
247  }
248 
250  {
251  OffsetValueType r = FindRoot(n);
252 
253  if ( r != p )
254  {
255  if ( Criterion(r, p) )
256  {
257  m_AuxData[p] += m_AuxData[r];
258  m_Parent[r] = p;
259  }
260  else
261  {
262  m_AuxData[p] = m_Lambda;
263  }
264  }
265  }
266 
267 #else
268  // version from ISMM paper
269  void MakeSet(OffsetValueType x)
270  {
271  m_Parent[x] = ACTIVE;
273  }
274 
275  void Link(OffsetValueType x, OffsetValueType y)
276  {
277  if ( ( m_Parent[y] == ACTIVE ) && ( m_Parent[x] == ACTIVE ) )
278  {
279  // should be a call to MergeAuxData
280  m_AuxData[y] = m_AuxData[x] + m_AuxData[y];
281  m_AuxData[x] = -m_AttributeValuePerPixel;
282  }
283  else if ( m_Parent[x] == ACTIVE )
284  {
286  }
287  else
288  {
290  m_Parent[y] = INACTIVE;
291  }
292  m_Parent[x] = y;
293  }
294 
296  {
297  if ( m_Parent[x] >= 0 )
298  {
299  m_Parent[x] = FindRoot(m_Parent[x]);
300  return ( m_Parent[x] );
301  }
302  else
303  {
304  return ( x );
305  }
306  }
307 
308  bool Equiv(OffsetValueType x, OffsetValueType y)
309  {
310  return ( ( m_Raw[x] == m_Raw[y] ) || ( m_Parent[x] == ACTIVE ) );
311  }
312 
314  {
315  OffsetValueType r = FindRoot(n);
316 
317  if ( r != p )
318  {
319  if ( Equiv(r, p) )
320  {
321  Link(r, p);
322  }
323  else if ( m_Parent[p] == ACTIVE )
324  {
325  m_Parent[p] = INACTIVE;
327  }
328  }
329  }
330 
331 #endif
332 };
333 } // end namespace itk
334 
335 #ifndef ITK_MANUAL_INSTANTIATION
336 #include "itkAttributeMorphologyBaseImageFilter.hxx"
337 #endif
338 
339 #endif
void EnlargeOutputRequestedRegion(DataObject *) override
bool Criterion(OffsetValueType x, OffsetValueType y)
bool operator()(GreyAndPos const &l, GreyAndPos const &r) const
signed long OffsetValueType
Definition: itkIntTypes.h:154
ImageToImageFilter< TInputImage, TOutputImage > Superclass
InputImageType::Pointer InputImagePointer
Base class for all process objects that output image data.
void Union(OffsetValueType n, OffsetValueType p)
void SetupOffsetVec(OffsetDirectVecType &PosOffsets, OffsetVecType &Offsets)
void PrintSelf(std::ostream &os, Indent indent) const override
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
Base class for all data objects in ITK.