SimpleITK  
ImageRegistrationMethod2/ImageRegistrationMethod2.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 with histogram
21  mutual information as a similarity measure. """
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 process is running. """
30  print(
31  f"{method.GetOptimizerIteration():3} "
32  + f" = {method.GetMetricValue():7.5f} "
33  + f" : {method.GetOptimizerPosition()}"
34  )
35 
36 
37 def main(args):
38  """ A SimpleITK example demonstrating image registration with histogram
39  mutual information as a similarity measure. """
40  if len(args) < 3:
41  print(
42  "Usage:",
43  "ImageRegistrationMethod2",
44  "<fixedImageFilter> <movingImageFile> <outputTransformFile>",
45  )
46  sys.exit(1)
47 
48  fixed = sitk.ReadImage(args[1], sitk.sitkFloat32)
49  fixed = sitk.Normalize(fixed)
50  fixed = sitk.DiscreteGaussian(fixed, 2.0)
51 
52  moving = sitk.ReadImage(args[2], sitk.sitkFloat32)
53  moving = sitk.Normalize(moving)
54  moving = sitk.DiscreteGaussian(moving, 2.0)
55 
57 
58  R.SetMetricAsJointHistogramMutualInformation()
59 
60  R.SetOptimizerAsGradientDescentLineSearch(
61  learningRate=1.0,
62  numberOfIterations=200,
63  convergenceMinimumValue=1e-5,
64  convergenceWindowSize=5,
65  )
66 
67  R.SetInitialTransform(sitk.TranslationTransform(fixed.GetDimension()))
68 
69  R.SetInterpolator(sitk.sitkLinear)
70 
71  R.AddCommand(sitk.sitkIterationEvent, lambda: command_iteration(R))
72 
73  outTx = R.Execute(fixed, moving)
74 
75  print("-------")
76  print(outTx)
77  print(f"Optimizer stop condition: {R.GetOptimizerStopConditionDescription()}")
78  print(f" Iteration: {R.GetOptimizerIteration()}")
79  print(f" Metric value: {R.GetMetricValue()}")
80 
81  sitk.WriteTransform(outTx, args[3])
82 
83  resampler = sitk.ResampleImageFilter()
84  resampler.SetReferenceImage(fixed)
85  resampler.SetInterpolator(sitk.sitkLinear)
86  resampler.SetDefaultPixelValue(1)
87  resampler.SetTransform(outTx)
88 
89  out = resampler.Execute(moving)
90 
91  simg1 = sitk.Cast(sitk.RescaleIntensity(fixed), sitk.sitkUInt8)
92  simg2 = sitk.Cast(sitk.RescaleIntensity(out), sitk.sitkUInt8)
93  cimg = sitk.Compose(simg1, simg2, simg1 // 2.0 + simg2 // 2.0)
94 
95  return {"fixed": fixed, "moving": moving, "composition": cimg}
96 
97 
98 if __name__ == "__main__":
99  return_dict = main(sys.argv)
100  if "SITK_NOSHOW" not in os.environ:
101  sitk.Show(return_dict["composition"], "ImageRegistration2 Composition")
itk::simple::Normalize
Image Normalize(const Image &image1)
Normalize an image by setting its mean to zero and variance to one.
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::DiscreteGaussian
Image DiscreteGaussian(const Image &image1, std::vector< double > variance=std::vector< double >(3, 1.0), unsigned int maximumKernelWidth=32u, std::vector< double > maximumError=std::vector< double >(3, 0.01), bool useImageSpacing=true)
Blurs an image by separable convolution with discrete gaussian kernels. This filter performs Gaussian...
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...