ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkLabelSetMorphBaseImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
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  *=========================================================================*/
18 #ifndef itkLabelSetMorphBaseImageFilter_h
19 #define itkLabelSetMorphBaseImageFilter_h
20 
21 #include "itkNumericTraits.h"
22 #include "itkImageToImageFilter.h"
23 
24 namespace itk
25 {
26 #if ITK_VERSION_MAJOR < 4
27 typedef int ThreadIdType;
28 typedef int RegionIndexType;
29 #else
30 typedef unsigned int RegionIndexType;
31 #endif
32 
45 template< typename TInputImage, bool doDilate,
46  typename TOutputImage = TInputImage >
48  public ImageToImageFilter< TInputImage, TOutputImage >
49 {
50 public:
51 
57 
59  itkNewMacro(Self);
60 
63 
65  typedef TInputImage InputImageType;
66  typedef TOutputImage OutputImageType;
67  typedef typename TInputImage::PixelType PixelType;
69  typedef typename TOutputImage::PixelType OutputPixelType;
71 
74 
76  typedef typename TInputImage::Pointer InputImagePointer;
77  typedef typename TInputImage::ConstPointer InputImageConstPointer;
80 
83 
86  typedef typename OutputImageType::RegionType OutputImageRegionType;
87 
88  // set all of the scales the same
89  void SetRadius(ScalarRealType scale);
90 
91  itkSetMacro(Radius, RadiusType);
92  itkGetConstReferenceMacro(Radius, RadiusType);
93 
98  itkSetMacro(UseImageSpacing, bool);
99  itkGetConstReferenceMacro(UseImageSpacing, bool);
100  itkBooleanMacro(UseImageSpacing);
102 
104  itkStaticConstMacro(ImageDimension, unsigned int,
105  TInputImage::ImageDimension);
106  itkStaticConstMacro(OutputImageDimension, unsigned int,
107  TOutputImage::ImageDimension);
108  itkStaticConstMacro(InputImageDimension, unsigned int,
109  TInputImage::ImageDimension);
111 
116  void writeDist(std::string fname);
117 
118 protected:
121 
122  RegionIndexType SplitRequestedRegion(RegionIndexType i, RegionIndexType num,
123  OutputImageRegionType & splitRegion) ITK_OVERRIDE;
124 
125  virtual void ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread,
126  ThreadIdType threadId) ITK_OVERRIDE;
127 
128  void GenerateData(void) ITK_OVERRIDE;
129 
130  // Override since the filter produces the entire dataset.
131  void EnlargeOutputRequestedRegion(DataObject *output) ITK_OVERRIDE;
132 
133  bool m_UseImageSpacing;
134  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
135 
136  RadiusType m_Radius;
137  RadiusType m_Scale;
138  typedef typename itk::Image< RealType, TInputImage::ImageDimension > DistanceImageType;
139  typename TInputImage::PixelType m_Extreme;
140 
141  typename DistanceImageType::Pointer m_DistanceImage;
142 
143  int m_MagnitudeSign;
144  int m_CurrentDimension;
145  bool m_FirstPassDone;
146 
147  // this is the first non-zero entry in the radius. Needed to
148  // support elliptical operations
149  RealType m_BaseSigma;
150 private:
151  LabelSetMorphBaseImageFilter(const Self &); //purposely not implemented
152  void operator=(const Self &); //purposely not implemented
153 };
154 } // end namespace itk
155 
156 #ifndef ITK_MANUAL_INSTANTIATION
157 #include "itkLabelSetMorphBaseImageFilter.hxx"
158 #endif
159 
160 #endif
ImageToImageFilter< TInputImage, TOutputImage > Superclass
signed long IndexValueType
Definition: itkIntTypes.h:150
OutputImageType::IndexValueType OutputIndexValueType
Base class for all process objects that output image data.
Base class for binary morphological erosion of label images.
NumericTraits< PixelType >::FloatType RealType
unsigned int ThreadIdType
Definition: itkIntTypes.h:159
itk::FixedArray< ScalarRealType, TInputImage::ImageDimension > RadiusType
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Define additional traits for native types such as int or float.
NumericTraits< PixelType >::ScalarRealType ScalarRealType
Base class for all data objects in ITK.
Templated n-dimensional image class.
Definition: itkImage.h:75