ITK  6.0.0
Insight Toolkit
itkCentralDifferenceImageFunction.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 itkCentralDifferenceImageFunction_h
19 #define itkCentralDifferenceImageFunction_h
20 
21 #include "itkImageFunction.h"
22 #include "itkCovariantVector.h"
25 
26 namespace itk
27 {
73 template <typename TInputImage,
74  typename TCoordRep = float,
75  typename TOutputType = CovariantVector<double, TInputImage::ImageDimension>>
76 class ITK_TEMPLATE_EXPORT CentralDifferenceImageFunction : public ImageFunction<TInputImage, TOutputType, TCoordRep>
77 {
78 public:
79  ITK_DISALLOW_COPY_AND_MOVE(CentralDifferenceImageFunction);
80 
82  static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
83 
89 
91  itkOverrideGetNameOfClassMacro(CentralDifferenceImageFunction);
92 
94  itkNewMacro(Self);
95 
97  using InputImageType = TInputImage;
98 
100  using InputPixelType = typename InputImageType::PixelType;
101 
104 
106  using typename Superclass::OutputType;
107 
110 
113 
116 
118  using typename Superclass::IndexType;
119 
121  using typename Superclass::ContinuousIndexType;
122 
124  using typename Superclass::PointType;
125 
127  using SpacingType = typename TInputImage::SpacingType;
128 
132 
134  void
135  SetInputImage(const TInputImage * inputData) override;
136 
139  virtual void
140  SetInterpolator(InterpolatorType * interpolator);
141 
143  itkGetModifiableObjectMacro(Interpolator, InterpolatorType);
144 
155  OutputType
156  EvaluateAtIndex(const IndexType & index) const override;
157 
171  OutputType
172  Evaluate(const PointType & point) const override;
173 
185  OutputType
186  EvaluateAtContinuousIndex(const ContinuousIndexType & cindex) const override;
187 
202  itkSetMacro(UseImageDirection, bool);
203  itkGetConstMacro(UseImageDirection, bool);
204  itkBooleanMacro(UseImageDirection);
207 protected:
209  ~CentralDifferenceImageFunction() override = default;
210  void
211  PrintSelf(std::ostream & os, Indent indent) const override;
212 
213 private:
215  template <typename T>
217  {
218  using Type = T;
219  };
220 
222  template <typename Type>
223  inline void
224  EvaluateAtIndexSpecialized(const IndexType & index,
225  OutputType & orientedDerivative,
227 
229  template <typename Type>
230  inline void
231  EvaluateAtIndexSpecialized(const IndexType & index,
232  OutputType & derivative,
234 
236  template <typename Type>
237  inline void
238  EvaluateAtContinuousIndexSpecialized(const ContinuousIndexType & cindex,
239  OutputType & orientedDerivative,
241 
243  template <typename Type>
244  inline void
245  EvaluateAtContinuousIndexSpecialized(const ContinuousIndexType & cindex,
246  OutputType & derivative,
248 
250  // NOTE: for some unknown reason, making these methods inline (as those above are inlined) makes them run *slower*.
251  template <typename Type>
252  void
253  EvaluateSpecialized(const PointType & point,
254  OutputType & orientedDerivative,
256 
258  template <typename Type>
259  void
260  EvaluateSpecialized(const PointType & point, OutputType & derivative, OutputTypeSpecializationStructType<Type>) const;
261 
262  // flag to take or not the image direction into account
263  // when computing the derivatives.
264  bool m_UseImageDirection{ true };
265 
266  // interpolator
267  InterpolatorPointer m_Interpolator{};
268 };
269 } // end namespace itk
270 
271 #ifndef ITK_MANUAL_INSTANTIATION
272 # include "itkCentralDifferenceImageFunction.hxx"
273 #endif
274 
275 #endif
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itkImageFunction.h
itk::CentralDifferenceImageFunction::OutputTypeSpecializationStructType::Type
T Type
Definition: itkCentralDifferenceImageFunction.h:218
itkCovariantVector.h
itk::ImageFunction< TInputImage, TOutputType, TCoordRep >::IndexType
typename InputImageType::IndexType IndexType
Definition: itkImageFunction.h:90
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::CentralDifferenceImageFunction
Calculate the derivative by central differencing.
Definition: itkCentralDifferenceImageFunction.h:76
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::FunctionBase< Point< TCoordRep, TInputImage::ImageDimension >, TOutputType >::OutputType
TOutputType OutputType
Definition: itkFunctionBase.h:62
itk::CentralDifferenceImageFunction::OutputTypeSpecializationStructType
Definition: itkCentralDifferenceImageFunction.h:216
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itkDefaultConvertPixelTraits.h
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::CentralDifferenceImageFunction::InputPixelType
typename InputImageType::PixelType InputPixelType
Definition: itkCentralDifferenceImageFunction.h:100
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::DefaultConvertPixelTraits
Traits class used to by ConvertPixels to convert blocks of pixels.
Definition: itkDefaultConvertPixelTraits.h:41
itk::CentralDifferenceImageFunction::InputImageType
TInputImage InputImageType
Definition: itkCentralDifferenceImageFunction.h:97
itk::CentralDifferenceImageFunction::OutputValueType
typename OutputConvertType::ComponentType OutputValueType
Definition: itkCentralDifferenceImageFunction.h:112
itk::CovariantVector
A templated class holding a n-Dimensional covariant vector.
Definition: itkCovariantVector.h:70
itk::CentralDifferenceImageFunction::SpacingType
typename TInputImage::SpacingType SpacingType
Definition: itkCentralDifferenceImageFunction.h:127
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::ContinuousIndex< TCoordRep, Self::ImageDimension >
itk::Point
A templated class holding a geometric point in n-Dimensional space.
Definition: itkPoint.h:53
itkInterpolateImageFunction.h
itk::InterpolateImageFunction
Base class for all image interpolators.
Definition: itkInterpolateImageFunction.h:45
itk::CentralDifferenceImageFunction::InterpolatorPointer
typename InterpolatorType::Pointer InterpolatorPointer
Definition: itkCentralDifferenceImageFunction.h:131
itk::DefaultConvertPixelTraits::ComponentType
typename PixelType::ComponentType ComponentType
Definition: itkDefaultConvertPixelTraits.h:45