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((original_size[1] - 1) * original_spacing[1] / new_spacing[1]),
51  ]
52  image = sitk.Resample(
53  image1=image,
54  size=new_size,
55  transform=sitk.Transform(),
56  interpolator=sitk.sitkLinear,
57  outputOrigin=image.GetOrigin(),
58  outputSpacing=new_spacing,
59  outputDirection=image.GetDirection(),
60  defaultPixelValue=0,
61  outputPixelType=image.GetPixelID(),
62  )
63  # If a single channel image, rescale to [0,255]. Also modify the
64  # intensity values based on the photometric interpretation. If
65  # MONOCHROME2 (minimum should be displayed as black) we don't need to
66  # do anything, if image has MONOCRHOME1 (minimum should be displayed as
67  # white) we flip # the intensities. This is a constraint imposed by ITK
68  # which always assumes MONOCHROME2.
69  if image.GetNumberOfComponentsPerPixel() == 1:
70  image = sitk.RescaleIntensity(image, 0, 255)
71  if image_file_reader.GetMetaData("0028|0004").strip() == "MONOCHROME1":
72  image = sitk.InvertIntensity(image, maximum=255)
73  image = sitk.Cast(image, sitk.sitkUInt8)
74  sitk.WriteImage(image, output_file_name)
75  return True
76  except BaseException:
77  return False
78 
79 
80 def convert_images(input_file_names, output_file_names, new_width):
81  MAX_PROCESSES = 15
82  with multiprocessing.Pool(processes=MAX_PROCESSES) as pool:
83  return pool.starmap(
84  functools.partial(convert_image, new_width=new_width),
85  zip(input_file_names, output_file_names),
86  )
87 
88 
89 def positive_int(int_str):
90  value = int(int_str)
91  if value <= 0:
92  raise argparse.ArgumentTypeError(int_str + " is not a positive integer value")
93  return value
94 
95 
96 def directory(dir_name):
97  if not os.path.isdir(dir_name):
98  raise argparse.ArgumentTypeError(dir_name + " is not a valid directory name")
99  return dir_name
100 
101 
102 def main(argv=None):
103  parser = argparse.ArgumentParser(
104  description="Convert and resize DICOM files to common image types."
105  )
106  parser.add_argument(
107  "root_of_data_directory",
108  type=directory,
109  help="Path to the topmost directory containing data.",
110  )
111  parser.add_argument(
112  "output_file_extension",
113  help="Image file extension, this determines output file type " "(e.g. png) .",
114  )
115  parser.add_argument("--w", type=positive_int, help="Width of converted images.")
116  parser.add_argument("--od", type=directory, help="Output directory.")
117  args = parser.parse_args(argv)
118 
119  input_file_names = []
120  for dir_name, subdir_names, file_names in os.walk(args.root_of_data_directory):
121  input_file_names += [
122  os.path.join(os.path.abspath(dir_name), fname) for fname in file_names
123  ]
124  if args.od:
125  # if all output files are written to the same directory we need them
126  # to have a unique name, so use an index.
127  file_names = [
128  os.path.join(os.path.abspath(args.od), str(i))
129  for i in range(len(input_file_names))
130  ]
131  else:
132  file_names = input_file_names
133  output_file_names = [
134  file_name + "." + args.output_file_extension for file_name in file_names
135  ]
136 
137  res = convert_images(input_file_names, output_file_names, args.w)
138  input_file_names = list(itertools.compress(input_file_names, res))
139  output_file_names = list(itertools.compress(output_file_names, res))
140 
141  # save csv file mapping input and output file names.
142  # using csv module and not pandas so as not to create more dependencies
143  # for the examples. pandas based code is more elegant/shorter.
144  dir_name = args.od if args.od else os.getcwd()
145  with open(os.path.join(dir_name, "file_names.csv"), mode="w") as fp:
146  fp_writer = csv.writer(
147  fp, delimiter=",", quotechar='"', quoting=csv.QUOTE_MINIMAL
148  )
149  fp_writer.writerow(["input file name", "output file name"])
150  for data in zip(input_file_names, output_file_names):
151  fp_writer.writerow(data)
152 
153 
154 if __name__ == "__main__":
155  sys.exit(main())
itk::simple::Transform
A simplified wrapper around a variety of ITK transforms.
Definition: sitkTransform.h:86
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::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:74
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
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...