ITK  5.4.0
Insight Toolkit
itkBSplineInterpolateImageFunction.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright NumFOCUS
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  * https://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  * Portions of this file are subject to the VTK Toolkit Version 3 copyright.
21  *
22  * Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
23  *
24  * For complete copyright, license and disclaimer of warranty information
25  * please refer to the NOTICE file at the top of the ITK source tree.
26  *
27  *=========================================================================*/
28 #ifndef itkBSplineInterpolateImageFunction_h
29 #define itkBSplineInterpolateImageFunction_h
30 
32 #include "vnl/vnl_matrix.h"
33 
35 #include "itkConceptChecking.h"
36 #include "itkCovariantVector.h"
37 
38 #include <memory> // For unique_ptr.
39 #include <vector>
40 
41 namespace itk
42 {
82 template <typename TImageType, typename TCoordRep = double, typename TCoefficientType = double>
83 class ITK_TEMPLATE_EXPORT BSplineInterpolateImageFunction : public InterpolateImageFunction<TImageType, TCoordRep>
84 {
85 public:
86  ITK_DISALLOW_COPY_AND_MOVE(BSplineInterpolateImageFunction);
87 
93 
95  itkOverrideGetNameOfClassMacro(BSplineInterpolateImageFunction);
96 
98  itkNewMacro(Self);
99 
101  using typename Superclass::OutputType;
102 
104  using typename Superclass::InputImageType;
105 
107  static constexpr unsigned int ImageDimension = Superclass::ImageDimension;
108 
110  using typename Superclass::IndexType;
111 
113  using typename Superclass::SizeType;
114 
116  using typename Superclass::ContinuousIndexType;
117 
119  using typename Superclass::PointType;
120 
123 
125  using CoefficientDataType = TCoefficientType;
127 
131 
134 
143  OutputType
144  Evaluate(const PointType & point) const override
145  {
146  const ContinuousIndexType index =
147  this->GetInputImage()->template TransformPhysicalPointToContinuousIndex<TCoordRep>(point);
148  // No thread info passed in, so call method that doesn't need thread ID.
149  return (this->EvaluateAtContinuousIndex(index));
150  }
153  virtual OutputType
154  Evaluate(const PointType & point, ThreadIdType threadId) const
155  {
156  const ContinuousIndexType index =
157  this->GetInputImage()->template TransformPhysicalPointToContinuousIndex<TCoordRep>(point);
158  return (this->EvaluateAtContinuousIndex(index, threadId));
159  }
160 
161  OutputType
162  EvaluateAtContinuousIndex(const ContinuousIndexType & index) const override
163  {
164  // Don't know thread information, make evaluateIndex, weights on the stack.
165  // Slower, but safer.
166  vnl_matrix<long> evaluateIndex(ImageDimension, (m_SplineOrder + 1));
167  vnl_matrix<double> weights(ImageDimension, (m_SplineOrder + 1));
168 
169  // Pass evaluateIndex, weights by reference. They're only good as long
170  // as this method is in scope.
171  return this->EvaluateAtContinuousIndexInternal(index, evaluateIndex, weights);
172  }
173 
174  virtual OutputType
176  {
177  // Pass evaluateIndex, weights by reference. Different threadIDs get different instances.
178  return this->EvaluateAtContinuousIndexInternal(x, m_ThreadedEvaluateIndex[threadId], m_ThreadedWeights[threadId]);
179  }
180 
181  CovariantVectorType
183  {
184  const ContinuousIndexType index =
185  this->GetInputImage()->template TransformPhysicalPointToContinuousIndex<TCoordRep>(point);
186 
187  // No thread info passed in, so call method that doesn't need thread ID.
188  return (this->EvaluateDerivativeAtContinuousIndex(index));
189  }
190 
191  CovariantVectorType
193  {
194  const ContinuousIndexType index =
195  this->GetInputImage()->template TransformPhysicalPointToContinuousIndex<TCoordRep>(point);
196  return (this->EvaluateDerivativeAtContinuousIndex(index, threadId));
197  }
198 
199  CovariantVectorType
201  {
202  // Don't know thread information, make evaluateIndex, weights,
203  // weightsDerivative
204  // on the stack.
205  // Slower, but safer.
206  vnl_matrix<long> evaluateIndex(ImageDimension, (m_SplineOrder + 1));
207  vnl_matrix<double> weights(ImageDimension, (m_SplineOrder + 1));
208  vnl_matrix<double> weightsDerivative(ImageDimension, (m_SplineOrder + 1));
209 
210  // Pass evaluateIndex, weights, weightsDerivative by reference. They're only
211  // good
212  // as long as this method is in scope.
213  return this->EvaluateDerivativeAtContinuousIndexInternal(x, evaluateIndex, weights, weightsDerivative);
214  }
215 
216  CovariantVectorType
218  {
219  return this->EvaluateDerivativeAtContinuousIndexInternal(
220  x, m_ThreadedEvaluateIndex[threadId], m_ThreadedWeights[threadId], m_ThreadedWeightsDerivative[threadId]);
221  }
222 
223  void
225  {
226  const ContinuousIndexType index =
227  this->GetInputImage()->template TransformPhysicalPointToContinuousIndex<TCoordRep>(point);
228 
229  // No thread info passed in, so call method that doesn't need thread ID.
230  this->EvaluateValueAndDerivativeAtContinuousIndex(index, value, deriv);
231  }
232 
233  void
235  OutputType & value,
236  CovariantVectorType & deriv,
237  ThreadIdType threadId) const
238  {
239  const ContinuousIndexType index =
240  this->GetInputImage()->template TransformPhysicalPointToContinuousIndex<TCoordRep>(point);
241  this->EvaluateValueAndDerivativeAtContinuousIndex(index, value, deriv, threadId);
242  }
243 
244  void
246  OutputType & value,
247  CovariantVectorType & deriv) const
248  {
249  // Don't know thread information, make evaluateIndex, weights,
250  // weightsDerivative
251  // on the stack.
252  // Slower, but safer.
253  vnl_matrix<long> evaluateIndex(ImageDimension, (m_SplineOrder + 1));
254  vnl_matrix<double> weights(ImageDimension, (m_SplineOrder + 1));
255  vnl_matrix<double> weightsDerivative(ImageDimension, (m_SplineOrder + 1));
256 
257  // Pass evaluateIndex, weights, weightsDerivative by reference. They're only
258  // good
259  // as long as this method is in scope.
260  this->EvaluateValueAndDerivativeAtContinuousIndexInternal(
261  x, value, deriv, evaluateIndex, weights, weightsDerivative);
262  }
263 
264  void
266  OutputType & value,
267  CovariantVectorType & derivativeValue,
268  ThreadIdType threadId) const
269  {
270  this->EvaluateValueAndDerivativeAtContinuousIndexInternal(x,
271  value,
272  derivativeValue,
273  m_ThreadedEvaluateIndex[threadId],
274  m_ThreadedWeights[threadId],
275  m_ThreadedWeightsDerivative[threadId]);
276  }
277 
280  void
281  SetSplineOrder(unsigned int SplineOrder);
282 
283  itkGetConstMacro(SplineOrder, unsigned int);
284 
285  void
286  SetNumberOfWorkUnits(ThreadIdType numThreads);
287 
288  itkGetConstMacro(NumberOfWorkUnits, ThreadIdType);
289 
291  void
292  SetInputImage(const TImageType * inputData) override;
293 
304  itkSetMacro(UseImageDirection, bool);
305  itkGetConstMacro(UseImageDirection, bool);
306  itkBooleanMacro(UseImageDirection);
309  SizeType
310  GetRadius() const override
311  {
312  return SizeType::Filled(m_SplineOrder + 1);
313  }
314 
315 protected:
334  virtual OutputType
335  EvaluateAtContinuousIndexInternal(const ContinuousIndexType & x,
336  vnl_matrix<long> & evaluateIndex,
337  vnl_matrix<double> & weights) const;
338 
339  virtual void
340  EvaluateValueAndDerivativeAtContinuousIndexInternal(const ContinuousIndexType & x,
341  OutputType & value,
342  CovariantVectorType & derivativeValue,
343  vnl_matrix<long> & evaluateIndex,
344  vnl_matrix<double> & weights,
345  vnl_matrix<double> & weightsDerivative) const;
346 
347  virtual CovariantVectorType
348  EvaluateDerivativeAtContinuousIndexInternal(const ContinuousIndexType & x,
349  vnl_matrix<long> & evaluateIndex,
350  vnl_matrix<double> & weights,
351  vnl_matrix<double> & weightsDerivative) const;
352 
354  ~BSplineInterpolateImageFunction() override = default;
355  void
356  PrintSelf(std::ostream & os, Indent indent) const override;
357 
358  // These are needed by the smoothing spline routine.
359  // temp storage for processing of Coefficients
360  std::vector<CoefficientDataType> m_Scratch{};
361  // Image size
362  typename TImageType::SizeType m_DataLength{};
363  // User specified spline order (3rd or cubic is the default)
364  unsigned int m_SplineOrder{};
365 
366  // Spline coefficients
367  typename CoefficientImageType::ConstPointer m_Coefficients{};
368 
369 private:
371  void
372  SetInterpolationWeights(const ContinuousIndexType & x,
373  const vnl_matrix<long> & EvaluateIndex,
374  vnl_matrix<double> & weights,
375  unsigned int splineOrder) const;
376 
378  void
379  SetDerivativeWeights(const ContinuousIndexType & x,
380  const vnl_matrix<long> & EvaluateIndex,
381  vnl_matrix<double> & weights,
382  unsigned int splineOrder) const;
383 
386  void
387  GeneratePointsToIndex();
388 
390  void
391  DetermineRegionOfSupport(vnl_matrix<long> & evaluateIndex,
392  const ContinuousIndexType & x,
393  unsigned int splineOrder) const;
394 
397  void
398  ApplyMirrorBoundaryConditions(vnl_matrix<long> & evaluateIndex, unsigned int splineOrder) const;
399 
400  Iterator m_CIterator{}; // Iterator for
401  // traversing spline
402  // coefficients.
403  unsigned long m_MaxNumberInterpolationPoints{}; // number of
404  // neighborhood
405  // points used for
406  // interpolation
407  std::vector<IndexType> m_PointsToIndex{}; // Preallocation of
408  // interpolation
409  // neighborhood
410  // indices
411 
412  CoefficientFilterPointer m_CoefficientFilter{};
413 
414  // flag to take or not the image direction into account when computing the
415  // derivatives.
416  bool m_UseImageDirection{ true };
417 
418  ThreadIdType m_NumberOfWorkUnits{};
419  std::unique_ptr<vnl_matrix<long>[]> m_ThreadedEvaluateIndex;
420  std::unique_ptr<vnl_matrix<double>[]> m_ThreadedWeights;
421  std::unique_ptr<vnl_matrix<double>[]> m_ThreadedWeightsDerivative;
422 };
423 } // namespace itk
424 
425 #ifndef ITK_MANUAL_INSTANTIATION
426 # include "itkBSplineInterpolateImageFunction.hxx"
427 #endif
428 
429 #endif
itk::BSplineInterpolateImageFunction::m_ThreadedWeights
std::unique_ptr< vnl_matrix< double >[]> m_ThreadedWeights
Definition: itkBSplineInterpolateImageFunction.h:420
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::BSplineInterpolateImageFunction::EvaluateDerivativeAtContinuousIndex
CovariantVectorType EvaluateDerivativeAtContinuousIndex(const ContinuousIndexType &x, ThreadIdType threadId) const
Definition: itkBSplineInterpolateImageFunction.h:217
itkCovariantVector.h
itk::BSplineInterpolateImageFunction::Evaluate
virtual OutputType Evaluate(const PointType &point, ThreadIdType threadId) const
Definition: itkBSplineInterpolateImageFunction.h:154
itk::BSplineInterpolateImageFunction
Evaluates the B-Spline interpolation of an image. Spline order may be from 0 to 5.
Definition: itkBSplineInterpolateImageFunction.h:83
itk::BSplineDecompositionImageFilter
Calculates the B-Spline coefficients of an image. Spline order may be from 0 to 5.
Definition: itkBSplineDecompositionImageFilter.h:74
itk::BSplineInterpolateImageFunction::GetRadius
SizeType GetRadius() const override
Definition: itkBSplineInterpolateImageFunction.h:310
itk::ImageLinearIteratorWithIndex< TImageType >
itk::BSplineInterpolateImageFunction::EvaluateAtContinuousIndex
OutputType EvaluateAtContinuousIndex(const ContinuousIndexType &index) const override
Definition: itkBSplineInterpolateImageFunction.h:162
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itkConceptChecking.h
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::BSplineInterpolateImageFunction::EvaluateValueAndDerivativeAtContinuousIndex
void EvaluateValueAndDerivativeAtContinuousIndex(const ContinuousIndexType &x, OutputType &value, CovariantVectorType &deriv) const
Definition: itkBSplineInterpolateImageFunction.h:245
itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TImageType::PixelType >::CoefficientFilterPointer
typename CoefficientFilter::Pointer CoefficientFilterPointer
Definition: itkBSplineInterpolateImageFunction.h:130
itk::BSplineInterpolateImageFunction::m_ThreadedEvaluateIndex
std::unique_ptr< vnl_matrix< long >[]> m_ThreadedEvaluateIndex
Definition: itkBSplineInterpolateImageFunction.h:419
itk::BSplineInterpolateImageFunction::EvaluateValueAndDerivative
void EvaluateValueAndDerivative(const PointType &point, OutputType &value, CovariantVectorType &deriv) const
Definition: itkBSplineInterpolateImageFunction.h:224
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::FunctionBase< Point< TCoordRep, TImageType ::ImageDimension >, NumericTraits< TImageType ::PixelType >::RealType >::OutputType
NumericTraits< TImageType ::PixelType >::RealType OutputType
Definition: itkFunctionBase.h:62
itk::BSplineInterpolateImageFunction::EvaluateDerivativeAtContinuousIndex
CovariantVectorType EvaluateDerivativeAtContinuousIndex(const ContinuousIndexType &x) const
Definition: itkBSplineInterpolateImageFunction.h:200
itkBSplineDecompositionImageFilter.h
itk::ThreadIdType
unsigned int ThreadIdType
Definition: itkIntTypes.h:99
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::BSplineInterpolateImageFunction::EvaluateValueAndDerivative
void EvaluateValueAndDerivative(const PointType &point, OutputType &value, CovariantVectorType &deriv, ThreadIdType threadId) const
Definition: itkBSplineInterpolateImageFunction.h:234
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:55
itk::point
*par Constraints *The filter requires an image with at least two dimensions and a vector *length of at least The theory supports extension to scalar but *the implementation of the itk vector classes do not **The template parameter TRealType must be floating point(float or double) or *a user-defined "real" numerical type with arithmetic operations defined *sufficient to compute derivatives. **\par Performance *This filter will automatically multithread if run with *SetUsePrincipleComponents
itk::BSplineInterpolateImageFunction::m_ThreadedWeightsDerivative
std::unique_ptr< vnl_matrix< double >[]> m_ThreadedWeightsDerivative
Definition: itkBSplineInterpolateImageFunction.h:421
itk::BSplineInterpolateImageFunction::EvaluateDerivative
CovariantVectorType EvaluateDerivative(const PointType &point) const
Definition: itkBSplineInterpolateImageFunction.h:182
itk::BSplineInterpolateImageFunction::EvaluateValueAndDerivativeAtContinuousIndex
void EvaluateValueAndDerivativeAtContinuousIndex(const ContinuousIndexType &x, OutputType &value, CovariantVectorType &derivativeValue, ThreadIdType threadId) const
Definition: itkBSplineInterpolateImageFunction.h:265
itk::CovariantVector
A templated class holding a n-Dimensional covariant vector.
Definition: itkCovariantVector.h:70
itk::BSplineInterpolateImageFunction::EvaluateDerivative
CovariantVectorType EvaluateDerivative(const PointType &point, ThreadIdType threadId) const
Definition: itkBSplineInterpolateImageFunction.h:192
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ContinuousIndex< TCoordRep, Self::ImageDimension >
itk::Point
A templated class holding a geometric point in n-Dimensional space.
Definition: itkPoint.h:53
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
itk::BSplineInterpolateImageFunction::Evaluate
OutputType Evaluate(const PointType &point) const override
Definition: itkBSplineInterpolateImageFunction.h:144
itkInterpolateImageFunction.h
itk::BSplineInterpolateImageFunction::EvaluateAtContinuousIndex
virtual OutputType EvaluateAtContinuousIndex(const ContinuousIndexType &x, ThreadIdType threadId) const
Definition: itkBSplineInterpolateImageFunction.h:175
itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TImageType::PixelType >::CoefficientDataType
TImageType::PixelType CoefficientDataType
Definition: itkBSplineInterpolateImageFunction.h:125
itk::InterpolateImageFunction
Base class for all image interpolators.
Definition: itkInterpolateImageFunction.h:45