ITK  4.8.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 #include "itkEnableIf.h"
26 #include "itkIsSame.h"
27 
28 namespace itk
29 {
75 template<
76  typename TInputImage,
77  typename TCoordRep = float,
78  typename TOutputType = CovariantVector<double, TInputImage::ImageDimension >
79  >
81  public ImageFunction< TInputImage,
82  TOutputType,
83  TCoordRep >
84 {
85 public:
86 
88  itkStaticConstMacro(ImageDimension, unsigned int,
89  TInputImage::ImageDimension);
90 
93  typedef ImageFunction< TInputImage,
94  TOutputType,
95  TCoordRep > Superclass;
98 
101 
103  itkNewMacro(Self);
104 
106  typedef TInputImage InputImageType;
107 
109  typedef typename InputImageType::PixelType InputPixelType;
110 
113 
116 
119 
122 
125 
128 
131 
134 
136  typedef typename TInputImage::SpacingType SpacingType;
137 
141 
143  virtual void SetInputImage(const TInputImage *inputData) ITK_OVERRIDE;
144 
147  virtual void SetInterpolator(InterpolatorType *interpolator);
148 
150  itkGetModifiableObjectMacro(Interpolator, InterpolatorType );
151 
162  virtual OutputType EvaluateAtIndex(const IndexType & index) const ITK_OVERRIDE;
163 
177  virtual OutputType Evaluate(const PointType & point) const ITK_OVERRIDE;
178 
190  virtual OutputType EvaluateAtContinuousIndex( const ContinuousIndexType & cindex) const ITK_OVERRIDE;
191 
206  itkSetMacro(UseImageDirection, bool);
207  itkGetConstMacro(UseImageDirection, bool);
208  itkBooleanMacro(UseImageDirection);
210 
211 protected:
214  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
215 
216 private:
217  CentralDifferenceImageFunction(const Self &); //purposely not implemented
218  void operator=(const Self &); //purposely not implemented
219 
220 
222  template<typename T>
224  {
225  typedef T Type;
226  };
227 
229  template< typename Type >
231  template< typename Type >
232  inline void EvaluateAtIndexSpecialized( const IndexType & index, OutputType & derivative, OutputTypeSpecializationStructType<Type>) const;
234 
236  template< typename Type >
238  template< typename Type >
241 
243  // NOTE: for some unknown reason, making these methods inline (as those above are inlined) makes them run *slower*.
244  template< typename Type >
246  template< typename Type >
247  void EvaluateSpecialized( const PointType & point, OutputType & derivative, OutputTypeSpecializationStructType<Type>) const;
249 
250  // flag to take or not the image direction into account
251  // when computing the derivatives.
253 
254  // interpolator
256 };
257 } // end namespace itk
258 
259 #ifndef ITK_MANUAL_INSTANTIATION
260 #include "itkCentralDifferenceImageFunction.hxx"
261 #endif
262 
263 #endif
Light weight base class for most itk classes.
Point< TCoordRep, itkGetStaticConstMacro(ImageDimension) > PointType
CovariantVector< OutputValueType, itkGetStaticConstMacro(ImageDimension) > ScalarDerivativeType
InterpolateImageFunction< TInputImage, TCoordRep > InterpolatorType
void EvaluateAtIndexSpecialized(const IndexType &index, OutputType &derivative, OutputTypeSpecializationStructType< OutputType >) const
Traits class used to by ConvertPixels to convert blocks of pixels.
Calculate the derivative by central differencing.
void PrintSelf(std::ostream &os, Indent indent) const override
ImageFunction< TInputImage, TOutputType, TCoordRep > Superclass
virtual OutputType EvaluateAtIndex(const IndexType &index) const override
virtual void SetInputImage(const TInputImage *inputData) override
DefaultConvertPixelTraits< InputPixelType > InputPixelConvertType
virtual OutputType Evaluate(const PointType &point) const override
ContinuousIndex< TCoordRep, itkGetStaticConstMacro(ImageDimension) > ContinuousIndexType
virtual void SetInterpolator(InterpolatorType *interpolator)
Base class for all image interpolaters.
void EvaluateAtContinuousIndexSpecialized(const ContinuousIndexType &index, OutputType &derivative, OutputTypeSpecializationStructType< OutputType >) const
Control indentation during Print() invocation.
Definition: itkIndent.h:49
virtual OutputType EvaluateAtContinuousIndex(const ContinuousIndexType &cindex) const override
void EvaluateSpecialized(const PointType &point, OutputType &derivative, OutputTypeSpecializationStructType< OutputType >) const
A templated class holding a n-Dimensional covariant vector.
Evaluates a function of an image at specified position.
DefaultConvertPixelTraits< OutputType > OutputConvertType