OTB  9.0.0
Orfeo Toolbox
otbBSplineInterpolateImageFunction.hxx
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2005-2022 Centre National d'Etudes Spatiales (CNES)
3  *
4  * This file is part of Orfeo Toolbox
5  *
6  * https://www.orfeo-toolbox.org/
7  *
8  * Licensed under the Apache License, Version 2.0 (the "License");
9  * you may not use this file except in compliance with the License.
10  * You may obtain a copy of the License at
11  *
12  * http://www.apache.org/licenses/LICENSE-2.0
13  *
14  * Unless required by applicable law or agreed to in writing, software
15  * distributed under the License is distributed on an "AS IS" BASIS,
16  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
17  * See the License for the specific language governing permissions and
18  * limitations under the License.
19  */
20 
21 #ifndef otbBSplineInterpolateImageFunction_hxx
22 #define otbBSplineInterpolateImageFunction_hxx
24 #include "itkImageLinearIteratorWithIndex.h"
25 #include "itkImageRegionConstIteratorWithIndex.h"
26 #include "itkImageRegionIterator.h"
27 
28 #include "itkVector.h"
29 
30 #include "itkMatrix.h"
31 
32 namespace otb
33 {
34 
38 template <class TImageType, class TCoordRep, class TCoefficientType>
40  :
41  m_SplineOrder(0),
42  m_Coefficients(CoefficientImageType::New()), // TODO: Should we store coefficients in a variable or retrieve from filter?
43  m_CoefficientFilter(CoefficientFilter::New())
44 {
45  this->SetSplineOrder(3);
46 }
47 
51 template <class TImageType, class TCoordRep, class TCoefficientType>
53 {
54  Superclass::PrintSelf(os, indent);
55  os << indent << "Spline Order: " << m_SplineOrder << std::endl;
56 }
58 
59 template <class TImageType, class TCoordRep, class TCoefficientType>
61 {
62  m_CoefficientFilter->GetOutput()->UpdateOutputInformation();
63  m_CoefficientFilter->GetOutput()->SetRequestedRegion(m_CoefficientFilter->GetInput()->GetBufferedRegion());
64  m_CoefficientFilter->GetOutput()->PropagateRequestedRegion();
65  m_CoefficientFilter->GetOutput()->UpdateOutputData();
66  m_Coefficients = m_CoefficientFilter->GetOutput();
67  m_CurrentBufferedRegion = m_CoefficientFilter->GetInput()->GetBufferedRegion();
68 }
69 template <class TImageType, class TCoordRep, class TCoefficientType>
71 {
72  if (inputData)
73  {
74  m_CoefficientFilter->SetInput(inputData);
75 
76  // the Coefficient Filter requires that the spline order and the input data be set.
77  // TODO: We need to ensure that this is only run once and only after both input and
78  // spline order have been set. Should we force an update after the
79  // splineOrder has been set also?
80 
81  UpdateCoefficientsFilter();
82 
83  // Call the Superclass implementation after, in case the filter
84  // pulls in more of the input image
85  Superclass::SetInputImage(inputData);
86 
87  m_DataLength = inputData->GetBufferedRegion().GetSize();
88  }
89  else
90  {
91  m_Coefficients = nullptr;
92  }
93 }
94 
95 template <class TImageType, class TCoordRep, class TCoefficientType>
97 {
98  if (SplineOrder == m_SplineOrder)
99  {
100  return;
101  }
102  m_SplineOrder = SplineOrder;
103  m_CoefficientFilter->SetSplineOrder(SplineOrder);
104 
105  // this->SetPoles();
106  m_MaxNumberInterpolationPoints = 1;
107  for (unsigned int n = 0; n < ImageDimension; ++n)
108  {
109  m_MaxNumberInterpolationPoints *= (m_SplineOrder + 1);
110  }
111  this->GeneratePointsToIndex();
112 }
113 
114 template <class TImageType, class TCoordRep, class TCoefficientType>
117 {
118  // UpdateCoefficientsFilter();
119  vnl_matrix<long> EvaluateIndex(ImageDimension, (m_SplineOrder + 1));
120 
121  // compute the interpolation indexes
122  this->DetermineRegionOfSupport(EvaluateIndex, x, m_SplineOrder);
123 
124  // Determine weights
125  vnl_matrix<double> weights(ImageDimension, (m_SplineOrder + 1));
126  SetInterpolationWeights(x, EvaluateIndex, weights, m_SplineOrder);
127 
128  // std::cout<<"EvaluateIndex: "<<std::endl;
129  // std::cout<<EvaluateIndex[0][0]<<"\t"<<EvaluateIndex[0][1]<<"\t"
130  // <<EvaluateIndex[0][2]<<"\t"<<EvaluateIndex[0][3]<<std::endl;
131  // std::cout<<EvaluateIndex[1][0]<<"\t"<<EvaluateIndex[1][1]<<"\t"
132  // <<EvaluateIndex[1][2]<<"\t"<<EvaluateIndex[1][3]<<std::endl;
133 
134  // Modify EvaluateIndex at the boundaries using mirror boundary conditions
135  this->ApplyMirrorBoundaryConditions(EvaluateIndex, m_SplineOrder);
136 
137  // perform interpolation
138  double interpolated = 0.0;
139  IndexType coefficientIndex;
140  // Step through eachpoint in the N-dimensional interpolation cube.
141  for (unsigned int p = 0; p < m_MaxNumberInterpolationPoints; ++p)
142  {
143  // translate each step into the N-dimensional index.
144  // IndexType pointIndex = PointToIndex( p );
145 
146  double w = 1.0;
147  for (unsigned int n = 0; n < ImageDimension; ++n)
148  {
149  w *= weights[n][m_PointsToIndex[p][n]];
150  coefficientIndex[n] = EvaluateIndex[n][m_PointsToIndex[p][n]]; // Build up ND index for coefficients.
151  // std::cout<<"From inside: "<<n<<" "<<p<<" "<<m_PointsToIndex[p][n]<<" "<< EvaluateIndex[n][m_PointsToIndex[p][n]]<<std::endl;
152  }
153  // std::cout<<"CoefficientIndex: "<<coefficientIndex<<std::endl;
154  // Convert our step p to the appropriate point in ND space in the
155  // m_Coefficients cube.
156  interpolated += w * m_Coefficients->GetPixel(coefficientIndex);
157  }
158 
159  /* double interpolated = 0.0;
160  IndexType coefficientIndex;
161  // Step through eachpoint in the N-dimensional interpolation cube.
162  for (unsigned int sp = 0; sp <= m_SplineOrder; ++sp)
163  {
164  for (unsigned int sp1=0; sp1 <= m_SplineOrder; sp1++)
165  {
166 
167  double w = 1.0;
168  for (unsigned int n1 = 0; n1 < ImageDimension; n1++ )
169  {
170  w *= weights[n1][ sp1 ];
171  coefficientIndex[n1] = EvaluateIndex[n1][sp]; // Build up ND index for coefficients.
172  }
173 
174  interpolated += w * m_Coefficients->GetPixel(coefficientIndex);
175  }
176  }
177  */
178  return (interpolated);
179 }
180 
181 template <class TImageType, class TCoordRep, class TCoefficientType>
184 {
185  UpdateCoefficientsFilter();
186  vnl_matrix<long> EvaluateIndex(ImageDimension, (m_SplineOrder + 1));
187 
188  // compute the interpolation indexes
189  // TODO: Do we need to revisit region of support for the derivatives?
190  this->DetermineRegionOfSupport(EvaluateIndex, x, m_SplineOrder);
191 
192  // Determine weights
193  vnl_matrix<double> weights(ImageDimension, (m_SplineOrder + 1));
194  SetInterpolationWeights(x, EvaluateIndex, weights, m_SplineOrder);
195 
196  vnl_matrix<double> weightsDerivative(ImageDimension, (m_SplineOrder + 1));
197  SetDerivativeWeights(x, EvaluateIndex, weightsDerivative, (m_SplineOrder));
198 
199  // Modify EvaluateIndex at the boundaries using mirror boundary conditions
200  this->ApplyMirrorBoundaryConditions(EvaluateIndex, m_SplineOrder);
201 
202  // Calculate derivative
203  CovariantVectorType derivativeValue;
204  double tempValue;
205  IndexType coefficientIndex;
206  for (unsigned int n = 0; n < ImageDimension; ++n)
207  {
208  derivativeValue[n] = 0.0;
209  for (unsigned int p = 0; p < m_MaxNumberInterpolationPoints; ++p)
210  {
211  tempValue = 1.0;
212  for (unsigned int n1 = 0; n1 < ImageDimension; n1++)
213  {
214  // coefficientIndex[n1] = EvaluateIndex[n1][sp];
215  coefficientIndex[n1] = EvaluateIndex[n1][m_PointsToIndex[p][n1]];
216 
217  if (n1 == n)
218  {
219  // w *= weights[n][ m_PointsToIndex[p][n] ];
220  tempValue *= weightsDerivative[n1][m_PointsToIndex[p][n1]];
221  }
222  else
223  {
224  tempValue *= weights[n1][m_PointsToIndex[p][n1]];
225  }
226  }
227  derivativeValue[n] += m_Coefficients->GetPixel(coefficientIndex) * tempValue;
228  }
229  derivativeValue[n] /= this->GetInputImage()->GetSignedSpacing()[n]; // take spacing into account
230  }
231 
232  return (derivativeValue);
233 }
234 
235 template <class TImageType, class TCoordRep, class TCoefficientType>
237  const vnl_matrix<long>& EvaluateIndex,
238  vnl_matrix<double>& weights,
239  unsigned int splineOrder) const
240 {
241  // For speed improvements we could make each case a separate function and use
242  // function pointers to reference the correct weight order.
243  // Left as is for now for readability.
244  double w, w2, w4, t, t0, t1;
245 
246  switch (splineOrder)
247  {
248  case 3:
249  for (unsigned int n = 0; n < ImageDimension; ++n)
250  {
251  w = x[n] - (double)EvaluateIndex[n][1];
252  weights[n][3] = (1.0 / 6.0) * w * w * w;
253  weights[n][0] = (1.0 / 6.0) + 0.5 * w * (w - 1.0) - weights[n][3];
254  weights[n][2] = w + weights[n][0] - 2.0 * weights[n][3];
255  weights[n][1] = 1.0 - weights[n][0] - weights[n][2] - weights[n][3];
256  }
257  break;
258  case 0:
259  for (unsigned int n = 0; n < ImageDimension; ++n)
260  {
261  weights[n][0] = 1; // implements nearest neighbor
262  }
263  break;
264  case 1:
265  for (unsigned int n = 0; n < ImageDimension; ++n)
266  {
267  w = x[n] - (double)EvaluateIndex[n][0];
268  weights[n][1] = w;
269  weights[n][0] = 1.0 - w;
270  }
271  break;
272  case 2:
273  for (unsigned int n = 0; n < ImageDimension; ++n)
274  {
275  /* x */
276  w = x[n] - (double)EvaluateIndex[n][1];
277  weights[n][1] = 0.75 - w * w;
278  weights[n][2] = 0.5 * (w - weights[n][1] + 1.0);
279  weights[n][0] = 1.0 - weights[n][1] - weights[n][2];
280  }
281  break;
282  case 4:
283  for (unsigned int n = 0; n < ImageDimension; ++n)
284  {
285  /* x */
286  w = x[n] - (double)EvaluateIndex[n][2];
287  w2 = w * w;
288  t = (1.0 / 6.0) * w2;
289  weights[n][0] = 0.5 - w;
290  weights[n][0] *= weights[n][0];
291  weights[n][0] *= (1.0 / 24.0) * weights[n][0];
292  t0 = w * (t - 11.0 / 24.0);
293  t1 = 19.0 / 96.0 + w2 * (0.25 - t);
294  weights[n][1] = t1 + t0;
295  weights[n][3] = t1 - t0;
296  weights[n][4] = weights[n][0] + t0 + 0.5 * w;
297  weights[n][2] = 1.0 - weights[n][0] - weights[n][1] - weights[n][3] - weights[n][4];
298  }
299  break;
300  case 5:
301  for (unsigned int n = 0; n < ImageDimension; ++n)
302  {
303  /* x */
304  w = x[n] - (double)EvaluateIndex[n][2];
305  w2 = w * w;
306  weights[n][5] = (1.0 / 120.0) * w * w2 * w2;
307  w2 -= w;
308  w4 = w2 * w2;
309  w -= 0.5;
310  t = w2 * (w2 - 3.0);
311  weights[n][0] = (1.0 / 24.0) * (1.0 / 5.0 + w2 + w4) - weights[n][5];
312  t0 = (1.0 / 24.0) * (w2 * (w2 - 5.0) + 46.0 / 5.0);
313  t1 = (-1.0 / 12.0) * w * (t + 4.0);
314  weights[n][2] = t0 + t1;
315  weights[n][3] = t0 - t1;
316  t0 = (1.0 / 16.0) * (9.0 / 5.0 - t);
317  t1 = (1.0 / 24.0) * w * (w4 - w2 - 5.0);
318  weights[n][1] = t0 + t1;
319  weights[n][4] = t0 - t1;
320  }
321  break;
322  default:
323  // SplineOrder not implemented yet.
324  itk::ExceptionObject err(__FILE__, __LINE__);
325  err.SetLocation(ITK_LOCATION);
326  err.SetDescription("SplineOrder must be between 0 and 5. Requested spline order has not been implemented yet.");
327  throw err;
328  break;
329  }
330 }
331 
332 template <class TImageType, class TCoordRep, class TCoefficientType>
334  const vnl_matrix<long>& EvaluateIndex,
335  vnl_matrix<double>& weights, unsigned int splineOrder) const
336 {
337  // For speed improvements we could make each case a separate function and use
338  // function pointers to reference the correct weight order.
339  // Another possibility would be to loop inside the case statement (reducing the number
340  // of switch statement executions to one per routine call.
341  // Left as is for now for readability.
342  double w, w1, w2, w3, w4, w5, t, t0, t1, t2;
343  int derivativeSplineOrder = (int)splineOrder - 1;
344 
345  switch (derivativeSplineOrder)
346  {
347 
348  // Calculates B(splineOrder) ( (x + 1/2) - xi) - B(splineOrder -1) ( (x - 1/2) - xi)
349  case -1:
350  // Why would we want to do this?
351  for (unsigned int n = 0; n < ImageDimension; ++n)
352  {
353  weights[n][0] = 0.0;
354  }
355  break;
356  case 0:
357  for (unsigned int n = 0; n < ImageDimension; ++n)
358  {
359  weights[n][0] = -1.0;
360  weights[n][1] = 1.0;
361  }
362  break;
363  case 1:
364  for (unsigned int n = 0; n < ImageDimension; ++n)
365  {
366  w = x[n] + 0.5 - (double)EvaluateIndex[n][1];
367  // w2 = w;
368  w1 = 1.0 - w;
369 
370  weights[n][0] = 0.0 - w1;
371  weights[n][1] = w1 - w;
372  weights[n][2] = w;
373  }
374  break;
375  case 2:
376 
377  for (unsigned int n = 0; n < ImageDimension; ++n)
378  {
379  w = x[n] + .5 - (double)EvaluateIndex[n][2];
380  w2 = 0.75 - w * w;
381  w3 = 0.5 * (w - w2 + 1.0);
382  w1 = 1.0 - w2 - w3;
383 
384  weights[n][0] = 0.0 - w1;
385  weights[n][1] = w1 - w2;
386  weights[n][2] = w2 - w3;
387  weights[n][3] = w3;
388  }
389  break;
390  case 3:
391 
392  for (unsigned int n = 0; n < ImageDimension; ++n)
393  {
394  w = x[n] + 0.5 - (double)EvaluateIndex[n][2];
395  w4 = (1.0 / 6.0) * w * w * w;
396  w1 = (1.0 / 6.0) + 0.5 * w * (w - 1.0) - w4;
397  w3 = w + w1 - 2.0 * w4;
398  w2 = 1.0 - w1 - w3 - w4;
399 
400  weights[n][0] = 0.0 - w1;
401  weights[n][1] = w1 - w2;
402  weights[n][2] = w2 - w3;
403  weights[n][3] = w3 - w4;
404  weights[n][4] = w4;
405  }
406  break;
407  case 4:
408  for (unsigned int n = 0; n < ImageDimension; ++n)
409  {
410  w = x[n] + .5 - (double)EvaluateIndex[n][3];
411  t2 = w * w;
412  t = (1.0 / 6.0) * t2;
413  w1 = 0.5 - w;
414  w1 *= w1;
415  w1 *= (1.0 / 24.0) * w1;
416  t0 = w * (t - 11.0 / 24.0);
417  t1 = 19.0 / 96.0 + t2 * (0.25 - t);
418  w2 = t1 + t0;
419  w4 = t1 - t0;
420  w5 = w1 + t0 + 0.5 * w;
421  w3 = 1.0 - w1 - w2 - w4 - w5;
422 
423  weights[n][0] = 0.0 - w1;
424  weights[n][1] = w1 - w2;
425  weights[n][2] = w2 - w3;
426  weights[n][3] = w3 - w4;
427  weights[n][4] = w4 - w5;
428  weights[n][5] = w5;
429  }
430  break;
431 
432  default:
433  // SplineOrder not implemented yet.
434  itk::ExceptionObject err(__FILE__, __LINE__);
435  err.SetLocation(ITK_LOCATION);
436  err.SetDescription("SplineOrder (for derivatives) must be between 1 and 5. Requested spline order has not been implemented yet.");
437  throw err;
438  break;
439  }
440 }
441 
442 // Generates m_PointsToIndex;
443 template <class TImageType, class TCoordRep, class TCoefficientType>
445 {
446  // m_PointsToIndex is used to convert a sequential location to an N-dimension
447  // index vector. This is precomputed to save time during the interpolation routine.
448  m_PointsToIndex.resize(m_MaxNumberInterpolationPoints);
449  for (unsigned int p = 0; p < m_MaxNumberInterpolationPoints; ++p)
450  {
451  int pp = p;
452  unsigned long indexFactor[ImageDimension];
453  indexFactor[0] = 1;
454  for (int j = 1; j < static_cast<int>(ImageDimension); ++j)
455  {
456  indexFactor[j] = indexFactor[j - 1] * (m_SplineOrder + 1);
457  }
458  for (int j = (static_cast<int>(ImageDimension) - 1); j >= 0; j--)
459  {
460  m_PointsToIndex[p][j] = pp / indexFactor[j];
461  pp = pp % indexFactor[j];
462  }
463  }
464 }
465 
466 template <class TImageType, class TCoordRep, class TCoefficientType>
468  const ContinuousIndexType& x,
469  unsigned int splineOrder) const
470 {
471  long indx;
472 
473  // compute the interpolation indexes
474  for (unsigned int n = 0; n < ImageDimension; ++n)
475  {
476  if (splineOrder & 1) // Use this index calculation for odd splineOrder
477  {
478  indx = (long)std::floor((float)x[n]) - splineOrder / 2;
479  // std::cout<<"x: "<<x<<std::endl;
480  // std::cout<<"splineOrder: "<<splineOrder<<std::endl;
481  // std::cout<<"indx: "<<indx<<std::endl;
482  for (unsigned int k = 0; k <= splineOrder; ++k)
483  {
484  evaluateIndex[n][k] = indx++;
485  }
486  }
487  else // Use this index calculation for even splineOrder
488  {
489 
490  indx = (long)std::floor((float)(x[n] + 0.5)) - splineOrder / 2;
491  // std::cout<<"x: "<<x<<std::endl;
492  // std::cout<<"splineOrder: "<<splineOrder<<std::endl;
493  // std::cout<<"indx: "<<indx<<std::endl;
494  for (unsigned int k = 0; k <= splineOrder; ++k)
495  {
496  evaluateIndex[n][k] = indx++;
497  }
498  }
499  }
500 }
501 
502 template <class TImageType, class TCoordRep, class TCoefficientType>
504  unsigned int splineOrder) const
505 {
506 
507  for (unsigned int n = 0; n < ImageDimension; ++n)
508  {
509  long dataLength = m_DataLength[n];
510  long dataOffset = m_CurrentBufferedRegion.GetIndex()[n];
511 
512  // apply the mirror boundary conditions
513  // TODO: We could implement other boundary options beside mirror
514  if (m_DataLength[n] == 1)
515  {
516  for (unsigned int k = 0; k <= splineOrder; ++k)
517  {
518  evaluateIndex[n][k] = 0;
519  }
520  }
521  else
522  {
523  for (unsigned int k = 0; k <= splineOrder; ++k)
524  {
525  // btw - Think about this couldn't this be replaced with a more elagent modulus method?
526  evaluateIndex[n][k] = (evaluateIndex[n][k] < dataOffset) ? (dataOffset + (dataOffset - evaluateIndex[n][k]) % dataLength) : (evaluateIndex[n][k]);
527  if ((long)dataLength + dataOffset <= evaluateIndex[n][k])
528  {
529  evaluateIndex[n][k] = dataOffset + dataLength - (evaluateIndex[n][k] - dataOffset - dataLength) % dataLength;
530  }
531  }
532  }
533  }
534 }
535 
536 } // namespace otb
537 
538 #endif
otbBSplineInterpolateImageFunction.h
otb::BSplineInterpolateImageFunction::ApplyMirrorBoundaryConditions
void ApplyMirrorBoundaryConditions(vnl_matrix< long > &evaluateIndex, unsigned int splineOrder) const
Definition: otbBSplineInterpolateImageFunction.hxx:503
otb::BSplineInterpolateImageFunction::GeneratePointsToIndex
void GeneratePointsToIndex()
Definition: otbBSplineInterpolateImageFunction.hxx:444
otb::BSplineInterpolateImageFunction::CoefficientImageType
itk::Image< CoefficientDataType, itkGetStaticConstMacro(ImageDimension)> CoefficientImageType
Definition: otbBSplineInterpolateImageFunction.h:89
otb
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
Definition: otbJoinContainer.h:32
otb::BSplineDecompositionImageFilter
This class is an evolution of the itk::BSplineDecompositionImageFilter to handle huge images with thi...
Definition: otbBSplineDecompositionImageFilter.h:43
otb::BSplineInterpolateImageFunction::SetInputImage
void SetInputImage(const TImageType *inputData) override
Definition: otbBSplineInterpolateImageFunction.hxx:70
otb::BSplineInterpolateImageFunction::SetSplineOrder
void SetSplineOrder(unsigned int SplineOrder)
Definition: otbBSplineInterpolateImageFunction.hxx:96
otb::BSplineInterpolateImageFunction::PrintSelf
void PrintSelf(std::ostream &os, itk::Indent indent) const override
Definition: otbBSplineInterpolateImageFunction.hxx:52
otb::BSplineInterpolateImageFunction::UpdateCoefficientsFilter
virtual void UpdateCoefficientsFilter(void)
Definition: otbBSplineInterpolateImageFunction.hxx:60
otb::BSplineInterpolateImageFunction::DetermineRegionOfSupport
void DetermineRegionOfSupport(vnl_matrix< long > &evaluateIndex, const ContinuousIndexType &x, unsigned int splineOrder) const
Definition: otbBSplineInterpolateImageFunction.hxx:467
otb::BSplineInterpolateImageFunction::SetInterpolationWeights
void SetInterpolationWeights(const ContinuousIndexType &x, const vnl_matrix< long > &EvaluateIndex, vnl_matrix< double > &weights, unsigned int splineOrder) const
Definition: otbBSplineInterpolateImageFunction.hxx:236
otb::BSplineInterpolateImageFunction::IndexType
Superclass::IndexType IndexType
Definition: otbBSplineInterpolateImageFunction.h:73
otb::BSplineInterpolateImageFunction::EvaluateAtContinuousIndex
OutputType EvaluateAtContinuousIndex(const ContinuousIndexType &index) const override
Definition: otbBSplineInterpolateImageFunction.hxx:116
otb::BSplineInterpolateImageFunction::SetDerivativeWeights
void SetDerivativeWeights(const ContinuousIndexType &x, const vnl_matrix< long > &EvaluateIndex, vnl_matrix< double > &weights, unsigned int splineOrder) const
Definition: otbBSplineInterpolateImageFunction.hxx:333
otb::BSplineInterpolateImageFunction::EvaluateDerivativeAtContinuousIndex
CovariantVectorType EvaluateDerivativeAtContinuousIndex(const ContinuousIndexType &x) const
Definition: otbBSplineInterpolateImageFunction.hxx:183
otb::BSplineInterpolateImageFunction::OutputType
Superclass::OutputType OutputType
Definition: otbBSplineInterpolateImageFunction.h:61
otb::BSplineInterpolateImageFunction::CovariantVectorType
itk::CovariantVector< OutputType, itkGetStaticConstMacro(ImageDimension)> CovariantVectorType
Definition: otbBSplineInterpolateImageFunction.h:106
otb::BSplineInterpolateImageFunction::ContinuousIndexType
Superclass::ContinuousIndexType ContinuousIndexType
Definition: otbBSplineInterpolateImageFunction.h:79
otb::BSplineInterpolateImageFunction::BSplineInterpolateImageFunction
BSplineInterpolateImageFunction()
Definition: otbBSplineInterpolateImageFunction.hxx:39