ITK  4.6.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 >
122  public ImageToImageFilter< TInputImage, TOutputImage >
123 {
124 public:
130 
132  itkNewMacro(Self);
133 
136 
139  typedef typename TOutputImage::PixelType OutputPixelType;
140  typedef typename TInputImage::PixelType InputPixelType;
142 
145  itkStaticConstMacro(ImageDimension, unsigned int,
146  TOutputImage::ImageDimension);
147 
149  typedef TInputImage InputImageType;
150  typedef typename InputImageType::Pointer InputImagePointer;
151  typedef TOutputImage OutputImageType;
152  typedef typename OutputImageType::Pointer OutputImagePointer;
153 
156 
158  itkSetMacro(ForegroundValue, InputPixelType);
159  itkGetConstMacro(ForegroundValue, InputPixelType);
161 
165  const std::vector< double > & GetSpecificity() const
166  {
167  return m_Specificity;
168  }
169 
173  const std::vector< double > & GetSensitivity() const
174  {
175  return m_Sensitivity;
176  }
177 
180  double GetSensitivity(unsigned int i)
181  {
182  if ( i > this->GetNumberOfIndexedInputs() )
183  {
184  itkExceptionMacro(<< "Array reference out of bounds.");
185  }
186  return m_Sensitivity[i];
187  }
189 
192  double GetSpecificity(unsigned int i)
193  {
194  if ( i > this->GetNumberOfIndexedInputs() )
195  {
196  itkExceptionMacro(<< "Array reference out of bounds.");
197  }
198  return m_Specificity[i];
199  }
201 
205  itkSetMacro(MaximumIterations, unsigned int);
206  itkGetConstMacro(MaximumIterations, unsigned int);
208 
216  itkSetMacro(ConfidenceWeight, double);
217  itkGetConstMacro(ConfidenceWeight, double);
219 
221  itkGetConstMacro(ElapsedIterations, unsigned int);
222 
223 #ifdef ITK_USE_CONCEPT_CHECKING
224  // Begin concept checking
225  itkConceptMacro( InputHasNumericTraitsCheck,
227  // End concept checking
228 #endif
229 
230 protected:
232  {
236  m_ConfidenceWeight = 1.0;
237  }
238 
239  virtual ~STAPLEImageFilter() {}
240  void GenerateData();
241 
242  void PrintSelf(std::ostream &, Indent) const;
243 
244 private:
245  STAPLEImageFilter(const Self &); //purposely not implemented
246  void operator=(const Self &); //purposely not implemented
247 
249  unsigned int m_ElapsedIterations;
250  unsigned int m_MaximumIterations;
251 
253 
254  std::vector< double > m_Sensitivity;
255  std::vector< double > m_Specificity;
256 };
257 } // end namespace itk
258 
259 #ifndef ITK_MANUAL_INSTANTIATION
260 #include "itkSTAPLEImageFilter.hxx"
261 #endif
262 
263 #endif
std::vector< double > m_Sensitivity
NumericTraits< InputPixelType >::RealType RealType
The STAPLE filter implements the Simultaneous Truth and Performance Level Estimation algorithm for ge...
double GetSensitivity(unsigned int i)
InputImageType::Pointer InputImagePointer
void operator=(const Self &)
double GetSpecificity(unsigned int i)
const std::vector< double > & GetSpecificity() const
void PrintSelf(std::ostream &, Indent) const
Base class for all process objects that output image data.
std::vector< double > m_Specificity
TOutputImage::PixelType OutputPixelType
ImageToImageFilter< TInputImage, TOutputImage > Superclass
SmartPointer< const Self > ConstPointer
TInputImage::PixelType InputPixelType
static T max(const T &)
static const unsigned int ImageDimension
SmartPointer< Self > Pointer
DataObjectPointerArraySizeType GetNumberOfIndexedInputs() const
const std::vector< double > & GetSensitivity() const
Superclass::OutputImageRegionType OutputImageRegionType
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
Superclass::OutputImageRegionType OutputImageRegionType
Define additional traits for native types such as int or float.
OutputImageType::Pointer OutputImagePointer
#define itkConceptMacro(name, concept)