ITK  6.0.0
Insight Toolkit
itkCovarianceImageFunction.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 itkCovarianceImageFunction_h
19 #define itkCovarianceImageFunction_h
20 
21 #include "itkImageFunction.h"
22 #include "itkNumericTraits.h"
23 
24 namespace itk
25 {
44 template <typename TInputImage, typename TCoordinate = float>
45 class ITK_TEMPLATE_EXPORT CovarianceImageFunction
46  : public ImageFunction<TInputImage,
47  vnl_matrix<typename NumericTraits<typename TInputImage::PixelType::ValueType>::RealType>,
48  TCoordinate>
49 {
50 public:
51  ITK_DISALLOW_COPY_AND_MOVE(CovarianceImageFunction);
52 
55  using Superclass =
56  ImageFunction<TInputImage,
57  vnl_matrix<typename NumericTraits<typename TInputImage::PixelType::ValueType>::RealType>,
58  TCoordinate>;
59 
62 
64  itkOverrideGetNameOfClassMacro(CovarianceImageFunction);
65 
67  itkNewMacro(Self);
68 
70  using InputImageType = TInputImage;
71 
73  using typename Superclass::OutputType;
74 
76  using typename Superclass::IndexType;
77 
79  using typename Superclass::ContinuousIndexType;
80 
82  using typename Superclass::PointType;
83 
85  static constexpr unsigned int ImageDimension = InputImageType::ImageDimension;
86 
88  using RealType = vnl_matrix<typename NumericTraits<typename InputImageType::PixelType::ValueType>::RealType>;
89 
91  RealType
92  EvaluateAtIndex(const IndexType & index) const override;
93 
95  RealType
96  Evaluate(const PointType & point) const override
97  {
98  IndexType index;
99 
100  this->ConvertPointToNearestIndex(point, index);
101  return this->EvaluateAtIndex(index);
102  }
103 
104  RealType
105  EvaluateAtContinuousIndex(const ContinuousIndexType & cindex) const override
106  {
107  IndexType index;
108 
109  this->ConvertContinuousIndexToNearestIndex(cindex, index);
110  return this->EvaluateAtIndex(index);
111  }
112 
115  itkSetMacro(NeighborhoodRadius, unsigned int);
116  itkGetConstReferenceMacro(NeighborhoodRadius, unsigned int);
119 protected:
121  ~CovarianceImageFunction() override = default;
122  void
123  PrintSelf(std::ostream & os, Indent indent) const override;
124 
125 private:
126  unsigned int m_NeighborhoodRadius{ 1 };
127 };
128 } // end namespace itk
129 
130 #ifndef ITK_MANUAL_INSTANTIATION
131 # include "itkCovarianceImageFunction.hxx"
132 #endif
133 
134 #endif
itkImageFunction.h
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::CovarianceImageFunction
Calculate the covariance matrix in the neighborhood of a pixel in a Vector image.
Definition: itkCovarianceImageFunction.h:45
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::CovarianceImageFunction::Evaluate
RealType Evaluate(const PointType &point) const override
Definition: itkCovarianceImageFunction.h:96
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::CovarianceImageFunction::InputImageType
TInputImage InputImageType
Definition: itkCovarianceImageFunction.h:70
itk::CovarianceImageFunction::RealType
vnl_matrix< typename NumericTraits< typename InputImageType::PixelType::ValueType >::RealType > RealType
Definition: itkCovarianceImageFunction.h:88
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::CovarianceImageFunction::EvaluateAtContinuousIndex
RealType EvaluateAtContinuousIndex(const ContinuousIndexType &cindex) const override
Definition: itkCovarianceImageFunction.h:105
itk::ImageFunction< TInputImage, vnl_matrix< NumericTraits< TInputImage::PixelType::ValueType >::RealType >, TCoordinate >::IndexType
typename InputImageType::IndexType IndexType
Definition: itkImageFunction.h:94