SimpleITK  
ImageRegistrationMethod4/ImageRegistrationMethod4.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 """ A SimpleITK example demonstrating image registration using Mattes mutual
21  information as the metric. """
22 
23 import sys
24 import os
25 import SimpleITK as sitk
26 
27 
28 def command_iteration(method):
29  """ Callback invoked when the optimization has an iteration. """
30  print(
31  f"{method.GetOptimizerIteration():3} "
32  + f"= {method.GetMetricValue():10.5f} "
33  + f": {method.GetOptimizerPosition()}"
34  )
35 
36 
37 def main(args):
38  """ A SimpleITK example demonstrating image registration using Mattes mutual
39  information as the metric. """
40 
41  if len(args) < 3:
42  print(
43  "Usage:",
44  "ImageRegistrationMethod4",
45  "<fixedImageFilter> <movingImageFile>",
46  "<outputTransformFile> <numberOfBins> <samplingPercentage>",
47  )
48  sys.exit(1)
49 
50  fixed = sitk.ReadImage(args[1], sitk.sitkFloat32)
51  moving = sitk.ReadImage(args[2], sitk.sitkFloat32)
52 
53  numberOfBins = 24
54  samplingPercentage = 0.10
55 
56  if len(args) > 4:
57  numberOfBins = int(args[4])
58  if len(args) > 5:
59  samplingPercentage = float(args[5])
60 
62  R.SetMetricAsMattesMutualInformation(numberOfBins)
63  R.SetMetricSamplingPercentage(samplingPercentage, sitk.sitkWallClock)
64  R.SetMetricSamplingStrategy(R.RANDOM)
65  R.SetOptimizerAsRegularStepGradientDescent(1.0, 0.001, 200)
66  R.SetInitialTransform(sitk.TranslationTransform(fixed.GetDimension()))
67  R.SetInterpolator(sitk.sitkLinear)
68 
69  R.AddCommand(sitk.sitkIterationEvent, lambda: command_iteration(R))
70 
71  outTx = R.Execute(fixed, moving)
72 
73  print("-------")
74  print(outTx)
75  print(f"Optimizer stop condition: {R.GetOptimizerStopConditionDescription()}")
76  print(f" Iteration: {R.GetOptimizerIteration()}")
77  print(f" Metric value: {R.GetMetricValue()}")
78 
79  sitk.WriteTransform(outTx, args[3])
80 
81  resampler = sitk.ResampleImageFilter()
82  resampler.SetReferenceImage(fixed)
83  resampler.SetInterpolator(sitk.sitkLinear)
84  resampler.SetDefaultPixelValue(100)
85  resampler.SetTransform(outTx)
86 
87  out = resampler.Execute(moving)
88  simg1 = sitk.Cast(sitk.RescaleIntensity(fixed), sitk.sitkUInt8)
89  simg2 = sitk.Cast(sitk.RescaleIntensity(out), sitk.sitkUInt8)
90  cimg = sitk.Compose(simg1, simg2, simg1 // 2.0 + simg2 // 2.0)
91 
92  return {"fixed": fixed, "moving": moving, "composition": cimg}
93 
94 
95 if __name__ == "__main__":
96  return_dict = main(sys.argv)
97  if "SITK_NOSHOW" not in os.environ:
98  sitk.Show(return_dict["composition"], "ImageRegistration4 Composition")
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::TranslationTransform
Translation of a 2D or 3D coordinate space.
Definition: sitkTranslationTransform.h:33
itk::simple::Cast
Image Cast(const Image &image, PixelIDValueEnum pixelID)
itk::simple::ResampleImageFilter
Resample an image via a coordinate transform.
Definition: sitkResampleImageFilter.h:52
itk::simple::WriteTransform
SITKCommon_EXPORT void WriteTransform(const Transform &transform, const PathType &filename)
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:87
itk::simple::ReadImage
SITKIO_EXPORT Image ReadImage(const std::vector< PathType > &fileNames, PixelIDValueEnum outputPixelType=sitkUnknown, const std::string &imageIO="")
ReadImage is a procedural interface to the ImageSeriesReader class which is convenient for most image...