ITK  6.0.0
Insight Toolkit
Examples/RegistrationITKv4/ImageRegistration4.py
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 # https://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 
19 import itk
20 from sys import argv
21 
22 
23 #
24 # Check input parameters
25 # INPUTS(fixedImage): {BrainProtonDensitySliceBorder20.png}
26 # INPUTS(movingImage): {BrainProtonDensitySliceShifted13x17y.png}
27 #
28 if len(argv) < 4:
29  print("Missing Parameters")
30  print(
31  "Usage: ImageRegistration4.py fixedImageFile movingImageFile outputImagefile"
32  )
33  exit()
34 
35 
36 #
37 # Define data types
38 #
39 FixedImageType = itk.Image[itk.F, 2]
40 MovingImageType = itk.Image[itk.F, 2]
41 TransformType = itk.TranslationTransform[itk.D, 2]
42 OptimizerType = itk.RegularStepGradientDescentOptimizerv4[itk.D]
43 RegistrationType = itk.ImageRegistrationMethodv4[FixedImageType, MovingImageType]
45  FixedImageType, MovingImageType
46 ]
47 
48 
49 #
50 # Read the fixed and moving images using filenames
51 # from the command line arguments
52 #
53 fixedImageReader = itk.ImageFileReader[FixedImageType].New()
54 movingImageReader = itk.ImageFileReader[MovingImageType].New()
55 
56 fixedImageReader.SetFileName(argv[1])
57 movingImageReader.SetFileName(argv[2])
58 
59 fixedImageReader.Update()
60 movingImageReader.Update()
61 
62 fixedImage = fixedImageReader.GetOutput()
63 movingImage = movingImageReader.GetOutput()
64 
65 
66 #
67 # Instantiate the classes for the registration framework
68 #
69 registration = RegistrationType.New()
70 imageMetric = MetricType.New()
71 transform = TransformType.New()
72 optimizer = OptimizerType.New()
73 
74 registration.SetOptimizer(optimizer)
75 registration.SetMetric(imageMetric)
76 
77 numberOfBins = 24
78 
79 imageMetric.SetNumberOfHistogramBins(numberOfBins)
80 imageMetric.SetUseMovingImageGradientFilter(False)
81 imageMetric.SetUseFixedImageGradientFilter(False)
82 
83 registration.SetFixedImage(fixedImage)
84 registration.SetMovingImage(movingImage)
85 
86 registration.SetInitialTransform(transform)
87 
88 
89 #
90 # Define optimizer parameters
91 #
92 optimizer.SetLearningRate(8.00)
93 optimizer.SetMinimumStepLength(0.001)
94 optimizer.SetNumberOfIterations(100)
95 optimizer.ReturnBestParametersAndValueOn()
96 optimizer.SetRelaxationFactor(0.8)
97 
98 
99 #
100 # One level registration process without shrinking and smoothing.
101 #
102 registration.SetNumberOfLevels(1)
103 registration.SetSmoothingSigmasPerLevel([0])
104 registration.SetShrinkFactorsPerLevel([1])
105 
106 registration.SetMetricSamplingStrategy(RegistrationType.RANDOM)
107 registration.SetMetricSamplingPercentage(0.20)
108 
109 
110 #
111 # Iteration Observer
112 #
113 def iterationUpdate():
114  currentParameter = registration.GetOutput().Get().GetParameters()
115  print(
116  "M: %f P: %f %f "
117  % (
118  optimizer.GetValue(),
119  currentParameter.GetElement(0),
120  currentParameter.GetElement(1),
121  )
122  )
123 
124 
125 iterationCommand = itk.PyCommand.New()
126 iterationCommand.SetCommandCallable(iterationUpdate)
127 optimizer.AddObserver(itk.IterationEvent(), iterationCommand)
128 
129 print("Starting registration")
130 
131 
132 #
133 # Start the registration process
134 #
135 registration.Update()
136 
137 
138 #
139 # Get the final parameters of the transformation
140 #
141 finalParameters = registration.GetOutput().Get().GetParameters()
142 
143 print("Final Registration Parameters ")
144 print(f"Translation X = {finalParameters.GetElement(0):f}")
145 print(f"Translation Y = {finalParameters.GetElement(1):f}")
146 
147 
148 #
149 # Now, we use the final transform for resampling the
150 # moving image.
151 #
152 resampler = itk.ResampleImageFilter[MovingImageType, FixedImageType].New()
153 resampler.SetTransform(registration.GetTransform())
154 resampler.SetInput(movingImageReader.GetOutput())
155 
156 region = fixedImage.GetLargestPossibleRegion()
157 
158 resampler.SetSize(region.GetSize())
159 resampler.SetOutputOrigin(fixedImage.GetOrigin())
160 resampler.SetOutputSpacing(fixedImage.GetSpacing())
161 resampler.SetOutputDirection(fixedImage.GetDirection())
162 resampler.SetDefaultPixelValue(100)
163 
164 OutputImageType = itk.Image[itk.UC, 2]
165 outputCast = itk.CastImageFilter[FixedImageType, OutputImageType].New()
166 outputCast.SetInput(resampler.GetOutput())
167 
168 
169 #
170 # Write the resampled image
171 #
172 writer = itk.ImageFileWriter[OutputImageType].New()
173 writer.SetFileName(argv[3])
174 writer.SetInput(outputCast.GetOutput())
175 writer.Update()
itk::CastImageFilter
Casts input pixels to output pixel type.
Definition: itkCastImageFilter.h:100
itk::ImageFileReader
Data source that reads image data from a single file.
Definition: itkImageFileReader.h:75
itk::RegularStepGradientDescentOptimizerv4
Regular Step Gradient descent optimizer.
Definition: itkRegularStepGradientDescentOptimizerv4.h:47
itk::PyCommand::New
static Pointer New()
! Method for creation through the object factory.
itk::ImageFileWriter
Writes image data to a single file.
Definition: itkImageFileWriter.h:90
itk::TranslationTransform
Translation transformation of a vector space (e.g. space coordinates)
Definition: itkTranslationTransform.h:43
itk::ImageRegistrationMethodv4
Interface method for the current registration framework.
Definition: itkImageRegistrationMethodv4.h:117
itk::ResampleImageFilter
Resample an image via a coordinate transform.
Definition: itkResampleImageFilter.h:90
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
New
static Pointer New()
itk::MattesMutualInformationImageToImageMetricv4
Computes the mutual information between two images to be registered using the method of Mattes et al.
Definition: itkMattesMutualInformationImageToImageMetricv4.h:103