ITK  6.0.0
Insight Toolkit
itkSTAPLEImageFilter.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 itkSTAPLEImageFilter_h
19 #define itkSTAPLEImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 #include <vector>
23 
24 namespace itk
25 {
121 template <typename TInputImage, typename TOutputImage>
122 class ITK_TEMPLATE_EXPORT STAPLEImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
123 {
124 public:
125  ITK_DISALLOW_COPY_AND_MOVE(STAPLEImageFilter);
126 
132 
134  itkNewMacro(Self);
135 
137  itkOverrideGetNameOfClassMacro(STAPLEImageFilter);
138 
141  using OutputPixelType = typename TOutputImage::PixelType;
142  using InputPixelType = typename TInputImage::PixelType;
144 
147  static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
148 
150  using InputImageType = TInputImage;
152  using OutputImageType = TOutputImage;
154 
156  using typename Superclass::OutputImageRegionType;
157 
159  itkSetMacro(ForegroundValue, InputPixelType);
160  itkGetConstMacro(ForegroundValue, InputPixelType);
166  const std::vector<double> &
168  {
169  return m_Specificity;
170  }
171 
175  const std::vector<double> &
177  {
178  return m_Sensitivity;
179  }
180 
183  double
184  GetSensitivity(unsigned int i)
185  {
186  if (i > this->GetNumberOfIndexedInputs())
187  {
188  itkExceptionMacro("Array reference out of bounds.");
189  }
190  return m_Sensitivity[i];
191  }
196  double
197  GetSpecificity(unsigned int i)
198  {
199  if (i > this->GetNumberOfIndexedInputs())
200  {
201  itkExceptionMacro("Array reference out of bounds.");
202  }
203  return m_Specificity[i];
204  }
210  itkSetClampMacro(MaximumIterations, unsigned int, 1, NumericTraits<unsigned int>::max());
211  itkGetConstMacro(MaximumIterations, unsigned int);
221  itkSetMacro(ConfidenceWeight, double);
222  itkGetConstMacro(ConfidenceWeight, double);
226  itkGetConstMacro(ElapsedIterations, unsigned int);
227 
228 #ifdef ITK_USE_CONCEPT_CHECKING
229  // Begin concept checking
230  itkConceptMacro(InputHasNumericTraitsCheck, (Concept::HasNumericTraits<InputPixelType>));
231  // End concept checking
232 #endif
233 
234 protected:
236  {
237  m_ForegroundValue = NumericTraits<InputPixelType>::OneValue();
238  m_MaximumIterations = NumericTraits<unsigned int>::max();
239  m_ElapsedIterations = 0;
240  m_ConfidenceWeight = 1.0;
241  }
242 
243  ~STAPLEImageFilter() override = default;
244  void
245  GenerateData() override;
246 
247  void
248  PrintSelf(std::ostream &, Indent) const override;
249 
250 private:
251  InputPixelType m_ForegroundValue{};
252  unsigned int m_ElapsedIterations{};
253  unsigned int m_MaximumIterations{};
254 
255  double m_ConfidenceWeight{};
256 
257  std::vector<double> m_Sensitivity{};
258  std::vector<double> m_Specificity{};
259 };
260 } // end namespace itk
261 
262 #ifndef ITK_MANUAL_INSTANTIATION
263 # include "itkSTAPLEImageFilter.hxx"
264 #endif
265 
266 #endif
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::STAPLEImageFilter::OutputPixelType
typename TOutputImage::PixelType OutputPixelType
Definition: itkSTAPLEImageFilter.h:141
itk::STAPLEImageFilter::InputPixelType
typename TInputImage::PixelType InputPixelType
Definition: itkSTAPLEImageFilter.h:142
itk::STAPLEImageFilter::GetSpecificity
const std::vector< double > & GetSpecificity() const
Definition: itkSTAPLEImageFilter.h:167
itk::Concept::HasNumericTraits
Definition: itkConceptChecking.h:716
itk::ImageSource::OutputImagePointer
typename OutputImageType::Pointer OutputImagePointer
Definition: itkImageSource.h:91
itk::STAPLEImageFilter::STAPLEImageFilter
STAPLEImageFilter()
Definition: itkSTAPLEImageFilter.h:235
itk::STAPLEImageFilter::GetSensitivity
const std::vector< double > & GetSensitivity() const
Definition: itkSTAPLEImageFilter.h:176
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::ImageToImageFilter
Base class for filters that take an image as input and produce an image as output.
Definition: itkImageToImageFilter.h:108
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::ImageToImageFilter::InputImagePointer
typename InputImageType::Pointer InputImagePointer
Definition: itkImageToImageFilter.h:130
itk::NumericTraits::OneValue
static T OneValue()
Definition: itkNumericTraits.h:158
itk::ImageToImageFilter::InputImageType
TInputImage InputImageType
Definition: itkImageToImageFilter.h:129
itkImageToImageFilter.h
itk::NumericTraits
Define additional traits for native types such as int or float.
Definition: itkNumericTraits.h:60
itk::NumericTraits::max
static constexpr T max(const T &)
Definition: itkNumericTraits.h:169
itkConceptMacro
#define itkConceptMacro(name, concept)
Definition: itkConceptChecking.h:65
itk::STAPLEImageFilter
The STAPLE filter implements the Simultaneous Truth and Performance Level Estimation algorithm for ge...
Definition: itkSTAPLEImageFilter.h:122
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::STAPLEImageFilter::RealType
typename NumericTraits< InputPixelType >::RealType RealType
Definition: itkSTAPLEImageFilter.h:143
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:139
itk::STAPLEImageFilter::GetSpecificity
double GetSpecificity(unsigned int i)
Definition: itkSTAPLEImageFilter.h:197
itk::ImageSource::OutputImageType
TOutputImage OutputImageType
Definition: itkImageSource.h:90
itk::STAPLEImageFilter::GetSensitivity
double GetSensitivity(unsigned int i)
Definition: itkSTAPLEImageFilter.h:184