PCL
BicubicFilterInterpolation.h
Go to the documentation of this file.
1 // ____ ______ __
2 // / __ \ / ____// /
3 // / /_/ // / / /
4 // / ____// /___ / /___ PixInsight Class Library
5 // /_/ \____//_____/ PCL 2.8.5
6 // ----------------------------------------------------------------------------
7 // pcl/BicubicFilterInterpolation.h - Released 2024-12-28T16:53:48Z
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_BicubicFilterInterpolation_h
53 #define __PCL_BicubicFilterInterpolation_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 #include <pcl/Vector.h>
64 
65 namespace pcl
66 {
67 
68 // ----------------------------------------------------------------------------
69 
86 class PCL_CLASS CubicFilter
87 {
88 public:
89 
94  CubicFilter( double B, double C )
95  : m_B( B )
96  , m_C( C )
97  {
98  m_k31 = ( 12 - 9*m_B - 6*m_C)/6;
99  m_k21 = (-18 + 12*m_B + 6*m_C)/6;
100  m_k01 = ( 6 - 2*m_B )/6;
101  m_k32 = ( -m_B - 6*m_C)/6;
102  m_k22 = ( 6*m_B + 30*m_C)/6;
103  m_k12 = ( -12*m_B - 48*m_C)/6;
104  m_k02 = ( 8*m_B + 24*m_C)/6;
105  }
106 
110  CubicFilter( const CubicFilter& ) = default;
111 
115  virtual ~CubicFilter()
116  {
117  }
118 
122  CubicFilter& operator =( const CubicFilter& ) = default;
123 
134  PCL_HOT_FUNCTION
135  double operator()( double x ) const noexcept
136  {
137  if ( x < 0 )
138  x = -x;
139  return (x < 1) ? m_k01 + x*x*(m_k21 + x*m_k31) :
140  m_k02 + x*(m_k12 + x*(m_k22 + x*m_k32));
141  }
142 
151  double Width() const noexcept
152  {
153  return 2.0;
154  }
155 
159  virtual String Description() const
160  {
161  return String().Format( "Cubic filter (B=%.6f, C=%.6f)", m_B, m_C );
162  }
163 
167  virtual CubicFilter* Clone() const
168  {
169  return new CubicFilter( *this );
170  }
171 
172 protected:
173 
174  double m_k31, m_k21, m_k01, m_k32, m_k22, m_k12, m_k02;
175 
176 private:
177 
178  double m_B, m_C; // for reference only; not used in calculations
179 };
180 
181 // ----------------------------------------------------------------------------
182 
193 {
194 public:
195 
200  : CubicFilter( 1.0/3, 1.0/3 )
201  {
202  }
203 
208 
213  {
214  }
215 
218  String Description() const override
219  {
220  return "Mitchell-Netravali cubic filter";
221  }
222 
225  CubicFilter* Clone() const override
226  {
227  return new MitchellNetravaliCubicFilter( *this );
228  }
229 };
230 
231 // ----------------------------------------------------------------------------
232 
242 class PCL_CLASS CatmullRomSplineFilter : public CubicFilter
243 {
244 public:
245 
250  : CubicFilter( 0, 0.5 )
251  {
252  }
253 
258 
263  {
264  }
265 
268  String Description() const override
269  {
270  return "Catmull-Rom spline filter";
271  }
272 
275  CubicFilter* Clone() const override
276  {
277  return new CatmullRomSplineFilter( *this );
278  }
279 };
280 
281 // ----------------------------------------------------------------------------
282 
292 class PCL_CLASS CubicBSplineFilter : public CubicFilter
293 {
294 public:
295 
300  : CubicFilter( 1, 0 )
301  {
302  }
303 
308 
313  {
314  }
315 
318  String Description() const override
319  {
320  return "Cubic B-spline filter";
321  }
322 
325  CubicFilter* Clone() const override
326  {
327  return new CubicBSplineFilter( *this );
328  }
329 };
330 
331 // ----------------------------------------------------------------------------
332 
333 #define m_width this->m_width
334 #define m_height this->m_height
335 #define m_fillBorder this->m_fillBorder
336 #define m_fillValue this->m_fillValue
337 #define m_data this->m_data
338 
339 // ----------------------------------------------------------------------------
340 
356 template <typename T>
358 {
359 public:
360 
372  BicubicFilterInterpolation( int rx, int ry, const CubicFilter& filter )
373  : m_rx( Max( 1, rx ) )
374  , m_ry( Max( 1, ry ) )
375  , m_filter( filter )
376  {
377  PCL_PRECONDITION( rx >= 1 )
378  PCL_PRECONDITION( ry >= 1 )
379  Initialize();
380  }
381 
386 
391  {
392  }
393 
396  void Initialize( const T* data, int dataWidth, int dataHeight ) override
397  {
398  BidimensionalInterpolation<T>::Initialize( data, dataWidth, dataHeight );
399 
400  if ( m_rx > m_width || m_ry > m_height )
401  {
402  m_rx = Min( m_rx, m_width );
403  m_ry = Min( m_ry, m_height );
404  Initialize();
405  }
406  }
407 
418  PCL_HOT_FUNCTION
419  double operator()( double x, double y ) const override
420  {
421  PCL_PRECONDITION( m_data != nullptr )
422  PCL_PRECONDITION( m_width > 0 && m_height > 0 )
423  PCL_PRECONDITION( x >= 0 && x < m_width )
424  PCL_PRECONDITION( y >= 0 && y < m_height )
425  PCL_CHECK( m_rx >= 1 && m_rx <= m_width )
426  PCL_CHECK( m_ry >= 1 && m_ry <= m_height )
427  PCL_CHECK( !m_k.IsEmpty() )
428 
429  int x0 = Range( RoundInt( x ), 0, m_width-1 );
430  int y0 = Range( RoundInt( y ), 0, m_height-1 );
431 
432  double dx = x - x0;
433  double dy = y - y0;
434 
435  int64 d = int64( y0 - m_ry )*m_width + x0 - m_rx;
436 
437  if ( m_rx >= m_ry )
438  {
439  double wx = 0;
440  for ( int xi = -m_rx, c = 0; xi <= m_rx; ++xi, ++c )
441  wx += (m_k[c] = m_filter( m_sx*(xi - dx) ));
442 
443  double sy = 0, wy = 0;
444  for ( int yi = -m_ry; yi <= m_ry; ++yi )
445  {
446  double ky = m_filter( m_sy*(yi - dy) );
447  wy += ky;
448 
449  const T* p;
450  int y = y0 + yi;
451 
452  if ( y < 0 )
453  {
454  d += m_width;
455 
456  if ( m_fillBorder )
457  {
458  sy += m_fillValue * ky;
459  continue;
460  }
461 
462  p = m_data + (d - 2*int64( m_width )*(y + 1));
463  }
464  else if ( y >= m_height )
465  {
466  if ( m_fillBorder )
467  {
468  sy += m_fillValue * ky;
469  continue;
470  }
471 
472  p = m_data + (d -= m_width);
473  }
474  else
475  {
476  p = m_data + d;
477  d += m_width;
478  }
479 
480  int x = x0 - m_rx;
481  int x2 = x0 + m_rx;
482  int x1 = Min( x2, m_width-1 );
483  double sx = 0;
484  const double* kx = m_k.Begin();
485 
486  while ( x < 0 )
487  {
488  ++x;
489  ++p;
490  sx += (m_fillBorder ? m_fillValue : double( *(p - (x + x)) )) * *kx++;
491  }
492 
493  for ( ;; )
494  {
495  sx += *p++ * *kx++;
496  if ( ++x > x1 )
497  break;
498  }
499 
500  while ( x <= x2 )
501  {
502  sx += (m_fillBorder ? m_fillValue : double( *--p )) * *kx++;
503  ++x;
504  }
505 
506  sy += sx/wx * ky;
507  }
508 
509  return sy/wy;
510  }
511  else
512  {
513  double wy = 0;
514  for ( int yi = -m_ry, r = 0; yi <= m_ry; ++yi, ++r )
515  wy += (m_k[r] = m_filter( m_sy*(yi - dy) ));
516 
517  double sx = 0, wx = 0;
518  for ( int xi = -m_rx; xi <= m_rx; ++xi )
519  {
520  double kx = m_filter( m_sx*(xi - dx) );
521  wx += kx;
522 
523  const T* p;
524  int x = x0 + xi;
525 
526  if ( x < 0 )
527  {
528  ++d;
529 
530  if ( m_fillBorder )
531  {
532  sx += m_fillValue * kx;
533  continue;
534  }
535 
536  p = m_data + (d - 2*(x + 1));
537  }
538  else if ( x >= m_width )
539  {
540  if ( m_fillBorder )
541  {
542  sx += m_fillValue * kx;
543  continue;
544  }
545 
546  p = m_data + --d;
547  }
548  else
549  p = m_data + d++;
550 
551  int y = y0 - m_ry;
552  int y2 = y0 + m_ry;
553  int y1 = Min( y2, m_height-1 );
554  double sy = 0;
555  const double* ky = m_k.Begin();
556 
557  while ( y < 0 )
558  {
559  ++y;
560  p += m_width;
561  sy += (m_fillBorder ? m_fillValue : double( *(p - (y + y)*m_width) )) * *ky++;
562  }
563 
564  for ( ;; )
565  {
566  sy += *p * *ky++;
567  p += m_width;
568  if ( ++y > y1 )
569  break;
570  }
571 
572  while ( y <= y2 )
573  {
574  sy += (m_fillBorder ? m_fillValue : double( *(p -= m_width) )) * *ky++;
575  ++y;
576  }
577 
578  sx += sy/wy * kx;
579  }
580 
581  return sx/wx;
582  }
583  }
584 
588  int HorizontalRadius() const noexcept
589  {
590  return m_rx;
591  }
592 
596  int VerticalRadius() const noexcept
597  {
598  return m_ry;
599  }
600 
609  void SetRadii( int rx, int ry )
610  {
611  if ( rx != m_rx || ry != m_ry )
612  {
613  m_rx = Max( 1, rx );
614  m_ry = Max( 1, ry );
615 
616  if ( this->data != nullptr )
617  {
618  if ( m_rx > m_width )
619  m_rx = m_width;
620  if ( m_ry > m_height )
621  m_ry = m_height;
622  }
623 
624  Initialize();
625  }
626  }
627 
632  const CubicFilter& Filter() const noexcept
633  {
634  return m_filter;
635  }
636 
640  void SetFilter( const CubicFilter& filter )
641  {
642  m_filter = filter;
643  Initialize();
644  }
645 
646 protected:
647 
648  int m_rx, m_ry; // horizontal and vertical radii
649  double m_sx, m_sy; // filter scaling ratios
650  mutable DVector m_k; // workspace for interpolated values
651  CubicFilter m_filter;
652 
653  void Initialize()
654  {
655  m_sx = m_filter.Width()/(m_rx + 0.5);
656  m_sy = m_filter.Width()/(m_ry + 0.5);
657  m_k = DVector( 1 + (Max( m_rx, m_ry ) << 1) );
658  }
659 };
660 
661 // ----------------------------------------------------------------------------
662 
663 #undef m_width
664 #undef m_height
665 #undef m_fillBorder
666 #undef m_fillValue
667 #undef m_data
668 
669 // ----------------------------------------------------------------------------
670 
671 } // pcl
672 
673 #endif // __PCL_BicubicFilterInterpolation_h
674 
675 // ----------------------------------------------------------------------------
676 // EOF pcl/BicubicFilterInterpolation.h - Released 2024-12-28T16:53:48Z
Bicubic filter interpolation algorithms.
const CubicFilter & Filter() const noexcept
void Initialize(const T *data, int dataWidth, int dataHeight) override
BicubicFilterInterpolation(const BicubicFilterInterpolation &)=default
PCL_HOT_FUNCTION double operator()(double x, double y) const override
BicubicFilterInterpolation(int rx, int ry, const CubicFilter &filter)
void SetFilter(const CubicFilter &filter)
A generic interface to two-dimensional interpolation algorithms.
virtual void Initialize(const T *data, int width, int height)
CubicFilter * Clone() const override
CatmullRomSplineFilter(const CatmullRomSplineFilter &)=default
CubicBSplineFilter(const CubicBSplineFilter &)=default
CubicFilter * Clone() const override
String Description() const override
Mitchell-Netravali parameterized cubic filters.
virtual CubicFilter * Clone() const
CubicFilter(double B, double C)
PCL_HOT_FUNCTION double operator()(double x) const noexcept
CubicFilter(const CubicFilter &)=default
virtual String Description() const
double Width() const noexcept
Generic vector of arbitrary length.
Definition: Vector.h:107
Mitchell-Netravali cubic filter with B=C=1/3.
CubicFilter * Clone() const override
MitchellNetravaliCubicFilter(const MitchellNetravaliCubicFilter &)=default
Unicode (UTF-16) string.
Definition: String.h:8146
String & Format(const_c_string8 fmt,...)
Definition: String.h:11129
int RoundInt(T x) noexcept
Definition: Math.h:1550
signed long long int64
Definition: Defs.h:676
constexpr const T & Min(const T &a, const T &b) noexcept
Definition: Utility.h:90
constexpr const T & Range(const T &x, const T &a, const T &b) noexcept
Definition: Utility.h:190
constexpr const T & Max(const T &a, const T &b) noexcept
Definition: Utility.h:119
PCL root namespace.
Definition: AbstractImage.h:77