PCL
BicubicInterpolation.h
Go to the documentation of this file.
1 // ____ ______ __
2 // / __ \ / ____// /
3 // / /_/ // / / /
4 // / ____// /___ / /___ PixInsight Class Library
5 // /_/ \____//_____/ PCL 2.7.0
6 // ----------------------------------------------------------------------------
7 // pcl/BicubicInterpolation.h - Released 2024-06-18T15:48:54Z
8 // ----------------------------------------------------------------------------
9 // This file is part of the PixInsight Class Library (PCL).
10 // PCL is a multiplatform C++ framework for development of PixInsight modules.
11 //
12 // Copyright (c) 2003-2024 Pleiades Astrophoto S.L. All Rights Reserved.
13 //
14 // Redistribution and use in both source and binary forms, with or without
15 // modification, is permitted provided that the following conditions are met:
16 //
17 // 1. All redistributions of source code must retain the above copyright
18 // notice, this list of conditions and the following disclaimer.
19 //
20 // 2. All redistributions in binary form must reproduce the above copyright
21 // notice, this list of conditions and the following disclaimer in the
22 // documentation and/or other materials provided with the distribution.
23 //
24 // 3. Neither the names "PixInsight" and "Pleiades Astrophoto", nor the names
25 // of their contributors, may be used to endorse or promote products derived
26 // from this software without specific prior written permission. For written
27 // permission, please contact info@pixinsight.com.
28 //
29 // 4. All products derived from this software, in any form whatsoever, must
30 // reproduce the following acknowledgment in the end-user documentation
31 // and/or other materials provided with the product:
32 //
33 // "This product is based on software from the PixInsight project, developed
34 // by Pleiades Astrophoto and its contributors (https://pixinsight.com/)."
35 //
36 // Alternatively, if that is where third-party acknowledgments normally
37 // appear, this acknowledgment must be reproduced in the product itself.
38 //
39 // THIS SOFTWARE IS PROVIDED BY PLEIADES ASTROPHOTO AND ITS CONTRIBUTORS
40 // "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED
41 // TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
42 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL PLEIADES ASTROPHOTO OR ITS
43 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
44 // EXEMPLARY OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, BUSINESS
45 // INTERRUPTION; PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; AND LOSS OF USE,
46 // DATA OR PROFITS) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
47 // CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
48 // ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
49 // POSSIBILITY OF SUCH DAMAGE.
50 // ----------------------------------------------------------------------------
51 
52 #ifndef __PCL_BicubicInterpolation_h
53 #define __PCL_BicubicInterpolation_h
54 
56 
57 #include <pcl/Defs.h>
58 #include <pcl/Diagnostics.h>
59 
61 #include <pcl/Math.h>
62 #include <pcl/Utility.h>
63 
64 namespace pcl
65 {
66 
67 #define m_width this->m_width
68 #define m_height this->m_height
69 #define m_fillBorder this->m_fillBorder
70 #define m_fillValue this->m_fillValue
71 #define m_data this->m_data
72 
73 // ----------------------------------------------------------------------------
74 
84 template <typename T>
86 {
87 public:
88 
93 
98 
99 protected:
100 
101  void InitXY( int& i1, int& j1, double* p0, double* p1, double* p2, double* p3, double x, double y ) const noexcept
102  {
103  PCL_PRECONDITION( int( x ) >= 0 )
104  PCL_PRECONDITION( int( x ) < m_width )
105  PCL_PRECONDITION( int( y ) >= 0 )
106  PCL_PRECONDITION( int( y ) < m_height )
107 
108  // Central grid coordinates
109  i1 = pcl::Range( TruncInt( y ), 0, m_height-1 );
110  j1 = pcl::Range( TruncInt( x ), 0, m_width-1 );
111 
112  // Set up source matrix
113 
114  int j0 = j1 - 1;
115  int i0 = i1 - 1;
116 
117  int j2 = j1 + 1;
118  int i2 = i1 + 1;
119 
120  int j3 = j1 + 2;
121  int i3 = i1 + 2;
122 
123  const T* fp = m_data + (int64( i0 )*m_width + j0);
124 
125  // Row 0
126 
127  if ( i0 < 0 )
128  {
129  fp += m_width;
130 
131  if ( m_fillBorder )
132  {
133  p0[0] = p0[1] = p0[2] = p0[3] = m_fillValue;
134  goto __row1;
135  }
136  }
137 
138  GetRow( p0, fp, j0, j2, j3 );
139 
140  if ( i0 >= 0 )
141  fp += m_width;
142 
143 __row1: // Row 1
144 
145  GetRow( p1, fp, j0, j2, j3 );
146 
147  // Row 2
148 
149  if ( i2 < m_height )
150  fp += m_width;
151 
152  GetRow( p2, fp, j0, j2, j3 );
153 
154  // Row 3
155 
156  if ( i3 < m_height )
157  fp += m_width;
158  else if ( m_fillBorder )
159  {
160  p3[0] = p3[1] = p3[2] = p3[3] = m_fillValue;
161  goto __done;
162  }
163 
164  GetRow( p3, fp, j0, j2, j3 );
165 
166 __done:
167  ;
168  }
169 
170  void InitYX( int& i1, int& j1, double* p0, double* p1, double* p2, double* p3, double x, double y ) const noexcept
171  {
172  PCL_PRECONDITION( int( x ) >= 0 )
173  PCL_PRECONDITION( int( x ) < m_width )
174  PCL_PRECONDITION( int( y ) >= 0 )
175  PCL_PRECONDITION( int( y ) < m_height )
176 
177  // Central grid coordinates
178  i1 = pcl::Range( TruncInt( y ), 0, m_height-1 );
179  j1 = pcl::Range( TruncInt( x ), 0, m_width-1 );
180 
181  // Set up source matrix
182 
183  int j0 = j1 - 1;
184  int i0 = i1 - 1;
185 
186  int j2 = j1 + 1;
187  int i2 = i1 + 1;
188 
189  int j3 = j1 + 2;
190  int i3 = i1 + 2;
191 
192  const T* fp = m_data + (int64( i0 )*m_width + j0);
193 
194  // Column 0
195 
196  if ( j0 < 0 )
197  {
198  ++fp;
199 
200  if ( m_fillBorder )
201  {
202  p0[0] = p0[1] = p0[2] = p0[3] = m_fillValue;
203  goto __col1;
204  }
205  }
206 
207  GetColumn( p0, fp, i0, i2, i3 );
208 
209  if ( j0 >= 0 )
210  ++fp;
211 
212 __col1: // Column 1
213 
214  GetColumn( p1, fp, i0, i2, i3 );
215 
216  // Column 2
217 
218  if ( j2 < m_width )
219  ++fp;
220 
221  GetColumn( p2, fp, i0, i2, i3 );
222 
223  // Column 3
224 
225  if ( j3 < m_width )
226  ++fp;
227  else if ( m_fillBorder )
228  {
229  p3[0] = p3[1] = p3[2] = p3[3] = m_fillValue;
230  goto __done;
231  }
232 
233  GetColumn( p3, fp, i0, i2, i3 );
234 
235 __done:
236  ;
237  }
238 
239 private:
240 
241  void GetRow( double* p, const T* fp, int j0, int j2, int j3 ) const noexcept
242  {
243  if ( m_fillBorder )
244  {
245  *p = (j0 >= 0) ? *fp : m_fillValue;
246  *++p = *++fp;
247 
248  if ( j2 < m_width )
249  {
250  *++p = *++fp;
251  *++p = (j3 < m_width) ? *++fp : m_fillValue;
252  }
253  else
254  p[1] = p[2] = m_fillValue;
255  }
256  else
257  {
258  if ( j0 < 0 )
259  ++fp;
260  *p = *fp;
261  if ( j0 >= 0 )
262  ++fp;
263  *++p = *fp;
264 
265  if ( j2 < m_width )
266  {
267  *++p = *++fp;
268  if ( j3 < m_width )
269  ++fp;
270  *++p = *fp;
271  }
272  else
273  {
274  *++p = *fp;
275  *++p = *(fp - 1);
276  }
277  }
278  }
279 
280  void GetColumn( double* p, const T* fp, int i0, int i2, int i3 ) const noexcept
281  {
282  if ( m_fillBorder )
283  {
284  *p = (i0 >= 0) ? *fp : m_fillValue;
285  *++p = *(fp += m_width);
286 
287  if ( i2 < m_height )
288  {
289  *++p = *(fp += m_width);
290  *++p = (i3 < m_height) ? *(fp += m_width) : m_fillValue;
291  }
292  else
293  p[1] = p[2] = m_fillValue;
294  }
295  else
296  {
297  if ( i0 < 0 )
298  fp += m_width;
299  *p = *fp;
300  if ( i0 >= 0 )
301  fp += m_width;
302  *++p = *fp;
303 
304  if ( i2 < m_height )
305  {
306  *++p = *(fp += m_width);
307  if ( i3 < m_height )
308  fp += m_width;
309  *++p = *fp;
310  }
311  else
312  {
313  *++p = *fp;
314  *++p = *(fp - m_width);
315  }
316  }
317  }
318 };
319 
320 // ----------------------------------------------------------------------------
321 
322 #define InitXY this->InitXY
323 #define InitYX this->InitYX
324 
325 /*
326  * Undefine the following to use any free constant value a != -1/2
327  */
328 #define __PCL_BICUBIC_SPLINE_A_IS_MINUS_ONE_HALF 1
329 
330 /*
331  * The "a" constant controls the depth of the negative lobes in the bicubic
332  * spline interpolation function, -1 <= a < 0. Larger values of a (in absolute
333  * value) lead to more ringing (sharpness). The default is -1/2, as recommended
334  * in [Keys 1981].
335  */
336 #define __PCL_BICUBIC_SPLINE_A -0.5
337 
338 /*
339  * Default clamping threshold for linear images. This value has been optimized
340  * for a = -1/2. If you use another value for a, this must also be fine tuned.
341  */
342 #ifndef __PCL_BICUBIC_SPLINE_CLAMPING_THRESHOLD
343 #define __PCL_BICUBIC_SPLINE_CLAMPING_THRESHOLD 0.3F
344 #endif
345 
376 template <typename T>
378 {
379 public:
380 
390  BicubicSplineInterpolation( float clamp = __PCL_BICUBIC_SPLINE_CLAMPING_THRESHOLD )
391  : m_clamp( Range( clamp, 0.0F, 1.0F ) )
392  {
393  PCL_PRECONDITION( 0 <= clamp && clamp <= 1 )
394  }
395 
400 
405  {
406  }
407 
417  double operator()( double x, double y ) const override
418  {
419  PCL_PRECONDITION( m_data != nullptr )
420  PCL_PRECONDITION( m_width > 0 && m_height > 0 )
421 
422  // Initialize grid coordinates and source matrix
423  int i1, j1;
424  double p0[ 4 ], p1[ 4 ], p2[ 4 ], p3[ 4 ];
425  InitXY( i1, j1, p0, p1, p2, p3, x, y );
426 
427  // Cubic spline coefficients
428  double C[ 4 ];
429  GetSplineCoefficients( C, x-j1 );
430 
431  // Interpolate neighbor rows
432  double c[ 4 ];
433  c[0] = Interpolate( p0, C );
434  c[1] = Interpolate( p1, C );
435  c[2] = Interpolate( p2, C );
436  c[3] = Interpolate( p3, C );
437 
438  // Interpolate result vertically
439  GetSplineCoefficients( C, y-i1 );
440  return Interpolate( c, C );
441  }
442 
457  void InterpolateX( double fx[], double x, double y ) const noexcept
458  {
459  PCL_PRECONDITION( m_data != nullptr )
460  PCL_PRECONDITION( m_width > 0 && m_height > 0 )
461 
462  // Initialize grid coordinates and source matrix
463  int i1, j1;
464  double p0[ 4 ], p1[ 4 ], p2[ 4 ], p3[ 4 ];
465  InitXY( i1, j1, p0, p1, p2, p3, x, y );
466 
467  // Cubic spline coefficients
468  double C[ 4 ];
469  GetSplineCoefficients( C, x-j1 );
470 
471  // Interpolate neighbor rows
472  fx[0] = Interpolate( p0, C );
473  fx[1] = Interpolate( p1, C );
474  fx[2] = Interpolate( p2, C );
475  fx[3] = Interpolate( p3, C );
476  }
477 
492  void InterpolateY( double fy[], double x, double y ) const noexcept
493  {
494  PCL_PRECONDITION( m_data != nullptr )
495  PCL_PRECONDITION( m_width > 0 && m_height > 0 )
496 
497  // Initialize grid coordinates and source matrix
498  int i1, j1;
499  double p0[ 4 ], p1[ 4 ], p2[ 4 ], p3[ 4 ];
500  InitYX( i1, j1, p0, p1, p2, p3, x, y );
501 
502  // Cubic spline coefficients
503  double C[ 4 ];
504  GetSplineCoefficients( C, y-i1 );
505 
506  // Interpolate neighbor columns
507  fy[0] = Interpolate( p0, C );
508  fy[1] = Interpolate( p1, C );
509  fy[2] = Interpolate( p2, C );
510  fy[3] = Interpolate( p3, C );
511  }
512 
532  double InterpolateVector( const double c[], double dx ) const noexcept
533  {
534  double C[ 4 ];
535  GetSplineCoefficients( C, dx );
536  return Interpolate( c, C );
537  }
538 
546  float ClampingThreshold() const noexcept
547  {
548  return m_clamp;
549  }
550 
574  void SetClampingThreshold( float c ) noexcept
575  {
576  PCL_PRECONDITION( 0 <= c && c <= 1 )
577  m_clamp = Range( c, 0.0F, 1.0F );
578  }
579 
580 private:
581 
582  double m_clamp;
583 
584  PCL_HOT_FUNCTION
585  double Interpolate( const double* __restrict__ p, const double* __restrict__ C ) const noexcept
586  {
587  // Unclamped code:
588  //return p[0]*C[0] + p[1]*C[1] + p[2]*C[2] + p[3]*C[3];
589  double f12 = p[1]*C[1] + p[2]*C[2];
590  double f03 = p[0]*C[0] + p[3]*C[3];
591  return (-f03 < f12*m_clamp) ? f12 + f03 : f12/(C[1] + C[2]);
592  }
593 
594  PCL_HOT_FUNCTION
595  void GetSplineCoefficients( double* __restrict__ C, double dx ) const noexcept
596  {
597  double dx2 = dx*dx;
598  double dx3 = dx2*dx;
599 
600 #ifdef __PCL_BICUBIC_SPLINE_A_IS_MINUS_ONE_HALF
601  // Optimized code for a = -1/2
602  // We *really* need optimization here since this routine is called twice
603  // for each interpolated pixel.
604  double dx1_2 = dx/2;
605  double dx2_2 = dx2/2;
606  double dx3_2 = dx3/2;
607  double dx22 = dx2 + dx2;
608  double dx315 = dx3 + dx3_2;
609  C[0] = dx2 - dx3_2 - dx1_2;
610  C[1] = dx315 - dx22 - dx2_2 + 1;
611  C[2] = dx22 - dx315 + dx1_2;
612  C[3] = dx3_2 - dx2_2;
613 #else
614 # define a (__PCL_BICUBIC_SPLINE_A)
615  C[0] = a*dx3 - 2*a*dx2 + a*dx;
616  C[1] = (a + 2)*dx3 - (a + 3)*dx2 + 1;
617  C[2] = -(a + 2)*dx3 + (2*a + 3)*dx2 - a*dx;
618  C[3] = -a*dx3 + a*dx2;
619 # undef a
620 #endif
621  }
622 };
623 
635 template <typename T>
637 {
638 public:
639 
643  BicubicInterpolation() = default;
644 
649 
654  {
655  }
656 };
657 
658 // ----------------------------------------------------------------------------
659 
676 template <typename T>
678 {
679 public:
680 
685 
690 
695  {
696  }
697 
704  double operator()( double x, double y ) const override
705  {
706  PCL_PRECONDITION( f != 0 )
707  PCL_PRECONDITION( m_width > 0 && m_height > 0 )
708 
709  // Initialize grid coordinates and source matrix
710  int i1, j1;
711  double p0[ 4 ], p1[ 4 ], p2[ 4 ], p3[ 4 ];
712  InitXY( i1, j1, p0, p1, p2, p3, x, y );
713 
714  // Bicubic B-Spline interpolation functions
715 
716  double dx = x - j1;
717  double dx0 = -1 - dx;
718  double dx1 = -dx;
719  double dx2 = 1 - dx;
720  double dx3 = 2 - dx;
721 
722  double dy = y - i1;
723  double dy0 = -1 - dy;
724  double dy1 = -dy;
725  double dy2 = 1 - dy;
726  double dy3 = 2 - dy;
727 
728  return
729  BXY( p0, 0, dx0, dy0 ) +
730  BXY( p0, 1, dx1, dy0 ) +
731  BXY( p0, 2, dx2, dy0 ) +
732  BXY( p0, 3, dx3, dy0 ) +
733  BXY( p1, 0, dx0, dy1 ) +
734  BXY( p1, 1, dx1, dy1 ) +
735  BXY( p1, 2, dx2, dy1 ) +
736  BXY( p1, 3, dx3, dy1 ) +
737  BXY( p2, 0, dx0, dy2 ) +
738  BXY( p2, 1, dx1, dy2 ) +
739  BXY( p2, 2, dx2, dy2 ) +
740  BXY( p2, 3, dx3, dy2 ) +
741  BXY( p3, 0, dx0, dy3 ) +
742  BXY( p3, 1, dx1, dy3 ) +
743  BXY( p3, 2, dx2, dy3 ) +
744  BXY( p3, 3, dx3, dy3 );
745  }
746 
747 private:
748 
749  double BXY( const double* __restrict__ p, int j, double x, double y ) const noexcept
750  {
751  return p[j]*B( x )*B( y );
752  }
753 
754  double B( double x ) const noexcept
755  {
756  double fx = (x > 0) ? x*x*x : 0;
757 
758  double fxp1 = x + 1;
759  fxp1 = (fxp1 > 0) ? fxp1*fxp1*fxp1 : 0;
760 
761  double fxp2 = x + 2;
762  fxp2 = (fxp2 > 0) ? fxp2*fxp2*fxp2 : 0;
763 
764  double fxm1 = x - 1;
765  fxm1 = (fxm1 > 0) ? fxm1*fxm1*fxm1 : 0;
766 
767  return (fxp2 - 4*fxp1 + 6*fx - 4*fxm1)/6;
768  }
769 };
770 
771 // ----------------------------------------------------------------------------
772 
773 #undef InitXY
774 #undef InitYX
775 
776 // ----------------------------------------------------------------------------
777 
778 #undef m_width
779 #undef m_height
780 #undef m_fillBorder
781 #undef m_fillValue
782 #undef m_data
783 
784 // ----------------------------------------------------------------------------
785 
786 } // pcl
787 
788 #endif // __PCL_BicubicInterpolation_h
789 
790 // ----------------------------------------------------------------------------
791 // EOF pcl/BicubicInterpolation.h - Released 2024-06-18T15:48:54Z
Bicubic B-Spline Interpolation Algorithm.
double operator()(double x, double y) const override
BicubicBSplineInterpolation(const BicubicBSplineInterpolation &)=default
Base class for bicubic interpolation algorithm implementations.
BicubicInterpolationBase(const BicubicInterpolationBase &)=default
Bicubic interpolation - an alias for BicubicSplineInterpolation.
BicubicInterpolation(const BicubicInterpolation &)=default
Bicubic spline interpolation algorithm.
BicubicSplineInterpolation(const BicubicSplineInterpolation &)=default
void InterpolateY(double fy[], double x, double y) const noexcept
float ClampingThreshold() const noexcept
BicubicSplineInterpolation(float clamp=__PCL_BICUBIC_SPLINE_CLAMPING_THRESHOLD)
double InterpolateVector(const double c[], double dx) const noexcept
void SetClampingThreshold(float c) noexcept
double operator()(double x, double y) const override
void InterpolateX(double fx[], double x, double y) const noexcept
A generic interface to two-dimensional interpolation algorithms.
int TruncInt(T x) noexcept
Definition: Math.h:1132
signed long long int64
Definition: Defs.h:676
constexpr const T & Range(const T &x, const T &a, const T &b) noexcept
Definition: Utility.h:190
PCL root namespace.
Definition: AbstractImage.h:77