ITK  5.0.0
Insight Segmentation and Registration Toolkit
SphinxExamples/src/Segmentation/Watersheds/SegmentWithWatershedImageFilter/Code.py
1 #!/usr/bin/env python
2 
3 # Copyright Insight Software Consortium
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 # Run with:
18 # ./WatershedImageFilter.py <InputFileName> <OutputFileName> <Threshold> <Level>
19 # e.g.
20 # ./WatershedImageFilter.py BrainProtonDensitySlice.png OutBrainWatershed.png 0.005 .5
21 # (A rule of thumb is to set the Threshold to be about 1 / 100 of the Level.)
22 #
23 # threshold: absolute minimum height value used during processing.
24 # Raising this threshold percentage effectively decreases the number of local minima in the input,
25 # resulting in an initial segmentation with fewer regions.
26 # The assumption is that the shallow regions that thresholding removes are of of less interest.
27 # level: parameter controls the depth of metaphorical flooding of the image.
28 # That is, it sets the maximum saliency value of interest in the result.
29 # Raising and lowering the Level influences the number of segments
30 # in the basic segmentation that are merged to produce the final output.
31 # A level of 1.0 is analogous to flooding the image up to a
32 # depth that is 100 percent of the maximum value in the image.
33 # A level of 0.0 produces the basic segmentation, which will typically be very oversegmented.
34 # Level values of interest are typically low (i.e. less than about 0.40 or 40%),
35 # since higher values quickly start to undersegment the image.
36 
37 import sys
38 import itk
39 
40 if len(sys.argv) != 5:
41  print('Usage: ' + sys.argv[0] +
42  ' <InputFileName> <OutputFileName> <Threshold> <Level>')
43  sys.exit(1)
44 
45 inputFileName = sys.argv[1]
46 outputFileName = sys.argv[2]
47 
48 Dimension = 2
49 
50 FloatPixelType = itk.ctype('float')
51 FloatImageType = itk.Image[FloatPixelType, Dimension]
52 
53 reader = itk.ImageFileReader[FloatImageType].New()
54 reader.SetFileName(inputFileName)
55 
56 gradientMagnitude = \
57  itk.GradientMagnitudeImageFilter.New(Input=reader.GetOutput())
58 
59 watershed = \
60  itk.WatershedImageFilter.New(Input=gradientMagnitude.GetOutput())
61 
62 threshold = float(sys.argv[3])
63 level = float(sys.argv[4])
64 watershed.SetThreshold(threshold)
65 watershed.SetLevel(level)
66 
67 LabeledImageType = type(watershed.GetOutput())
68 
69 PixelType = itk.ctype('unsigned char')
70 RGBPixelType = itk.RGBPixel[PixelType]
71 RGBImageType = itk.Image[RGBPixelType, Dimension]
72 
73 ScalarToRGBColormapFilterType = \
74  itk.ScalarToRGBColormapImageFilter[LabeledImageType, RGBImageType]
75 colormapImageFilter = ScalarToRGBColormapFilterType.New()
76 colormapImageFilter.SetColormap(ScalarToRGBColormapFilterType.Jet)
77 colormapImageFilter.SetInput(watershed.GetOutput())
78 colormapImageFilter.Update()
79 
80 WriterType = itk.ImageFileWriter[RGBImageType]
81 writer = WriterType.New()
82 writer.SetFileName(outputFileName)
83 writer.SetInput(colormapImageFilter.GetOutput())
84 writer.Update()