SimpleITK  
ImageRegistrationMethodBSpline2/ImageRegistrationMethodBSpline2.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():10.5f}"
29  )
30  print("\t#: ", len(method.GetOptimizerPosition()))
31 
32 
33 def command_multi_iteration(method):
34  print("--------- Resolution Changing ---------")
35 
36 
37 if len(sys.argv) < 4:
38  print(
39  "Usage:",
40  sys.argv[0],
41  "<fixedImageFilter> <movingImageFile>",
42  "<outputTransformFile>",
43  )
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 transformDomainMeshSize = [10] * moving.GetDimension()
51 tx = sitk.BSplineTransformInitializer(fixed, transformDomainMeshSize)
52 
53 print("Initial Parameters:")
54 print(tx.GetParameters())
55 
57 R.SetMetricAsMattesMutualInformation(50)
58 R.SetOptimizerAsGradientDescentLineSearch(
59  5.0, 100, convergenceMinimumValue=1e-4, convergenceWindowSize=5
60 )
61 R.SetOptimizerScalesFromPhysicalShift()
62 R.SetInitialTransform(tx)
63 R.SetInterpolator(sitk.sitkLinear)
64 
65 R.SetShrinkFactorsPerLevel([6, 2, 1])
66 R.SetSmoothingSigmasPerLevel([6, 2, 1])
67 
68 R.AddCommand(sitk.sitkIterationEvent, lambda: command_iteration(R))
69 R.AddCommand(
70  sitk.sitkMultiResolutionIterationEvent, lambda: command_multi_iteration(R)
71 )
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, sys.argv[3])
82 
83 if "SITK_NOSHOW" not in os.environ:
84  resampler = sitk.ResampleImageFilter()
85  resampler.SetReferenceImage(fixed)
86  resampler.SetInterpolator(sitk.sitkLinear)
87  resampler.SetDefaultPixelValue(100)
88  resampler.SetTransform(outTx)
89 
90  out = resampler.Execute(moving)
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  sitk.Show(cimg, "ImageRegistration1 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::BSplineTransformInitializer
BSplineTransform BSplineTransformInitializer(const Image &image1, const std::vector< uint32_t > &transformDomainMeshSize=std::vector< uint32_t >(3, 1u), unsigned int order=3u)
BSplineTransformInitializerFilter is a helper class intended to initialize the control point grid suc...
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:52
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