SimpleITK  
DicomConvert/DicomConvert.py
1 # =========================================================================
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 # http://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 # =========================================================================
18 
19 import argparse
20 import csv
21 import functools
22 import itertools
23 import multiprocessing
24 import os
25 import sys
26 
27 import SimpleITK as sitk
28 
29 
30 def convert_image(input_file_name, output_file_name, new_width=None):
31  try:
32  image_file_reader = sitk.ImageFileReader()
33  # only read DICOM images
34  image_file_reader.SetImageIO("GDCMImageIO")
35  image_file_reader.SetFileName(input_file_name)
36  image_file_reader.ReadImageInformation()
37  image_size = list(image_file_reader.GetSize())
38  if len(image_size) == 3 and image_size[2] == 1:
39  image_size[2] = 0
40  image_file_reader.SetExtractSize(image_size)
41  image = image_file_reader.Execute()
42  if new_width:
43  original_size = image.GetSize()
44  original_spacing = image.GetSpacing()
45  new_spacing = [
46  (original_size[0] - 1) * original_spacing[0] / (new_width - 1)
47  ] * 2
48  new_size = [
49  new_width,
50  int(
51  (original_size[1] - 1)
52  * original_spacing[1]
53  / new_spacing[1]
54  ),
55  ]
56  image = sitk.Resample(
57  image1=image,
58  size=new_size,
59  transform=sitk.Transform(),
60  interpolator=sitk.sitkLinear,
61  outputOrigin=image.GetOrigin(),
62  outputSpacing=new_spacing,
63  outputDirection=image.GetDirection(),
64  defaultPixelValue=0,
65  outputPixelType=image.GetPixelID(),
66  )
67  # If a single channel image, rescale to [0,255]. Also modify the
68  # intensity values based on the photometric interpretation. If
69  # MONOCHROME2 (minimum should be displayed as black) we don't need to
70  # do anything, if image has MONOCRHOME1 (minimum should be displayed as
71  # white) we flip # the intensities. This is a constraint imposed by ITK
72  # which always assumes MONOCHROME2.
73  if image.GetNumberOfComponentsPerPixel() == 1:
74  image = sitk.RescaleIntensity(image, 0, 255)
75  if (
76  image_file_reader.GetMetaData("0028|0004").strip()
77  == "MONOCHROME1"
78  ):
79  image = sitk.InvertIntensity(image, maximum=255)
80  image = sitk.Cast(image, sitk.sitkUInt8)
81  sitk.WriteImage(image, output_file_name)
82  return True
83  except BaseException:
84  return False
85 
86 
87 def convert_images(input_file_names, output_file_names, new_width):
88  MAX_PROCESSES = 15
89  with multiprocessing.Pool(processes=MAX_PROCESSES) as pool:
90  return pool.starmap(
91  functools.partial(convert_image, new_width=new_width),
92  zip(input_file_names, output_file_names),
93  )
94 
95 
96 def positive_int(int_str):
97  value = int(int_str)
98  if value <= 0:
99  raise argparse.ArgumentTypeError(
100  int_str + " is not a positive integer value"
101  )
102  return value
103 
104 
105 def directory(dir_name):
106  if not os.path.isdir(dir_name):
107  raise argparse.ArgumentTypeError(
108  dir_name + " is not a valid directory name"
109  )
110  return dir_name
111 
112 
113 def main(argv=None):
114  parser = argparse.ArgumentParser(
115  description="Convert and resize DICOM files to common image types."
116  )
117  parser.add_argument(
118  "root_of_data_directory",
119  type=directory,
120  help="Path to the topmost directory containing data.",
121  )
122  parser.add_argument(
123  "output_file_extension",
124  help="Image file extension, this determines output file type "
125  "(e.g. png) .",
126  )
127  parser.add_argument(
128  "--w", type=positive_int, help="Width of converted images."
129  )
130  parser.add_argument("--od", type=directory, help="Output directory.")
131  args = parser.parse_args(argv)
132 
133  input_file_names = []
134  for dir_name, subdir_names, file_names in os.walk(
135  args.root_of_data_directory
136  ):
137  input_file_names += [
138  os.path.join(os.path.abspath(dir_name), fname)
139  for fname in file_names
140  ]
141  if args.od:
142  # if all output files are written to the same directory we need them
143  # to have a unique name, so use an index.
144  file_names = [
145  os.path.join(os.path.abspath(args.od), str(i))
146  for i in range(len(input_file_names))
147  ]
148  else:
149  file_names = input_file_names
150  output_file_names = [
151  file_name + "." + args.output_file_extension
152  for file_name in file_names
153  ]
154 
155  res = convert_images(input_file_names, output_file_names, args.w)
156  input_file_names = list(itertools.compress(input_file_names, res))
157  output_file_names = list(itertools.compress(output_file_names, res))
158 
159  # save csv file mapping input and output file names.
160  # using csv module and not pandas so as not to create more dependencies
161  # for the examples. pandas based code is more elegant/shorter.
162  dir_name = args.od if args.od else os.getcwd()
163  with open(os.path.join(dir_name, "file_names.csv"), mode="w") as fp:
164  fp_writer = csv.writer(
165  fp, delimiter=",", quotechar='"', quoting=csv.QUOTE_MINIMAL
166  )
167  fp_writer.writerow(["input file name", "output file name"])
168  for data in zip(input_file_names, output_file_names):
169  fp_writer.writerow(data)
170 
171 
172 if __name__ == "__main__":
173  sys.exit(main())
itk::simple::Transform
A simplified wrapper around a variety of ITK transforms.
Definition: sitkTransform.h:81
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::WriteImage
SITKIO_EXPORT void WriteImage(const Image &image, const std::vector< std::string > &fileNames, bool useCompression=false, int compressionLevel=-1)
WriteImage is a procedural interface to the ImageSeriesWriter. class which is convenient for many ima...
itk::simple::InvertIntensity
Image InvertIntensity(const Image &image1, double maximum=255)
Invert the intensity of an image.
itk::simple::ImageFileReader
Read an image file and return a SimpleITK Image.
Definition: sitkImageFileReader.h:60
itk::simple::Cast
Image Cast(const Image &image, PixelIDValueEnum pixelID)
itk::simple::Resample
Image Resample(const Image &image1, const std::vector< uint32_t > &size, Transform transform=itk::simple::Transform(), InterpolatorEnum interpolator=itk::simple::sitkLinear, const std::vector< double > &outputOrigin=std::vector< double >(3, 0.0), const std::vector< double > &outputSpacing=std::vector< double >(3, 1.0), const std::vector< double > &outputDirection=std::vector< double >(), double defaultPixelValue=0.0, PixelIDValueEnum outputPixelType=sitkUnknown, bool useNearestNeighborExtrapolator=false)
itk::simple::ResampleImageFilter Procedural Interface