ITK  4.13.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 
86  itkStaticConstMacro(ImageDimension, unsigned int,
87  TInputImage::ImageDimension);
88 
91  typedef ImageFunction< TInputImage,
92  TOutputType,
93  TCoordRep > Superclass;
96 
99 
101  itkNewMacro(Self);
102 
104  typedef TInputImage InputImageType;
105 
107  typedef typename InputImageType::PixelType InputPixelType;
108 
111 
113  typedef typename Superclass::OutputType OutputType;
114 
117 
120 
123 
126 
128  typedef typename Superclass::ContinuousIndexType ContinuousIndexType;
129 
132 
134  typedef typename TInputImage::SpacingType SpacingType;
135 
139 
141  virtual void SetInputImage(const TInputImage *inputData) ITK_OVERRIDE;
142 
145  virtual void SetInterpolator(InterpolatorType *interpolator);
146 
148  itkGetModifiableObjectMacro(Interpolator, InterpolatorType );
149 
160  virtual OutputType EvaluateAtIndex(const IndexType & index) const ITK_OVERRIDE;
161 
175  virtual OutputType Evaluate(const PointType & point) const ITK_OVERRIDE;
176 
188  virtual OutputType EvaluateAtContinuousIndex( const ContinuousIndexType & cindex) const ITK_OVERRIDE;
189 
204  itkSetMacro(UseImageDirection, bool);
205  itkGetConstMacro(UseImageDirection, bool);
206  itkBooleanMacro(UseImageDirection);
208 
209 protected:
212  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
213 
214 private:
215  ITK_DISALLOW_COPY_AND_ASSIGN(CentralDifferenceImageFunction);
216 
217 
219  template<typename T>
221  {
222  typedef T Type;
223  };
224 
226  template< typename Type >
227  inline void EvaluateAtIndexSpecialized( const IndexType & index, OutputType & derivative, OutputTypeSpecializationStructType<OutputType>) const;
228  template< typename Type >
229  inline void EvaluateAtIndexSpecialized( const IndexType & index, OutputType & derivative, OutputTypeSpecializationStructType<Type>) const;
231 
233  template< typename Type >
234  inline void EvaluateAtContinuousIndexSpecialized( const ContinuousIndexType & index, OutputType & derivative, OutputTypeSpecializationStructType<OutputType>) const;
235  template< typename Type >
236  inline void EvaluateAtContinuousIndexSpecialized( const ContinuousIndexType & index, OutputType & derivative, OutputTypeSpecializationStructType<Type>) const;
238 
240  // NOTE: for some unknown reason, making these methods inline (as those above are inlined) makes them run *slower*.
241  template< typename Type >
242  void EvaluateSpecialized( const PointType & point, OutputType & derivative, OutputTypeSpecializationStructType<OutputType>) const;
243  template< typename Type >
244  void EvaluateSpecialized( const PointType & point, OutputType & derivative, OutputTypeSpecializationStructType<Type>) const;
246 
247  // flag to take or not the image direction into account
248  // when computing the derivatives.
250 
251  // interpolator
253 };
254 } // end namespace itk
255 
256 #ifndef ITK_MANUAL_INSTANTIATION
257 #include "itkCentralDifferenceImageFunction.hxx"
258 #endif
259 
260 #endif
Light weight base class for most itk classes.
CovariantVector< OutputValueType, itkGetStaticConstMacro(ImageDimension) > ScalarDerivativeType
InterpolateImageFunction< TInputImage, TCoordRep > InterpolatorType
Traits class used to by ConvertPixels to convert blocks of pixels.
Calculate the derivative by central differencing.
ImageFunction< TInputImage, TOutputType, TCoordRep > Superclass
DefaultConvertPixelTraits< InputPixelType > InputPixelConvertType
Base class for all image interpolaters.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
A templated class holding a n-Dimensional covariant vector.
Evaluates a function of an image at specified position.
DefaultConvertPixelTraits< OutputType > OutputConvertType