ITK  4.12.0
Insight Segmentation and Registration Toolkit
itkTestingComparisonImageFilter.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 itkTestingComparisonImageFilter_h
19 #define itkTestingComparisonImageFilter_h
20 
21 #include "itkArray.h"
22 #include "itkNumericTraits.h"
23 #include "itkImageSource.h"
24 
25 namespace itk
26 {
27 namespace Testing
28 {
41 template< typename TInputImage, typename TOutputImage >
42 class ITK_TEMPLATE_EXPORT ComparisonImageFilter:
43  public ImageSource< TOutputImage >
44 {
45 public:
51 
53  itkNewMacro(Self);
54 
57 
59  typedef TInputImage InputImageType;
60  typedef typename InputImageType::PixelType InputPixelType;
61  typedef TOutputImage OutputImageType;
62  typedef typename OutputImageType::PixelType OutputPixelType;
63  typedef typename OutputImageType::RegionType OutputImageRegionType;
66 
68  virtual void SetValidInput(const InputImageType *validImage);
69 
71  virtual void SetTestInput(const InputImageType *testImage);
72 
75  itkSetMacro(ToleranceRadius, int);
76  itkGetConstMacro(ToleranceRadius, int);
78 
81  itkSetMacro(DifferenceThreshold, OutputPixelType);
82  itkGetConstMacro(DifferenceThreshold, OutputPixelType);
84 
88  itkSetMacro(IgnoreBoundaryPixels, bool);
89  itkGetConstMacro(IgnoreBoundaryPixels, bool);
91 
94  itkGetConstMacro(MinimumDifference, OutputPixelType);
95  itkGetConstMacro(MaximumDifference, OutputPixelType);
96  itkGetConstMacro(MeanDifference, RealType);
97  itkGetConstMacro(TotalDifference, AccumulateType);
98  itkGetConstMacro(NumberOfPixelsWithDifferences, SizeValueType);
100 
102  using Superclass::SetInput;
103  virtual void SetInput(const TInputImage *image);
104  virtual void SetInput(unsigned int, const TInputImage *image);
105  const TInputImage * GetInput() const;
106  const TInputImage * GetInput(unsigned int idx) const;
108 
109 protected:
112 
113  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
114 
124  void ThreadedGenerateData(const OutputImageRegionType & threadRegion,
125  ThreadIdType threadId) ITK_OVERRIDE;
126 
127  void BeforeThreadedGenerateData() ITK_OVERRIDE;
128 
129  void AfterThreadedGenerateData() ITK_OVERRIDE;
130 
131  OutputPixelType m_DifferenceThreshold;
132 
133  RealType m_MeanDifference;
134  OutputPixelType m_MinimumDifference;
135  OutputPixelType m_MaximumDifference;
136 
137  AccumulateType m_TotalDifference;
138 
139  SizeValueType m_NumberOfPixelsWithDifferences;
140 
141  int m_ToleranceRadius;
142 
143  Array< AccumulateType > m_ThreadDifferenceSum;
144  Array< SizeValueType > m_ThreadNumberOfPixels;
145 
146  Array< OutputPixelType > m_ThreadMinimumDifference;
147  Array< OutputPixelType > m_ThreadMaximumDifference;
148 
149 private:
150  ITK_DISALLOW_COPY_AND_ASSIGN(ComparisonImageFilter);
151 
152  bool m_IgnoreBoundaryPixels;
153 };
154 } // end namespace Testing
155 } // end namespace itk
156 
157 #ifndef ITK_MANUAL_INSTANTIATION
158 #include "itkTestingComparisonImageFilter.hxx"
159 #endif
160 
161 
162 #endif
Array class with size defined at construction time.
Definition: itkArray.h:50
virtual void PrintSelf(std::ostream &os, Indent indent) const override
Base class for all process objects that output image data.
unsigned long SizeValueType
Definition: itkIntTypes.h:143
NumericTraits< RealType >::AccumulateType AccumulateType
NumericTraits< OutputPixelType >::RealType RealType
unsigned int ThreadIdType
Definition: itkIntTypes.h:159
Implements comparison between two images.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Define additional traits for native types such as int or float.