ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkCentralDifferenceImageFunction.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 itkCentralDifferenceImageFunction_h
19 #define itkCentralDifferenceImageFunction_h
20 
21 #include "itkImageFunction.h"
22 #include "itkCovariantVector.h"
25 
26 namespace itk
27 {
73 template<
74  typename TInputImage,
75  typename TCoordRep = float,
76  typename TOutputType = CovariantVector<double, TInputImage::ImageDimension >
77  >
78 class ITK_TEMPLATE_EXPORT CentralDifferenceImageFunction:
79  public ImageFunction< TInputImage,
80  TOutputType,
81  TCoordRep >
82 {
83 public:
84  ITK_DISALLOW_COPY_AND_ASSIGN(CentralDifferenceImageFunction);
85 
87  static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
88 
91  using Superclass = ImageFunction< TInputImage,
92  TOutputType,
93  TCoordRep >;
96 
99 
101  itkNewMacro(Self);
102 
104  using InputImageType = TInputImage;
105 
107  using InputPixelType = typename InputImageType::PixelType;
108 
111 
113  using OutputType = typename Superclass::OutputType;
114 
117 
120 
123 
126 
128  using ContinuousIndexType = typename Superclass::ContinuousIndexType;
129 
132 
134  using SpacingType = typename TInputImage::SpacingType;
135 
139 
141  void SetInputImage(const TInputImage *inputData) override;
142 
145  virtual void SetInterpolator(InterpolatorType *interpolator);
146 
148  itkGetModifiableObjectMacro(Interpolator, InterpolatorType );
149 
160  OutputType EvaluateAtIndex(const IndexType & index) const override;
161 
175  OutputType Evaluate(const PointType & point) const override;
176 
188  OutputType EvaluateAtContinuousIndex( const ContinuousIndexType & cindex) const override;
189 
204  itkSetMacro(UseImageDirection, bool);
205  itkGetConstMacro(UseImageDirection, bool);
206  itkBooleanMacro(UseImageDirection);
208 
209 protected:
211  ~CentralDifferenceImageFunction() override = default;
212  void PrintSelf(std::ostream & os, Indent indent) const override;
213 
214 private:
215 
217  template<typename T>
219  {
220  using Type = T;
221  };
222 
224  template< typename Type >
225  inline void EvaluateAtIndexSpecialized( const IndexType & index, OutputType & derivative, OutputTypeSpecializationStructType<OutputType>) const;
226  template< typename Type >
227  inline void EvaluateAtIndexSpecialized( const IndexType & index, OutputType & derivative, OutputTypeSpecializationStructType<Type>) const;
229 
231  template< typename Type >
232  inline void EvaluateAtContinuousIndexSpecialized( const ContinuousIndexType & index, OutputType & derivative, OutputTypeSpecializationStructType<OutputType>) const;
233  template< typename Type >
234  inline void EvaluateAtContinuousIndexSpecialized( const ContinuousIndexType & index, OutputType & derivative, OutputTypeSpecializationStructType<Type>) const;
236 
238  // NOTE: for some unknown reason, making these methods inline (as those above are inlined) makes them run *slower*.
239  template< typename Type >
240  void EvaluateSpecialized( const PointType & point, OutputType & derivative, OutputTypeSpecializationStructType<OutputType>) const;
241  template< typename Type >
242  void EvaluateSpecialized( const PointType & point, OutputType & derivative, OutputTypeSpecializationStructType<Type>) const;
244 
245  // flag to take or not the image direction into account
246  // when computing the derivatives.
248 
249  // interpolator
251 };
252 } // end namespace itk
253 
254 #ifndef ITK_MANUAL_INSTANTIATION
255 #include "itkCentralDifferenceImageFunction.hxx"
256 #endif
257 
258 #endif
typename Superclass::ContinuousIndexType ContinuousIndexType
Light weight base class for most itk classes.
typename OutputConvertType::ComponentType OutputValueType
Traits class used to by ConvertPixels to convert blocks of pixels.
Calculate the derivative by central differencing.
typename InterpolatorType::Pointer InterpolatorPointer
Base class for all image interpolaters.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
typename PixelType::ComponentType ComponentType
A templated class holding a n-Dimensional covariant vector.
Evaluates a function of an image at specified position.
typename InputImageType::PixelType InputPixelType