SimpleITK  
ImageRegistrationMethodBSpline1/ImageRegistrationMethodBSpline1.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(f"{method.GetOptimizerIteration():3} " + f"= {method.GetMetricValue():10.5f}")
27 
28 
29 def main(args):
30  if len(args) < 4:
31  print(
32  "Usage:",
33  sys.args[0],
34  "<fixedImageFilter> <movingImageFile>",
35  "<outputTransformFile>",
36  )
37  sys.exit(1)
38 
39  fixed = sitk.ReadImage(args[1], sitk.sitkFloat32)
40 
41  moving = sitk.ReadImage(args[2], sitk.sitkFloat32)
42 
43  transformDomainMeshSize = [8] * moving.GetDimension()
44  tx = sitk.BSplineTransformInitializer(fixed, transformDomainMeshSize)
45 
46  print("Initial Parameters:")
47  print(tx.GetParameters())
48 
50  R.SetMetricAsCorrelation()
51 
52  R.SetOptimizerAsLBFGSB(
53  gradientConvergenceTolerance=1e-5,
54  numberOfIterations=100,
55  maximumNumberOfCorrections=5,
56  maximumNumberOfFunctionEvaluations=1000,
57  costFunctionConvergenceFactor=1e7,
58  )
59  R.SetInitialTransform(tx, True)
60  R.SetInterpolator(sitk.sitkLinear)
61 
62  R.AddCommand(sitk.sitkIterationEvent, lambda: command_iteration(R))
63 
64  outTx = R.Execute(fixed, moving)
65 
66  print("-------")
67  print(outTx)
68  print(f"Optimizer stop condition: {R.GetOptimizerStopConditionDescription()}")
69  print(f" Iteration: {R.GetOptimizerIteration()}")
70  print(f" Metric value: {R.GetMetricValue()}")
71 
72  sitk.WriteTransform(outTx, args[3])
73 
74  resampler = sitk.ResampleImageFilter()
75  resampler.SetReferenceImage(fixed)
76  resampler.SetInterpolator(sitk.sitkLinear)
77  resampler.SetDefaultPixelValue(100)
78  resampler.SetTransform(outTx)
79 
80  out = resampler.Execute(moving)
81  simg1 = sitk.Cast(sitk.RescaleIntensity(fixed), sitk.sitkUInt8)
82  simg2 = sitk.Cast(sitk.RescaleIntensity(out), sitk.sitkUInt8)
83  cimg = sitk.Compose(simg1, simg2, simg1 // 2.0 + simg2 // 2.0)
84 
85  return_images = {"fixed": fixed, "moving": moving, "composition": cimg}
86  return return_images
87 
88 
89 if __name__ == "__main__":
90  return_dict = main(sys.argv)
91  if "SITK_NOSHOW" not in os.environ:
92  sitk.Show(
93  return_dict["composition"], "ImageRegistrationMethodBSpline1 Composition"
94  )
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::Cast
Image Cast(const Image &image, PixelIDValueEnum pixelID)
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...