PCL
ShepardInterpolation.h
Go to the documentation of this file.
1 // ____ ______ __
2 // / __ \ / ____// /
3 // / /_/ // / / /
4 // / ____// /___ / /___ PixInsight Class Library
5 // /_/ \____//_____/ PCL 2.8.5
6 // ----------------------------------------------------------------------------
7 // pcl/ShepardInterpolation.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_ShepardInterpolation_h
53 #define __PCL_ShepardInterpolation_h
54 
56 
57 #include <pcl/Defs.h>
58 #include <pcl/Diagnostics.h>
59 
60 #include <pcl/QuadTree.h>
61 #include <pcl/Vector.h>
62 
63 #ifndef __PCL_BUILDING_PIXINSIGHT_APPLICATION
64 # include <pcl/Console.h>
65 #endif
66 
67 namespace pcl
68 {
69 
70 // ----------------------------------------------------------------------------
71 
72 /*
73  * Default normalized search radius for local Shepard interpolation. This is an
74  * initial search radius relative to the unit circle for the adaptive quadtree
75  * search algorithm.
76  */
77 #define __PCL_SHEPARD_DEFAULT_SEARCH_RADIUS 0.10
78 
79 /*
80  * Default power parameter for local Shepard interpolation. Larger values tend
81  * to yield more accurate interpolation devices. Small powers lead to more
82  * approximating (smoothing) devices. The chosen default value is intermediate.
83  */
84 #define __PCL_SHEPARD_DEFAULT_POWER 4
85 
86 /*
87  * Default regularization (smoothing) factor for local Shepard interpolation,
88  * in the range [0,1). This is a clipping fraction for Winsorization of nearby
89  * function values in the point interpolation routine.
90  */
91 #define __PCL_SHEPARD_DEFAULT_REGULARIZATION 0
92 
93 // ----------------------------------------------------------------------------
94 
124 template <typename T>
125 class PCL_CLASS ShepardInterpolation
126 {
127 public:
128 
133 
137  using scalar = typename vector_type::scalar;
138 
143 
148 
152  constexpr static int BucketCapacity = 16;
153 
157  ShepardInterpolation() = default;
158 
165 
170 
175  {
176  }
177 
183  ShepardInterpolation& operator =( const ShepardInterpolation& ) = delete;
184 
188  ShepardInterpolation& operator =( ShepardInterpolation&& ) = default;
189 
194  bool IsValid() const
195  {
196  return !m_T.IsEmpty();
197  }
198 
220  void SetPower( int m )
221  {
222  PCL_PRECONDITION( m > 0 )
223  m_mu = (m > 0) ? m : __PCL_SHEPARD_DEFAULT_POWER;
224  }
225 
231  int Power() const
232  {
233  return m_mu;
234  }
235 
262  void SetRadius( double R )
263  {
264  PCL_PRECONDITION( R > 0 )
265  m_R = (R < 0 || 1 + R == 1) ? __PCL_SHEPARD_DEFAULT_SEARCH_RADIUS : R;
266  }
267 
272  double Radius() const
273  {
274  return m_R;
275  }
276 
294  void SetSmoothing( float r )
295  {
296  PCL_PRECONDITION( r >= 0 && r < 1 )
297  m_r = (r < 0 || r >= 1) ? __PCL_SHEPARD_DEFAULT_REGULARIZATION : r;
298  }
299 
304  float Smoothing() const
305  {
306  return m_r;
307  }
308 
332  void Initialize( const T* x, const T* y, const T* z, int n )
333  {
334  DoInitialize( nullptr/*rect*/, x, y, z, n );
335  }
336 
372  void Initialize( const search_rect& rect, const T* x, const T* y, const T* z, int n )
373  {
374  DoInitialize( &rect, x, y, z, n );
375  }
376 
381  void Clear()
382  {
383  m_T.Clear();
384  }
385 
403  T operator ()( double x, double y ) const
404  {
405  PCL_PRECONDITION( !m_T.IsEmpty() )
406 
407  const double dx = m_r0*(x - m_x0);
408  const double dy = m_r0*(y - m_y0);
409 
410  for ( double R = m_R; ; R += m_R )
411  {
412  int m = 0;
413  double R2 = R*R;
414  Array<DPoint> V;
415  m_T.Search( DRect( dx-R, dy-R, dx+R, dy+R ),
416  [&]( const vector_type& v, void* )
417  {
418  double x = dx - v[0];
419  double y = dy - v[1];
420  double r2 = x*x + y*y;
421  if ( r2 < R2 )
422  {
423  ++m;
424  double w = PowI( 1 - Sqrt( r2 )/R, m_mu );
425  /*
426  * N.B. The equivalent code below is about 400 times
427  * slower than the above call to PowI() for m_mu=16.
428  * Measured on a TR 2990WX.
429  *
430  double w = 1 - Sqrt( r2 )/R;
431  for ( int i = 1; i < m_mu; ++i )
432  w *= w;
433  */
434  V << DPoint( w, w*v[2] );
435  }
436  },
437  nullptr );
438  if ( m >= 3 )
439  {
440  if ( m_r > 0 )
441  {
442  /*
443  * Regularization by Winsorization of the weighted sample.
444  */
445  int r = Min( TruncInt( m_r * m ), (m >> 1) - (m & 1)^1 );
446  if ( r > 0 )
447  {
448  V.Sort( []( const DPoint& v1, const DPoint& v2 )
449  {
450  return v1.y < v2.y;
451  } );
452  for ( int i = 0; i < r; ++i )
453  {
454  V[i].y = V[r].y;
455  V[m-i-1].y = V[m-r-1].y;
456  }
457  }
458  }
459  DPoint Wz( 0 );
460  for ( const DPoint& v : V )
461  Wz += v;
462  if ( 1 + Wz.x != 1 )
463  return T( Wz.y/Wz.x );
464  if ( R >= 1 )
465  break; // degenerate!
466  }
467  }
468 
469  return 0; // empty!?
470  }
471 
477  template <typename Tp>
478  T operator ()( const GenericPoint<Tp>& p ) const
479  {
480  return operator ()( double( p.x ), double( p.y ) );
481  }
482 
494  template <typename T1>
495  void Evaluate( T1* Z, const T1* X, const T1* Y, size_type n ) const
496  {
497  PCL_PRECONDITION( !m_T.IsEmpty() )
498  PCL_PRECONDITION( X != nullptr && Y != nullptr && Z != nullptr )
499  PCL_PRECONDITION( n > 0 )
500  for ( size_type i = 0; i < n; ++i )
501  Z[i] = T1( operator()( double( X[i] ), double( Y[i] ) ) );
502  }
503 
514  {
515  return false;
516  }
517 
518 protected:
519 
520  double m_r0 = 1; // scaling factor for normalization of node coordinates
521  double m_x0 = 0; // zero offset for normalization of X node coordinates
522  double m_y0 = 0; // zero offset for normalization of Y node coordinates
523  int m_mu = __PCL_SHEPARD_DEFAULT_POWER; // power parameter (leveling factor)
524  double m_R = __PCL_SHEPARD_DEFAULT_SEARCH_RADIUS; // initial search radius
525  float m_r = __PCL_SHEPARD_DEFAULT_REGULARIZATION; // regularization (clipping fraction)
526  search_tree m_T; // tree points store input coordinates and function values
527 
533  void DoInitialize( const search_rect* rect, const T* x, const T* y, const T* z, int n )
534  {
535  PCL_PRECONDITION( x != nullptr && y != nullptr && z != nullptr )
536  PCL_PRECONDITION( n > 2 )
537 
538  if ( n < 3 )
539  throw Error( "ShepardInterpolation::Initialize(): At least three input nodes must be specified." );
540 
541  Clear();
542 
543  try
544  {
545  if ( rect == nullptr )
546  {
547  /*
548  * Find mean coordinate values.
549  */
550  m_x0 = m_y0 = 0;
551  for ( int i = 0; i < n; ++i )
552  {
553  m_x0 += x[i];
554  m_y0 += y[i];
555  }
556  m_x0 /= n;
557  m_y0 /= n;
558 
559  /*
560  * Find radius of unit circle.
561  */
562  m_r0 = 0;
563  for ( int i = 0; i < n; ++i )
564  {
565  double dx = x[i] - m_x0;
566  double dy = y[i] - m_y0;
567  double r = Sqrt( dx*dx + dy*dy );
568  if ( r > m_r0 )
569  m_r0 = r;
570  }
571  }
572  else
573  {
574  m_x0 = rect->CenterX();
575  m_y0 = rect->CenterY();
576  m_r0 = rect->Diagonal()/2;
577  }
578 
579  if ( 1 + m_r0 == 1 )
580  throw Error( "ShepardInterpolation::Initialize(): Empty or insignificant interpolation space." );
581  m_r0 = 1/m_r0;
582 
583  /*
584  * Build working vector. Transform coordinates to the unit circle.
585  */
587  for ( int i = 0; i < n; ++i )
588  V << vector_type( m_r0*(x[i] - m_x0), m_r0*(y[i] - m_y0), z[i] );
589 
590  /*
591  * Find and remove duplicate input nodes. Two nodes are equal if their
592  * coordinates don't differ more than the machine epsilon for the
593  * floating point type T.
594  */
595  V.Sort( []( const vector_type& p, const vector_type& q )
596  {
597  return (p[0] != q[0]) ? p[0] < q[0] : p[1] < q[1];
598  } );
599  Array<int> remove;
600  for ( int i = 0, j = 1; j < n; ++i, ++j )
601  if ( Abs( V[i][0] - V[j][0] ) <= std::numeric_limits<T>::epsilon() )
602  if ( Abs( V[i][1] - V[j][1] ) <= std::numeric_limits<T>::epsilon() )
603  remove << i;
604  if ( !remove.IsEmpty() )
605  {
607  int i = 0;
608  for ( int j : remove )
609  {
610  for ( ; i < j; ++i )
611  U << V[i];
612  ++i;
613  }
614  for ( ; i < n; ++i )
615  U << V[i];
616  if ( U.Length() < 3 )
617  throw Error( "ShepardInterpolation::Initialize(): Less than three input nodes left after sanitization." );
618  V = U;
619  }
620 
621  /*
622  * Build the point search tree.
623  */
624  if ( rect == nullptr )
625  m_T.Build( V, BucketCapacity );
626  else
627  {
628  m_T.Build( *rect, V, BucketCapacity );
629  if ( m_T.Length() < 3 )
630  throw Error( "ShepardInterpolation::Initialize(): Less than three input nodes in the specified search region." );
631  }
632  }
633  catch ( ... )
634  {
635  Clear();
636  throw;
637  }
638  }
639 };
640 
641 // ----------------------------------------------------------------------------
642 
654 template <class P = DPoint>
656 {
657 public:
658 
662  using point = P;
663 
668 
673 
679 
684 
689 
698  int power = __PCL_SHEPARD_DEFAULT_POWER,
699  double radius = __PCL_SHEPARD_DEFAULT_SEARCH_RADIUS )
700  {
701  Initialize( P1, P2, power, radius );
702  }
703 
711  PointShepardInterpolation( const surface& Sx, const surface& Sy )
712  {
713  Initialize( Sx, Sy );
714  }
715 
721  PointShepardInterpolation& operator =( const PointShepardInterpolation& ) = delete;
722 
727 
760  void Initialize( const point_list& P1, const point_list& P2,
761  int power = __PCL_SHEPARD_DEFAULT_POWER,
762  double radius = __PCL_SHEPARD_DEFAULT_SEARCH_RADIUS,
763  float smoothing = __PCL_SHEPARD_DEFAULT_REGULARIZATION )
764  {
765  PCL_PRECONDITION( P1.Length() >= 3 )
766  PCL_PRECONDITION( P1.Length() <= P2.Length() )
767  PCL_PRECONDITION( power > 0 )
768  PCL_PRECONDITION( radius > 0 )
769  PCL_PRECONDITION( smoothing >= 0 && smoothing < 1 )
770 
771  m_Sx.Clear();
772  m_Sy.Clear();
773 
774  m_Sx.SetPower( power );
775  m_Sy.SetPower( power );
776 
777  m_Sx.SetRadius( radius );
778  m_Sy.SetRadius( radius );
779 
780  m_Sx.SetSmoothing( smoothing );
781  m_Sy.SetSmoothing( smoothing );
782 
783  if ( P1.Length() < 3 || P2.Length() < 3 )
784  throw Error( "PointShepardInterpolation::Initialize(): At least three input nodes must be specified." );
785 
786  if ( P2.Length() < P1.Length() )
787  throw Error( "PointShepardInterpolation::Initialize(): Incompatible point array lengths." );
788 
789  DVector X( int( P1.Length() ) ),
790  Y( int( P1.Length() ) ),
791  Zx( int( P1.Length() ) ),
792  Zy( int( P1.Length() ) );
793  for ( int i = 0; i < X.Length(); ++i )
794  {
795  X[i] = P1[i].x;
796  Y[i] = P1[i].y;
797  Zx[i] = P2[i].x;
798  Zy[i] = P2[i].y;
799  }
800  m_Sx.Initialize( X.Begin(), Y.Begin(), Zx.Begin(), X.Length() );
801  m_Sy.Initialize( X.Begin(), Y.Begin(), Zy.Begin(), X.Length() );
802  }
803 
808  void Clear()
809  {
810  m_Sx.Clear();
811  m_Sy.Clear();
812  }
813 
818  bool IsValid() const
819  {
820  return m_Sx.IsValid() && m_Sy.IsValid();
821  }
822 
827  const surface& SurfaceX() const
828  {
829  return m_Sx;
830  }
831 
836  const surface& SurfaceY() const
837  {
838  return m_Sy;
839  }
840 
844  template <typename T>
845  DPoint operator ()( T x, T y ) const
846  {
847  return DPoint( m_Sx( x, y ), m_Sy( x, y ) );
848  }
849 
853  template <typename T>
854  DPoint operator ()( const GenericPoint<T>& p ) const
855  {
856  return operator ()( p.x, p.y );
857  }
858 
870  template <typename T>
871  void Evaluate( T* ZX, T* ZY, const T* X, const T* Y, size_type n ) const
872  {
873  PCL_PRECONDITION( ISValid() )
874  m_Sx.Evaluate( ZX, X, Y, n );
875  m_Sy.Evaluate( ZY, X, Y, n );
876  }
877 
888  {
889  return false;
890  }
891 
892 private:
893 
894  /*
895  * The surface interpolations in the X and Y plane directions.
896  */
897  surface m_Sx, m_Sy;
898 
899  friend class DrizzleData;
900  friend class DrizzleDataDecoder;
901 };
902 
903 // ----------------------------------------------------------------------------
904 
905 } // pcl
906 
907 #endif // __PCL_ShepardInterpolation_h
908 
909 // ----------------------------------------------------------------------------
910 // EOF pcl/ShepardInterpolation.h - Released 2024-12-28T16:53:48Z
bool IsEmpty() const noexcept
Definition: Array.h:324
size_type Length() const noexcept
Definition: Array.h:278
Drizzle integration data parser and generator.
Definition: DrizzleData.h:143
A simple exception with an associated error message.
Definition: Exception.h:239
A generic point in the two-dimensional space.
Definition: Point.h:100
component x
Abscissa (horizontal, or X-axis coordinate).
Definition: Point.h:111
component y
Ordinate (vertical, or Y-axis coordinate).
Definition: Point.h:112
A generic rectangle in the two-dimensional space.
Definition: Rectangle.h:316
Generic vector of arbitrary length.
Definition: Vector.h:107
iterator Begin()
Definition: Vector.h:1918
Vector Shepard interpolation/approximation in two dimensions.
const surface & SurfaceX() const
PointShepardInterpolation(const surface &Sx, const surface &Sy)
PointShepardInterpolation(const point_list &P1, const point_list &P2, int power=__PCL_SHEPARD_DEFAULT_POWER, double radius=__PCL_SHEPARD_DEFAULT_SEARCH_RADIUS)
PointShepardInterpolation(const PointShepardInterpolation &)=default
const surface & SurfaceY() const
void Evaluate(T *ZX, T *ZY, const T *X, const T *Y, size_type n) const
void Initialize(const point_list &P1, const point_list &P2, int power=__PCL_SHEPARD_DEFAULT_POWER, double radius=__PCL_SHEPARD_DEFAULT_SEARCH_RADIUS, float smoothing=__PCL_SHEPARD_DEFAULT_REGULARIZATION)
PointShepardInterpolation(PointShepardInterpolation &&)=default
Two-dimensional surface interpolation with the local Shepard method.
void Evaluate(T1 *Z, const T1 *X, const T1 *Y, size_type n) const
ShepardInterpolation(const ShepardInterpolation &)=delete
ShepardInterpolation(ShepardInterpolation &&)=default
void DoInitialize(const search_rect *rect, const T *x, const T *y, const T *z, int n)
typename vector_type::scalar scalar
typename search_tree::rectangle search_rect
void Initialize(const search_rect &rect, const T *x, const T *y, const T *z, int n)
void Initialize(const T *x, const T *y, const T *z, int n)
Complex< T > Sqrt(const Complex< T > &c) noexcept
Definition: Complex.h:688
T Abs(const Complex< T > &c) noexcept
Definition: Complex.h:443
T PowI(T x, int n) noexcept
Definition: Math.h:1786
int TruncInt(T x) noexcept
Definition: Math.h:1203
size_t size_type
Definition: Defs.h:609
constexpr const T & Min(const T &a, const T &b) noexcept
Definition: Utility.h:90
PCL root namespace.
Definition: AbstractImage.h:77