ITK  5.2.0
Insight Toolkit
itkMeanSquaresImageToImageMetric.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  * 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 itkMeanSquaresImageToImageMetric_h
19 #define itkMeanSquaresImageToImageMetric_h
20 
21 #include "itkImageToImageMetric.h"
22 #include "itkPoint.h"
23 #include "itkIndex.h"
24 
25 
26 namespace itk
27 {
28 
38 template <typename TFixedImage, typename TMovingImage>
39 class ITK_TEMPLATE_EXPORT MeanSquaresImageToImageMetric : public ImageToImageMetric<TFixedImage, TMovingImage>
40 {
41 public:
42  ITK_DISALLOW_COPY_AND_MOVE(MeanSquaresImageToImageMetric);
43 
49 
51  itkNewMacro(Self);
52 
55 
57  using TransformType = typename Superclass::TransformType;
58  using TransformPointer = typename Superclass::TransformPointer;
59  using TransformJacobianType = typename Superclass::TransformJacobianType;
60  using InterpolatorType = typename Superclass::InterpolatorType;
61  using MeasureType = typename Superclass::MeasureType;
62  using DerivativeType = typename Superclass::DerivativeType;
63  using ParametersType = typename Superclass::ParametersType;
64  using FixedImageType = typename Superclass::FixedImageType;
65  using MovingImageType = typename Superclass::MovingImageType;
66  using MovingImagePointType = typename Superclass::MovingImagePointType;
67  using FixedImageConstPointer = typename Superclass::FixedImageConstPointer;
68  using MovingImageConstPointer = typename Superclass::MovingImageConstPointer;
69  using CoordinateRepresentationType = typename Superclass::CoordinateRepresentationType;
70  using FixedImageSampleContainer = typename Superclass::FixedImageSampleContainer;
71  using ImageDerivativesType = typename Superclass::ImageDerivativesType;
72  using WeightsValueType = typename Superclass::WeightsValueType;
74 
75  // Needed for evaluation of Jacobian.
76  using FixedImagePointType = typename Superclass::FixedImagePointType;
77 
79  static constexpr unsigned int MovingImageDimension = MovingImageType::ImageDimension;
80 
88  void
89  Initialize() override;
90 
93  GetValue(const ParametersType & parameters) const override;
94 
96  void
97  GetDerivative(const ParametersType & parameters, DerivativeType & derivative) const override;
98 
100  void
101  GetValueAndDerivative(const ParametersType & parameters,
102  MeasureType & value,
103  DerivativeType & derivative) const override;
104 
105 protected:
107  ~MeanSquaresImageToImageMetric() override;
108  void
109  PrintSelf(std::ostream & os, Indent indent) const override;
110 
111 private:
112  bool
113  GetValueThreadProcessSample(ThreadIdType threadId,
114  SizeValueType fixedImageSample,
115  const MovingImagePointType & mappedPoint,
116  double movingImageValue) const override;
117 
118  bool
119  GetValueAndDerivativeThreadProcessSample(ThreadIdType threadId,
120  SizeValueType fixedImageSample,
121  const MovingImagePointType & mappedPoint,
122  double movingImageValue,
123  const ImageDerivativesType & movingImageGradientValue) const override;
124 
125  struct PerThreadS
126  {
130  };
131 
132  itkAlignedTypedef(64, PerThreadS, AlignedPerThreadType);
133  AlignedPerThreadType * m_PerThread;
134 };
135 } // end namespace itk
136 
137 #ifndef ITK_MANUAL_INSTANTIATION
138 # include "itkMeanSquaresImageToImageMetric.hxx"
139 #endif
140 
141 #endif
itk::MeanSquaresImageToImageMetric::PerThreadS::m_MSEDerivative
DerivativeType m_MSEDerivative
Definition: itkMeanSquaresImageToImageMetric.h:129
itk::ImageToImageMetric::FixedImageConstPointer
typename FixedImageType::ConstPointer FixedImageConstPointer
Definition: itkImageToImageMetric.h:77
itk::ImageToImageMetric
Computes similarity between regions of two images.
Definition: itkImageToImageMetric.h:52
itk::SingleValuedCostFunction::MeasureType
double MeasureType
Definition: itkSingleValuedCostFunction.h:50
itk::OptimizerParameters< TInternalComputationValueType >
itk::ImageToImageMetric::WeightsValueType
typename BSplineTransformWeightsType::ValueType WeightsValueType
Definition: itkImageToImageMetric.h:450
itk::MeanSquaresImageToImageMetric::PerThreadS
Definition: itkMeanSquaresImageToImageMetric.h:125
itk::ImageToImageMetric::FixedImagePointType
typename TransformType::InputPointType FixedImagePointType
Definition: itkImageToImageMetric.h:97
itkPoint.h
itk::ImageToImageMetric::TransformJacobianType
typename TransformType::JacobianType TransformJacobianType
Definition: itkImageToImageMetric.h:91
itk::ImageToImageMetric::MovingImagePointType
typename TransformType::OutputPointType MovingImagePointType
Definition: itkImageToImageMetric.h:98
itk::ImageToImageMetric::FixedImageType
TFixedImage FixedImageType
Definition: itkImageToImageMetric.h:75
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::ImageToImageMetric::IndexValueType
typename BSplineTransformIndexArrayType::ValueType IndexValueType
Definition: itkImageToImageMetric.h:454
itk::ThreadIdType
unsigned int ThreadIdType
Definition: itkIntTypes.h:99
itk::MeanSquaresImageToImageMetric::m_PerThread
AlignedPerThreadType * m_PerThread
Definition: itkMeanSquaresImageToImageMetric.h:133
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:59
itkImageToImageMetric.h
itk::ImageToImageMetric::MovingImageConstPointer
typename MovingImageType::ConstPointer MovingImageConstPointer
Definition: itkImageToImageMetric.h:72
itk::MeanSquaresImageToImageMetric::PerThreadS::m_Jacobian
TransformJacobianType m_Jacobian
Definition: itkMeanSquaresImageToImageMetric.h:127
itk::ImageToImageMetric::FixedImageSampleContainer
std::vector< FixedImageSamplePoint > FixedImageSampleContainer
Definition: itkImageToImageMetric.h:380
itkIndex.h
itk::MeanSquaresImageToImageMetric
TODO.
Definition: itkMeanSquaresImageToImageMetric.h:39
itk::MeanSquaresImageToImageMetric::PerThreadS::m_MSE
MeasureType m_MSE
Definition: itkMeanSquaresImageToImageMetric.h:128
itk::CovariantVector
A templated class holding a n-Dimensional covariant vector.
Definition: itkCovariantVector.h:70
itk::ImageToImageMetric::TransformPointer
typename TransformType::Pointer TransformPointer
Definition: itkImageToImageMetric.h:87
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ImageToImageMetric::MovingImageType
TMovingImage MovingImageType
Definition: itkImageToImageMetric.h:70
itk::IndexValueType
signed long IndexValueType
Definition: itkIntTypes.h:90
itk::Array
Array class with size defined at construction time.
Definition: itkArray.h:47
itk::Transform
Transform points and vectors from an input space to an output space.
Definition: itkTransform.h:82
itk::ImageToImageMetric::CoordinateRepresentationType
typename Superclass::ParametersValueType CoordinateRepresentationType
Definition: itkImageToImageMetric.h:64
itk::InterpolateImageFunction
Base class for all image interpolators.
Definition: itkInterpolateImageFunction.h:45
itk::SizeValueType
unsigned long SizeValueType
Definition: itkIntTypes.h:83