ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkMutualInformationImageToImageMetric.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 itkMutualInformationImageToImageMetric_h
19 #define itkMutualInformationImageToImageMetric_h
20 
21 #include "itkImageToImageMetric.h"
22 #include "itkPoint.h"
23 
24 #include "itkIndex.h"
25 #include "itkKernelFunctionBase.h"
26 
27 namespace itk
28 {
93 template< typename TFixedImage, typename TMovingImage >
94 class ITK_TEMPLATE_EXPORT MutualInformationImageToImageMetric:
95  public ImageToImageMetric< TFixedImage, TMovingImage >
96 {
97 public:
98  ITK_DISALLOW_COPY_AND_ASSIGN(MutualInformationImageToImageMetric);
99 
105 
107  itkNewMacro(Self);
108 
111 
113  using TransformType = typename Superclass::TransformType;
114  using TransformPointer = typename Superclass::TransformPointer;
115  using TransformJacobianType = typename Superclass::TransformJacobianType;
116  using InterpolatorType = typename Superclass::InterpolatorType;
117  using MeasureType = typename Superclass::MeasureType;
118  using DerivativeType = typename Superclass::DerivativeType;
119  using ParametersType = typename Superclass::ParametersType;
120  using FixedImageType = typename Superclass::FixedImageType;
121  using MovingImageType = typename Superclass::MovingImageType;
122  using FixedImageConstPointer = typename Superclass::FixedImageConstPointer;
123  using MovingImageCosntPointer = typename Superclass::MovingImageConstPointer;
124 
131 
133 
135  static constexpr unsigned int MovingImageDimension = MovingImageType::ImageDimension;
136 
138  void GetDerivative(
139  const ParametersType & parameters,
140  DerivativeType & Derivative) const override;
141 
143  MeasureType GetValue(const ParametersType & parameters) const override;
144 
146  void GetValueAndDerivative(const ParametersType & parameters,
147  MeasureType & Value, DerivativeType & Derivative) const override;
148 
153  void SetNumberOfSpatialSamples(unsigned int num);
154 
156  itkGetConstReferenceMacro(NumberOfSpatialSamples, unsigned int);
157 
163  itkSetClampMacro( MovingImageStandardDeviation, double,
165  itkGetConstReferenceMacro(MovingImageStandardDeviation, double);
167 
173  itkSetClampMacro( FixedImageStandardDeviation, double,
175  itkGetConstMacro(FixedImageStandardDeviation, double);
177 
180  itkSetObjectMacro(KernelFunction, KernelFunctionType);
181  itkGetModifiableObjectMacro(KernelFunction, KernelFunctionType);
183 
184 
185 protected:
187  ~MutualInformationImageToImageMetric() override = default;
188  void PrintSelf(std::ostream & os, Indent indent) const override;
189 
190 private:
197  {
198 public:
200  { FixedImagePointValue.Fill(0.0); }
201  ~SpatialSample()= default;
202 
204  double FixedImageValue{0.0};
205  double MovingImageValue{0.0};
206  };
207 
209  using SpatialSampleContainer = std::vector< SpatialSample >;
210 
214 
218 
223 
225 
231  virtual void SampleFixedImageDomain(SpatialSampleContainer & samples) const;
232 
236  void CalculateDerivatives(const FixedImagePointType &, DerivativeType &, TransformJacobianType &) const;
237 
238  using CoordinateRepresentationType = typename Superclass::CoordinateRepresentationType;
241 
243 
244 };
245 } // end namespace itk
246 
247 #ifndef ITK_MANUAL_INSTANTIATION
248 #include "itkMutualInformationImageToImageMetric.hxx"
249 #endif
250 
251 #endif
Array class with size defined at construction time.
Definition: itkArray.h:46
Light weight base class for most itk classes.
typename TransformType::InputPointType FixedImagePointType
Define numeric traits for std::vector.
typename FixedImageIndexType::IndexValueType FixedImageIndexValueType
typename TransformType::JacobianType TransformJacobianType
Calculate the derivative by central differencing.
typename FixedImageType::IndexType FixedImageIndexType
Transform points and vectors from an input space to an output space.
Definition: itkTransform.h:83
Computes the mutual information between two images to be registered.
signed long IndexValueType
Definition: itkIntTypes.h:90
typename Superclass::ParametersValueType CoordinateRepresentationType
typename FixedImageType::ConstPointer FixedImageConstPointer
Base class for all image interpolaters.
typename TransformType::Pointer TransformPointer
typename MovingImageType::IndexType MovingImageIndexType
Control indentation during Print() invocation.
Definition: itkIndent.h:49
typename Superclass::MovingImageConstPointer MovingImageCosntPointer
Computes similarity between regions of two images.
A templated class holding a geometric point in n-Dimensional space.
Definition: itkPoint.h:52
typename TransformType::OutputPointType MovingImagePointType