ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkCovarianceImageFunction.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 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 TCoordRep = float >
45 class ITK_TEMPLATE_EXPORT CovarianceImageFunction:
46  public ImageFunction< TInputImage,
47  vnl_matrix<
48  typename NumericTraits< typename TInputImage::PixelType::ValueType >::RealType >,
49  TCoordRep >
50 {
51 public:
52  ITK_DISALLOW_COPY_AND_ASSIGN(CovarianceImageFunction);
53 
56  using Superclass = ImageFunction< TInputImage,
57  vnl_matrix<
59  TCoordRep >;
60 
63 
66 
68  itkNewMacro(Self);
69 
71  using InputImageType = TInputImage;
72 
74  using OutputType = typename Superclass::OutputType;
75 
77  using IndexType = typename Superclass::IndexType;
78 
80  using ContinuousIndexType = typename Superclass::ContinuousIndexType;
81 
83  using PointType = typename Superclass::PointType;
84 
86  static constexpr unsigned int ImageDimension = InputImageType::ImageDimension;
87 
89  using RealType = vnl_matrix<typename NumericTraits<typename InputImageType::PixelType::ValueType>::RealType>;
90 
92  RealType EvaluateAtIndex(const IndexType & index) const override;
93 
95  RealType Evaluate(const PointType & point) const override
96  {
97  IndexType index;
98 
99  this->ConvertPointToNearestIndex(point, index);
100  return this->EvaluateAtIndex(index);
101  }
102 
104  const ContinuousIndexType & cindex) const override
105  {
106  IndexType index;
107 
108  this->ConvertContinuousIndexToNearestIndex(cindex, index);
109  return this->EvaluateAtIndex(index);
110  }
111 
114  itkSetMacro(NeighborhoodRadius, unsigned int);
115  itkGetConstReferenceMacro(NeighborhoodRadius, unsigned int);
117 
118 protected:
120  ~CovarianceImageFunction() override = default;
121  void PrintSelf(std::ostream & os, Indent indent) const override;
122 
123 private:
124  unsigned int m_NeighborhoodRadius{ 1 };
125 };
126 } // end namespace itk
127 
128 #ifndef ITK_MANUAL_INSTANTIATION
129 #include "itkCovarianceImageFunction.hxx"
130 #endif
131 
132 #endif
Light weight base class for most itk classes.
Define numeric traits for std::vector.
Calculate the covariance matrix in the neighborhood of a pixel in a Vector image. ...
RealType Evaluate(const PointType &point) const override
typename Superclass::OutputType OutputType
vnl_matrix< typename NumericTraits< typename InputImageType::PixelType::ValueType >::RealType > RealType
RealType EvaluateAtContinuousIndex(const ContinuousIndexType &cindex) const override
typename Superclass::ContinuousIndexType ContinuousIndexType
typename Superclass::IndexType IndexType
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Evaluates a function of an image at specified position.
typename Superclass::PointType PointType