ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkTestingComparisonImageFilter.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 itkTestingComparisonImageFilter_h
19#define itkTestingComparisonImageFilter_h
20
21#include "itkArray.h"
22#include "itkNumericTraits.h"
24#include <mutex>
25
26namespace itk
27{
28namespace Testing
29{
43template <typename TInputImage, typename TOutputImage>
44class ITK_TEMPLATE_EXPORT ComparisonImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
45{
46public:
47 ITK_DISALLOW_COPY_AND_MOVE(ComparisonImageFilter);
48
54
56 itkNewMacro(Self);
57
59 itkOverrideGetNameOfClassMacro(ComparisonImageFilter);
60
62 using InputImageType = TInputImage;
63 using InputPixelType = typename InputImageType::PixelType;
64 using OutputImageType = TOutputImage;
65 using OutputPixelType = typename OutputImageType::PixelType;
66 using OutputImageRegionType = typename OutputImageType::RegionType;
69
72
75
78 itkSetMacro(VerifyInputInformation, bool);
79 itkGetConstMacro(VerifyInputInformation, bool);
80 itkBooleanMacro(VerifyInputInformation);
81
84 itkSetMacro(ToleranceRadius, int);
85 itkGetConstMacro(ToleranceRadius, int);
86
89 itkSetMacro(DifferenceThreshold, OutputPixelType);
90 itkGetConstMacro(DifferenceThreshold, OutputPixelType);
91
95 itkSetMacro(IgnoreBoundaryPixels, bool);
96 itkGetConstMacro(IgnoreBoundaryPixels, bool);
97 itkBooleanMacro(IgnoreBoundaryPixels);
98
101 itkGetConstMacro(MinimumDifference, OutputPixelType);
102 itkGetConstMacro(MaximumDifference, OutputPixelType);
103 itkGetConstMacro(MeanDifference, RealType);
104 itkGetConstMacro(TotalDifference, AccumulateType);
105 itkGetConstMacro(NumberOfPixelsWithDifferences, SizeValueType);
106
107protected:
109 ~ComparisonImageFilter() override = default;
110
111 void
112 PrintSelf(std::ostream & os, Indent indent) const override;
113
117 void
119
120 void
122
123 void
125
126 void
127 VerifyInputInformation() const override;
128
130
135
137
139
141
142private:
144
145 std::mutex m_Mutex;
146};
147} // end namespace Testing
148} // end namespace itk
149
150#ifndef ITK_MANUAL_INSTANTIATION
151# include "itkTestingComparisonImageFilter.hxx"
152#endif
153
154
155#endif
Control indentation during Print() invocation.
Definition itkIndent.h:50
static constexpr T NonpositiveMin()
static constexpr T max(const T &)
Implements transparent reference counting.
typename OutputImageType::PixelType OutputPixelType
void VerifyInputInformation() const override
Verifies that the inputs meta-data is consistent and valid for continued execution of the pipeline,...
typename InputImageType::PixelType InputPixelType
void DynamicThreadedGenerateData(const OutputImageRegionType &) override
void BeforeThreadedGenerateData() override
ImageToImageFilter< TInputImage, TOutputImage > Superclass
typename NumericTraits< RealType >::AccumulateType AccumulateType
~ComparisonImageFilter() override=default
typename OutputImageType::RegionType OutputImageRegionType
itkSetInputMacro(ValidInput, InputImageType)
itkSetInputMacro(TestInput, InputImageType)
void PrintSelf(std::ostream &os, Indent indent) const override
typename NumericTraits< OutputPixelType >::RealType RealType
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned long SizeValueType
Definition itkIntTypes.h:86