ITK  5.2.0
Insight Toolkit
itkVectorMeanImageFunction.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  * 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 itkVectorMeanImageFunction_h
19 #define itkVectorMeanImageFunction_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 VectorMeanImageFunction
46  : public ImageFunction<TInputImage, typename NumericTraits<typename TInputImage::PixelType>::RealType, TCoordRep>
47 {
48 public:
49  ITK_DISALLOW_COPY_AND_MOVE(VectorMeanImageFunction);
50 
53  using Superclass =
55 
58 
61 
63  itkNewMacro(Self);
64 
66  using InputImageType = TInputImage;
67 
69  using OutputType = typename Superclass::OutputType;
70 
72  using IndexType = typename Superclass::IndexType;
73 
75  using ContinuousIndexType = typename Superclass::ContinuousIndexType;
76 
78  using PointType = typename Superclass::PointType;
79 
81  static constexpr unsigned int ImageDimension = InputImageType::ImageDimension;
82 
85 
87  RealType
88  EvaluateAtIndex(const IndexType & index) const override;
89 
91  RealType
92  Evaluate(const PointType & point) const override
93  {
94  IndexType index;
95 
96  this->ConvertPointToNearestIndex(point, index);
97  return this->EvaluateAtIndex(index);
98  }
99 
100  RealType
101  EvaluateAtContinuousIndex(const ContinuousIndexType & cindex) const override
102  {
103  IndexType index;
104 
105  this->ConvertContinuousIndexToNearestIndex(cindex, index);
106  return this->EvaluateAtIndex(index);
107  }
108 
111  itkSetMacro(NeighborhoodRadius, unsigned int);
112  itkGetConstReferenceMacro(NeighborhoodRadius, unsigned int);
114 
115 protected:
117  ~VectorMeanImageFunction() override = default;
118  void
119  PrintSelf(std::ostream & os, Indent indent) const override;
120 
121 private:
122  unsigned int m_NeighborhoodRadius;
123 };
124 } // end namespace itk
125 
126 #ifndef ITK_MANUAL_INSTANTIATION
127 # include "itkVectorMeanImageFunction.hxx"
128 #endif
129 
130 #endif
itk::VectorMeanImageFunction::IndexType
typename Superclass::IndexType IndexType
Definition: itkVectorMeanImageFunction.h:72
itkImageFunction.h
itk::VectorMeanImageFunction::Evaluate
RealType Evaluate(const PointType &point) const override
Definition: itkVectorMeanImageFunction.h:92
itk::VectorMeanImageFunction::EvaluateAtContinuousIndex
RealType EvaluateAtContinuousIndex(const ContinuousIndexType &cindex) const override
Definition: itkVectorMeanImageFunction.h:101
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::VectorMeanImageFunction
Calculate the mean value in the neighborhood of a pixel in a Vector image.
Definition: itkVectorMeanImageFunction.h:45
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::VectorMeanImageFunction::m_NeighborhoodRadius
unsigned int m_NeighborhoodRadius
Definition: itkVectorMeanImageFunction.h:122
itk::VectorMeanImageFunction::RealType
typename NumericTraits< typename TInputImage::PixelType >::RealType RealType
Definition: itkVectorMeanImageFunction.h:84
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:59
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::VectorMeanImageFunction::InputImageType
TInputImage InputImageType
Definition: itkVectorMeanImageFunction.h:66
itkNumericTraits.h
itk::NumericTraits::RealType
double RealType
Definition: itkNumericTraits.h:84
itk::VectorMeanImageFunction::ContinuousIndexType
typename Superclass::ContinuousIndexType ContinuousIndexType
Definition: itkVectorMeanImageFunction.h:75
itk::VectorMeanImageFunction::OutputType
typename Superclass::OutputType OutputType
Definition: itkVectorMeanImageFunction.h:69
itk::VectorMeanImageFunction::PointType
typename Superclass::PointType PointType
Definition: itkVectorMeanImageFunction.h:78