ITK  4.6.0
Insight Segmentation and Registration Toolkit
itkLabelGeometryImageFilter.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 __itkLabelGeometryImageFilter_h
19 #define __itkLabelGeometryImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 #include "itkNumericTraits.h"
23 #include "itkArray.h"
25 #include "itksys/hash_map.hxx"
26 #include "itkFastMutexLock.h"
27 #include <vector>
28 #include "vnl/algo/vnl_symmetric_eigensystem.h"
29 #include "vnl/vnl_det.h"
30 #include "vnl/vnl_math.h"
31 
32 namespace itk
33 {
77 template< typename TLabelImage, typename TIntensityImage = TLabelImage >
79  public ImageToImageFilter< TLabelImage, TIntensityImage >
80 {
81 public:
87 
89  itkNewMacro(Self);
90 
93 
95  typedef TIntensityImage IntensityImageType;
96  typedef typename TIntensityImage::Pointer InputImagePointer;
97  typedef typename TIntensityImage::RegionType RegionType;
98  typedef typename TIntensityImage::SizeType SizeType;
99  typedef typename TIntensityImage::IndexType IndexType;
100  typedef typename TIntensityImage::PixelType PixelType;
101 
103  typedef TLabelImage LabelImageType;
104  typedef typename TLabelImage::Pointer LabelImagePointer;
105  typedef typename TLabelImage::RegionType LabelRegionType;
106  typedef typename TLabelImage::SizeType LabelSizeType;
107  typedef typename TLabelImage::IndexType LabelIndexType;
108  typedef typename TLabelImage::PixelType LabelPixelType;
109  typedef typename TLabelImage::PointType LabelPointType;
110 
112  itkStaticConstMacro(ImageDimension, unsigned int,
113  TLabelImage::ImageDimension);
114 
117 
120 
123 
126  itkGetStaticConstMacro(ImageDimension) *2 > BoundingBoxType;
127  typedef itk::FixedArray< float,
128  itkGetStaticConstMacro(ImageDimension) *2 > BoundingBoxFloatType;
130 
131  //typedef itk::FixedArray<
132  // LabelPointType,std::pow(2.0,itkGetStaticConstMacro(ImageDimension))>
133  // BoundingBoxVerticesType;
134  typedef std::vector< LabelPointType > BoundingBoxVerticesType;
135 
138 
141  itkGetStaticConstMacro(ImageDimension) > IndexArrayType;
142 
144  typedef std::vector< LabelPixelType > LabelsType;
145 
147  typedef std::vector< LabelIndexType > LabelIndicesType;
148 
150  typedef std::vector< double > VectorType;
151 
153  typedef vnl_matrix< double > MatrixType;
154 
160  {
161 public:
162  // default constructor
164  {
165  // initialized to the default values
166  this->m_Label = 0;
168 
169  const unsigned int imageDimension = itkGetStaticConstMacro(ImageDimension);
170 
171  //m_BoundingBox.resize(imageDimension*2);
172  for ( unsigned int i = 0; i < imageDimension * 2; i += 2 )
173  {
176  }
177 
179  m_BoundingBoxSize.Fill(0);
180  m_PixelIndices.clear();
181  m_Centroid.Fill(0);
182  m_WeightedCentroid.Fill(0);
183  m_ZeroOrderMoment = 0;
187  m_Eigenvalues.clear();
189  m_Eigenvectors.fill(0);
190  m_AxesLength.Fill(0);
191  m_Eccentricity = 1;
192  m_Elongation = 1;
193  m_Orientation = 0;
194  LabelPointType emptyPoint;
195  emptyPoint.Fill(0);
196  unsigned int numberOfVertices = (unsigned int)std::pow( (double)2, (int)ImageDimension );
197  m_OrientedBoundingBoxVertices.resize(numberOfVertices, emptyPoint);
200  m_OrientedLabelImage = LabelImageType::New();
201  m_OrientedIntensityImage = IntensityImageType::New();
204  m_RotationMatrix.fill(0.0);
205 
208  for ( unsigned int i = 0; i < ImageDimension; i++ )
209  {
210  for ( unsigned int j = 0; j < ImageDimension; j++ )
211  {
212  m_SecondOrderRawMoments(i, j) = 0;
213  m_SecondOrderCentralMoments(i, j) = 0;
214  }
215  }
216  }
217 
242  typename LabelImageType::Pointer m_OrientedLabelImage;
243  typename IntensityImageType::Pointer m_OrientedIntensityImage;
246  };
247 
249  // Map from the label to the class storing all of the geometry information.
250  typedef itksys::hash_map< LabelPixelType, LabelGeometry > MapType;
251  typedef typename itksys::hash_map< LabelPixelType, LabelGeometry >::iterator MapIterator;
252  typedef typename itksys::hash_map< LabelPixelType, LabelGeometry >::const_iterator MapConstIterator;
253 
254  // Macros for enabling the calculation of additional features.
255  itkGetMacro(CalculatePixelIndices, bool);
256  itkBooleanMacro(CalculatePixelIndices);
257  void SetCalculatePixelIndices(const bool value)
258  {
259  // CalculateOrientedBoundingBox, CalculateOrientedLabelImage, and
260  // CalculateOrientedIntensityImage all need CalculatePixelIndices to be
261  // turned
262  // on if they are turned on. So, CalculatePixelIndices cannot be
263  // turned off if any of these flags are turned on.
264  if ( value == false )
265  {
266  if ( ( this->m_CalculateOrientedBoundingBox == true )
267  || ( this->m_CalculateOrientedLabelRegions == true )
268  || ( this->m_CalculateOrientedIntensityRegions == true ) )
269  {
270  // We cannot change the value, so return.
271  return;
272  }
273  }
274 
275  if ( this->m_CalculatePixelIndices != value )
276  {
277  this->m_CalculatePixelIndices = value;
278  this->Modified();
279  }
280  }
281 
282  itkGetMacro(CalculateOrientedBoundingBox, bool);
283  itkBooleanMacro(CalculateOrientedBoundingBox);
284  void SetCalculateOrientedBoundingBox(const bool value)
285  {
286  if ( this->m_CalculateOrientedBoundingBox != value )
287  {
288  this->m_CalculateOrientedBoundingBox = value;
289  this->Modified();
290  }
291 
292  // CalculateOrientedBoundingBox needs
293  // CalculatePixelIndices to be turned on.
294  if ( value == true )
295  {
296  this->SetCalculatePixelIndices(true);
297  }
298  }
299 
300  itkGetMacro(CalculateOrientedLabelRegions, bool);
301  itkBooleanMacro(CalculateOrientedLabelRegions);
302  void SetCalculateOrientedLabelRegions(const bool value)
303  {
304  if ( this->m_CalculateOrientedLabelRegions != value )
305  {
306  this->m_CalculateOrientedLabelRegions = value;
307  this->Modified();
308 
309  // CalculateOrientedLabelImage needs
310  // CalculateOrientedBoundingBox to be turned on.
311  if ( value == true )
312  {
314  }
315  }
316  }
317 
318  itkGetMacro(CalculateOrientedIntensityRegions, bool);
319  itkBooleanMacro(CalculateOrientedIntensityRegions);
321  {
322  if ( this->m_CalculateOrientedIntensityRegions != value )
323  {
325  this->Modified();
326 
327  // CalculateOrientedIntensityImage needs
328  // CalculateOrientedBoundingBox to be turned on.
329  if ( value == true )
330  {
332  }
333  }
334  }
335 
337  void SetIntensityInput(const TIntensityImage *input)
338  {
339  // Process object is not const-correct so the const casting is required.
340  this->SetNthInput( 1, const_cast< TIntensityImage * >( input ) );
341  }
342 
344  const TIntensityImage * GetIntensityInput() const
345  {
346  return static_cast< TIntensityImage * >( const_cast< DataObject * >( this->ProcessObject::GetInput(1) ) );
347  }
348 
351  bool HasLabel(LabelPixelType label) const
352  {
353  return m_LabelGeometryMapper.find(label) != m_LabelGeometryMapper.end();
354  }
355 
358  {
359  return m_LabelGeometryMapper.size();
360  }
361 
363  {
364  return this->GetNumberOfObjects();
365  }
366 
368  std::vector< LabelPixelType > GetLabels() const
369  {
370  return m_AllLabels;
371  }
372 
375 
379 
381  //std::vector< SizeValueType > GetAllCounts() const;
382 
385 
388 
391 
394 
397 
400 
404 
408 
411 
415 
418 
423 
426 
429 
437 
440 
443 
446 
450 
452  RegionType GetRegion(LabelPixelType label) const;
453 
455  TLabelImage * GetOrientedLabelImage(LabelPixelType label) const;
456 
459  TIntensityImage * GetOrientedIntensityImage(LabelPixelType label) const;
460 
461 #ifdef ITK_USE_CONCEPT_CHECKING
462  // Begin concept checking
463  itkConceptMacro( InputHasNumericTraitsCheck,
465  // End concept checking
466 #endif
467 
468 protected:
471  void PrintSelf(std::ostream & os, Indent indent) const;
472 
473  void GenerateData();
474 
475 private:
476  LabelGeometryImageFilter(const Self &); //purposely not implemented
477  void operator=(const Self &); //purposely not implemented
478 
479  bool CalculateOrientedBoundingBoxVertices(vnl_symmetric_eigensystem< double > eig, LabelGeometry & m_LabelGeometry);
480 
485 
489 
491 }; // end of class
492 
493 template< typename TLabelImage, typename TIntensityImage >
495  vnl_symmetric_eigensystem< double > eig);
496 
497 template< typename TLabelImage, typename TIntensityImage, typename TGenericImage >
500  vnl_symmetric_eigensystem< double > eig,
502  bool useLabelImage = true);
503 } // end namespace itk
504 
505 #ifndef ITK_MANUAL_INSTANTIATION
506 #include "itkLabelGeometryImageFilter.hxx"
507 #endif
508 
509 #endif
LabelIndicesType GetPixelIndices(LabelPixelType label) const
itk::FixedArray< float, itkGetStaticConstMacro(ImageDimension)*2 > BoundingBoxFloatType
Critical section locking class that can be allocated on the stack.
bool HasLabel(LabelPixelType label) const
LabelPointType GetOrientedBoundingBoxOrigin(LabelPixelType label) const
bool CalculateOrientedBoundingBoxVertices(vnl_symmetric_eigensystem< double > eig, LabelGeometry &m_LabelGeometry)
RealType GetMinorAxisLength(LabelPixelType label) const
Light weight base class for most itk classes.
itksys::hash_map< LabelPixelType, LabelGeometry >::const_iterator MapConstIterator
itk::FixedArray< RealType, itkGetStaticConstMacro(ImageDimension) > AxesLengthType
SimpleDataObjectDecorator< RealType > RealObjectType
void SetCalculateOrientedIntensityRegions(const bool value)
RealType GetEccentricity(LabelPixelType label) const
VectorType GetEigenvalues(LabelPixelType label) const
itksys::hash_map< LabelPixelType, LabelGeometry >::iterator MapIterator
signed long IndexValueType
Definition: itkIntTypes.h:150
RealType GetOrientation(LabelPixelType label) const
TIntensityImage::IndexType IndexType
RealType GetElongation(LabelPixelType label) const
LabelGeometryImageFilter< TLabelImage, TIntensityImage >::MatrixType CalculateRotationMatrix(vnl_symmetric_eigensystem< double > eig)
static const unsigned int ImageDimension
std::vector< LabelIndexType > LabelIndicesType
void SetCalculateOrientedBoundingBox(const bool value)
itksys::hash_map< LabelPixelType, LabelGeometry > MapType
unsigned long SizeValueType
Definition: itkIntTypes.h:143
RealType GetOrientedBoundingBoxVolume(LabelPixelType label) const
bool CalculateOrientedImage(LabelGeometryImageFilter< TLabelImage, TIntensityImage > *filter, vnl_symmetric_eigensystem< double > eig, typename LabelGeometryImageFilter< TLabelImage, TIntensityImage >::LabelGeometry &labelGeometry, bool useLabelImage=true)
MatrixType GetRotationMatrix(LabelPixelType label) const
void Fill(const ValueType &)
TLabelImage * GetOrientedLabelImage(LabelPixelType label) const
Simulate a standard C array with copy semnatics.
Definition: itkFixedArray.h:50
TIntensityImage::RegionType RegionType
std::vector< LabelPixelType > LabelsType
itk::FixedArray< typename LabelIndexType::IndexValueType, itkGetStaticConstMacro(ImageDimension) > IndexArrayType
void SetCalculateOrientedLabelRegions(const bool value)
MatrixType GetEigenvectors(LabelPixelType label) const
RealType GetIntegratedIntensity(LabelPixelType label) const
Decorates any &quot;simple&quot; data type (data types without smart pointers) with a DataObject API...
void SetCalculatePixelIndices(const bool value)
TIntensityImage::Pointer InputImagePointer
static T max(const T &)
LabelPointType GetOrientedBoundingBoxSize(LabelPixelType label) const
std::vector< LabelPixelType > GetLabels() const
RegionType GetRegion(LabelPixelType label) const
Given a label map and an optional intensity image, compute geometric features.
virtual void Modified() const
DataObject * GetInput(const DataObjectIdentifierType &key)
RealType GetBoundingBoxVolume(LabelPixelType label) const
LabelSizeType GetBoundingBoxSize(LabelPixelType label) const
void SetIntensityInput(const TIntensityImage *input)
SmartPointer< const Self > ConstPointer
std::vector< LabelPointType > BoundingBoxVerticesType
ImageToImageFilter< TLabelImage, TIntensityImage > Superclass
BoundingBoxVerticesType GetOrientedBoundingBoxVertices(LabelPixelType label) const
LabelPointType GetWeightedCentroid(LabelPixelType label) const
Base class for filters that take an image as input and produce an image as output.
void PrintSelf(std::ostream &os, Indent indent) const
FixedArray< float, itkGetStaticConstMacro(ImageDimension) > m_AxesLength
TIntensityImage * GetOrientedIntensityImage(LabelPixelType label) const
const TIntensityImage * GetIntensityInput() const
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Define additional traits for native types such as int or float.
virtual void SetNthInput(DataObjectPointerArraySizeType num, DataObject *input)
BoundingBoxType GetBoundingBox(LabelPixelType label) const
void operator=(const Self &)
#define itkConceptMacro(name, concept)
itk::FixedArray< typename LabelIndexType::IndexValueType, itkGetStaticConstMacro(ImageDimension)*2 > BoundingBoxType
LabelPointType GetCentroid(LabelPixelType label) const
NumericTraits< PixelType >::RealType RealType
AxesLengthType GetAxesLength(LabelPixelType label) const
Base class for all data objects in ITK.
RealType GetMajorAxisLength(LabelPixelType label) const
SizeValueType GetVolume(LabelPixelType label) const