ITK  6.0.0
Insight Toolkit
itkContourDirectedMeanDistanceImageFilter.h
Go to the documentation of this file.
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  * https://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 itkContourDirectedMeanDistanceImageFilter_h
19 #define itkContourDirectedMeanDistanceImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 #include "itkNumericTraits.h"
23 #include "itkArray.h"
24 #include "itkImage.h"
25 
26 namespace itk
27 {
63 template <typename TInputImage1, typename TInputImage2>
64 class ITK_TEMPLATE_EXPORT ContourDirectedMeanDistanceImageFilter : public ImageToImageFilter<TInputImage1, TInputImage1>
65 {
66 public:
67  ITK_DISALLOW_COPY_AND_MOVE(ContourDirectedMeanDistanceImageFilter);
68 
74 
76  itkNewMacro(Self);
77 
79  itkOverrideGetNameOfClassMacro(ContourDirectedMeanDistanceImageFilter);
80 
82  using InputImage1Type = TInputImage1;
83  using InputImage2Type = TInputImage2;
88 
90  using SizeType = typename TInputImage1::SizeType;
92 
93  using InputImage1PixelType = typename TInputImage1::PixelType;
94  using InputImage2PixelType = typename TInputImage2::PixelType;
95 
97  static constexpr unsigned int ImageDimension = TInputImage1::ImageDimension;
98 
101 
103  void
104  SetInput1(const InputImage1Type * image);
105 
107  void
108  SetInput2(const InputImage2Type * image);
109 
111  const InputImage1Type *
112  GetInput1();
113 
115  const InputImage2Type *
116  GetInput2();
117 
119  itkGetConstMacro(ContourDirectedMeanDistance, RealType);
120 
122  itkSetMacro(UseImageSpacing, bool);
123  itkGetConstMacro(UseImageSpacing, bool);
124  itkBooleanMacro(UseImageSpacing);
127 #ifdef ITK_USE_CONCEPT_CHECKING
128  // Begin concept checking
130  // End concept checking
131 #endif
132 
133 protected:
135  ~ContourDirectedMeanDistanceImageFilter() override = default;
136  void
137  PrintSelf(std::ostream & os, Indent indent) const override;
138 
141  void
142  AllocateOutputs() override;
143 
145  void
146  BeforeThreadedGenerateData() override;
147 
150  void
151  AfterThreadedGenerateData() override;
152 
154  void
155  ThreadedGenerateData(const RegionType & outputRegionForThread, ThreadIdType threadId) override;
156 
157  void
159  {
160  itkExceptionMacro("This class requires threadId so it must use classic multi-threading model");
161  }
162 
163  // Override since the filter needs all the data for the algorithm
164  void
165  GenerateInputRequestedRegion() override;
166 
167  // Override since the filter produces all of its output
168  void
169  EnlargeOutputRequestedRegion(DataObject * data) override;
170 
171 private:
173 
174  typename DistanceMapType::Pointer m_DistanceMap{};
175 
176  Array<RealType> m_MeanDistance{};
178  RealType m_ContourDirectedMeanDistance{};
179  bool m_UseImageSpacing{ true };
180 };
181 } // end namespace itk
182 
183 #ifndef ITK_MANUAL_INSTANTIATION
184 # include "itkContourDirectedMeanDistanceImageFilter.hxx"
185 #endif
186 
187 #endif
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAddImageFilter.h:94
itk::ContourDirectedMeanDistanceImageFilter::InputImage1Type
TInputImage1 InputImage1Type
Definition: itkContourDirectedMeanDistanceImageFilter.h:82
itk::Concept::HasNumericTraits
Definition: itkConceptChecking.h:716
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::ContourDirectedMeanDistanceImageFilter::InputImage1PixelType
typename TInputImage1::PixelType InputImage1PixelType
Definition: itkContourDirectedMeanDistanceImageFilter.h:93
itkImage.h
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::ContourDirectedMeanDistanceImageFilter::InputImage2PixelType
typename TInputImage2::PixelType InputImage2PixelType
Definition: itkContourDirectedMeanDistanceImageFilter.h:94
itk::ContourDirectedMeanDistanceImageFilter::InputImage2ConstPointer
typename TInputImage2::ConstPointer InputImage2ConstPointer
Definition: itkContourDirectedMeanDistanceImageFilter.h:87
itk::ContourDirectedMeanDistanceImageFilter::InputImage2Type
TInputImage2 InputImage2Type
Definition: itkContourDirectedMeanDistanceImageFilter.h:83
itk::ThreadIdType
unsigned int ThreadIdType
Definition: itkIntTypes.h:102
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::ImageToImageFilter
Base class for filters that take an image as input and produce an image as output.
Definition: itkImageToImageFilter.h:108
itk::ContourDirectedMeanDistanceImageFilter
Computes the directed Mean distance between the boundaries of non-zero pixel regions of two images.
Definition: itkContourDirectedMeanDistanceImageFilter.h:64
itk::ContourDirectedMeanDistanceImageFilter::RealType
typename NumericTraits< InputImage1PixelType >::RealType RealType
Definition: itkContourDirectedMeanDistanceImageFilter.h:100
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:55
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::ContourDirectedMeanDistanceImageFilter::RegionType
typename TInputImage1::RegionType RegionType
Definition: itkContourDirectedMeanDistanceImageFilter.h:89
itk::ContourDirectedMeanDistanceImageFilter::InputImage1ConstPointer
typename TInputImage1::ConstPointer InputImage1ConstPointer
Definition: itkContourDirectedMeanDistanceImageFilter.h:86
itkImageToImageFilter.h
itk::ContourDirectedMeanDistanceImageFilter::IndexType
typename TInputImage1::IndexType IndexType
Definition: itkContourDirectedMeanDistanceImageFilter.h:91
itk::NumericTraits
Define additional traits for native types such as int or float.
Definition: itkNumericTraits.h:60
itkArray.h
itkConceptMacro
#define itkConceptMacro(name, concept)
Definition: itkConceptChecking.h:65
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ContourDirectedMeanDistanceImageFilter::SizeType
typename TInputImage1::SizeType SizeType
Definition: itkContourDirectedMeanDistanceImageFilter.h:90
itk::Array< RealType >
itkNumericTraits.h
itk::ContourDirectedMeanDistanceImageFilter::DynamicThreadedGenerateData
void DynamicThreadedGenerateData(const RegionType &) override
Definition: itkContourDirectedMeanDistanceImageFilter.h:158
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
itk::ContourDirectedMeanDistanceImageFilter::InputImage2Pointer
typename TInputImage2::Pointer InputImage2Pointer
Definition: itkContourDirectedMeanDistanceImageFilter.h:85
itk::ContourDirectedMeanDistanceImageFilter::InputImage1Pointer
typename TInputImage1::Pointer InputImage1Pointer
Definition: itkContourDirectedMeanDistanceImageFilter.h:84
itk::DataObject
Base class for all data objects in ITK.
Definition: itkDataObject.h:293