ITK  5.0.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 <vector>
27 #include "vnl/algo/vnl_symmetric_eigensystem.h"
28 #include "vnl/vnl_det.h"
29 #include "itkMath.h"
30 
31 namespace itk
32 {
76 template< typename TLabelImage, typename TIntensityImage = TLabelImage >
77 class ITK_TEMPLATE_EXPORT LabelGeometryImageFilter:
78  public ImageToImageFilter< TLabelImage, TIntensityImage >
79 {
80 public:
81  ITK_DISALLOW_COPY_AND_ASSIGN(LabelGeometryImageFilter);
82 
88 
90  itkNewMacro(Self);
91 
94 
96  using IntensityImageType = TIntensityImage;
97  using InputImagePointer = typename TIntensityImage::Pointer;
101  using PixelType = typename TIntensityImage::PixelType;
102 
104  using LabelImageType = TLabelImage;
105  using LabelImagePointer = typename TLabelImage::Pointer;
109  using LabelPixelType = typename TLabelImage::PixelType;
111 
113  static constexpr unsigned int ImageDimension = TLabelImage::ImageDimension;
114 
117 
120 
123 
126  Self::ImageDimension *2 >;
127  using BoundingBoxFloatType = itk::FixedArray< float,
128  Self::ImageDimension *2 >;
129 
130  // using BoundingBoxVerticesType = itk::FixedArray<
131  // LabelPointType,std::pow(2.0,Self::ImageDimension)>;
132  using BoundingBoxVerticesType = std::vector< LabelPointType >;
133 
136 
139  Self::ImageDimension >;
140 
142  using LabelsType = std::vector< LabelPixelType >;
143 
145  using LabelIndicesType = std::vector< LabelIndexType >;
146 
148  using VectorType = std::vector< double >;
149 
151  using MatrixType = vnl_matrix< double >;
152 
158  {
159 public:
160  // default constructor
162  {
163  // initialized to the default values
164  this->m_Label = 0;
165  this->m_Sum = NumericTraits< RealType >::ZeroValue();
166 
167  const unsigned int imageDimension = Self::ImageDimension;
168 
169  //m_BoundingBox.resize(imageDimension*2);
170  for ( unsigned int i = 0; i < imageDimension * 2; i += 2 )
171  {
174  }
175 
176  m_BoundingBoxVolume = 0;
177  m_BoundingBoxSize.Fill(0);
178  m_PixelIndices.clear();
179  m_Centroid.Fill(0);
180  m_WeightedCentroid.Fill(0);
181  m_ZeroOrderMoment = 0;
182  m_FirstOrderRawMoments.Fill(0);
183  m_FirstOrderWeightedRawMoments.Fill(0);
184  m_Eigenvalues.resize(ImageDimension);
185  m_Eigenvalues.clear();
186  m_Eigenvectors.set_size(ImageDimension, ImageDimension);
187  m_Eigenvectors.fill(0);
188  m_AxesLength.Fill(0);
189  m_Eccentricity = 1;
190  m_Elongation = 1;
191  m_Orientation = 0;
192  LabelPointType emptyPoint;
193  emptyPoint.Fill(0);
194  unsigned int numberOfVertices = 1 << ImageDimension;
195  m_OrientedBoundingBoxVertices.resize(numberOfVertices, emptyPoint);
196  m_OrientedBoundingBoxVolume = 0;
197  m_OrientedBoundingBoxSize.Fill(0);
198  m_OrientedLabelImage = LabelImageType::New();
199  m_OrientedIntensityImage = IntensityImageType::New();
200  m_OrientedBoundingBoxOrigin.Fill(0);
201  m_RotationMatrix.set_size(ImageDimension, ImageDimension);
202  m_RotationMatrix.fill(0.0);
203 
204  m_SecondOrderRawMoments.set_size(ImageDimension, ImageDimension);
205  m_SecondOrderCentralMoments.set_size(ImageDimension, ImageDimension);
206  for ( unsigned int i = 0; i < ImageDimension; i++ )
207  {
208  for ( unsigned int j = 0; j < ImageDimension; j++ )
209  {
210  m_SecondOrderRawMoments(i, j) = 0;
211  m_SecondOrderCentralMoments(i, j) = 0;
212  }
213  }
214  }
215 
240  typename LabelImageType::Pointer m_OrientedLabelImage;
241  typename IntensityImageType::Pointer m_OrientedIntensityImage;
244  };
245 
247  // Map from the label to the class storing all of the geometry information.
248  using MapType = itksys::hash_map< LabelPixelType, LabelGeometry >;
249  using MapIterator = typename itksys::hash_map< LabelPixelType, LabelGeometry >::iterator;
250  using MapConstIterator = typename itksys::hash_map< LabelPixelType, LabelGeometry >::const_iterator;
251 
252  // Macros for enabling the calculation of additional features.
253  itkGetMacro(CalculatePixelIndices, bool);
254  itkBooleanMacro(CalculatePixelIndices);
255  void SetCalculatePixelIndices(const bool value)
256  {
257  // CalculateOrientedBoundingBox, CalculateOrientedLabelImage, and
258  // CalculateOrientedIntensityImage all need CalculatePixelIndices to be
259  // turned
260  // on if they are turned on. So, CalculatePixelIndices cannot be
261  // turned off if any of these flags are turned on.
262  if ( value == false )
263  {
264  if ( ( this->m_CalculateOrientedBoundingBox == true )
265  || ( this->m_CalculateOrientedLabelRegions == true )
266  || ( this->m_CalculateOrientedIntensityRegions == true ) )
267  {
268  // We cannot change the value, so return.
269  return;
270  }
271  }
272 
273  if ( this->m_CalculatePixelIndices != value )
274  {
275  this->m_CalculatePixelIndices = value;
276  this->Modified();
277  }
278  }
279 
280  itkGetMacro(CalculateOrientedBoundingBox, bool);
281  itkBooleanMacro(CalculateOrientedBoundingBox);
282  void SetCalculateOrientedBoundingBox(const bool value)
283  {
284  if ( this->m_CalculateOrientedBoundingBox != value )
285  {
286  this->m_CalculateOrientedBoundingBox = value;
287  this->Modified();
288  }
289 
290  // CalculateOrientedBoundingBox needs
291  // CalculatePixelIndices to be turned on.
292  if ( value == true )
293  {
294  this->SetCalculatePixelIndices(true);
295  }
296  }
297 
298  itkGetMacro(CalculateOrientedLabelRegions, bool);
299  itkBooleanMacro(CalculateOrientedLabelRegions);
300  void SetCalculateOrientedLabelRegions(const bool value)
301  {
302  if ( this->m_CalculateOrientedLabelRegions != value )
303  {
304  this->m_CalculateOrientedLabelRegions = value;
305  this->Modified();
306 
307  // CalculateOrientedLabelImage needs
308  // CalculateOrientedBoundingBox to be turned on.
309  if ( value == true )
310  {
311  SetCalculateOrientedBoundingBox(true);
312  }
313  }
314  }
315 
316  itkGetMacro(CalculateOrientedIntensityRegions, bool);
317  itkBooleanMacro(CalculateOrientedIntensityRegions);
319  {
320  if ( this->m_CalculateOrientedIntensityRegions != value )
321  {
322  this->m_CalculateOrientedIntensityRegions = value;
323  this->Modified();
324 
325  // CalculateOrientedIntensityImage needs
326  // CalculateOrientedBoundingBox to be turned on.
327  if ( value == true )
328  {
329  this->SetCalculateOrientedBoundingBox(true);
330  }
331  }
332  }
333 
335  void SetIntensityInput(const TIntensityImage *input)
336  {
337  // Process object is not const-correct so the const casting is required.
338  this->SetNthInput( 1, const_cast< TIntensityImage * >( input ) );
339  }
340 
342  const TIntensityImage * GetIntensityInput() const
343  {
344  return static_cast< TIntensityImage * >( const_cast< DataObject * >( this->ProcessObject::GetInput(1) ) );
345  }
346 
349  bool HasLabel(LabelPixelType label) const
350  {
351  return m_LabelGeometryMapper.find(label) != m_LabelGeometryMapper.end();
352  }
353 
356  {
357  return m_LabelGeometryMapper.size();
358  }
359 
361  {
362  return this->GetNumberOfObjects();
363  }
364 
366  std::vector< LabelPixelType > GetLabels() const
367  {
368  return m_AllLabels;
369  }
370 
372  LabelIndicesType GetPixelIndices(LabelPixelType label) const;
373 
376  SizeValueType GetVolume(LabelPixelType label) const;
377 
379  //std::vector< SizeValueType > GetAllCounts() const;
380 
382  RealType GetIntegratedIntensity(LabelPixelType label) const;
383 
385  LabelPointType GetCentroid(LabelPixelType label) const;
386 
388  LabelPointType GetWeightedCentroid(LabelPixelType label) const;
389 
391  VectorType GetEigenvalues(LabelPixelType label) const;
392 
394  MatrixType GetEigenvectors(LabelPixelType label) const;
395 
397  AxesLengthType GetAxesLength(LabelPixelType label) const;
398 
401  RealType GetMinorAxisLength(LabelPixelType label) const;
402 
405  RealType GetMajorAxisLength(LabelPixelType label) const;
406 
408  RealType GetEccentricity(LabelPixelType label) const;
409 
412  RealType GetElongation(LabelPixelType label) const;
413 
415  RealType GetOrientation(LabelPixelType label) const;
416 
420  BoundingBoxType GetBoundingBox(LabelPixelType label) const;
421 
423  RealType GetBoundingBoxVolume(LabelPixelType label) const;
424 
426  LabelSizeType GetBoundingBoxSize(LabelPixelType label) const;
427 
434  BoundingBoxVerticesType GetOrientedBoundingBoxVertices(LabelPixelType label) const;
435 
437  RealType GetOrientedBoundingBoxVolume(LabelPixelType label) const;
438 
440  LabelPointType GetOrientedBoundingBoxSize(LabelPixelType label) const;
441 
443  LabelPointType GetOrientedBoundingBoxOrigin(LabelPixelType label) const;
444 
447  MatrixType GetRotationMatrix(LabelPixelType label) const;
448 
450  RegionType GetRegion(LabelPixelType label) const;
451 
453  TLabelImage * GetOrientedLabelImage(LabelPixelType label) const;
454 
457  TIntensityImage * GetOrientedIntensityImage(LabelPixelType label) const;
458 
459 #ifdef ITK_USE_CONCEPT_CHECKING
460  // Begin concept checking
461  itkConceptMacro( InputHasNumericTraitsCheck,
463  // End concept checking
464 #endif
465 
466 protected:
468  ~LabelGeometryImageFilter() override = default;
469  void PrintSelf(std::ostream & os, Indent indent) const override;
470 
471  void GenerateData() override;
472 
473 private:
474  bool CalculateOrientedBoundingBoxVertices(vnl_symmetric_eigensystem< double > eig, LabelGeometry & m_LabelGeometry);
475 
480 
484 }; // end of class
485 
486 } // end namespace itk
487 
488 #ifndef ITK_MANUAL_INSTANTIATION
489 #include "itkLabelGeometryImageFilter.hxx"
490 #endif
491 
492 #endif
typename TLabelImage::RegionType LabelRegionType
bool HasLabel(LabelPixelType label) const
FixedArray< float, Self::ImageDimension > m_AxesLength
typename TLabelImage::PixelType LabelPixelType
typename itksys::hash_map< LabelPixelType, LabelGeometry >::const_iterator MapConstIterator
Light weight base class for most itk classes.
Define numeric traits for std::vector.
void SetCalculateOrientedIntensityRegions(const bool value)
unsigned long SizeValueType
Definition: itkIntTypes.h:83
std::vector< LabelPointType > BoundingBoxVerticesType
void SetCalculateOrientedBoundingBox(const bool value)
typename TIntensityImage::RegionType RegionType
typename TIntensityImage::IndexType IndexType
typename TLabelImage::IndexType LabelIndexType
std::vector< LabelPixelType > LabelsType
void SetCalculateOrientedLabelRegions(const bool value)
typename TIntensityImage::Pointer InputImagePointer
Decorates any &quot;simple&quot; data type (data types without smart pointers) with a DataObject API...
void SetCalculatePixelIndices(const bool value)
typename NumericTraits< PixelType >::RealType RealType
typename TLabelImage::SizeType LabelSizeType
signed long IndexValueType
Definition: itkIntTypes.h:90
std::vector< LabelPixelType > GetLabels() const
Given a label map and an optional intensity image, compute geometric features.
DataObject * GetInput(const DataObjectIdentifierType &key)
Return an input.
void SetIntensityInput(const TIntensityImage *input)
typename TLabelImage::PointType LabelPointType
std::vector< LabelIndexType > LabelIndicesType
Base class for filters that take an image as input and produce an image as output.
typename TIntensityImage::SizeType SizeType
typename TIntensityImage::PixelType PixelType
const TIntensityImage * GetIntensityInput() const
typename itksys::hash_map< LabelPixelType, LabelGeometry >::iterator MapIterator
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itksys::hash_map< LabelPixelType, LabelGeometry > MapType
SmartPointer< Self > Pointer
typename TLabelImage::Pointer LabelImagePointer
#define itkConceptMacro(name, concept)
Base class for all data objects in ITK.