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""" A SimpleITK example demonstrating how to convert and resize DICOM files
20 to common image types. """
21
22import argparse
23import csv
24import functools
25import itertools
26import multiprocessing
27import os
28import sys
29
30import SimpleITK as sitk
31
32
33def convert_image(input_file_name, output_file_name, new_width=None):
34 """ Convert a single DICOM image to a common image type. """
35 try:
36 image_file_reader = sitk.ImageFileReader()
37 # only read DICOM images
38 image_file_reader.SetImageIO("GDCMImageIO")
39 image_file_reader.SetFileName(input_file_name)
40 image_file_reader.ReadImageInformation()
41 image_size = list(image_file_reader.GetSize())
42 if len(image_size) == 3 and image_size[2] == 1:
43 image_size[2] = 0
44 image_file_reader.SetExtractSize(image_size)
45 image = image_file_reader.Execute()
46 if new_width:
47 original_size = image.GetSize()
48 original_spacing = image.GetSpacing()
49 new_spacing = [
50 (original_size[0] - 1) * original_spacing[0] / (new_width - 1)
51 ] * 2
52 new_size = [
53 new_width,
54 int((original_size[1] - 1) * original_spacing[1] / new_spacing[1]),
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 image_file_reader.GetMetaData("0028|0004").strip() == "MONOCHROME1":
76 image = sitk.InvertIntensity(image, maximum=255)
77 image = sitk.Cast(image, sitk.sitkUInt8)
78 sitk.WriteImage(image, output_file_name)
79 return True
80 except RuntimeError:
81 return False
82
83
84def convert_images(input_file_names, output_file_names, new_width):
85 """ Convert multiple DICOM images in parallel to a common image type. """
86 MAX_PROCESSES = 15
87 with multiprocessing.Pool(processes=MAX_PROCESSES) as pool:
88 return pool.starmap(
89 functools.partial(convert_image, new_width=new_width),
90 zip(input_file_names, output_file_names),
91 )
92
93
94def positive_int(int_str):
95 """ Custom argparse type for positive integers. """
96 value = int(int_str)
97 if value <= 0:
98 raise argparse.ArgumentTypeError(int_str + " is not a positive integer value")
99 return value
100
101
102def directory(dir_name):
103 """ Custom argparse type for directory. """
104 if not os.path.isdir(dir_name):
105 raise argparse.ArgumentTypeError(dir_name + " is not a valid directory name")
106 return dir_name
107
108
109def main(argv=None):
110 """ Main function. """
111 parser = argparse.ArgumentParser(
112 description="Convert and resize DICOM files to common image types."
113 )
114 parser.add_argument(
115 "root_of_data_directory",
116 type=directory,
117 help="Path to the topmost directory containing data.",
118 )
119 parser.add_argument(
120 "output_file_extension",
121 help="Image file extension, this determines output file type " "(e.g. png) .",
122 )
123 parser.add_argument("--w", type=positive_int, help="Width of converted images.")
124 parser.add_argument("--od", type=directory, help="Output directory.")
125 args = parser.parse_args(argv)
126
127 input_file_names = []
128 for dir_name, _, file_names in os.walk(args.root_of_data_directory):
129 input_file_names += [
130 os.path.join(os.path.abspath(dir_name), fname) for fname in file_names
131 ]
132 if args.od:
133 # if all output files are written to the same directory we need them
134 # to have a unique name, so use an index.
135 file_names = [
136 os.path.join(os.path.abspath(args.od), str(i))
137 for i in range(len(input_file_names))
138 ]
139 else:
140 file_names = input_file_names
141 output_file_names = [
142 file_name + "." + args.output_file_extension for file_name in file_names
143 ]
144
145 res = convert_images(input_file_names, output_file_names, args.w)
146 input_file_names = list(itertools.compress(input_file_names, res))
147 output_file_names = list(itertools.compress(output_file_names, res))
148
149 # save csv file mapping input and output file names.
150 # using csv module and not pandas so as not to create more dependencies
151 # for the examples. pandas based code is more elegant/shorter.
152 dir_name = args.od if args.od else os.getcwd()
153 with open(os.path.join(dir_name, "file_names.csv"), mode="w", encoding='utf-8') as fp:
154 fp_writer = csv.writer(
155 fp, delimiter=",", quotechar='"', quoting=csv.QUOTE_MINIMAL
156 )
157 fp_writer.writerow(["input file name", "output file name"])
158 for data in zip(input_file_names, output_file_names):
159 fp_writer.writerow(data)
160
161
162if __name__ == "__main__":
163 sys.exit(main())
Read an image file and return a SimpleITK Image.
A simplified wrapper around a variety of ITK transforms.
Image InvertIntensity(Image &&image1, double maximum=255)
Invert the intensity of an image.
SITKBasicFilters_EXPORT Image Resample(const Image &image1, Transform transform=itk::simple::Transform(), InterpolatorEnum interpolator=itk::simple::sitkLinear, double defaultPixelValue=0.0, PixelIDValueEnum outputPixelType=sitkUnknown, bool useNearestNeighborExtrapolator=false)
itk::simple::ResampleImageFilter Procedural Interface
SITKIO_EXPORT void WriteImage(const Image &image, const PathType &fileName, bool useCompression=false, int compressionLevel=-1)
WriteImage is a procedural interface to the ImageFileWriter. class which is convenient for many image...
Image RescaleIntensity(Image &&image1, double outputMinimum=0, double outputMaximum=255)
Applies a linear transformation to the intensity levels of the input Image .
Image Cast(const Image &image, PixelIDValueEnum pixelID)