ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkMeanSquaresImageToImageMetric.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 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:
40  public ImageToImageMetric< TFixedImage, TMovingImage >
41 {
42 public:
43  ITK_DISALLOW_COPY_AND_ASSIGN(MeanSquaresImageToImageMetric);
44 
50 
52  itkNewMacro(Self);
53 
56 
58  using TransformType = typename Superclass::TransformType;
59  using TransformPointer = typename Superclass::TransformPointer;
60  using TransformJacobianType = typename Superclass::TransformJacobianType;
61  using InterpolatorType = typename Superclass::InterpolatorType;
62  using MeasureType = typename Superclass::MeasureType;
63  using DerivativeType = typename Superclass::DerivativeType;
64  using ParametersType = typename Superclass::ParametersType;
65  using FixedImageType = typename Superclass::FixedImageType;
66  using MovingImageType = typename Superclass::MovingImageType;
67  using MovingImagePointType = typename Superclass::MovingImagePointType;
68  using FixedImageConstPointer = typename Superclass::FixedImageConstPointer;
69  using MovingImageConstPointer = typename Superclass::MovingImageConstPointer;
70  using CoordinateRepresentationType = typename Superclass::CoordinateRepresentationType;
71  using FixedImageSampleContainer = typename Superclass::FixedImageSampleContainer;
72  using ImageDerivativesType = typename Superclass::ImageDerivativesType;
73  using WeightsValueType = typename Superclass::WeightsValueType;
75 
76  // Needed for evaluation of Jacobian.
77  using FixedImagePointType = typename Superclass::FixedImagePointType;
78 
80  static constexpr unsigned int MovingImageDimension = MovingImageType::ImageDimension;
81 
89  void Initialize() override;
90 
92  MeasureType GetValue(const ParametersType & parameters) const override;
93 
95  void GetDerivative(const ParametersType & parameters,
96  DerivativeType & Derivative) const override;
97 
99  void GetValueAndDerivative(const ParametersType & parameters,
100  MeasureType & Value,
101  DerivativeType & Derivative) const override;
102 
103 protected:
104 
106  ~MeanSquaresImageToImageMetric() override;
107  void PrintSelf(std::ostream & os, Indent indent) const override;
108 
109 private:
110 
111  bool GetValueThreadProcessSample(ThreadIdType threadId,
112  SizeValueType fixedImageSample,
113  const MovingImagePointType & mappedPoint,
114  double movingImageValue) const override;
115 
116  bool GetValueAndDerivativeThreadProcessSample(ThreadIdType threadId,
117  SizeValueType fixedImageSample,
118  const MovingImagePointType & mappedPoint,
119  double movingImageValue,
120  const ImageDerivativesType &
121  movingImageGradientValue) const override;
122 
123  struct PerThreadS
124  {
128  };
129 
130  itkAlignedTypedef( 64, PerThreadS, AlignedPerThreadType );
131  AlignedPerThreadType *m_PerThread;
132 };
133 } // end namespace itk
134 
135 #ifndef ITK_MANUAL_INSTANTIATION
136 #include "itkMeanSquaresImageToImageMetric.hxx"
137 #endif
138 
139 #endif
Array class with size defined at construction time.
Definition: itkArray.h:46
typename BSplineTransformWeightsType::ValueType WeightsValueType
Light weight base class for most itk classes.
typename TransformType::InputPointType FixedImagePointType
typename TransformType::JacobianType TransformJacobianType
unsigned long SizeValueType
Definition: itkIntTypes.h:83
Transform points and vectors from an input space to an output space.
Definition: itkTransform.h:83
signed long IndexValueType
Definition: itkIntTypes.h:90
typename Superclass::ParametersValueType CoordinateRepresentationType
std::vector< FixedImageSamplePoint > FixedImageSampleContainer
typename FixedImageType::ConstPointer FixedImageConstPointer
unsigned int ThreadIdType
Definition: itkIntTypes.h:99
Base class for all image interpolaters.
typename TransformType::Pointer TransformPointer
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Computes similarity between regions of two images.
A templated class holding a n-Dimensional covariant vector.
typename MovingImageType::ConstPointer MovingImageConstPointer
typename BSplineTransformIndexArrayType::ValueType IndexValueType
typename TransformType::OutputPointType MovingImagePointType