ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkScatterMatrixImageFunction.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 itkScatterMatrixImageFunction_h
19 #define itkScatterMatrixImageFunction_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 ScatterMatrixImageFunction:
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(ScatterMatrixImageFunction);
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<
91 
93  RealType EvaluateAtIndex(const IndexType & index) const override;
94 
96  RealType Evaluate(const PointType & point) const override
97  {
98  IndexType index;
99 
100  this->ConvertPointToNearestIndex(point, index);
101  return this->EvaluateAtIndex(index);
102  }
103 
105  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);
118 
119 protected:
121  ~ScatterMatrixImageFunction() override = default;
122  void PrintSelf(std::ostream & os, Indent indent) const override;
123 
124 private:
125  unsigned int m_NeighborhoodRadius;
126 };
127 } // end namespace itk
128 
129 #ifndef ITK_MANUAL_INSTANTIATION
130 #include "itkScatterMatrixImageFunction.hxx"
131 #endif
132 
133 #endif
Calculate the scatter matrix in the neighborhood of a pixel in a Vector image.
Light weight base class for most itk classes.
Define numeric traits for std::vector.
typename Superclass::PointType PointType
typename Superclass::OutputType OutputType
RealType Evaluate(const PointType &point) const override
typename Superclass::IndexType IndexType
RealType EvaluateAtContinuousIndex(const ContinuousIndexType &cindex) const override
Control indentation during Print() invocation.
Definition: itkIndent.h:49
typename Superclass::ContinuousIndexType ContinuousIndexType
vnl_matrix< typename NumericTraits< typename InputImageType::PixelType::ValueType >::RealType > RealType
Evaluates a function of an image at specified position.