ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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
22#include <vector>
23
24namespace itk
25{
116template <typename TInputImage, typename TOutputImage>
117class ITK_TEMPLATE_EXPORT STAPLEImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
118{
119public:
120 ITK_DISALLOW_COPY_AND_MOVE(STAPLEImageFilter);
121
127
129 itkNewMacro(Self);
130
132 itkOverrideGetNameOfClassMacro(STAPLEImageFilter);
133
136 using OutputPixelType = typename TOutputImage::PixelType;
137 using InputPixelType = typename TInputImage::PixelType;
139
142 static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
143
145 using InputImageType = TInputImage;
146 using InputImagePointer = typename InputImageType::Pointer;
147 using OutputImageType = TOutputImage;
148 using OutputImagePointer = typename OutputImageType::Pointer;
149
152
154 itkSetMacro(ForegroundValue, InputPixelType);
155 itkGetConstMacro(ForegroundValue, InputPixelType);
157
161 const std::vector<double> &
163 {
164 return m_Specificity;
165 }
166
170 const std::vector<double> &
172 {
173 return m_Sensitivity;
174 }
175
178 double
179 GetSensitivity(unsigned int i)
180 {
181 if (i > this->GetNumberOfIndexedInputs())
182 {
183 itkExceptionMacro("Array reference out of bounds.");
184 }
185 return m_Sensitivity[i];
186 }
187
188
191 double
192 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 }
200
201
205 itkSetClampMacro(MaximumIterations, unsigned int, 1, NumericTraits<unsigned int>::max());
206 itkGetConstMacro(MaximumIterations, unsigned int);
208
216 itkSetMacro(ConfidenceWeight, double);
217 itkGetConstMacro(ConfidenceWeight, double);
219
221 itkGetConstMacro(ElapsedIterations, unsigned int);
222
223 itkConceptMacro(InputHasNumericTraitsCheck, (Concept::HasNumericTraits<InputPixelType>));
224
225protected:
233
234 ~STAPLEImageFilter() override = default;
235 void
236 GenerateData() override;
237
238 void
239 PrintSelf(std::ostream &, Indent) const override;
240
241private:
243 unsigned int m_ElapsedIterations{};
244 unsigned int m_MaximumIterations{};
245
247
248 std::vector<double> m_Sensitivity{};
249 std::vector<double> m_Specificity{};
250};
251} // end namespace itk
252
253#ifndef ITK_MANUAL_INSTANTIATION
254# include "itkSTAPLEImageFilter.hxx"
255#endif
256
257#endif
typename OutputImageType::RegionType OutputImageRegionType
Control indentation during Print() invocation.
Definition itkIndent.h:50
static constexpr T max(const T &)
DataObjectPointerArraySizeType GetNumberOfIndexedInputs() const
Get the number of defined Indexed inputs.
const std::vector< double > & GetSensitivity() const
typename TOutputImage::PixelType OutputPixelType
double GetSpecificity(unsigned int i)
~STAPLEImageFilter() override=default
typename InputImageType::Pointer InputImagePointer
std::vector< double > m_Specificity
static constexpr unsigned int ImageDimension
void GenerateData() override
ImageToImageFilter< TInputImage, TOutputImage > Superclass
void PrintSelf(std::ostream &, Indent) const override
typename OutputImageType::Pointer OutputImagePointer
std::vector< double > m_Sensitivity
SmartPointer< Self > Pointer
typename NumericTraits< InputPixelType >::RealType RealType
const std::vector< double > & GetSpecificity() const
SmartPointer< const Self > ConstPointer
double GetSensitivity(unsigned int i)
typename TInputImage::PixelType InputPixelType
Implements transparent reference counting.
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....