SimpleITK  2.1.0
ImageRegistrationMethodDisplacement1/ImageRegistrationMethodDisplacement1.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  if (method.GetOptimizerIteration() == 0):
27  print(f"\tLevel: {method.GetCurrentLevel()}")
28  print(f"\tScales: {method.GetOptimizerScales()}")
29  print(f"#{method.GetOptimizerIteration()}")
30  print(f"\tMetric Value: {method.GetMetricValue():10.5f}")
31  print(f"\tLearningRate: {method.GetOptimizerLearningRate():10.5f}")
32  if (method.GetOptimizerConvergenceValue() != sys.float_info.max):
33  print(f"\tConvergence Value: {method.GetOptimizerConvergenceValue():.5e}")
34 
35 
36 def command_multiresolution_iteration(method):
37  print(f"\tStop Condition: {method.GetOptimizerStopConditionDescription()}")
38  print("============= Resolution Change =============")
39 
40 
41 if len(sys.argv) < 4:
42  print("Usage:", sys.argv[0], "<fixedImageFilter> <movingImageFile>",
43  "<outputTransformFile>")
44  sys.exit(1)
45 
46 fixed = sitk.ReadImage(sys.argv[1], sitk.sitkFloat32)
47 
48 moving = sitk.ReadImage(sys.argv[2], sitk.sitkFloat32)
49 
50 initialTx = sitk.CenteredTransformInitializer(fixed, moving,
52  fixed.GetDimension()))
53 
55 
56 R.SetShrinkFactorsPerLevel([3, 2, 1])
57 R.SetSmoothingSigmasPerLevel([2, 1, 1])
58 
59 R.SetMetricAsJointHistogramMutualInformation(20)
60 R.MetricUseFixedImageGradientFilterOff()
61 
62 R.SetOptimizerAsGradientDescent(learningRate=1.0,
63  numberOfIterations=100,
64  estimateLearningRate=R.EachIteration)
65 R.SetOptimizerScalesFromPhysicalShift()
66 
67 R.SetInitialTransform(initialTx)
68 
69 R.SetInterpolator(sitk.sitkLinear)
70 
71 R.AddCommand(sitk.sitkIterationEvent, lambda: command_iteration(R))
72 R.AddCommand(sitk.sitkMultiResolutionIterationEvent,
73  lambda: command_multiresolution_iteration(R))
74 
75 outTx1 = R.Execute(fixed, moving)
76 
77 print("-------")
78 print(outTx1)
79 print(f"Optimizer stop condition: {R.GetOptimizerStopConditionDescription()}")
80 print(f" Iteration: {R.GetOptimizerIteration()}")
81 print(f" Metric value: {R.GetMetricValue()}")
82 
83 displacementField = sitk.Image(fixed.GetSize(), sitk.sitkVectorFloat64)
84 displacementField.CopyInformation(fixed)
85 displacementTx = sitk.DisplacementFieldTransform(displacementField)
86 del displacementField
87 displacementTx.SetSmoothingGaussianOnUpdate(varianceForUpdateField=0.0,
88  varianceForTotalField=1.5)
89 
90 R.SetMovingInitialTransform(outTx1)
91 R.SetInitialTransform(displacementTx, inPlace=True)
92 
93 R.SetMetricAsANTSNeighborhoodCorrelation(4)
94 R.MetricUseFixedImageGradientFilterOff()
95 
96 R.SetShrinkFactorsPerLevel([3, 2, 1])
97 R.SetSmoothingSigmasPerLevel([2, 1, 1])
98 
99 R.SetOptimizerScalesFromPhysicalShift()
100 R.SetOptimizerAsGradientDescent(learningRate=1,
101  numberOfIterations=300,
102  estimateLearningRate=R.EachIteration)
103 
104 R.Execute(fixed, moving)
105 
106 print("-------")
107 print(displacementTx)
108 print(f"Optimizer stop condition: {R.GetOptimizerStopConditionDescription()}")
109 print(f" Iteration: {R.GetOptimizerIteration()}")
110 print(f" Metric value: {R.GetMetricValue()}")
111 
112 compositeTx = sitk.CompositeTransform([outTx1, displacementTx])
113 sitk.WriteTransform(compositeTx, sys.argv[3])
114 
115 if ("SITK_NOSHOW" not in os.environ):
116  sitk.Show(displacementTx.GetDisplacementField(), "Displacement Field")
117 
118  resampler = sitk.ResampleImageFilter()
119  resampler.SetReferenceImage(fixed)
120  resampler.SetInterpolator(sitk.sitkLinear)
121  resampler.SetDefaultPixelValue(100)
122  resampler.SetTransform(compositeTx)
123 
124  out = resampler.Execute(moving)
125  simg1 = sitk.Cast(sitk.RescaleIntensity(fixed), sitk.sitkUInt8)
126  simg2 = sitk.Cast(sitk.RescaleIntensity(out), sitk.sitkUInt8)
127  cimg = sitk.Compose(simg1, simg2, simg1 // 2. + simg2 // 2.)
128  sitk.Show(cimg, "ImageRegistration1 Composition")
itk::simple::Image
The Image class for SimpleITK.
Definition: sitkImage.h:75
itk::simple::WriteTransform
SITKCommon_EXPORT void WriteTransform(const Transform &transform, const std::string &filename)
itk::simple::DisplacementFieldTransform
A dense deformable transform over a bounded spatial domain for 2D or 3D coordinates space.
Definition: sitkDisplacementFieldTransform.h:36
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::Cast
Image Cast(const Image &image, PixelIDValueEnum pixelID)
itk::simple::ResampleImageFilter
Resample an image via a coordinate transform.
Definition: sitkResampleImageFilter.h:53
itk::simple::CompositeTransform
This class contains a stack of transforms and concatenates them by composition.
Definition: sitkCompositeTransform.h:58
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
itk::simple::AffineTransform
An affine transformation about a fixed center with translation for a 2D or 3D coordinate.
Definition: sitkAffineTransform.h:35
itk::simple::CenteredTransformInitializer
Transform CenteredTransformInitializer(const Image &fixedImage, const Image &movingImage, const Transform &transform, CenteredTransformInitializerFilter::OperationModeType operationMode=itk::simple::CenteredTransformInitializerFilter::MOMENTS)
CenteredTransformInitializer is a helper class intended to initialize the center of rotation and the ...