go home Home | Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Data Structures | File List | Namespace Members | Data Fields | Globals | Related Pages
itkReducedDimensionBSplineInterpolateImageFunction.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright UMC Utrecht and contributors
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 
19 /*=========================================================================
20 
21  Program: Insight Segmentation & Registration Toolkit
22  Module: $RCSfile: itkReducedDimBSplineInterpolateImageFunction.h,v $
23  Language: C++
24  Date: $Date: 2009-04-25 12:27:05 $
25  Version: $Revision: 1.24 $
26 
27  Copyright (c) Insight Software Consortium. All rights reserved.
28  See ITKCopyright.txt or http://www.itk.org/HTML/Copyright.htm for details.
29 
30  Portions of this code are covered under the VTK copyright.
31  See VTKCopyright.txt or http://www.kitware.com/VTKCopyright.htm for details.
32 
33  This software is distributed WITHOUT ANY WARRANTY; without even
34  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
35  PURPOSE. See the above copyright notices for more information.
36 
37 =========================================================================*/
38 #ifndef __itkReducedDimensionBSplineInterpolateImageFunction_h
39 #define __itkReducedDimensionBSplineInterpolateImageFunction_h
40 
41 #include <vector>
42 
43 #include "itkImageLinearIteratorWithIndex.h"
44 #include "itkInterpolateImageFunction.h"
45 #include "vnl/vnl_matrix.h"
46 
48 #include "itkConceptChecking.h"
49 #include "itkCovariantVector.h"
50 
51 namespace itk
52 {
87 template<
88 class TImageType,
89 class TCoordRep = double,
90 class TCoefficientType = double >
92  public InterpolateImageFunction< TImageType, TCoordRep >
93 {
94 public:
95 
98  typedef InterpolateImageFunction< TImageType, TCoordRep > Superclass;
100  typedef SmartPointer< const Self > ConstPointer;
101 
103  itkTypeMacro( ReducedDimensionBSplineInterpolateImageFunction, InterpolateImageFunction );
104 
106  itkNewMacro( Self );
107 
109  typedef typename Superclass::OutputType OutputType;
110 
112  typedef typename Superclass::InputImageType InputImageType;
113 
115  itkStaticConstMacro( ImageDimension, unsigned int, Superclass::ImageDimension );
116 
118  typedef typename Superclass::IndexType IndexType;
119 
121  typedef typename Superclass::ContinuousIndexType ContinuousIndexType;
122 
124  typedef typename Superclass::PointType PointType;
125 
127  typedef ImageLinearIteratorWithIndex< TImageType > Iterator;
128 
130  typedef TCoefficientType CoefficientDataType;
131  typedef Image< CoefficientDataType,
132  itkGetStaticConstMacro( ImageDimension )
134 
138 
140 
150  const ContinuousIndexType & index ) const;
151 
153  typedef CovariantVector< OutputType,
154  itkGetStaticConstMacro( ImageDimension )
156 
158  {
159  ContinuousIndexType index;
160  this->GetInputImage()->TransformPhysicalPointToContinuousIndex( point, index );
161  return ( this->EvaluateDerivativeAtContinuousIndex( index ) );
162  }
163 
164 
166  const ContinuousIndexType & x ) const;
167 
170  void SetSplineOrder( unsigned int SplineOrder );
171 
172  itkGetConstMacro( SplineOrder, int );
173 
175  virtual void SetInputImage( const TImageType * inputData );
176 
189  itkSetMacro( UseImageDirection, bool );
190  itkGetConstMacro( UseImageDirection, bool );
191  itkBooleanMacro( UseImageDirection );
192 
193 protected:
194 
197  void PrintSelf( std::ostream & os, Indent indent ) const;
198 
199  // These are needed by the smoothing spline routine.
200  std::vector< CoefficientDataType > m_Scratch; // temp storage for processing of Coefficients
201  typename TImageType::SizeType m_DataLength; // Image size
202  unsigned int m_SplineOrder; // User specified spline order (3rd or cubic is the default)
203 
204  typename CoefficientImageType::ConstPointer m_Coefficients; // Spline coefficients
205 
206 private:
207 
208  ReducedDimensionBSplineInterpolateImageFunction( const Self & ); //purposely not implemented
209  void operator=( const Self & ); //purposely not implemented
210 
213  const vnl_matrix< long > & EvaluateIndex,
214  vnl_matrix< double > & weights,
215  unsigned int splineOrder ) const;
216 
219  const vnl_matrix< long > & EvaluateIndex,
220  vnl_matrix< double > & weights,
221  unsigned int splineOrder ) const;
222 
226 
228  void DetermineRegionOfSupport( vnl_matrix< long > & evaluateIndex,
229  const ContinuousIndexType & x,
230  unsigned int splineOrder ) const;
231 
234  void ApplyMirrorBoundaryConditions( vnl_matrix< long > & evaluateIndex,
235  unsigned int splineOrder ) const;
236 
237  Iterator m_CIterator; // Iterator for traversing spline coefficients.
238  unsigned long m_MaxNumberInterpolationPoints; // number of neighborhood points used for interpolation
239  std::vector< IndexType > m_PointsToIndex; // Preallocation of interpolation neighborhood indicies
240 
242 
243  // flag to take or not the image direction into account when computing the
244  // derivatives.
246 
247 };
248 
249 } // namespace itk
250 
251 #ifndef ITK_MANUAL_INSTANTIATION
252 #include "itkReducedDimensionBSplineInterpolateImageFunction.hxx"
253 #endif
254 
255 #endif
itk::ReducedDimensionBSplineInterpolateImageFunction::PrintSelf
void PrintSelf(std::ostream &os, Indent indent) const
itk::ReducedDimensionBSplineInterpolateImageFunction::m_DataLength
TImageType::SizeType m_DataLength
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:201
itk::ReducedDimensionBSplineInterpolateImageFunction::SetDerivativeWeights
void SetDerivativeWeights(const ContinuousIndexType &x, const vnl_matrix< long > &EvaluateIndex, vnl_matrix< double > &weights, unsigned int splineOrder) const
itk::ReducedDimensionBSplineInterpolateImageFunction::CoefficientDataType
TCoefficientType CoefficientDataType
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:130
itk::ReducedDimensionBSplineInterpolateImageFunction::~ReducedDimensionBSplineInterpolateImageFunction
virtual ~ReducedDimensionBSplineInterpolateImageFunction()
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:196
SmartPointer< Self >
itk::ReducedDimensionBSplineInterpolateImageFunction::OutputType
Superclass::OutputType OutputType
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:106
itk::ReducedDimensionBSplineInterpolateImageFunction::ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:100
itk::ReducedDimensionBSplineInterpolateImageFunction::Self
ReducedDimensionBSplineInterpolateImageFunction Self
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:97
itk::ReducedDimensionBSplineInterpolateImageFunction::Superclass
InterpolateImageFunction< TImageType, TCoordRep > Superclass
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:98
itk::MultiOrderBSplineDecompositionImageFilter
Calculates the B-Spline coefficients of an image. Spline order may be per dimension from 0 to 5 per.
Definition: itkMultiOrderBSplineDecompositionImageFilter.h:84
itk::ReducedDimensionBSplineInterpolateImageFunction
Evaluates the B-Spline interpolation of an image. Spline order may be from 0 to 5.
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:93
itk::ReducedDimensionBSplineInterpolateImageFunction::InputImageType
Superclass::InputImageType InputImageType
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:112
Image
ImageLinearIteratorWithIndex< InterpolatorBase< TElastix >::InputImageType >
itk::ReducedDimensionBSplineInterpolateImageFunction::GeneratePointsToIndex
void GeneratePointsToIndex()
itk::ReducedDimensionBSplineInterpolateImageFunction::ReducedDimensionBSplineInterpolateImageFunction
ReducedDimensionBSplineInterpolateImageFunction()
itk::ReducedDimensionBSplineInterpolateImageFunction::SetSplineOrder
void SetSplineOrder(unsigned int SplineOrder)
itk::ReducedDimensionBSplineInterpolateImageFunction::ReducedDimensionBSplineInterpolateImageFunction
ReducedDimensionBSplineInterpolateImageFunction(const Self &)
itk::ReducedDimensionBSplineInterpolateImageFunction::itkStaticConstMacro
itkStaticConstMacro(ImageDimension, unsigned int, Superclass::ImageDimension)
itk::ReducedDimensionBSplineInterpolateImageFunction::CoefficientFilter
MultiOrderBSplineDecompositionImageFilter< TImageType, CoefficientImageType > CoefficientFilter
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:137
itk::ReducedDimensionBSplineInterpolateImageFunction::CoefficientFilterPointer
CoefficientFilter::Pointer CoefficientFilterPointer
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:139
itk::ReducedDimensionBSplineInterpolateImageFunction::m_PointsToIndex
std::vector< IndexType > m_PointsToIndex
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:239
double
itkMultiOrderBSplineDecompositionImageFilter.h
itk::ReducedDimensionBSplineInterpolateImageFunction::Iterator
ImageLinearIteratorWithIndex< TImageType > Iterator
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:127
itk::ReducedDimensionBSplineInterpolateImageFunction::m_SplineOrder
unsigned int m_SplineOrder
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:202
itk::ReducedDimensionBSplineInterpolateImageFunction::ApplyMirrorBoundaryConditions
void ApplyMirrorBoundaryConditions(vnl_matrix< long > &evaluateIndex, unsigned int splineOrder) const
itk::ReducedDimensionBSplineInterpolateImageFunction::EvaluateDerivativeAtContinuousIndex
CovariantVectorType EvaluateDerivativeAtContinuousIndex(const ContinuousIndexType &x) const
itk::ReducedDimensionBSplineInterpolateImageFunction::ContinuousIndexType
Superclass::ContinuousIndexType ContinuousIndexType
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:121
itk::ReducedDimensionBSplineInterpolateImageFunction::CoefficientImageType
Image< CoefficientDataType, itkGetStaticConstMacro(ImageDimension) > CoefficientImageType
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:133
itk::ReducedDimensionBSplineInterpolateImageFunction::DetermineRegionOfSupport
void DetermineRegionOfSupport(vnl_matrix< long > &evaluateIndex, const ContinuousIndexType &x, unsigned int splineOrder) const
itk::ReducedDimensionBSplineInterpolateImageFunction::m_Coefficients
CoefficientImageType::ConstPointer m_Coefficients
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:204
itk::ReducedDimensionBSplineInterpolateImageFunction::CovariantVectorType
CovariantVector< OutputType, itkGetStaticConstMacro(ImageDimension) > CovariantVectorType
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:155
itk::ReducedDimensionBSplineInterpolateImageFunction::m_UseImageDirection
bool m_UseImageDirection
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:245
itk::ReducedDimensionBSplineInterpolateImageFunction::EvaluateDerivative
CovariantVectorType EvaluateDerivative(const PointType &point) const
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:157
itk::ReducedDimensionBSplineInterpolateImageFunction::SetInterpolationWeights
void SetInterpolationWeights(const ContinuousIndexType &x, const vnl_matrix< long > &EvaluateIndex, vnl_matrix< double > &weights, unsigned int splineOrder) const
itk::ReducedDimensionBSplineInterpolateImageFunction::m_Scratch
std::vector< CoefficientDataType > m_Scratch
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:200
vnl_matrix< double >
itk
Definition: itkAdvancedImageToImageMetric.h:40
itk::ReducedDimensionBSplineInterpolateImageFunction::Pointer
SmartPointer< Self > Pointer
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:99
itk::ReducedDimensionBSplineInterpolateImageFunction::m_CIterator
Iterator m_CIterator
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:237
itk::ReducedDimensionBSplineInterpolateImageFunction::SetInputImage
virtual void SetInputImage(const TImageType *inputData)
itk::ReducedDimensionBSplineInterpolateImageFunction::m_CoefficientFilter
CoefficientFilterPointer m_CoefficientFilter
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:241
itk::ReducedDimensionBSplineInterpolateImageFunction::IndexType
Superclass::IndexType IndexType
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:118
itk::ReducedDimensionBSplineInterpolateImageFunction::m_MaxNumberInterpolationPoints
unsigned long m_MaxNumberInterpolationPoints
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:238
itk::ReducedDimensionBSplineInterpolateImageFunction::PointType
Superclass::PointType PointType
Definition: itkReducedDimensionBSplineInterpolateImageFunction.h:124
itk::ReducedDimensionBSplineInterpolateImageFunction::operator=
void operator=(const Self &)
itk::ReducedDimensionBSplineInterpolateImageFunction::EvaluateAtContinuousIndex
virtual OutputType EvaluateAtContinuousIndex(const ContinuousIndexType &index) const


Generated on OURCE_DATE_EPOCH for elastix by doxygen 1.8.18 elastix logo