SimpleITK  
Python/FFTConvolution.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 """ Example script showing Fast Fourier Transform convolution in SimpleITK """
21 
22 import sys
23 import os
24 import math
25 import SimpleITK as sitk
26 
27 if len(sys.argv) < 4:
28  print("Usage: FFTConvolution <input> <kernel> <output>")
29  sys.exit(1)
30 
31 inputFileName = sys.argv[1]
32 kernelFileName = sys.argv[2]
33 outputFileName = sys.argv[3]
34 
35 # -- Input Image --
36 # read the input image
37 img = sitk.ReadImage(inputFileName)
38 
39 # save the type of pixel the input is, so that we can cast the result
40 # back out to the same type
41 pixelID = img.GetPixelID()
42 
43 # pad the image
44 img = sitk.MirrorPad(img, [128] * 2, [128] * 2)
45 size = img.GetSize()
46 
47 # perform the FFT
48 fftimg = sitk.ForwardFFT(sitk.Cast(img, sitk.sitkFloat32))
49 
50 # -- Kernel Image --
51 # Read the kernel image file
52 kernel = sitk.ReadImage(kernelFileName)
53 
54 # flip kernel about all axis
55 kernel = sitk.Flip(kernel, [True] * 2)
56 
57 # normalize the kernel to sum to ~1
59 stats.Execute(kernel)
60 kernel = sitk.Cast(kernel / stats.GetSum(), sitk.sitkFloat32)
61 
62 upadding = [0] * 2
63 upadding[0] = int(math.floor((size[0] - kernel.GetSize()[0]) / 2.0))
64 upadding[1] = int(math.floor((size[1] - kernel.GetSize()[1]) / 2.0))
65 
66 lpadding = [0] * 2
67 lpadding[0] = int(math.ceil((size[0] - kernel.GetSize()[0]) / 2.0))
68 lpadding[1] = int(math.ceil((size[1] - kernel.GetSize()[1]) / 2.0))
69 
70 # pad the kernel to prevent edge artifacts
71 kernel = sitk.ConstantPad(kernel, upadding, lpadding, 0.0)
72 
73 # perform FFT on kernel
74 fftkernel = sitk.ForwardFFT(sitk.FFTShift(kernel))
75 
76 # meta-data must match for multiplication
77 fftkernel.SetSpacing(fftimg.GetSpacing())
78 fftkernel.SetOrigin(fftimg.GetOrigin())
79 fftkernel.SetDirection(fftimg.GetDirection())
80 
81 # -- Convolution --
82 # Finally perform the convolution in Fourier space by multiplication
83 img = sitk.InverseFFT(fftimg * fftkernel)
84 
85 # remove the padding
86 img = sitk.Crop(img, [128] * 2, [128] * 2)
87 
88 # -- Writing --
89 # write the output image the same type as the input
90 sitk.WriteImage(sitk.Cast(img, pixelID), outputFileName)
91 
92 if "SITK_NOSHOW" not in os.environ:
93  sitk.Show(sitk.ReadImage(inputFileName), "original")
94  sitk.Show(sitk.ReadImage(kernelFileName), "kernel")
95  sitk.Show(sitk.Cast(img, pixelID), "FFT_Convolution")
itk::simple::StatisticsImageFilter
Compute min, max, variance and mean of an Image .
Definition: sitkStatisticsImageFilter.h:45
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::Crop
Image Crop(const Image &image1, std::vector< unsigned int > lowerBoundaryCropSize=std::vector< unsigned int >(3, 0), std::vector< unsigned int > upperBoundaryCropSize=std::vector< unsigned int >(3, 0))
Decrease the image size by cropping the image by an itk::Size at both the upper and lower bounds of t...
itk::simple::Flip
Image Flip(const Image &image1, std::vector< bool > flipAxes=std::vector< bool >(3, false), bool flipAboutOrigin=false)
Flips an image across user specified axes.
itk::simple::ForwardFFT
Image ForwardFFT(const Image &image1)
Base class for forward Fast Fourier Transform .
itk::simple::InverseFFT
Image InverseFFT(const Image &image1)
Base class for inverse Fast Fourier Transform .
itk::simple::WriteImage
SITKIO_EXPORT void WriteImage(const Image &image, const std::vector< PathType > &fileNames, bool useCompression=false, int compressionLevel=-1)
WriteImage is a procedural interface to the ImageSeriesWriter. class which is convenient for many ima...
itk::simple::MirrorPad
Image MirrorPad(const Image &image1, std::vector< unsigned int > padLowerBound=std::vector< unsigned int >(3, 0), std::vector< unsigned int > padUpperBound=std::vector< unsigned int >(3, 0), double decayBase=1.0)
Increase the image size by padding with replicants of the input image value.
itk::simple::ConstantPad
Image ConstantPad(const Image &image1, std::vector< unsigned int > padLowerBound=std::vector< unsigned int >(3, 0), std::vector< unsigned int > padUpperBound=std::vector< unsigned int >(3, 0), double constant=0.0)
Increase the image size by padding with a constant value.
itk::simple::FFTShift
Image FFTShift(const Image &image1, bool inverse=false)
Shift the zero-frequency components of a Fourier transform to the center of the image.
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...