SimpleITK  2.0.0
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 from __future__ import print_function
21 
22 import SimpleITK as sitk
23 import sys
24 import os
25 
26 if len(sys.argv) < 4:
27  print("Usage:", sys.argv[0], "<fixedImageFilter> <movingImageFile>",
28  "<outputTransformFile> <numberOfBins> <samplingPercentage>")
29  sys.exit(1)
30 
31 
32 def command_iteration(method):
33  print("{0:3} = {1:10.5f} : {2}".format(method.GetOptimizerIteration(),
34  method.GetMetricValue(),
35  method.GetOptimizerPosition()))
36 
37 
38 fixed = sitk.ReadImage(sys.argv[1], sitk.sitkFloat32)
39 moving = sitk.ReadImage(sys.argv[2], sitk.sitkFloat32)
40 
41 numberOfBins = 24
42 samplingPercentage = 0.10
43 
44 if len(sys.argv) > 4:
45  numberOfBins = int(sys.argv[4])
46 if len(sys.argv) > 5:
47  samplingPercentage = float(sys.argv[5])
48 
50 R.SetMetricAsMattesMutualInformation(numberOfBins)
51 R.SetMetricSamplingPercentage(samplingPercentage, sitk.sitkWallClock)
52 R.SetMetricSamplingStrategy(R.RANDOM)
53 R.SetOptimizerAsRegularStepGradientDescent(1.0, .001, 200)
54 R.SetInitialTransform(sitk.TranslationTransform(fixed.GetDimension()))
55 R.SetInterpolator(sitk.sitkLinear)
56 
57 R.AddCommand(sitk.sitkIterationEvent, lambda: command_iteration(R))
58 
59 outTx = R.Execute(fixed, moving)
60 
61 print("-------")
62 print(outTx)
63 print("Optimizer stop condition: {0}"
64  .format(R.GetOptimizerStopConditionDescription()))
65 print(" Iteration: {0}".format(R.GetOptimizerIteration()))
66 print(" Metric value: {0}".format(R.GetMetricValue()))
67 
68 sitk.WriteTransform(outTx, sys.argv[3])
69 
70 if ("SITK_NOSHOW" not in os.environ):
71  resampler = sitk.ResampleImageFilter()
72  resampler.SetReferenceImage(fixed)
73  resampler.SetInterpolator(sitk.sitkLinear)
74  resampler.SetDefaultPixelValue(100)
75  resampler.SetTransform(outTx)
76 
77  out = resampler.Execute(moving)
78  simg1 = sitk.Cast(sitk.RescaleIntensity(fixed), sitk.sitkUInt8)
79  simg2 = sitk.Cast(sitk.RescaleIntensity(out), sitk.sitkUInt8)
80  cimg = sitk.Compose(simg1, simg2, simg1 // 2. + simg2 // 2.)
81  sitk.Show(cimg, "ImageRegistration4 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::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::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