ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkSTAPLEImageFilter.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 itkSTAPLEImageFilter_h
19 #define itkSTAPLEImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 #include <vector>
23 
24 namespace itk
25 {
120 template< typename TInputImage, typename TOutputImage >
121 class ITK_TEMPLATE_EXPORT STAPLEImageFilter:
122  public ImageToImageFilter< TInputImage, TOutputImage >
123 {
124 public:
125  ITK_DISALLOW_COPY_AND_ASSIGN(STAPLEImageFilter);
126 
132 
134  itkNewMacro(Self);
135 
137  itkTypeMacro(STAPLEImageFilter, ImageToImageFilter);
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;
151  using InputImagePointer = typename InputImageType::Pointer;
152  using OutputImageType = TOutputImage;
153  using OutputImagePointer = typename OutputImageType::Pointer;
154 
156  using OutputImageRegionType = typename Superclass::OutputImageRegionType;
157 
159  itkSetMacro(ForegroundValue, InputPixelType);
160  itkGetConstMacro(ForegroundValue, InputPixelType);
162 
166  const std::vector< double > & GetSpecificity() const
167  {
168  return m_Specificity;
169  }
170 
174  const std::vector< double > & GetSensitivity() const
175  {
176  return m_Sensitivity;
177  }
178 
181  double GetSensitivity(unsigned int i)
182  {
183  if ( i > this->GetNumberOfIndexedInputs() )
184  {
185  itkExceptionMacro(<< "Array reference out of bounds.");
186  }
187  return m_Sensitivity[i];
188  }
190 
193  double GetSpecificity(unsigned int i)
194  {
195  if ( i > this->GetNumberOfIndexedInputs() )
196  {
197  itkExceptionMacro(<< "Array reference out of bounds.");
198  }
199  return m_Specificity[i];
200  }
202 
206  itkSetMacro(MaximumIterations, unsigned int);
207  itkGetConstMacro(MaximumIterations, unsigned int);
209 
217  itkSetMacro(ConfidenceWeight, double);
218  itkGetConstMacro(ConfidenceWeight, double);
220 
222  itkGetConstMacro(ElapsedIterations, unsigned int);
223 
224 #ifdef ITK_USE_CONCEPT_CHECKING
225  // Begin concept checking
226  itkConceptMacro( InputHasNumericTraitsCheck,
228  // End concept checking
229 #endif
230 
231 protected:
233  {
234  m_ForegroundValue = NumericTraits< InputPixelType >::OneValue();
235  m_MaximumIterations = NumericTraits< unsigned int >::max();
236  m_ElapsedIterations = 0;
237  m_ConfidenceWeight = 1.0;
238  }
239 
240  ~STAPLEImageFilter() override = default;
241  void GenerateData() override;
242 
243  void PrintSelf(std::ostream &, Indent) const override;
244 
245 private:
247  unsigned int m_ElapsedIterations;
248  unsigned int m_MaximumIterations;
249 
251 
252  std::vector< double > m_Sensitivity;
253  std::vector< double > m_Specificity;
254 };
255 } // end namespace itk
256 
257 #ifndef ITK_MANUAL_INSTANTIATION
258 #include "itkSTAPLEImageFilter.hxx"
259 #endif
260 
261 #endif
std::vector< double > m_Sensitivity
The STAPLE filter implements the Simultaneous Truth and Performance Level Estimation algorithm for ge...
typename OutputImageType::Pointer OutputImagePointer
double GetSensitivity(unsigned int i)
Define numeric traits for std::vector.
double GetSpecificity(unsigned int i)
const std::vector< double > & GetSpecificity() const
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
typename TOutputImage::PixelType OutputPixelType
Base class for all process objects that output image data.
std::vector< double > m_Specificity
typename NumericTraits< InputPixelType >::RealType RealType
typename InputImageType::Pointer InputImagePointer
typename OutputImageType::RegionType OutputImageRegionType
TOutputImage OutputImageType
typename TInputImage::PixelType InputPixelType
const std::vector< double > & GetSensitivity() const
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
#define itkConceptMacro(name, concept)