Compute Mean Squares Metric Between Two Images#

Synopsis#

Compute the mean squares metric between two images.

Results#

input image

Input image 1.#

input image

Input image 2.#

Data from 2 images passed through Output:

[-10, -10]: 23101.7
[-10, -5]: 23205.7
[-10, 0]: 23260.4
[-10, 5]: 23064.5
[-10, 10]: 22914.5
[-5, -10]: 23271.1
[-5, -5]: 23351.3
[-5, 0]: 23401
[-5, 5]: 23185.1
[-5, 10]: 23026.5
[0, -10]: 23486.5
[0, -5]: 23538.2
[0, 0]: 23566.2
[0, 5]: 23352.1
[0, 10]: 23175.2
[5, -10]: 23590.7
[5, -5]: 23625.7
[5, 0]: 23633.4
[5, 5]: 23401.1
[5, 10]: 23196.7
[10, -10]: 23723.5
[10, -5]: 23762.9
[10, 0]: 23767.1
[10, 5]: 23504.9
[10, 10]: 23298.3

Code#

C++#

#include "itkImage.h"
#include "itkImageFileReader.h"
#include "itkMeanSquaresImageToImageMetric.h"
#include "itkLinearInterpolateImageFunction.h"
#include "itkTranslationTransform.h"

int
main(int argc, char * argv[])
{
  using ImageType = itk::Image<double, 2>;

  if (argc < 3)
  {
    std::cout << "Usage: " << argv[0] << " imageFile1 imageFile2" << std::endl;
    return EXIT_FAILURE;
  }

  ImageType::Pointer fixedImage = itk::ReadImage<ImageType>(argv[1]);
  ImageType::Pointer movingImage = itk::ReadImage<ImageType>(argv[2]);

  using MetricType = itk::MeanSquaresImageToImageMetric<ImageType, ImageType>;
  using InterpolatorType = itk::LinearInterpolateImageFunction<ImageType, double>;
  using TransformType = itk::TranslationTransform<double, 2>;

  auto metric = MetricType::New();
  auto transform = TransformType::New();

  auto interpolator = InterpolatorType::New();
  interpolator->SetInputImage(fixedImage);

  metric->SetFixedImage(fixedImage);
  metric->SetMovingImage(movingImage);
  metric->SetFixedImageRegion(fixedImage->GetLargestPossibleRegion());
  metric->SetTransform(transform);
  metric->SetInterpolator(interpolator);

  TransformType::ParametersType params(transform->GetNumberOfParameters());
  params.Fill(0.0);

  metric->Initialize();
  for (double x = -10.0; x <= 10.0; x += 5.0)
  {
    params(0) = x;
    for (double y = -10.0; y <= 10.0; y += 5.0)
    {
      params(1) = y;
      std::cout << params << ": " << metric->GetValue(params) << std::endl;
    }
  }

  return EXIT_SUCCESS;
}

Python …

import itk
import argparse

parser = argparse.ArgumentParser(description="Compute Mean Square Between Two Images.")
parser.add_argument("input_image_1")
parser.add_argument("input_image_2")
args = parser.parse_args()

fixed_image = itk.imread(args.input_image_1, itk.F)
moving_image = itk.imread(args.input_image_2, itk.F)

metric = itk.MeanSquaresImageToImageMetric[type(fixed_image), type(moving_image)].New()
transform = itk.TranslationTransform[itk.D, fixed_image.GetImageDimension()].New()
interpolator = itk.LinearInterpolateImageFunction[type(fixed_image), itk.D].New()

metric.SetFixedImage(fixed_image)
metric.SetMovingImage(moving_image)
metric.SetFixedImageRegion(fixed_image.GetLargestPossibleRegion())
metric.SetTransform(transform)
metric.SetInterpolator(interpolator)
metric.Initialize()

params = itk.OptimizerParameters[itk.D]()
params.SetSize(2)

for x in range(-10, 15, 5):
    params.SetElement(0, x)
    for y in range(-10, 15, 5):
        params.SetElement(1, y)
        print(f"{list(params)}: {metric.GetValue(params):.1f}")

Classes demonstrated#

template<typename TFixedImage, typename TMovingImage>
class MeanSquaresImageToImageMetric : public itk::ImageToImageMetric<TFixedImage, TMovingImage>

TODO.

ITK Sphinx Examples:

See itk::MeanSquaresImageToImageMetric for additional documentation.