ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkImageMomentsCalculator.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 itkImageMomentsCalculator_h
19 #define itkImageMomentsCalculator_h
20 
21 #include "itkAffineTransform.h"
22 #include "itkImage.h"
23 #include "itkSpatialObject.h"
24 
25 #include "vnl/vnl_vector_fixed.h"
26 #include "vnl/vnl_matrix_fixed.h"
27 #include "vnl/vnl_diag_matrix.h"
28 
29 namespace itk
30 {
59 template< typename TImage >
60 class ITK_TEMPLATE_EXPORT ImageMomentsCalculator:public Object
61 {
62 public:
65  typedef Object Superclass;
68 
70  itkNewMacro(Self);
71 
73  itkTypeMacro(ImageMomentsCalculator, Object);
74 
76  itkStaticConstMacro(ImageDimension, unsigned int,
77  TImage::ImageDimension);
78 
80  typedef double ScalarType;
81 
84 
87 
91 
93  typedef Matrix< ScalarType,
94  itkGetStaticConstMacro(ImageDimension),
95  itkGetStaticConstMacro(ImageDimension) > MatrixType;
96 
98  typedef TImage ImageType;
99 
101  typedef typename ImageType::Pointer ImagePointer;
102  typedef typename ImageType::ConstPointer ImageConstPointer;
103 
107 
109  virtual void SetImage(const ImageType *image)
110  {
111  if ( m_Image != image )
112  {
113  m_Image = image;
114  this->Modified();
115  m_Valid = false;
116  }
117  }
119 
121  virtual void SetSpatialObjectMask(const SpatialObject< itkGetStaticConstMacro(ImageDimension) > *so)
122  {
123  if ( m_SpatialObjectMask != so )
124  {
125  m_SpatialObjectMask = so;
126  this->Modified();
127  m_Valid = false;
128  }
129  }
131 
137  void Compute();
138 
143  ScalarType GetTotalMass() const;
144 
150  VectorType GetFirstMoments() const;
151 
157  MatrixType GetSecondMoments() const;
158 
163  VectorType GetCenterOfGravity() const;
164 
169  MatrixType GetCentralMoments() const;
170 
177  VectorType GetPrincipalMoments() const;
178 
191  MatrixType GetPrincipalAxes() const;
193 
197  AffineTransformPointer GetPrincipalAxesToPhysicalAxesTransform() const;
198 
203  AffineTransformPointer GetPhysicalAxesToPrincipalAxesTransform() const;
204 
205 protected:
207  virtual ~ImageMomentsCalculator() ITK_OVERRIDE;
208  virtual void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
209 
210 private:
211  ITK_DISALLOW_COPY_AND_ASSIGN(ImageMomentsCalculator);
212 
213  bool m_Valid; // Have moments been computed yet?
214  ScalarType m_M0; // Zeroth moment
215  VectorType m_M1; // First moments about origin
216  MatrixType m_M2; // Second moments about origin
217  VectorType m_Cg; // Center of gravity (physical units)
218  MatrixType m_Cm; // Second central moments (physical)
219  VectorType m_Pm; // Principal moments (physical)
220  MatrixType m_Pa; // Principal axes (physical)
221 
223  SpatialObjectConstPointer m_SpatialObjectMask;
224 }; // class ImageMomentsCalculator
225 } // end namespace itk
226 
227 #ifndef ITK_MANUAL_INSTANTIATION
228 #include "itkImageMomentsCalculator.hxx"
229 #endif
230 
231 #endif /* itkImageMomentsCalculator_h */
SpatialObjectType::ConstPointer SpatialObjectConstPointer
A templated class holding a M x N size Matrix.
Definition: itkMatrix.h:53
Light weight base class for most itk classes.
SmartPointer< const Self > ConstPointer
SpatialObjectType::Pointer SpatialObjectPointer
virtual void SetSpatialObjectMask(const SpatialObject< itkGetStaticConstMacro(ImageDimension) > *so)
Matrix< ScalarType, itkGetStaticConstMacro(ImageDimension), itkGetStaticConstMacro(ImageDimension) > MatrixType
Compute moments of an n-dimensional image.
Implementation of the composite pattern.
SpatialObject< itkGetStaticConstMacro(ImageDimension) > SpatialObjectType
AffineTransform< double, itkGetStaticConstMacro(ImageDimension) > AffineTransformType
ImageType::ConstPointer ImageConstPointer
ImageMomentsCalculator< TImage > Self
Control indentation during Print() invocation.
Definition: itkIndent.h:49
AffineTransformType::Pointer AffineTransformPointer
virtual void SetImage(const ImageType *image)
Base class for most ITK classes.
Definition: itkObject.h:59
Vector< ScalarType, itkGetStaticConstMacro(ImageDimension) > VectorType