ITK  5.4.0
Insight Toolkit
SphinxExamples/src/IO/ImageBase/ProcessImageChunks/Code.py
1 #!/usr/bin/env python
2 
3 # Copyright NumFOCUS
4 #
5 # Licensed under the Apache License, Version 2.0 (the "License");
6 # you may not use this file except in compliance with the License.
7 # You may obtain a copy of the License at
8 #
9 # https://www.apache.org/licenses/LICENSE-2.0.txt
10 #
11 # Unless required by applicable law or agreed to in writing, software
12 # distributed under the License is distributed on an "AS IS" BASIS,
13 # WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14 # See the License for the specific language governing permissions and
15 # limitations under the License.
16 
17 import itk
18 import argparse
19 
20 itk.auto_progress(2)
21 
22 parser = argparse.ArgumentParser(description="Process Image Chunks.")
23 parser.add_argument("input_image")
24 parser.add_argument("output_image")
25 args = parser.parse_args()
26 
27 xDiv = 6
28 yDiv = 4
29 zDiv = 1
30 # 1 for 2D input
31 
32 Dimension = 3
33 
34 PixelType = itk.UC
35 ImageType = itk.Image[PixelType, Dimension]
36 
37 ReaderType = itk.ImageFileReader[ImageType]
38 reader = ReaderType.New()
39 reader.SetFileName(args.input_image)
40 reader.UpdateOutputInformation()
41 fullRegion = reader.GetOutput().GetLargestPossibleRegion()
42 fullSize = fullRegion.GetSize()
43 # when reading image from file, start index is always 0
44 
45 # create variables of the proper types
46 start = itk.Index[Dimension]()
47 end = itk.Index[Dimension]()
48 size = itk.Size[Dimension]()
49 
50 MedianType = itk.MedianImageFilter[ImageType, ImageType]
51 median = MedianType.New()
52 median.SetInput(reader.GetOutput())
53 median.SetRadius(2)
54 
55 WriterType = itk.ImageFileWriter[ImageType]
56 writer = WriterType.New()
57 writer.SetFileName(args.output_image)
58 
59 # Use for loops to split the image into chunks.
60 # This way of splitting gives us easy control over it.
61 # For example, we could make the regions always be of equal size
62 # in order to create samples for a deep learning algorithm.
63 # ImageRegionSplitterMultidimensional has a similar
64 # functionality to what is implemented below.
65 for z in range(zDiv):
66  start[2] = int(fullSize[2] * float(z) / zDiv)
67  end[2] = int(fullSize[2] * (z + 1.0) / zDiv)
68  size[2] = end[2] - start[2]
69 
70  for y in range(yDiv):
71  start[1] = int(fullSize[1] * float(y) / yDiv)
72  end[1] = int(fullSize[1] * (y + 1.0) / yDiv)
73  size[1] = end[1] - start[1]
74 
75  for x in range(xDiv):
76  start[0] = int(fullSize[0] * float(x) / xDiv)
77  end[0] = int(fullSize[0] * (x + 1.0) / xDiv)
78  size[0] = end[0] - start[0]
79 
80  region = itk.ImageRegion[Dimension]()
81  region.SetIndex(start)
82  region.SetSize(size)
83 
84  # now that we have our chunk, request the filter to only process that
85  median.GetOutput().SetRequestedRegion(region)
86  median.Update()
87  result = median.GetOutput()
88 
89  # only needed in case of further manual processing
90  result.DisconnectPipeline()
91 
92  # possible additional non-ITK pipeline processing
93 
94  writer.SetInput(result)
95 
96  # convert region into IO region
97  ioRegion = itk.ImageIORegion(Dimension)
98  ioRegion.SetIndex(0, start[0])
99  ioRegion.SetIndex(1, start[1])
100  ioRegion.SetIndex(2, start[2])
101  ioRegion.SetSize(0, region.GetSize()[0])
102  ioRegion.SetSize(1, region.GetSize()[1])
103  ioRegion.SetSize(2, region.GetSize()[2])
104  # tell the writer this is only a chunk of a larger image
105  writer.SetIORegion(ioRegion)
106 
107  try:
108  writer.Update()
109  except Exception as e:
110  print("Exception for chunk:", x, y, z, "\n", e)
111  sys.exit(1)
itk::Index
Represent a n-dimensional index in a n-dimensional image.
Definition: itkIndex.h:70
itk::Size
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition: itkSize.h:71
itk::ImageRegion
An image region represents a structured region of data.
Definition: itkImageRegion.h:80
itk::ImageIORegion
An ImageIORegion represents a structured region of data.
Definition: itkImageIORegion.h:52
itk::ImageFileReader
Data source that reads image data from a single file.
Definition: itkImageFileReader.h:75
itk::ImageFileWriter
Writes image data to a single file.
Definition: itkImageFileWriter.h:88
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
itk::MedianImageFilter
Applies a median filter to an image.
Definition: itkMedianImageFilter.h:53