ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkImageMomentsCalculator.h
Go to the documentation of this file.
1/*=========================================================================
2 *
3 * Copyright NumFOCUS
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 * https://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
29namespace itk
30{
59template <typename TImage>
60class ITK_TEMPLATE_EXPORT ImageMomentsCalculator : public Object
61{
62public:
63 ITK_DISALLOW_COPY_AND_MOVE(ImageMomentsCalculator);
64
70
72 itkNewMacro(Self);
73
75 itkOverrideGetNameOfClassMacro(ImageMomentsCalculator);
76
78 static constexpr unsigned int ImageDimension = TImage::ImageDimension;
79
81 using ScalarType = double;
82
85
88
92
95
97 using ImageType = TImage;
98
100 using ImagePointer = typename ImageType::Pointer;
101 using ImageConstPointer = typename ImageType::ConstPointer;
102
106
108 virtual void
109 SetImage(const ImageType * image)
110 {
111 if (m_Image != image)
112 {
113 m_Image = image;
114 this->Modified();
115 m_Valid = false;
116 }
117 }
118
119
121 virtual void
123 {
124 if (m_SpatialObjectMask != so)
125 {
127 this->Modified();
128 m_Valid = false;
129 }
130 }
131
132
138 void
140
147
155
163
170
177
186
202
208
215
216protected:
218 ~ImageMomentsCalculator() override = default;
219 void
220 PrintSelf(std::ostream & os, Indent indent) const override;
221
222private:
223 bool m_Valid{}; // Have moments been computed yet?
224 ScalarType m_M0{}; // Zeroth moment
225 VectorType m_M1{}; // First moments about origin
226 MatrixType m_M2{}; // Second moments about origin
227 VectorType m_Cg{}; // Center of gravity (physical units)
228 MatrixType m_Cm{}; // Second central moments (physical)
229 VectorType m_Pm{}; // Principal moments (physical)
230 MatrixType m_Pa{}; // Principal axes (physical)
231
234}; // class ImageMomentsCalculator
235} // end namespace itk
236
237#ifndef ITK_MANUAL_INSTANTIATION
238# include "itkImageMomentsCalculator.hxx"
239#endif
240
241#endif /* itkImageMomentsCalculator_h */
VectorType GetCenterOfGravity() const
MatrixType GetCentralMoments() const
typename AffineTransformType::Pointer AffineTransformPointer
ImageMomentsCalculator< TImage > Self
virtual void SetSpatialObjectMask(const SpatialObject< Self::ImageDimension > *so)
Matrix< ScalarType, Self::ImageDimension, Self::ImageDimension > MatrixType
MatrixType GetPrincipalAxes() const
void PrintSelf(std::ostream &os, Indent indent) const override
VectorType GetFirstMoments() const
~ImageMomentsCalculator() override=default
SmartPointer< const Self > ConstPointer
typename SpatialObjectType::ConstPointer SpatialObjectConstPointer
typename SpatialObjectType::Pointer SpatialObjectPointer
AffineTransform< double, Self::ImageDimension > AffineTransformType
virtual void SetImage(const ImageType *image)
typename ImageType::Pointer ImagePointer
VectorType GetPrincipalMoments() const
AffineTransformPointer GetPhysicalAxesToPrincipalAxesTransform() const
typename ImageType::ConstPointer ImageConstPointer
AffineTransformPointer GetPrincipalAxesToPhysicalAxesTransform() const
ScalarType GetTotalMass() const
MatrixType GetSecondMoments() const
SpatialObject< Self::ImageDimension > SpatialObjectType
Vector< ScalarType, Self::ImageDimension > VectorType
Control indentation during Print() invocation.
Definition itkIndent.h:50
A templated class holding a M x N size Matrix.
Definition itkMatrix.h:53
virtual void Modified() const
Implements transparent reference counting.
Implementation of the composite pattern.
A templated class holding a n-Dimensional vector.
Definition itkVector.h:63
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....