ITK  5.4.0
Insight Toolkit
itkLabelOverlapMeasuresImageFilter.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 itkLabelOverlapMeasuresImageFilter_h
19#define itkLabelOverlapMeasuresImageFilter_h
20
21#include "itkImageSink.h"
22#include "itkNumericTraits.h"
23#include <mutex>
24#include <unordered_map>
25
26namespace itk
27{
28
44template <typename TLabelImage>
45class ITK_TEMPLATE_EXPORT LabelOverlapMeasuresImageFilter : public ImageSink<TLabelImage>
46{
47public:
48 ITK_DISALLOW_COPY_AND_MOVE(LabelOverlapMeasuresImageFilter);
49
55
57 itkNewMacro(Self);
58
60 itkOverrideGetNameOfClassMacro(LabelOverlapMeasuresImageFilter);
61
63 using LabelImageType = TLabelImage;
66
70
71 using LabelType = typename TLabelImage::PixelType;
72
75
81 {
82 public:
83 // default constructor/copy/move etc...
84
85 SizeValueType m_Source{ 0 };
86 SizeValueType m_Target{ 0 };
87 SizeValueType m_Union{ 0 };
88 SizeValueType m_Intersection{ 0 };
89 SizeValueType m_SourceComplement{ 0 };
90 SizeValueType m_TargetComplement{ 0 };
91 };
92
94 using MapType = std::unordered_map<LabelType, LabelSetMeasures>;
95 using MapIterator = typename MapType::iterator;
96 using MapConstIterator = typename MapType::const_iterator;
97
99 static constexpr unsigned int ImageDimension = TLabelImage::ImageDimension;
100
110 MapType
112 {
113 return this->m_LabelSetMeasures;
114 }
115
116 // Overlap agreement metrics
117
119 RealType
121
124
130 {
131 return this->GetUnionOverlap();
132 }
140 {
141 return this->GetUnionOverlap(label);
142 }
146 RealType
150 {
151 return this->GetMeanOverlap();
152 }
160 {
161 return this->GetMeanOverlap(label);
162 }
166 RealType
168
171
172 // Overlap error metrics
173
177
180
184
187
191
194
195#ifdef ITK_USE_CONCEPT_CHECKING
196 // Begin concept checking
197 itkConceptMacro(Input1HasNumericTraitsCheck, (Concept::HasNumericTraits<LabelType>));
198 // End concept checking
199#endif
200
201protected:
204
207
208 void
209 PrintSelf(std::ostream & os, Indent indent) const override;
210
211 void
213
214 void
216
217 void
218 MergeMap(MapType & m1, MapType & m2) const;
219
220private:
221 MapType m_LabelSetMeasures{};
222
223 std::mutex m_Mutex{};
224}; // end of class
225
226} // end namespace itk
227
228#ifndef ITK_MANUAL_INSTANTIATION
229# include "itkLabelOverlapMeasuresImageFilter.hxx"
230#endif
231
232#endif
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Computes overlap measures between the set same set of labels of pixels of two images....
typename TLabelImage::ConstPointer LabelImageConstPointer
typename NumericTraits< LabelType >::RealType RealType
itkSetInputMacro(SourceImage, LabelImageType)
RealType GetUnionOverlap(LabelType) const
itkSetInputMacro(TargetImage, LabelImageType)
RealType GetFalsePositiveError(LabelType) const
RealType GetFalseDiscoveryRate(LabelType) const
RealType GetVolumeSimilarity(LabelType) const
RealType GetMeanOverlap(LabelType) const
void ThreadedStreamedGenerateData(const RegionType &) override
RealType GetTargetOverlap(LabelType) const
~LabelOverlapMeasuresImageFilter() override=default
typename NumericTraits< LabelType >::PrintType PrintType
itkGetInputMacro(SourceImage, LabelImageType)
void PrintSelf(std::ostream &os, Indent indent) const override
RealType GetFalseNegativeError(LabelType) const
itkGetInputMacro(TargetImage, LabelImageType)
std::unordered_map< LabelType, LabelSetMeasures > MapType
void MergeMap(MapType &m1, MapType &m2) const
Light weight base class for most itk classes.
SmartPointer< const Self > ConstPointer
SmartPointer< Self > Pointer
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned long SizeValueType
Definition: itkIntTypes.h:83