ITK  6.0.0
Insight Toolkit
itkVarianceImageFunction.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 itkVarianceImageFunction_h
19 #define itkVarianceImageFunction_h
20 
21 #include "itkImageFunction.h"
22 #include "itkNumericTraits.h"
23 
24 namespace itk
25 {
42 template <typename TInputImage, typename TCoordinate = float>
43 class ITK_TEMPLATE_EXPORT VarianceImageFunction
44  : public ImageFunction<TInputImage, typename NumericTraits<typename TInputImage::PixelType>::RealType, TCoordinate>
45 {
46 public:
47  ITK_DISALLOW_COPY_AND_MOVE(VarianceImageFunction);
48 
51  using Superclass =
53 
56 
58  itkOverrideGetNameOfClassMacro(VarianceImageFunction);
59 
61  itkNewMacro(Self);
62 
64  using InputImageType = TInputImage;
65 
67  using typename Superclass::OutputType;
68 
70  using typename Superclass::IndexType;
71 
73  using typename Superclass::ContinuousIndexType;
74 
76  using typename Superclass::PointType;
77 
79  static constexpr unsigned int ImageDimension = InputImageType::ImageDimension;
80 
83 
85  RealType
86  EvaluateAtIndex(const IndexType & index) const override;
87 
89  RealType
90  Evaluate(const PointType & point) const override
91  {
92  IndexType index;
93 
94  this->ConvertPointToNearestIndex(point, index);
95  return this->EvaluateAtIndex(index);
96  }
97 
98  RealType
99  EvaluateAtContinuousIndex(const ContinuousIndexType & cindex) const override
100  {
101  IndexType index;
102 
103  this->ConvertContinuousIndexToNearestIndex(cindex, index);
104  return this->EvaluateAtIndex(index);
105  }
106 
109  itkSetMacro(NeighborhoodRadius, unsigned int);
110  itkGetConstReferenceMacro(NeighborhoodRadius, unsigned int);
113 protected:
115  ~VarianceImageFunction() override = default;
116  void
117  PrintSelf(std::ostream & os, Indent indent) const override;
118 
119 private:
120  unsigned int m_NeighborhoodRadius{};
121 };
122 } // end namespace itk
123 
124 #ifndef ITK_MANUAL_INSTANTIATION
125 # include "itkVarianceImageFunction.hxx"
126 #endif
127 
128 #endif
itkImageFunction.h
itk::VarianceImageFunction::InputImageType
TInputImage InputImageType
Definition: itkVarianceImageFunction.h:64
itk::VarianceImageFunction::EvaluateAtContinuousIndex
RealType EvaluateAtContinuousIndex(const ContinuousIndexType &cindex) const override
Definition: itkVarianceImageFunction.h:99
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::VarianceImageFunction
Calculate the variance in the neighborhood of a pixel.
Definition: itkVarianceImageFunction.h:43
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::ImageFunction
Evaluates a function of an image at specified position.
Definition: itkImageFunction.h:55
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:55
itk::point
*par Constraints *The filter requires an image with at least two dimensions and a vector *length of at least The theory supports extension to scalar but *the implementation of the itk vector classes do not **The template parameter TRealType must be floating point(float or double) or *a user-defined "real" numerical type with arithmetic operations defined *sufficient to compute derivatives. **\par Performance *This filter will automatically multithread if run with *SetUsePrincipleComponents
itk::VarianceImageFunction::RealType
typename NumericTraits< typename InputImageType::PixelType >::RealType RealType
Definition: itkVarianceImageFunction.h:82
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::ContinuousIndex< TCoordinate, Self::ImageDimension >
itkNumericTraits.h
itk::Point
A templated class holding a geometric point in n-Dimensional space.
Definition: itkPoint.h:53
itk::VarianceImageFunction::Evaluate
RealType Evaluate(const PointType &point) const override
Definition: itkVarianceImageFunction.h:90
itk::NumericTraits::RealType
double RealType
Definition: itkNumericTraits.h:86
itk::ImageFunction< TInputImage, NumericTraits< TInputImage::PixelType >::RealType, TCoordinate >::IndexType
typename InputImageType::IndexType IndexType
Definition: itkImageFunction.h:94