ITK  4.9.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  >
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 
114 
117 
120 
123 
126 
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  CentralDifferenceImageFunction(const Self &) ITK_DELETE_FUNCTION;
216  void operator=(const Self &) ITK_DELETE_FUNCTION;
217 
218 
220  template<typename T>
222  {
223  typedef T Type;
224  };
225 
227  template< typename Type >
229  template< typename Type >
230  inline void EvaluateAtIndexSpecialized( const IndexType & index, OutputType & derivative, OutputTypeSpecializationStructType<Type>) const;
232 
234  template< typename Type >
236  template< typename Type >
239 
241  // NOTE: for some unknown reason, making these methods inline (as those above are inlined) makes them run *slower*.
242  template< typename Type >
244  template< typename Type >
245  void EvaluateSpecialized( const PointType & point, OutputType & derivative, OutputTypeSpecializationStructType<Type>) const;
247 
248  // flag to take or not the image direction into account
249  // when computing the derivatives.
251 
252  // interpolator
254 };
255 } // end namespace itk
256 
257 #ifndef ITK_MANUAL_INSTANTIATION
258 #include "itkCentralDifferenceImageFunction.hxx"
259 #endif
260 
261 #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