SimpleITK  2.0.0
ImageRegistrationMethodExhaustive/ImageRegistrationMethodExhaustive.py
1 #!/usr/bin/env python
2 # =========================================================================
3 #
4 # Copyright NumFOCUS
5 #
6 # Licensed under the Apache License, Version 2.0 (the "License");
7 # you may not use this file except in compliance with the License.
8 # You may obtain a copy of the License at
9 #
10 # http://www.apache.org/licenses/LICENSE-2.0.txt
11 #
12 # Unless required by applicable law or agreed to in writing, software
13 # distributed under the License is distributed on an "AS IS" BASIS,
14 # WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
15 # See the License for the specific language governing permissions and
16 # limitations under the License.
17 #
18 # =========================================================================
19 
20 """
21 This script demonstrates the use of the Exhaustive optimizer in the
22 ImageRegistrationMethod to estimate a good initial rotation position.
23 
24 Because gradient descent base optimization can get stuck in local
25 minima, a good initial transform is critical for reasonable
26 results. Search a reasonable space on a grid with brute force may be a
27 reliable way to get a starting location for further optimization.
28 
29 The initial translation and center of rotation for the transform is
30 initialized based on the first principle moments of the intensities of
31 the image. Then in either 2D or 3D a Euler transform is used to
32 exhaustively search a grid of the rotation space at a certain step
33 size. The resulting transform is a reasonable guess where to start
34 further registration.
35 """
36 
37 from __future__ import print_function
38 from __future__ import division
39 
40 import SimpleITK as sitk
41 import sys
42 import os
43 from math import pi
44 
45 
46 def command_iteration(method):
47  if (method.GetOptimizerIteration() == 0):
48  print("Scales: ", method.GetOptimizerScales())
49  print("{0:3} = {1:7.5f} : {2}".format(method.GetOptimizerIteration(),
50  method.GetMetricValue(),
51  method.GetOptimizerPosition()))
52 
53 
54 if len(sys.argv) < 4:
55  print("Usage:", sys.argv[0], "<fixedImageFilter> <movingImageFile>",
56  "<outputTransformFile>")
57  sys.exit(1)
58 
59 fixed = sitk.ReadImage(sys.argv[1], sitk.sitkFloat32)
60 
61 moving = sitk.ReadImage(sys.argv[2], sitk.sitkFloat32)
62 
64 
65 R.SetMetricAsMattesMutualInformation(numberOfHistogramBins=50)
66 
67 sample_per_axis = 12
68 if fixed.GetDimension() == 2:
70  # Set the number of samples (radius) in each dimension, with a
71  # default step size of 1.0
72  R.SetOptimizerAsExhaustive([sample_per_axis // 2, 0, 0])
73  # Utilize the scale to set the step size for each dimension
74  R.SetOptimizerScales([2.0 * pi / sample_per_axis, 1.0, 1.0])
75 elif fixed.GetDimension() == 3:
77  R.SetOptimizerAsExhaustive([sample_per_axis // 2, sample_per_axis // 2,
78  sample_per_axis // 4, 0, 0, 0])
79  R.SetOptimizerScales(
80  [2.0 * pi / sample_per_axis, 2.0 * pi / sample_per_axis,
81  2.0 * pi / sample_per_axis, 1.0, 1.0, 1.0])
82 
83 # Initialize the transform with a translation and the center of
84 # rotation from the moments of intensity.
85 tx = sitk.CenteredTransformInitializer(fixed, moving, tx)
86 
87 R.SetInitialTransform(tx)
88 
89 R.SetInterpolator(sitk.sitkLinear)
90 
91 R.AddCommand(sitk.sitkIterationEvent, lambda: command_iteration(R))
92 
93 outTx = R.Execute(fixed, moving)
94 
95 print("-------")
96 print(outTx)
97 print("Optimizer stop condition: {0}"
98  .format(R.GetOptimizerStopConditionDescription()))
99 print(" Iteration: {0}".format(R.GetOptimizerIteration()))
100 print(" Metric value: {0}".format(R.GetMetricValue()))
101 
102 sitk.WriteTransform(outTx, sys.argv[3])
103 
104 if ("SITK_NOSHOW" not in os.environ):
105  resampler = sitk.ResampleImageFilter()
106  resampler.SetReferenceImage(fixed)
107  resampler.SetInterpolator(sitk.sitkLinear)
108  resampler.SetDefaultPixelValue(1)
109  resampler.SetTransform(outTx)
110 
111  out = resampler.Execute(moving)
112 
113  simg1 = sitk.Cast(sitk.RescaleIntensity(fixed), sitk.sitkUInt8)
114  simg2 = sitk.Cast(sitk.RescaleIntensity(out), sitk.sitkUInt8)
115  cimg = sitk.Compose(simg1, simg2, simg1 // 2. + simg2 // 2.)
116  sitk.Show(cimg, "ImageRegistrationExhaustive Composition")
itk::simple::WriteTransform
SITKCommon_EXPORT void WriteTransform(const Transform &transform, const std::string &filename)
itk::simple::RescaleIntensity
Image RescaleIntensity(const Image &image1, double outputMinimum=0, double outputMaximum=255)
Applies a linear transformation to the intensity levels of the input Image .
itk::simple::Show
void SITKIO_EXPORT Show(const Image &image, const std::string &title="", const bool debugOn=ProcessObject::GetGlobalDefaultDebug())
itk::simple::ReadImage
SITKIO_EXPORT Image ReadImage(const std::vector< std::string > &fileNames, PixelIDValueEnum outputPixelType=sitkUnknown, const std::string &imageIO="")
ReadImage is a procedural interface to the ImageSeriesReader class which is convenient for most image...
itk::simple::Cast
Image Cast(const Image &image, PixelIDValueEnum pixelID)
itk::simple::Euler2DTransform
A rigid 2D transform with rotation in radians around a fixed center with translation.
Definition: sitkEuler2DTransform.h:35
itk::simple::ResampleImageFilter
Resample an image via a coordinate transform.
Definition: sitkResampleImageFilter.h:53
itk::simple::Compose
Image Compose(const Image &image1, const Image &image2, const Image &image3, const Image &image4, const Image &image5)
ComposeImageFilter combine several scalar images into a multicomponent image.
itk::simple::ImageRegistrationMethod
An interface method to the modular ITKv4 registration framework.
Definition: sitkImageRegistrationMethod.h:89
itk::simple::Euler3DTransform
A rigid 3D transform with rotation in radians around a fixed center with translation.
Definition: sitkEuler3DTransform.h:35
itk::simple::CenteredTransformInitializer
Transform CenteredTransformInitializer(const Image &fixedImage, const Image &movingImage, const Transform &transform, CenteredTransformInitializerFilter::OperationModeType operationMode=itk::simple::CenteredTransformInitializerFilter::MOMENTS)
CenteredTransformInitializer is a helper class intended to initialize the center of rotation and the ...