ITK  4.6.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 >
73  public ImageToImageFilter< TInputImage1, TInputImage1 >
74 {
75 public:
81 
83  itkNewMacro(Self);
84 
87 
89  typedef TInputImage1 InputImage1Type;
90  typedef TInputImage2 InputImage2Type;
91  typedef typename TInputImage1::Pointer InputImage1Pointer;
92  typedef typename TInputImage2::Pointer InputImage2Pointer;
93  typedef typename TInputImage1::ConstPointer InputImage1ConstPointer;
94  typedef typename TInputImage2::ConstPointer InputImage2ConstPointer;
95 
96  typedef typename TInputImage1::RegionType RegionType;
97  typedef typename TInputImage1::SizeType SizeType;
98  typedef typename TInputImage1::IndexType IndexType;
99 
100  typedef typename TInputImage1::PixelType InputImage1PixelType;
101  typedef typename TInputImage2::PixelType InputImage2PixelType;
102 
104  itkStaticConstMacro(ImageDimension, unsigned int,
105  TInputImage1::ImageDimension);
106 
109 
111  void SetInput1(const InputImage1Type *image);
112 
114  void SetInput2(const InputImage2Type *image);
115 
117  const InputImage1Type * GetInput1(void);
118 
120  const InputImage2Type * GetInput2(void);
121 
123  itkSetMacro(UseImageSpacing, bool);
124  itkGetConstMacro( UseImageSpacing, bool );
126 
128  itkGetConstMacro(DirectedHausdorffDistance, RealType);
129  itkGetConstMacro(AverageHausdorffDistance, RealType);
131 
132 #ifdef ITK_USE_CONCEPT_CHECKING
133  // Begin concept checking
134  itkConceptMacro( InputHasNumericTraitsCheck,
136  // End concept checking
137 #endif
138 
139 protected:
142  void PrintSelf(std::ostream & os, Indent indent) const;
143 
146  void AllocateOutputs();
147 
150 
154 
156  void ThreadedGenerateData(const RegionType &
157  outputRegionForThread,
158  ThreadIdType threadId);
159 
160  // Override since the filter needs all the data for the algorithm
162 
163  // Override since the filter produces all of its output
165 
166 private:
167  DirectedHausdorffDistanceImageFilter(const Self &); //purposely not implemented
168  void operator=(const Self &); //purposely not implemented
169 
172 
173 
175 
178 
180  std::vector< CompensatedSummationType > m_Sum;
181 
185 }; // end of class
186 } // end namespace itk
187 
188 #ifndef ITK_MANUAL_INSTANTIATION
189 #include "itkDirectedHausdorffDistanceImageFilter.hxx"
190 #endif
191 
192 #endif
Light weight base class for most itk classes.
Image< RealType, itkGetStaticConstMacro(ImageDimension) > DistanceMapType
const InputImage1Type * GetInput1(void)
NumericTraits< InputImage1PixelType >::RealType RealType
void SetInput1(const InputImage1Type *image)
Perform more precise accumulation of floating point numbers.
void PrintSelf(std::ostream &os, Indent indent) const
void SetInput2(const InputImage2Type *image)
void ThreadedGenerateData(const RegionType &outputRegionForThread, ThreadIdType threadId)
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 additional traits for native types such as int or float.
const InputImage2Type * GetInput2(void)
#define itkConceptMacro(name, concept)
Base class for all data objects in ITK.
Templated n-dimensional image class.
Definition: itkImage.h:75
void EnlargeOutputRequestedRegion(DataObject *data)
unsigned int ThreadIdType
Definition: itkIntTypes.h:159
ImageToImageFilter< TInputImage1, TInputImage1 > Superclass