ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkDirectedHausdorffDistanceImageFilter.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 itkDirectedHausdorffDistanceImageFilter_h
19 #define itkDirectedHausdorffDistanceImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 #include "itkNumericTraits.h"
23 #include "itkArray.h"
25 
26 namespace itk
27 {
71 template< typename TInputImage1, typename TInputImage2 >
72 class ITK_TEMPLATE_EXPORT DirectedHausdorffDistanceImageFilter:
73  public ImageToImageFilter< TInputImage1, TInputImage1 >
74 {
75 public:
76  ITK_DISALLOW_COPY_AND_ASSIGN(DirectedHausdorffDistanceImageFilter);
77 
83 
85  itkNewMacro(Self);
86 
89 
91  using InputImage1Type = TInputImage1;
92  using InputImage2Type = TInputImage2;
93  using InputImage1Pointer = typename TInputImage1::Pointer;
94  using InputImage2Pointer = typename TInputImage2::Pointer;
95  using InputImage1ConstPointer = typename TInputImage1::ConstPointer;
96  using InputImage2ConstPointer = typename TInputImage2::ConstPointer;
97 
99  using SizeType = typename TInputImage1::SizeType;
101 
102  using InputImage1PixelType = typename TInputImage1::PixelType;
103  using InputImage2PixelType = typename TInputImage2::PixelType;
104 
106  static constexpr unsigned int ImageDimension = TInputImage1::ImageDimension;
107 
110 
112  void SetInput1(const InputImage1Type *image);
113 
115  void SetInput2(const InputImage2Type *image);
116 
118  const InputImage1Type * GetInput1();
119 
121  const InputImage2Type * GetInput2();
122 
124  itkSetMacro(UseImageSpacing, bool);
125  itkGetConstMacro( UseImageSpacing, bool );
127 
129  itkGetConstMacro(DirectedHausdorffDistance, RealType);
130  itkGetConstMacro(AverageHausdorffDistance, RealType);
132 
133 #ifdef ITK_USE_CONCEPT_CHECKING
134  // Begin concept checking
135  itkConceptMacro( InputHasNumericTraitsCheck,
137  // End concept checking
138 #endif
139 
140 protected:
142  ~DirectedHausdorffDistanceImageFilter() override = default;
143  void PrintSelf(std::ostream & os, Indent indent) const override;
144 
147  void AllocateOutputs() override;
148 
150  void BeforeThreadedGenerateData() override;
151 
154  void AfterThreadedGenerateData() override;
155 
157  void ThreadedGenerateData(const RegionType &
158  outputRegionForThread,
159  ThreadIdType threadId) override;
160 
161  void DynamicThreadedGenerateData( const RegionType & ) override
162  {
163  itkExceptionMacro("This class requires threadId so it must use classic multi-threading model");
164  }
165 
166  // Override since the filter needs all the data for the algorithm
167  void GenerateInputRequestedRegion() override;
168 
169  // Override since the filter produces all of its output
170  void EnlargeOutputRequestedRegion(DataObject *data) override;
171 
172 private:
175 
176 
178 
181 
183  std::vector< CompensatedSummationType > m_Sum;
184 
188 }; // end of class
189 } // end namespace itk
190 
191 #ifndef ITK_MANUAL_INSTANTIATION
192 #include "itkDirectedHausdorffDistanceImageFilter.hxx"
193 #endif
194 
195 #endif
Light weight base class for most itk classes.
Define numeric traits for std::vector.
typename NumericTraits< InputImage1PixelType >::RealType RealType
unsigned int ThreadIdType
Definition: itkIntTypes.h:99
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
Computes the directed Hausdorff distance between the set of non-zero pixels of two images...
#define itkConceptMacro(name, concept)
Base class for all data objects in ITK.
Templated n-dimensional image class.
Definition: itkImage.h:75