Loading [MathJax]/extensions/tex2jax.js
SimpleITK  2.2.0
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 import SimpleITK as sitk
21 import sys
22 import os
23 
24 
25 def command_iteration(method):
26  print(
27  f"{method.GetOptimizerIteration():3}"
28  + f" = {method.GetMetricValue():7.5f}"
29  + f" : {method.GetOptimizerPosition()}"
30  )
31 
32 
33 if len(sys.argv) < 4:
34  print(
35  "Usage:",
36  sys.argv[0],
37  "<fixedImageFilter> <movingImageFile> <outputTransformFile>",
38  )
39  sys.exit(1)
40 
41 pixelType = sitk.sitkFloat32
42 
43 fixed = sitk.ReadImage(sys.argv[1], sitk.sitkFloat32)
44 fixed = sitk.Normalize(fixed)
45 fixed = sitk.DiscreteGaussian(fixed, 2.0)
46 
47 moving = sitk.ReadImage(sys.argv[2], sitk.sitkFloat32)
48 moving = sitk.Normalize(moving)
49 moving = sitk.DiscreteGaussian(moving, 2.0)
50 
52 
53 R.SetMetricAsJointHistogramMutualInformation()
54 
55 R.SetOptimizerAsGradientDescentLineSearch(
56  learningRate=1.0,
57  numberOfIterations=200,
58  convergenceMinimumValue=1e-5,
59  convergenceWindowSize=5,
60 )
61 
62 R.SetInitialTransform(sitk.TranslationTransform(fixed.GetDimension()))
63 
64 R.SetInterpolator(sitk.sitkLinear)
65 
66 R.AddCommand(sitk.sitkIterationEvent, lambda: command_iteration(R))
67 
68 outTx = R.Execute(fixed, moving)
69 
70 print("-------")
71 print(outTx)
72 print(f"Optimizer stop condition: {R.GetOptimizerStopConditionDescription()}")
73 print(f" Iteration: {R.GetOptimizerIteration()}")
74 print(f" Metric value: {R.GetMetricValue()}")
75 
76 sitk.WriteTransform(outTx, sys.argv[3])
77 
78 if "SITK_NOSHOW" not in os.environ:
79  resampler = sitk.ResampleImageFilter()
80  resampler.SetReferenceImage(fixed)
81  resampler.SetInterpolator(sitk.sitkLinear)
82  resampler.SetDefaultPixelValue(1)
83  resampler.SetTransform(outTx)
84 
85  out = resampler.Execute(moving)
86 
87  simg1 = sitk.Cast(sitk.RescaleIntensity(fixed), sitk.sitkUInt8)
88  simg2 = sitk.Cast(sitk.RescaleIntensity(out), sitk.sitkUInt8)
89  cimg = sitk.Compose(simg1, simg2, simg1 // 2.0 + simg2 // 2.0)
90  sitk.Show(cimg, "ImageRegistration2 Composition")
itk::simple::WriteTransform
SITKCommon_EXPORT void WriteTransform(const Transform &transform, const std::string &filename)
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::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::TranslationTransform
Translation of a 2D or 3D coordinate space.
Definition: sitkTranslationTransform.h:35
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: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