PCL
SurfaceSpline.h
Go to the documentation of this file.
1 // ____ ______ __
2 // / __ \ / ____// /
3 // / /_/ // / / /
4 // / ____// /___ / /___ PixInsight Class Library
5 // /_/ \____//_____/ PCL 2.9.3
6 // ----------------------------------------------------------------------------
7 // pcl/SurfaceSpline.h - Released 2025-02-21T12:13:32Z
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-2025 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_SurfaceSpline_h
53 #define __PCL_SurfaceSpline_h
54 
56 
57 #include <pcl/Defs.h>
58 #include <pcl/Diagnostics.h>
59 
60 #include <pcl/AbstractImage.h>
61 #include <pcl/Array.h>
62 #include <pcl/AutoPointer.h>
63 #include <pcl/Homography.h>
64 #include <pcl/ParallelProcess.h>
65 #include <pcl/Point.h>
66 #include <pcl/QuadTree.h>
67 #include <pcl/SurfaceSimplifier.h>
68 #include <pcl/Thread.h>
69 #include <pcl/Vector.h>
70 
71 #ifndef __PCL_BUILDING_PIXINSIGHT_APPLICATION
72 # include <pcl/Console.h>
73 # include <pcl/StandardStatus.h>
74 #endif
75 
76 namespace pcl
77 {
78 
79 // ----------------------------------------------------------------------------
80 
81 /*
82  * The maximum allowed number of points in a point surface spline by default.
83  */
84 #define __PCL_PSSPLINE_DEFAULT_MAX_POINTS 2100
85 
86 /*
87  * Default quadtree bucket capacity for recursive surface subspline generation.
88  */
89 #define __PCL_RSSPLINE_DEFAULT_TREE_BUCKET_CAPACITY 64
90 
91 /*
92  * Default maximum spline length for a non-recursive quadtree node in a
93  * recursive surface spline.
94  */
95 #define __PCL_RSSPLINE_DEFAULT_SPLINE_MAX_LENGTH 1600
96 
97 /*
98  * Whether to allow extrapolation outside the interpolation region for
99  * recursive surface splines. Extrapolation is disabled by default because
100  * recursively defined subsplines are slightly more prone to oscillation than
101  * normal surface splines.
102  */
103 #define __PCL_RSSPLINE_DEFAULT_ALLOW_EXTRAPOLATION false
104 
105 // ----------------------------------------------------------------------------
106 
141 namespace RadialBasisFunction
142 {
143  enum value_type
144  {
145  Unknown = -1,
146  VariableOrder = 0, // phi(r) = (r^2)^(m-1) * Ln( r^2 )
147  ThinPlateSpline, // phi(r) = r^2 * Ln( r )
148  Gaussian, // phi(r) = Exp( -(eps*r)^2 )
149  Multiquadric, // phi(r) = Sqrt( 1 + (eps*r)^2 )
150  InverseMultiquadric, // phi(r) = 1/Sqrt( 1 + (eps*r)^2 )
151  InverseQuadratic, // phi(r) = 1/( 1 + (eps*r)^2 )
152  __number_of_pcl_implemented_rbfs__,
153  DDMVariableOrder = 100, // Variable-order polyharmonic DDM-RBF
154  DDMThinPlateSpline = 101, // Thin plate spline DDM-RBF
155  DDMMultiquadric = 102, // Multiquadric DDM-RBF
156  Default = ThinPlateSpline
157  };
158 
163  inline static bool Validate( int rbf )
164  {
165  return rbf >= 0 && rbf < __number_of_pcl_implemented_rbfs__
166  || rbf == DDMThinPlateSpline
167  || rbf == DDMVariableOrder
168  || rbf == DDMMultiquadric;
169  }
170 
175  inline static bool HasCoreImplementation( int rbf )
176  {
177  return rbf > __number_of_pcl_implemented_rbfs__;
178  }
179 
185  inline static bool HasDDMImplementation( int rbf )
186  {
187  return rbf == DDMThinPlateSpline
188  || rbf == DDMVariableOrder
189  || rbf == DDMMultiquadric;
190  }
191 
196  inline static bool IsVariableOrder( int rbf )
197  {
198  return rbf == DDMVariableOrder || rbf == VariableOrder;
199  }
200 
206  inline static bool HasShapeParameter( int rbf )
207  {
208  return rbf == Gaussian
209  || rbf == Multiquadric
210  || rbf == InverseMultiquadric
211  || rbf == InverseQuadratic;
212  }
213 }
214 
215 // ----------------------------------------------------------------------------
216 
221 class PCL_CLASS SurfaceSplineBase
222 {
223 public:
224 
229  using rbf_type = RadialBasisFunction::value_type;
230 
231 protected:
232 
236  SurfaceSplineBase() = default;
237 
242 
247 
252  {
253  }
254 
258  SurfaceSplineBase& operator =( const SurfaceSplineBase& ) = default;
259 
263  SurfaceSplineBase& operator =( SurfaceSplineBase&& ) = default;
264 
268  static void Generate( float* __restrict__ c, void** handle,
269  rbf_type, double e2, bool polynomial,
270  const float* __restrict__ x, const float* __restrict__ y, const float* __restrict__ z,
271  int n, int m, float r, const float* __restrict__ w );
272 
276  static void Generate( double* __restrict__ c, void** handle,
277  rbf_type, double e2, bool polynomial,
278  const double* __restrict__ x, const double* __restrict__ y, const double* __restrict__ z,
279  int n, int m, float r, const float* __restrict__ w );
280 
284  static double EvaluateHandle( const void* handle, double x, double y );
285 
290  static void EvaluateHandle( const void* handle, float* z, const float* x, const float* y, double x0, double y0, double r, size_type n );
291 
296  static void EvaluateHandle( const void* handle, double* z, const double* x, const double* y, double x0, double y0, double r, size_type n );
297 
301  static void SerializeHandle( IsoString& data, const void* handle );
302 
306  static void* DeserializeHandle( const IsoString& data );
307 
311  static void DestroyHandle( void* handle );
312 
316  static void* DuplicateHandle( const void* handle );
317 
318  friend class SplineWorldTransformation;
319 };
320 
321 // ----------------------------------------------------------------------------
322 
370 template <typename T>
371 class PCL_CLASS SurfaceSpline : private SurfaceSplineBase
372 {
373 public:
374 
379  using scalar = T;
380 
386 
391  using weight = float;
392 
397 
402  using rbf_type = SurfaceSplineBase::rbf_type;
403 
408  SurfaceSpline() = default;
409 
414  : m_rbf( S.m_rbf )
415  , m_havePolynomial( S.m_havePolynomial )
416  , m_eps( S.m_eps )
417  , m_eps2( S.m_eps2 )
418  , m_x( S.m_x )
419  , m_y( S.m_y )
420  , m_z( S.m_z )
421  , m_r0( S.m_r0 )
422  , m_x0( S.m_x0 )
423  , m_y0( S.m_y0 )
424  , m_order( S.m_order )
425  , m_smoothing( S.m_smoothing )
426  , m_weights( S.m_weights )
427  , m_spline( S.m_spline )
428  , m_serialization( S.m_serialization )
429  {
430  m_handle = DuplicateHandle( S.m_handle );
431  }
432 
437  : m_rbf( S.m_rbf )
438  , m_havePolynomial( S.m_havePolynomial )
439  , m_eps( S.m_eps )
440  , m_eps2( S.m_eps2 )
441  , m_x( std::move( S.m_x ) )
442  , m_y( std::move( S.m_y ) )
443  , m_z( std::move( S.m_z ) )
444  , m_r0( S.m_r0 )
445  , m_x0( S.m_x0 )
446  , m_y0( S.m_y0 )
447  , m_order( S.m_order )
448  , m_smoothing( S.m_smoothing )
449  , m_weights( std::move( S.m_weights ) )
450  , m_spline( std::move( S.m_spline ) )
451  , m_serialization( S.m_serialization )
452  {
453  m_handle = S.m_handle;
454  S.m_handle = 0;
455  }
456 
460  ~SurfaceSpline() override
461  {
462  DestroyHandle( m_handle );
463  m_handle = 0;
464  }
465 
469  SurfaceSpline& operator =( const SurfaceSpline& S )
470  {
471  if ( &S != this )
472  {
473  Clear();
474  m_rbf = S.m_rbf;
475  m_havePolynomial = S.m_havePolynomial;
476  m_eps = S.m_eps;
477  m_eps2 = S.m_eps2;
478  m_x = S.m_x;
479  m_y = S.m_y;
480  m_z = S.m_z;
481  m_r0 = S.m_r0;
482  m_x0 = S.m_x0;
483  m_y0 = S.m_y0;
484  m_order = S.m_order;
485  m_smoothing = S.m_smoothing;
486  m_weights = S.m_weights;
487  m_spline = S.m_spline;
488  m_serialization = S.m_serialization;
489  m_handle = DuplicateHandle( S.m_handle );
490  }
491  return *this;
492  }
493 
497  SurfaceSpline& operator =( SurfaceSpline&& S )
498  {
499  if ( &S != this )
500  {
501  Clear();
502  m_rbf = S.m_rbf;
503  m_havePolynomial = S.m_havePolynomial;
504  m_eps = S.m_eps;
505  m_eps2 = S.m_eps2;
506  m_x = std::move( S.m_x );
507  m_y = std::move( S.m_y );
508  m_z = std::move( S.m_z );
509  m_r0 = S.m_r0;
510  m_x0 = S.m_x0;
511  m_y0 = S.m_y0;
512  m_order = S.m_order;
513  m_smoothing = S.m_smoothing;
514  m_weights = std::move( S.m_weights );
515  m_spline = std::move( S.m_spline );
516  m_serialization = S.m_serialization;
517  m_handle = S.m_handle;
518  S.m_handle = 0;
519  }
520  return *this;
521  }
522 
527  bool IsValid() const
528  {
529  return m_handle != 0 || !m_spline.IsEmpty();
530  }
531 
535  int Length() const
536  {
537  return m_x.Length();
538  }
539 
545  vector X() const
546  {
547  vector x( m_x.Length() );
548  if ( IsValid() )
549  for ( int i = 0; i < m_x.Length(); ++i )
550  x[i] = m_x0 + m_x[i]/m_r0;
551  return x;
552  }
553 
559  vector Y() const
560  {
561  vector y( m_y.Length() );
562  if ( IsValid() )
563  for ( int i = 0; i < m_y.Length(); ++i )
564  y[i] = m_y0 + m_y[i]/m_r0;
565  return y;
566  }
567 
573  const vector& Z() const
574  {
575  return m_z;
576  }
577 
584  const weight_vector& Weights() const
585  {
586  return m_weights;
587  }
588 
594  rbf_type RBFType() const
595  {
596  return m_rbf;
597  }
598 
607  void SetRBFType( rbf_type rbf )
608  {
609  Clear();
610  m_rbf = rbf;
611  }
612 
621  double ShapeParameter() const
622  {
623  return Sqrt( m_eps2 );
624  }
625 
641  void SetShapeParameter( double eps )
642  {
643  Clear();
644  m_eps = Abs( eps );
645  m_eps2 = m_eps*m_eps;
646  }
647 
651  int Order() const
652  {
653  return m_order;
654  }
655 
679  void SetOrder( int order )
680  {
681  PCL_PRECONDITION( order > 1 )
682  Clear();
683  m_order = pcl::Max( 2, order );
684  }
685 
694  bool IsPolynomialEnabled() const
695  {
696  return m_havePolynomial;
697  }
698 
706  void EnablePolynomial( bool enable = true )
707  {
708  Clear();
709  m_havePolynomial = enable;
710  }
711 
719  void DisablePolynomial( bool disable = true )
720  {
721  EnablePolynomial( !disable );
722  }
723 
728  float Smoothing() const
729  {
730  return m_smoothing;
731  }
732 
750  void SetSmoothing( float s )
751  {
752  PCL_PRECONDITION( s >= 0 )
753  Clear();
754  m_smoothing = pcl::Max( 0.0F, s );
755  }
756 
807  void Initialize( const scalar* x, const scalar* y, const scalar* z, int n, const weight* w = nullptr )
808  {
809  PCL_PRECONDITION( x != nullptr && y != nullptr && z != nullptr )
810  PCL_PRECONDITION( n > 2 )
811 
812  if ( x == nullptr || y == nullptr || z == nullptr )
813  throw Error( "SurfaceSpline::Initialize(): Null vector pointer(s)." );
814 
815  if ( n < 3 )
816  throw Error( "SurfaceSpline::Initialize(): At least three input nodes must be specified." );
817 
818  Clear();
819 
820  if ( m_smoothing <= 0 )
821  w = nullptr;
822 
823  try
824  {
825  /*
826  * Find mean coordinates.
827  */
828  m_x0 = m_y0 = 0;
829  for ( int i = 0; i < n; ++i )
830  {
831  m_x0 += x[i];
832  m_y0 += y[i];
833  }
834  m_x0 /= n;
835  m_y0 /= n;
836 
837  /*
838  * Find radius of largest containing circle.
839  */
840  m_r0 = 0;
841  for ( int i = 0; i < n; ++i )
842  {
843  double dx = x[i] - m_x0;
844  double dy = y[i] - m_y0;
845  double r = Sqrt( dx*dx + dy*dy );
846  if ( r > m_r0 )
847  m_r0 = r;
848  }
849  if ( 1 + m_r0 == 1 )
850  throw Error( "SurfaceSpline::Initialize(): Empty or insignificant interpolation space." );
851  m_r0 = 1/m_r0;
852 
853  /*
854  * If needed and requested, compute an optimal shape parameter.
855  */
856  if ( unlikely( RadialBasisFunction::HasShapeParameter( m_rbf ) ) && m_eps == 0 )
857  {
858  Array<double> R2;
859  for ( int i = 0; i < n; ++i )
860  for ( int j = i; ++j < n; )
861  {
862  double dx = x[i] - x[j];
863  double dy = y[i] - y[j];
864  R2 << dx*dx + dy*dy;
865  }
866  m_eps2 = 1/(m_r0*4*Sqrt( Median( R2.Begin(), R2.End() ) ));
867  }
868  else
869  m_eps2 = m_r0*m_eps;
870  m_eps2 *= m_eps2;
871 
872  /*
873  * Build point list with normalized node coordinates.
874  */
875  Array<NodeData> P;
876  for ( int i = 0; i < n; ++i )
877  P << NodeData( m_r0*(x[i] - m_x0),
878  m_r0*(y[i] - m_y0),
879  z[i],
880  (w != nullptr && w[i] > 0) ? w[i] : 1.0F );
881 
882  /*
883  * Find duplicate input nodes. Two nodes are considered equal if their
884  * (normalized) coordinates don't differ more than the machine epsilon
885  * for the floating point type T.
886  */
887  P.Sort();
888  Array<int> remove;
889  for ( int i = 0, j = 1; j < n; ++i, ++j )
890  if ( Abs( P[i].x - P[j].x ) <= std::numeric_limits<scalar>::epsilon() )
891  if ( Abs( P[i].y - P[j].y ) <= std::numeric_limits<scalar>::epsilon() )
892  remove << i;
893 
894  /*
895  * Build working vectors, excluding duplicate input nodes.
896  */
897  int N = n - int( remove.Length() );
898  if ( N < 3 )
899  throw Error( "SurfaceSpline::Initialize(): Less than three input nodes left after sanitization." );
900  m_x = vector( N );
901  m_y = vector( N );
902  m_z = vector( N );
903  if ( w != nullptr )
904  m_weights = weight_vector( N );
905  int i = 0, k = 0;
906  for ( int j : remove )
907  {
908  for ( ; i < j; ++i, ++k )
909  {
910  m_x[k] = P[i].x;
911  m_y[k] = P[i].y;
912  m_z[k] = P[i].z;
913  if ( w != nullptr )
914  m_weights[k] = P[i].w;
915  }
916  ++i;
917  }
918  for ( ; i < n; ++i, ++k )
919  {
920  m_x[k] = P[i].x;
921  m_y[k] = P[i].y;
922  m_z[k] = P[i].z;
923  if ( w != nullptr )
924  m_weights[k] = P[i].w;
925  }
926 
927  if ( !RadialBasisFunction::HasCoreImplementation( m_rbf ) )
928  m_spline = vector( scalar( 0 ), N + (m_havePolynomial ? ((m_order*(m_order + 1)) >> 1) : 0) );
929 
930  Generate( m_spline.Begin(), &m_handle, m_rbf, m_eps2, m_havePolynomial,
931  m_x.Begin(), m_y.Begin(), m_z.Begin(), N, m_order,
932  m_smoothing, m_weights.Begin() );
933  }
934  catch ( ... )
935  {
936  Clear();
937  throw;
938  }
939  }
940 
945  void Clear()
946  {
947  m_x.Clear();
948  m_y.Clear();
949  m_z.Clear();
950  m_weights.Clear();
951  m_spline.Clear();
952  DestroyHandle( m_handle );
953  m_handle = 0;
954  }
955 
970  {
971  IsoString data;
972  SerializeHandle( data, m_handle );
973  return data;
974  }
975 
985  scalar operator ()( double x, double y ) const
986  {
987  PCL_PRECONDITION( !m_x.IsEmpty() && !m_y.IsEmpty() )
988  PCL_CHECK( m_order >= 2 )
989  PCL_CHECK( !m_spline.IsEmpty() || m_handle != 0 )
990 
991  /*
992  * Normalized interpolation coordinates.
993  */
994  x = m_r0*(x - m_x0);
995  y = m_r0*(y - m_y0);
996 
997  /*
998  * Core RBF implementations
999  */
1000  if ( m_handle != 0 )
1001  return scalar( EvaluateHandle( m_handle, x, y ) );
1002 
1003  /*
1004  * Add polynomial part of the surface spline.
1005  */
1006  int n = m_x.Length();
1007  double z = 0;
1008  if ( m_havePolynomial )
1009  {
1010  z += m_spline[n];
1011  switch ( m_order )
1012  {
1013  case 2:
1014  z += m_spline[n+1]*x + m_spline[n+2]*y;
1015  break;
1016  case 3:
1017  z += (m_spline[n+1] + m_spline[n+3]*x + m_spline[n+4]*y)*x + (m_spline[n+2] + m_spline[n+5]*y)*y;
1018  break;
1019  default:
1020  for ( int i = 1, j = n+1, i1 = (m_order*(m_order + 1))>>1, ix = 0, iy = 0; i < i1; ++i, ++j )
1021  if ( ix == 0 )
1022  {
1023  ix = iy + 1;
1024  iy = 0;
1025  z += m_spline[j] * PowI( x, ix );
1026  }
1027  else
1028  {
1029  --ix;
1030  ++iy;
1031  z += m_spline[j] * PowI( x, ix ) * PowI( y, iy );
1032  }
1033  break;
1034  }
1035  }
1036 
1037  /*
1038  * Add radial basis functions.
1039  */
1040  switch ( m_rbf )
1041  {
1042  case RadialBasisFunction::VariableOrder:
1043  for ( int i = 0; i < n; ++i )
1044  {
1045  double dx = m_x[i] - x;
1046  double dy = m_y[i] - y;
1047  double r2 = dx*dx + dy*dy;
1048  if ( r2 != 0 )
1049  {
1050  double r2m1 = r2;
1051  for ( int j = m_order; --j > 1; )
1052  r2m1 *= r2;
1053  z += m_spline[i] * r2m1 * pcl::Ln( r2 );
1054  }
1055  }
1056  break;
1057  case RadialBasisFunction::ThinPlateSpline:
1058  for ( int i = 0; i < n; ++i )
1059  {
1060  double dx = m_x[i] - x;
1061  double dy = m_y[i] - y;
1062  double r2 = dx*dx + dy*dy;
1063  if ( r2 != 0 )
1064  z += m_spline[i] * r2 * pcl::Ln( Sqrt( r2 ) );
1065  }
1066  break;
1067  case RadialBasisFunction::Gaussian:
1068  for ( int i = 0; i < n; ++i )
1069  {
1070  double dx = m_x[i] - x;
1071  double dy = m_y[i] - y;
1072  z += m_spline[i] * pcl::Exp( -m_eps2 * (dx*dx + dy*dy) );
1073  }
1074  break;
1075  case RadialBasisFunction::Multiquadric:
1076  for ( int i = 0; i < n; ++i )
1077  {
1078  double dx = m_x[i] - x;
1079  double dy = m_y[i] - y;
1080  z += m_spline[i] * pcl::Sqrt( 1 + m_eps2 * (dx*dx + dy*dy) );
1081  }
1082  break;
1083  case RadialBasisFunction::InverseMultiquadric:
1084  for ( int i = 0; i < n; ++i )
1085  {
1086  double dx = m_x[i] - x;
1087  double dy = m_y[i] - y;
1088  z += m_spline[i] / pcl::Sqrt( 1 + m_eps2 * (dx*dx + dy*dy) );
1089  }
1090  break;
1091  case RadialBasisFunction::InverseQuadratic:
1092  for ( int i = 0; i < n; ++i )
1093  {
1094  double dx = m_x[i] - x;
1095  double dy = m_y[i] - y;
1096  z += m_spline[i] / (1 + m_eps2 * (dx*dx + dy*dy));
1097  }
1098  break;
1099  default:
1100  break;
1101  }
1102 
1103  return scalar( z );
1104  }
1105 
1111  template <typename Tp>
1112  scalar operator ()( const GenericPoint<Tp>& p ) const
1113  {
1114  return operator ()( double( p.x ), double( p.y ) );
1115  }
1116 
1131  void Evaluate( scalar* Z, const scalar* X, const scalar* Y, size_type n ) const
1132  {
1133  PCL_PRECONDITION( !m_x.IsEmpty() && !m_y.IsEmpty() )
1134  PCL_PRECONDITION( X != nullptr && Y != nullptr && Z != nullptr )
1135  PCL_PRECONDITION( n > 0 )
1136  PCL_CHECK( m_order >= 2 )
1137  PCL_CHECK( !m_spline.IsEmpty() || m_handle != 0 )
1138 
1139  if ( m_handle != 0 )
1140  EvaluateHandle( m_handle, Z, X, Y, m_x0, m_y0, m_r0, n );
1141  else
1142  for ( size_type i = 0; i < n; ++i )
1143  Z[i] = operator()( X[i], Y[i] );
1144  }
1145 
1154  vector Evaluate( const vector& X, const vector& Y ) const
1155  {
1156  PCL_PRECONDITION( !m_x.IsEmpty() && !m_y.IsEmpty() )
1157  PCL_CHECK( m_order >= 2 )
1158  PCL_CHECK( !m_spline.IsEmpty() || m_handle != 0 )
1159 
1160  int n = Min( X.Length(), Y.Length() );
1161  vector Z( n );
1162  Evaluate( Z.Begin(), X.Begin(), Y.Begin(), size_type( n ) );
1163  return Z;
1164  }
1165 
1176  {
1177  return m_handle != 0;
1178  }
1179 
1180 protected:
1181 
1182  rbf_type m_rbf = RadialBasisFunction::Default;
1183  bool m_havePolynomial = true; // use 1st order polynomial for stabilization
1184  double m_eps = 0; // shape parameter, 0 -> find optimal value automatically
1185  double m_eps2 = 0; // squared shape parameter
1186  vector m_x; // vector of normalized X node coordinates
1187  vector m_y; // vector of normalized Y node coordinates
1188  vector m_z; // vector of function values - for reference only
1189  double m_r0 = 1; // scaling factor for normalization of node coordinates
1190  double m_x0 = 0; // zero offset for normalization of X node coordinates
1191  double m_y0 = 0; // zero offset for normalization of Y node coordinates
1192  int m_order = 2; // derivative order >= 2
1193  float m_smoothing = 0; // <= 0 -> interpolating spline, > 0 -> smoothing factor of approximating spline
1194  weight_vector m_weights; // optional node weights for approximating spline
1195  vector m_spline; // coefficients of the 2-D surface spline, polynomial coeffs. at tail
1196  IsoString m_serialization; // internal use: SplineWorldTransformation for delayed deserialization
1197  void* m_handle = 0; // handle to a core-implemented RBF interpolation/smoothing object
1198 
1204  struct NodeData
1205  {
1206  scalar x, y, z;
1207  weight w;
1208 
1209  NodeData( scalar a_x, scalar a_y, scalar a_z, weight a_w )
1210  : x( a_x ), y( a_y ), z( a_z ), w( a_w )
1211  {
1212  }
1213 
1214  bool operator ==( const NodeData& p ) const
1215  {
1216  return x == p.x && y == p.y;
1217  }
1218 
1219  bool operator <( const NodeData& p ) const
1220  {
1221  return (x != p.x) ? x < p.x : y < p.y;
1222  }
1223  };
1224 
1225  friend class DrizzleData;
1226  friend class DrizzleDataDecoder;
1227  friend class SplineWorldTransformation;
1228 };
1229 
1230 // ----------------------------------------------------------------------------
1231 
1245 class PCL_CLASS PointSurfaceSpline
1246 {
1247 public:
1248 
1253 
1258  using rbf_type = spline::rbf_type;
1259 
1264  PointSurfaceSpline() = default;
1265 
1270 
1275 
1283  template <class point_list1, class point_list2, class weight_vector = FVector>
1284  PointSurfaceSpline( const point_list1& P1, const point_list2& P2,
1285  float smoothness = 0, int order = 2,
1286  const weight_vector& W = weight_vector(),
1287  rbf_type rbf = RadialBasisFunction::Default,
1288  double eps = 0,
1289  bool polynomial = true )
1290  {
1291  Initialize( P1, P2, smoothness, W, order, rbf, eps, polynomial );
1292  }
1293 
1301  PointSurfaceSpline( const spline& Sx, const spline& Sy )
1302  {
1303  Initialize( Sx, Sy );
1304  }
1305 
1309  PointSurfaceSpline& operator =( const PointSurfaceSpline& other ) = default;
1310 
1314  PointSurfaceSpline& operator =( PointSurfaceSpline&& ) = default;
1315 
1393  template <class point_list1, class point_list2, class weight_vector = FVector>
1394  void Initialize( const point_list1& P1, const point_list2& P2,
1395  float smoothness = 0, const weight_vector& W = weight_vector(), int order = 2,
1396  rbf_type rbf = RadialBasisFunction::Default,
1397  double eps = 0,
1398  bool polynomial = true )
1399  {
1400  PCL_PRECONDITION( P1.Length() >= 3 )
1401  PCL_PRECONDITION( P1.Length() <= P2.Length() )
1402  PCL_PRECONDITION( order >= 2 )
1403  PCL_PRECONDITION( W.IsEmpty() || P1.Length() <= size_type( W.Length() ) )
1404 
1405  Clear();
1406 
1407  int N = int( P1.Length() );
1408  if ( N < 3 || P2.Length() < 3 )
1409  throw Error( "PointSurfaceSpline::Initialize(): At least three input nodes must be specified." );
1410  if ( int( P2.Length() ) < N || !W.IsEmpty() && int( W.Length() ) < N )
1411  throw Error( "PointSurfaceSpline::Initialize(): Insufficient data." );
1412 
1413  m_Sx.SetRBFType( rbf );
1414  m_Sx.SetOrder( order );
1415  m_Sx.SetShapeParameter( eps );
1416  m_Sx.EnablePolynomial( polynomial );
1417  m_Sx.SetSmoothing( smoothness );
1418 
1419  m_Sy.SetRBFType( rbf );
1420  m_Sy.SetOrder( order );
1421  m_Sy.SetShapeParameter( eps );
1422  m_Sy.EnablePolynomial( polynomial );
1423  m_Sy.SetSmoothing( smoothness );
1424 
1425  DVector X( N ), Y( N ), ZX( N ), ZY( N );
1426  double zxMin = std::numeric_limits<double>::max();
1427  double zxMax = -std::numeric_limits<double>::max();
1428  double zyMin = std::numeric_limits<double>::max();
1429  double zyMax = -std::numeric_limits<double>::max();
1430  for ( int i = 0; i < N; ++i )
1431  {
1432  X[i] = double( P1[i].x );
1433  Y[i] = double( P1[i].y );
1434 
1435  if ( m_incremental )
1436  {
1437  ZX[i] = X[i];
1438  ZY[i] = Y[i];
1439  m_H.Apply( ZX[i], ZY[i] );
1440  ZX[i] = double( P2[i].x ) - ZX[i];
1441  ZY[i] = double( P2[i].y ) - ZY[i];
1442  }
1443  else
1444  {
1445  ZX[i] = double( P2[i].x );
1446  ZY[i] = double( P2[i].y );
1447  }
1448 
1449  if ( ZX[i] < zxMin )
1450  zxMin = ZX[i];
1451  if ( ZX[i] > zxMax )
1452  zxMax = ZX[i];
1453  if ( ZY[i] < zyMin )
1454  zyMin = ZY[i];
1455  if ( ZY[i] > zyMax )
1456  zyMax = ZY[i];
1457  }
1458 
1459  if ( m_useSimplifiers )
1460  {
1461  SurfaceSimplifier SS;
1462  SS.EnableRejection();
1463  SS.SetRejectFraction( m_simplifierRejectFraction );
1465 
1466  // Simplified surface, X axis.
1467  DVector XXS, YXS, ZXS;
1468  m_epsX = (zxMax - zxMin)/100;
1469  double epsLow = 0, epsHigh = 10*m_epsX;
1470  for ( int i = 0; i < 200; ++i )
1471  {
1472  SS.SetTolerance( m_epsX );
1473  SS.Simplify( XXS, YXS, ZXS, X, Y, ZX );
1474  if ( XXS.Length() > m_maxSplinePoints )
1475  {
1476  epsLow = m_epsX;
1477  if ( epsHigh - epsLow <= 4*std::numeric_limits<double>::epsilon() )
1478  epsHigh *= 2;
1479  }
1480  else
1481  {
1482  if ( epsHigh - epsLow <= 2*std::numeric_limits<double>::epsilon() )
1483  break;
1484  epsHigh = m_epsX;
1485  }
1486  m_epsX = (epsLow + epsHigh)/2;
1487  }
1488  if ( XXS.Length() > m_maxSplinePoints )
1489  {
1490  m_truncatedX = true;
1491  XXS = DVector( XXS.Begin(), m_maxSplinePoints );
1492  YXS = DVector( YXS.Begin(), m_maxSplinePoints );
1493  ZXS = DVector( ZXS.Begin(), m_maxSplinePoints );
1494  }
1495 
1496  // Simplified surface, Y axis.
1497  DVector XYS, YYS, ZYS;
1498  m_epsY = (zyMax - zyMin)/100;
1499  epsLow = 0, epsHigh = 10*m_epsY;
1500  for ( int i = 0; i < 200; ++i )
1501  {
1502  SS.SetTolerance( m_epsY );
1503  SS.Simplify( XYS, YYS, ZYS, X, Y, ZY );
1504  if ( XYS.Length() > m_maxSplinePoints )
1505  {
1506  epsLow = m_epsY;
1507  if ( epsHigh - epsLow <= 4*std::numeric_limits<double>::epsilon() )
1508  epsHigh *= 2;
1509  }
1510  else
1511  {
1512  if ( epsHigh - epsLow <= 2*std::numeric_limits<double>::epsilon() )
1513  break;
1514  epsHigh = m_epsY;
1515  }
1516  m_epsY = (epsLow + epsHigh)/2;
1517  }
1518  if ( XYS.Length() > m_maxSplinePoints )
1519  {
1520  m_truncatedY = true;
1521  XYS = DVector( XYS.Begin(), m_maxSplinePoints );
1522  YYS = DVector( YYS.Begin(), m_maxSplinePoints );
1523  ZYS = DVector( ZYS.Begin(), m_maxSplinePoints );
1524  }
1525 
1526  SplineGenerationThread<weight_vector> Tx( m_Sx, XXS, YXS, ZXS, XXS.Length(), weight_vector() );
1527  SplineGenerationThread<weight_vector> Ty( m_Sy, XYS, YYS, ZYS, XYS.Length(), weight_vector() );
1528  Tx.Start( ThreadPriority::DefaultMax );
1529  Ty.Start( ThreadPriority::DefaultMax );
1530  Tx.Wait();
1531  Ty.Wait();
1532  Tx.ValidateOrThrow();
1533  Ty.ValidateOrThrow();
1534  }
1535  else
1536  {
1537  m_truncatedX = m_truncatedY = N > m_maxSplinePoints;
1538  SplineGenerationThread<weight_vector> Tx( m_Sx, X, Y, ZX, Min( N, m_maxSplinePoints ), W );
1539  SplineGenerationThread<weight_vector> Ty( m_Sy, X, Y, ZY, Min( N, m_maxSplinePoints ), W );
1540  Tx.Start( ThreadPriority::DefaultMax );
1541  Ty.Start( ThreadPriority::DefaultMax );
1542  Tx.Wait();
1543  Ty.Wait();
1544  Tx.ValidateOrThrow();
1545  Ty.ValidateOrThrow();
1546  }
1547  }
1548 
1565  void Initialize( const spline& Sx, const spline& Sy )
1566  {
1567  Clear();
1568  m_useSimplifiers = m_incremental = false;
1569  m_H = Homography();
1570  if ( Sx.IsValid() && Sy.IsValid() )
1571  {
1572  m_Sx = Sx;
1573  m_Sy = Sy;
1574  }
1575  }
1576 
1581  void Clear()
1582  {
1583  m_Sx.Clear();
1584  m_Sy.Clear();
1585  m_epsX = m_epsY = 0;
1586  m_truncatedX = m_truncatedY = false;
1587  }
1588 
1593  bool IsValid() const
1594  {
1595  return m_Sx.IsValid() && m_Sy.IsValid() && (!m_incremental || m_H.IsValid());
1596  }
1597 
1602  const spline& SplineX() const
1603  {
1604  return m_Sx;
1605  }
1606 
1611  const spline& SplineY() const
1612  {
1613  return m_Sy;
1614  }
1615 
1620  int MaxSplinePoints() const
1621  {
1622  return m_maxSplinePoints;
1623  }
1624 
1633  void SetMaxSplinePoints( int n )
1634  {
1635  PCL_PRECONDITION( n >= 3 )
1636  m_maxSplinePoints = Max( 3, n );
1637  }
1638 
1645  bool SimplifiersEnabled() const
1646  {
1647  return m_useSimplifiers;
1648  }
1649 
1657  void EnableSimplifiers( bool enable = true )
1658  {
1659  m_useSimplifiers = enable;
1660  }
1661 
1669  void DisableSimplifiers( bool disable = true )
1670  {
1671  EnableSimplifiers( !disable );
1672  }
1673 
1680  {
1681  return m_simplifierRejectFraction;
1682  }
1683 
1690  void SetSimplifierRejectFraction( float rejectFraction )
1691  {
1692  PCL_PRECONDITION( rejectFraction > 0 && rejectFraction < 1 )
1693  m_simplifierRejectFraction = Range( rejectFraction, 0.0F, 1.0F );
1694  }
1695 
1700  double ErrorX() const
1701  {
1702  return m_epsX;
1703  }
1704 
1709  double ErrorY() const
1710  {
1711  return m_epsY;
1712  }
1713 
1719  bool TruncatedX() const
1720  {
1721  return m_truncatedX;
1722  }
1723 
1729  bool TruncatedY() const
1730  {
1731  return m_truncatedY;
1732  }
1733 
1739  bool Truncated() const
1740  {
1741  return TruncatedX() || TruncatedY();
1742  }
1743 
1751  const Matrix& LinearFunction() const
1752  {
1753  return m_H;
1754  }
1755 
1763  void SetLinearFunction( const Matrix& H )
1764  {
1765  m_H = H;
1766  }
1767 
1782  {
1783  return m_incremental;
1784  }
1785 
1791  void EnableIncrementalFunction( bool enable = true )
1792  {
1793  m_incremental = enable;
1794  }
1795 
1801  void DisableIncrementalFunction( bool disable = true )
1802  {
1803  EnableIncrementalFunction( !disable );
1804  }
1805 
1809  DPoint operator ()( double x, double y ) const
1810  {
1811  PCL_PRECONDITION( ISValid() )
1812  DPoint p( m_Sx( x, y ), m_Sy( x, y ) );
1813  if ( m_incremental )
1814  p += m_H( x, y );
1815  return p;
1816  }
1817 
1821  template <typename T>
1822  DPoint operator ()( const GenericPoint<T>& p ) const
1823  {
1824  return operator ()( double( p.x ), double( p.y ) );
1825  }
1826 
1842  template <typename T>
1843  void Evaluate( T* ZX, T* ZY, const T* X, const T* Y, size_type n ) const
1844  {
1845  PCL_PRECONDITION( ISValid() )
1846  m_Sx.Evaluate( ZX, X, Y, n );
1847  m_Sy.Evaluate( ZY, X, Y, n );
1848  if ( m_incremental )
1849  for ( size_type i = 0; i < n; ++i )
1850  {
1851  DPoint dxy = m_H( X[i], Y[i] );
1852  ZX[i] += dxy.x;
1853  ZY[i] += dxy.y;
1854  }
1855  }
1856 
1868  template <class V>
1869  Array<DPoint> Evaluate( const V& X, const V& Y ) const
1870  {
1871  PCL_PRECONDITION( ISValid() )
1872  size_type n = Min( X.Length(), Y.Length() );
1873  Array<double> ZX( n ), ZY( n );
1874  Evaluate( ZX.Begin(), ZY.Begin(), X.Begin(), Y.Begin(), n );
1875  Array<DPoint> Z( n );
1876  for ( size_type i = 0; i < n; ++i )
1877  {
1878  Z[i].x = ZX[i];
1879  Z[i].y = ZY[i];
1880  }
1881  return Z;
1882  }
1883 
1896  template <class PV>
1897  Array<DPoint> Evaluate( const PV& P ) const
1898  {
1899  PCL_PRECONDITION( ISValid() )
1900  size_type n = P.Length();
1901  Array<double> X( n ), Y( n );
1902  for ( size_type i = 0; i < n; ++i )
1903  {
1904  X[i] = P[i].x;
1905  Y[i] = P[i].y;
1906  }
1907  return Evaluate( X, Y );
1908  }
1909 
1920  {
1921  return m_Sx.HasFastVectorEvaluation() && m_Sy.HasFastVectorEvaluation();
1922  }
1923 
1924 private:
1925 
1926  /*
1927  * The surface splines in the X and Y plane directions.
1928  */
1929  spline m_Sx, m_Sy;
1930  int m_maxSplinePoints = __PCL_PSSPLINE_DEFAULT_MAX_POINTS;
1931 
1932  /*
1933  * Incremental surface splines.
1934  */
1935  bool m_incremental = false; // true => fit differences w.r.t. a projective transformation
1936  Homography m_H; // base projective transformation when m_incremental = true
1937 
1938  /*
1939  * Surface simplification.
1940  */
1941  bool m_useSimplifiers = false;
1942  float m_simplifierRejectFraction = 0.1F;
1943  double m_epsX = 0; // simplification error on the X axis
1944  double m_epsY = 0; // simplification error on the Y axis
1945  bool m_truncatedX = false; // true => truncated vectors in the X axis
1946  bool m_truncatedY = false; // true => truncated vectors in the Y axis
1947 
1948  template <class weight_vector>
1949  class SplineGenerationThread : public Thread
1950  {
1951  public:
1952 
1953  SplineGenerationThread( spline& S, const DVector& X, const DVector& Y, const DVector& Z, int N, const weight_vector& W )
1954  : m_S( S ), m_X( X ), m_Y( Y ), m_Z( Z ), m_N( N ), m_W( W )
1955  {
1956  }
1957 
1958  PCL_HOT_FUNCTION void Run() override
1959  {
1960  try
1961  {
1962  m_S.Initialize( m_X.Begin(), m_Y.Begin(), m_Z.Begin(), m_N, m_W.Begin() );
1963  return;
1964  }
1965  catch ( const Exception& x )
1966  {
1967  m_exception = x;
1968  }
1969  catch ( ... )
1970  {
1971  m_exception = Error( "Unknown exception" );
1972  }
1973  m_S.Clear();
1974  }
1975 
1976  void ValidateOrThrow() const
1977  {
1978  if ( !m_S.IsValid() )
1979  throw m_exception;
1980  }
1981 
1982  private:
1983 
1984  spline& m_S;
1985  const DVector& m_X;
1986  const DVector& m_Y;
1987  const DVector& m_Z;
1988  int m_N;
1989  weight_vector m_W; // ### N.B. do not use a reference: the W ctor. argument can be a temporary object.
1990  Exception m_exception;
1991  };
1992 
1993  friend class DrizzleData;
1994  friend class DrizzleDataDecoder;
1995  friend class SplineWorldTransformation;
1996 };
1997 
1998 // ----------------------------------------------------------------------------
1999 
2018 {
2019 public:
2020 
2026 
2034 
2042 
2047 
2052 
2060  template <class point_list1, class point_list2, class weight_vector = FVector>
2061  RecursivePointSurfaceSpline( const point_list1& P1, const point_list2& P2,
2062  float smoothness = 0,
2063  int order = 2,
2064  const weight_vector& W = weight_vector(),
2065  bool allowExtrapolation = __PCL_RSSPLINE_DEFAULT_ALLOW_EXTRAPOLATION,
2066  int maxSplineLength = __PCL_RSSPLINE_DEFAULT_SPLINE_MAX_LENGTH,
2067  int bucketCapacity = __PCL_RSSPLINE_DEFAULT_TREE_BUCKET_CAPACITY,
2068  bool verbose = true )
2069  {
2070  Initialize( P1, P2, smoothness, W, order, allowExtrapolation, maxSplineLength, bucketCapacity, verbose );
2071  }
2072 
2077  {
2078  Clear();
2079  }
2080 
2157  template <class point_list1, class point_list2, class weight_vector = FVector>
2158  void Initialize( const point_list1& P1, const point_list2& P2,
2159  float smoothness = 0,
2160  const weight_vector& W = weight_vector(),
2161  int order = 2,
2162  bool allowExtrapolation = __PCL_RSSPLINE_DEFAULT_ALLOW_EXTRAPOLATION,
2163  int maxSplineLength = __PCL_RSSPLINE_DEFAULT_SPLINE_MAX_LENGTH,
2164  int bucketCapacity = __PCL_RSSPLINE_DEFAULT_TREE_BUCKET_CAPACITY,
2165  bool verbose = true )
2166  {
2167  PCL_PRECONDITION( P1.Length() >= 3 )
2168  PCL_PRECONDITION( P1.Length() <= P2.Length() )
2169  PCL_PRECONDITION( order >= 2 )
2170  PCL_PRECONDITION( W.IsEmpty() || P1.Length() <= W.Length() )
2171 
2172  Clear();
2173 
2174  if ( P1.Length() < 3 || P2.Length() < 3 )
2175  throw Error( "RecursivePointSurfaceSpline::Initialize(): At least three input nodes must be specified." );
2176 
2177  bool weighted = smoothness > 0 && !W.IsEmpty();
2178 
2179  if ( P1.Length() > P2.Length() || weighted && P1.Length() > size_type( W.Length() ) )
2180  throw Error( "RecursivePointSurfaceSpline::Initialize(): Insufficient data." );
2181 
2182  m_extrapolate = allowExtrapolation;
2183 
2184  if ( P1.Length() <= size_type( maxSplineLength ) )
2185  {
2186  StatusMonitor monitor;
2187 #ifndef __PCL_BUILDING_PIXINSIGHT_APPLICATION
2188  StandardStatus status;
2189  if ( verbose )
2190  {
2191  monitor.SetCallback( &status );
2192  monitor.Initialize( "Building surface subsplines", 100 );
2193  }
2194 #endif
2195  for ( const auto& p : P1 )
2196  {
2197  double x = double( p.x );
2198  double y = double( p.y );
2199  if ( x < m_rect.x0 )
2200  m_rect.x0 = x;
2201  else if ( x > m_rect.x1 )
2202  m_rect.x1 = x;
2203  if ( y < m_rect.y0 )
2204  m_rect.y0 = y;
2205  else if ( y > m_rect.y1 )
2206  m_rect.y1 = y;
2207  }
2208 
2209  m_spline.Initialize( P1, P2, smoothness, W, order,
2210  (order == 2) ? RadialBasisFunction::ThinPlateSpline
2211  : RadialBasisFunction::VariableOrder );
2212 
2213 #ifndef __PCL_BUILDING_PIXINSIGHT_APPLICATION
2214  if ( verbose )
2215  monitor.Complete();
2216 #endif
2217  }
2218  else
2219  {
2220  /*
2221  * Find the total interpolation region.
2222  */
2223  search_rectangle rect = search_coordinate( 0 );
2224  node_list data;
2225  for ( size_type i = 0; i < P1.Length(); ++i )
2226  {
2227  const auto& p1 = P1[i];
2228  const auto& p2 = P2[i];
2229  data << (weighted ? Node( p1, p2, W[i] ) : Node( p1, p2 ));
2230  double x = p1.x;
2231  double y = p1.y;
2232  if ( x < rect.x0 )
2233  rect.x0 = x;
2234  else if ( x > rect.x1 )
2235  rect.x1 = x;
2236  if ( y < rect.y0 )
2237  rect.y0 = y;
2238  else if ( y > rect.y1 )
2239  rect.y1 = y;
2240  }
2241  // Force a square interpolation region for optimal quadtree behavior.
2242  if ( rect.Width() < rect.Height() )
2243  rect.InflateBy( (rect.Height() - rect.Width())/2, search_coordinate( 0 ) );
2244  else
2245  rect.InflateBy( search_coordinate( 0 ), (rect.Width() - rect.Height())/2 );
2246 
2247  /*
2248  * Build the interpolation quadtree.
2249  */
2250  m_tree.Build( rect, data, bucketCapacity );
2251 
2252  /*
2253  * Build subspline interpolation data.
2254  */
2255  Array<SubsplineData> subsplineData;
2256  m_tree.Traverse(
2257  [&]( const search_rectangle& rect, const node_list& points, void*& D )
2258  {
2259  /*
2260  * Ensure redundancy for all subsplines by gathering a set of
2261  * interpolation points in an extended rectangular region around
2262  * each quadtree node.
2263  */
2264  search_rectangle r = rect.InflatedBy( 1.5*Max( rect.Width(), rect.Height() ) );
2265  node_list N = m_tree.Search( r );
2266 
2267  /*
2268  * Sort the cloud of interpolation points by proximity to a
2269  * circle centered at the argument of the minimum RBF value
2270  * (approximately 0.61 in normalized coordinates).
2271  */
2272  if ( N.Length() > size_type( maxSplineLength ) )
2273  {
2274  DPoint c = rect.Center();
2275  double d = 0.61 * (Max( r.Width(), r.Height() )/2 + Max( rect.Width(), rect.Height() )/4);
2276  N.Sort(
2277  [&]( const auto& a, const auto& b )
2278  {
2279  return Abs( a.position.DistanceTo( c ) - d ) < Abs( b.position.DistanceTo( c ) - d );
2280  } );
2281  }
2282 
2283  /*
2284  * Get the optimal subset of at most maxSplineLength redundant
2285  * interpolation points for this quadtree node.
2286  */
2287  Array<DPoint> P1, P2;
2288  Array<float> PW;
2289  for ( int j = 0, l = Min( int( N.Length() ), maxSplineLength ); j < l; ++j )
2290  {
2291  P1 << N[j].position;
2292  P2 << N[j].value;
2293  if ( weighted )
2294  PW << N[j].weight;
2295  }
2296 
2297  subsplineData << SubsplineData( P1, P2, PW, D );
2298  } );
2299 
2300  /*
2301  * Generate all subsplines.
2302  */
2303  StatusMonitor monitor;
2304 #ifndef __PCL_BUILDING_PIXINSIGHT_APPLICATION
2305  StandardStatus status;
2306  if ( verbose )
2307  {
2308  monitor.SetCallback( &status );
2309  monitor.Initialize( "Building recursive surface subsplines", subsplineData.Length() );
2310  }
2311 #endif
2312  Array<size_type> L = Thread::OptimalThreadLoads( subsplineData.Length(),
2313  1/*overheadLimit*/,
2314  m_parallel ? m_maxProcessors : 1 );
2315  AbstractImage::ThreadData threadData( monitor, subsplineData.Length() );
2317  for ( int i = 0, n = 0; i < int( L.Length() ); n += int( L[i++] ) )
2318  threads.Add( new SubsplineGenerationThread( threadData,
2319  subsplineData,
2320  smoothness,
2321  order,
2322  allowExtrapolation,
2323  maxSplineLength,
2324  bucketCapacity,
2325  n,
2326  n + int( L[i] ) ) );
2327  AbstractImage::RunThreads( threads, threadData );
2328  threads.Destroy();
2329  }
2330  }
2331 
2336  void Clear()
2337  {
2338  m_tree.Traverse(
2339  []( const search_rectangle&, node_list&, void*& D )
2340  {
2341  delete reinterpret_cast<RecursivePointSurfaceSpline*&>( D );
2342  } );
2343  m_tree.Clear();
2344  m_spline.Clear();
2345  m_rect = search_coordinate( 0 );
2346  }
2347 
2358  bool IsRecursive() const
2359  {
2360  return !m_tree.IsEmpty();
2361  }
2362 
2367  bool IsValid() const
2368  {
2369  return IsRecursive() || m_spline.IsValid();
2370  }
2371 
2381  DPoint operator ()( double x, double y ) const
2382  {
2383  if ( m_spline.IsValid() )
2384  {
2385  if ( m_extrapolate || m_rect.IncludesFast( x, y ) )
2386  return m_spline( x, y );
2387  }
2388  else
2389  {
2390  const typename tree::Node* node = m_tree.NodeAt( search_point( x, y ) );
2391  if ( node == nullptr )
2392  {
2393  if ( !m_extrapolate )
2394  return 0.0;
2395 
2396  search_rectangle r0 = m_tree.Root()->rect;
2397  if ( x <= r0.x0 )
2398  {
2399  if ( y <= r0.y0 )
2400  node = m_tree.NodeAt( search_point( r0.x0 + SearchDelta, r0.y0 + SearchDelta ) );
2401  else if ( y >= r0.y1 )
2402  node = m_tree.NodeAt( search_point( r0.x0 + SearchDelta, r0.y1 - SearchDelta ) );
2403  else
2404  node = m_tree.NodeAt( search_point( r0.x0 + SearchDelta, search_coordinate( y ) ) );
2405  }
2406  else if ( x >= r0.x1 )
2407  {
2408  if ( y <= r0.y0 )
2409  node = m_tree.NodeAt( search_point( r0.x1 - SearchDelta, r0.y0 + SearchDelta ) );
2410  else if ( y >= r0.y1 )
2411  node = m_tree.NodeAt( search_point( r0.x1 - SearchDelta, r0.y1 - SearchDelta ) );
2412  else
2413  node = m_tree.NodeAt( search_point( r0.x1 - SearchDelta, search_coordinate( y ) ) );
2414  }
2415  else if ( y <= r0.y0 )
2416  node = m_tree.NodeAt( search_point( search_coordinate( x ), r0.y0 + SearchDelta ) );
2417  else // y >= r0.y1
2418  node = m_tree.NodeAt( search_point( search_coordinate( x ), r0.y1 - SearchDelta ) );
2419 
2420  if ( node == nullptr ) // ?!
2421  return 0.0;
2422  }
2423 
2424  if ( node->IsLeaf() )
2425  {
2426  const typename tree::LeafNode* leaf = static_cast<const typename tree::LeafNode*>( node );
2427  if ( leaf->data == nullptr ) // ?!
2428  return 0.0;
2429  return reinterpret_cast<RecursivePointSurfaceSpline*>( leaf->data )->operator()( x, y );
2430  }
2431 
2432  DPoint p = 0.0;
2433  int n = 0;
2434  if ( node->nw != nullptr )
2435  {
2436  const typename tree::LeafNode* leaf;
2437  if ( y <= node->nw->rect.y1 )
2438  leaf = m_tree.LeafNodeAt( search_point( node->nw->rect.x1 - SearchDelta, search_coordinate( y ) ) );
2439  else if ( x <= node->nw->rect.x1 )
2440  leaf = m_tree.LeafNodeAt( search_point( search_coordinate( x ), node->nw->rect.y1 - SearchDelta ) );
2441  else
2442  leaf = m_tree.LeafNodeAt( search_point( node->nw->rect.x1 - SearchDelta, node->nw->rect.y1 - SearchDelta ) );
2443 
2444  if ( leaf != nullptr )
2445  {
2446  p += reinterpret_cast<RecursivePointSurfaceSpline*>( leaf->data )->operator()( x, y );
2447  ++n;
2448  }
2449  }
2450  if ( node->ne != nullptr )
2451  {
2452  const typename tree::LeafNode* leaf;
2453  if ( y <= node->ne->rect.y1 )
2454  leaf = m_tree.LeafNodeAt( search_point( node->ne->rect.x0 + SearchDelta, search_coordinate( y ) ) );
2455  else if ( x <= node->ne->rect.x0 )
2456  leaf = m_tree.LeafNodeAt( search_point( node->ne->rect.x0 + SearchDelta, node->ne->rect.y1 - SearchDelta ) );
2457  else
2458  leaf = m_tree.LeafNodeAt( search_point( search_coordinate( x ), node->ne->rect.y1 - SearchDelta ) );
2459 
2460  if ( leaf != nullptr )
2461  {
2462  p += reinterpret_cast<RecursivePointSurfaceSpline*>( leaf->data )->operator()( x, y );
2463  ++n;
2464  }
2465  }
2466  if ( node->sw != nullptr )
2467  {
2468  const typename tree::LeafNode* leaf;
2469  if ( y >= node->sw->rect.y0 )
2470  leaf = m_tree.LeafNodeAt( search_point( node->sw->rect.x1 - SearchDelta, search_coordinate( y ) ) );
2471  else if ( x <= node->sw->rect.x1 )
2472  leaf = m_tree.LeafNodeAt( search_point( search_coordinate( x ), node->sw->rect.y0 + SearchDelta ) );
2473  else
2474  leaf = m_tree.LeafNodeAt( search_point( node->sw->rect.x1 - SearchDelta, node->sw->rect.y0 + SearchDelta ) );
2475 
2476  if ( leaf != nullptr )
2477  {
2478  p += reinterpret_cast<RecursivePointSurfaceSpline*>( leaf->data )->operator()( x, y );
2479  ++n;
2480  }
2481  }
2482  if ( node->se != nullptr )
2483  {
2484  const typename tree::LeafNode* leaf;
2485  if ( y >= node->se->rect.y0 )
2486  leaf = m_tree.LeafNodeAt( search_point( node->se->rect.x0 + SearchDelta, search_coordinate( y ) ) );
2487  else if ( x <= node->se->rect.x0 )
2488  leaf = m_tree.LeafNodeAt( search_point( node->se->rect.x0 + SearchDelta, node->se->rect.y0 + SearchDelta ) );
2489  else
2490  leaf = m_tree.LeafNodeAt( search_point( search_coordinate( x ), node->se->rect.y0 + SearchDelta ) );
2491 
2492  if ( leaf != nullptr )
2493  {
2494  p += reinterpret_cast<RecursivePointSurfaceSpline*>( leaf->data )->operator()( x, y );
2495  ++n;
2496  }
2497  }
2498 
2499  if ( n > 0 )
2500  return p/double( n );
2501  }
2502 
2503  return 0.0;
2504  }
2505 
2515  template <typename T>
2516  DPoint operator ()( const GenericPoint<T>& p ) const
2517  {
2518  return operator ()( double( p.x ), double( p.y ) );
2519  }
2520 
2521 private:
2522 
2523  /*
2524  * Interpolation data point, QuadTree-compatible.
2525  */
2526  struct Node
2527  {
2528  DPoint position, value;
2529  float weight;
2530 
2531  using component = DPoint::component;
2532 
2533  template <class point1, class point2>
2534  Node( const point1& p, const point2& v )
2535  : position( double( p.x ), double( p.y ) )
2536  , value( double( v.x ), double( v.y ) )
2537  {
2538  }
2539 
2540  template <class point1, class point2>
2541  Node( const point1& p, const point2& v, float w )
2542  : position( double( p.x ), double( p.y ) )
2543  , value( double( v.x ), double( v.y ) )
2544  , weight( w )
2545  {
2546  }
2547 
2548  Node( const Node& ) = default;
2549 
2550  double& operator []( int i )
2551  {
2552  return position[i];
2553  }
2554 
2555  double operator []( int i ) const
2556  {
2557  return position[i];
2558  }
2559  };
2560 
2561  using tree = QuadTree<Node>;
2562  using node_list = typename tree::point_list;
2563  using search_rectangle = typename tree::rectangle;
2564  using search_coordinate = typename tree::coordinate;
2565  using search_point = GenericPoint<search_coordinate>;
2566 
2567  tree m_tree; // the tree of subsplines
2568  PointSurfaceSpline m_spline; // final point spline if there is no further recursion
2569  search_rectangle m_rect = search_coordinate( 0 ); // the interpolation region for this subspline
2570  bool m_extrapolate = __PCL_RSSPLINE_DEFAULT_ALLOW_EXTRAPOLATION;
2571 
2572  static constexpr search_coordinate SearchDelta =
2573  2 * std::numeric_limits<search_coordinate>::epsilon();
2574 
2575  /*
2576  * Parallel subspline generation.
2577  */
2578  struct SubsplineData
2579  {
2580  Array<DPoint> P1, P2;
2581  Array<float> PW;
2582  mutable void** nodeData;
2583 
2584  SubsplineData( const Array<DPoint>& p1, const Array<DPoint>& p2, const Array<float>& pw, void*& nd )
2585  : P1( p1 )
2586  , P2( p2 )
2587  , PW( pw )
2588  , nodeData( &nd )
2589  {
2590  }
2591  };
2592 
2593  class SubsplineGenerationThread : public Thread
2594  {
2595  public:
2596 
2597  SubsplineGenerationThread( const AbstractImage::ThreadData& data,
2598  const Array<SubsplineData>& subsplineData,
2599  float smoothness,
2600  int order,
2601  bool allowExtrapolation,
2602  int maxSplineLength,
2603  int bucketCapacity,
2604  int startIndex, int endIndex )
2605  : m_data( data )
2606  , m_subsplineData( subsplineData )
2607  , m_smoothness( smoothness )
2608  , m_order( order )
2609  , m_allowExtrapolation( allowExtrapolation )
2610  , m_maxSplineLength( maxSplineLength )
2611  , m_bucketCapacity( bucketCapacity )
2612  , m_startIndex( startIndex )
2613  , m_endIndex( endIndex )
2614  {
2615  }
2616 
2617  PCL_HOT_FUNCTION void Run() override
2618  {
2620 
2621  for ( int i = m_startIndex; i < m_endIndex; ++i )
2622  {
2623  const SubsplineData& d = m_subsplineData[i];
2624  AutoPointer<RecursivePointSurfaceSpline> s(
2625  new RecursivePointSurfaceSpline( d.P1, d.P2, m_smoothness, m_order, d.PW,
2626  m_allowExtrapolation,
2627  m_maxSplineLength,
2628  m_bucketCapacity,
2629  false/*verbose*/ ) );
2630  *reinterpret_cast<RecursivePointSurfaceSpline**>( d.nodeData ) = s.Release();
2631 
2633  }
2634 
2635  m_success = true;
2636  }
2637 
2638  operator bool() const
2639  {
2640  return m_success;
2641  }
2642 
2643  private:
2644 
2645  const AbstractImage::ThreadData& m_data;
2646  const Array<SubsplineData>& m_subsplineData;
2647  float m_smoothness;
2648  int m_order;
2649  bool m_allowExtrapolation;
2650  int m_maxSplineLength;
2651  int m_bucketCapacity;
2652  int m_startIndex, m_endIndex;
2653  bool m_success = false;
2654  };
2655 
2656  friend class DrizzleData;
2657  friend class DrizzleDataDecoder;
2658 };
2659 
2660 // ----------------------------------------------------------------------------
2661 
2662 } // pcl
2663 
2664 #endif // __PCL_SurfaceSpline_h
2665 
2666 // ----------------------------------------------------------------------------
2667 // EOF pcl/SurfaceSpline.h - Released 2025-02-21T12:13:32Z
static void RunThreads(ReferenceArray< thread > &threads, ThreadData &data, bool useAffinity=true)
iterator Begin()
Definition: Array.h:438
size_type Length() const noexcept
Definition: Array.h:278
iterator End()
Definition: Array.h:463
64-bit floating point real vector.
Drizzle integration data parser and generator.
Definition: DrizzleData.h:143
A simple exception with an associated error message.
Definition: Exception.h:256
Generic dynamic matrix of arbitrary dimensions.
Definition: Matrix.h:123
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
iterator Begin()
Definition: Vector.h:1918
int Length() const noexcept
Definition: Vector.h:1802
Homography geometric transformation.
Definition: Homography.h:85
Eight-bit string (ISO/IEC-8859-1 or UTF-8 string)
Definition: String.h:5445
A process using multiple concurrent execution threads.
Vector surface spline interpolation/approximation in two dimensions.
const Matrix & LinearFunction() const
void SetLinearFunction(const Matrix &H)
void EnableSimplifiers(bool enable=true)
bool HasFastVectorEvaluation() const
PointSurfaceSpline(PointSurfaceSpline &&)=default
void Initialize(const point_list1 &P1, const point_list2 &P2, float smoothness=0, const weight_vector &W=weight_vector(), int order=2, rbf_type rbf=RadialBasisFunction::Default, double eps=0, bool polynomial=true)
bool SimplifiersEnabled() const
float SimplifierRejectFraction() const
const spline & SplineY() const
void DisableIncrementalFunction(bool disable=true)
const spline & SplineX() const
Array< DPoint > Evaluate(const PV &P) const
PointSurfaceSpline(const PointSurfaceSpline &)=default
void DisableSimplifiers(bool disable=true)
PointSurfaceSpline(const spline &Sx, const spline &Sy)
void Evaluate(T *ZX, T *ZY, const T *X, const T *Y, size_type n) const
void EnableIncrementalFunction(bool enable=true)
void SetSimplifierRejectFraction(float rejectFraction)
Array< DPoint > Evaluate(const V &X, const V &Y) const
PointSurfaceSpline(const point_list1 &P1, const point_list2 &P2, float smoothness=0, int order=2, const weight_vector &W=weight_vector(), rbf_type rbf=RadialBasisFunction::Default, double eps=0, bool polynomial=true)
bool IncrementalFunctionEnabled() const
void SetMaxSplinePoints(int n)
void Initialize(const spline &Sx, const spline &Sy)
Vector surface spline interpolation/approximation in two dimensions with recursive subspline generati...
RecursivePointSurfaceSpline(RecursivePointSurfaceSpline &&)=default
RecursivePointSurfaceSpline(const RecursivePointSurfaceSpline &)=delete
void Initialize(const point_list1 &P1, const point_list2 &P2, float smoothness=0, const weight_vector &W=weight_vector(), int order=2, bool allowExtrapolation=__PCL_RSSPLINE_DEFAULT_ALLOW_EXTRAPOLATION, int maxSplineLength=__PCL_RSSPLINE_DEFAULT_SPLINE_MAX_LENGTH, int bucketCapacity=__PCL_RSSPLINE_DEFAULT_TREE_BUCKET_CAPACITY, bool verbose=true)
RecursivePointSurfaceSpline(const point_list1 &P1, const point_list2 &P2, float smoothness=0, int order=2, const weight_vector &W=weight_vector(), bool allowExtrapolation=__PCL_RSSPLINE_DEFAULT_ALLOW_EXTRAPOLATION, int maxSplineLength=__PCL_RSSPLINE_DEFAULT_SPLINE_MAX_LENGTH, int bucketCapacity=__PCL_RSSPLINE_DEFAULT_TREE_BUCKET_CAPACITY, bool verbose=true)
Dynamic array of pointers to objects providing direct iteration and element access by reference.
void Destroy(iterator i, size_type n=1)
void Add(const ReferenceArray &x)
Surface spline world coordinate transformation.
A status monitoring callback that sends progress information to the process console.
An asynchronous status monitoring system.
void SetCallback(StatusCallback *callback)
void Initialize(const String &info, size_type count=0)
Shape-preserving simplification of 2-D surfaces.
void Simplify(C &xs, C &ys, C &zs, const C &x, const C &y, const C &z) const
void EnableCentroidInclusion(bool enable=true)
void EnableRejection(bool enabled=true)
void SetRejectFraction(float rejectFraction)
void SetTolerance(double tolerance)
Base class of two-dimensional surface splines.
SurfaceSplineBase(SurfaceSplineBase &&)=default
static void DestroyHandle(void *handle)
static void SerializeHandle(IsoString &data, const void *handle)
static void * DeserializeHandle(const IsoString &data)
static double EvaluateHandle(const void *handle, double x, double y)
static void * DuplicateHandle(const void *handle)
RadialBasisFunction::value_type rbf_type
static void EvaluateHandle(const void *handle, float *z, const float *x, const float *y, double x0, double y0, double r, size_type n)
static void Generate(double *__restrict__ c, void **handle, rbf_type, double e2, bool polynomial, const double *__restrict__ x, const double *__restrict__ y, const double *__restrict__ z, int n, int m, float r, const float *__restrict__ w)
static void EvaluateHandle(const void *handle, double *z, const double *x, const double *y, double x0, double y0, double r, size_type n)
static void Generate(float *__restrict__ c, void **handle, rbf_type, double e2, bool polynomial, const float *__restrict__ x, const float *__restrict__ y, const float *__restrict__ z, int n, int m, float r, const float *__restrict__ w)
SurfaceSplineBase(const SurfaceSplineBase &)=default
Two-dimensional interpolating/approximating surface spline.
SurfaceSpline()=default
void DisablePolynomial(bool disable=true)
int Length() const
bool HasFastVectorEvaluation() const
vector X() const
bool IsPolynomialEnabled() const
IsoString CoreSerialization() const
~SurfaceSpline() override
rbf_type RBFType() const
void SetShapeParameter(double eps)
void EnablePolynomial(bool enable=true)
const weight_vector & Weights() const
SurfaceSpline(SurfaceSpline &&S)
void SetRBFType(rbf_type rbf)
vector Y() const
const vector & Z() const
float Smoothing() const
void SetOrder(int order)
bool IsValid() const
SurfaceSpline(const SurfaceSpline &S)
void Initialize(const scalar *x, const scalar *y, const scalar *z, int n, const weight *w=nullptr)
void Evaluate(scalar *Z, const scalar *X, const scalar *Y, size_type n) const
vector Evaluate(const vector &X, const vector &Y) const
void SetSmoothing(float s)
double ShapeParameter() const
Client-side interface to a PixInsight thread.
Definition: Thread.h:233
static Array< size_type > OptimalThreadLoads(size_type count, size_type overheadLimit=1u, int maxThreads=PCL_MAX_PROCESSORS)
bool operator==(const Array< T, A > &x1, const Array< T, A > &x2) noexcept
Definition: Array.h:2285
bool operator<(const Array< T, A > &x1, const Array< T, A > &x2) noexcept
Definition: Array.h:2296
Complex< T > Sqrt(const Complex< T > &c) noexcept
Definition: Complex.h:688
Complex< T > Exp(const Complex< T > &c) noexcept
Definition: Complex.h:728
Complex< T > Ln(const Complex< T > &c) noexcept
Definition: Complex.h:739
T Abs(const Complex< T > &c) noexcept
Definition: Complex.h:443
T PowI(T x, int n) noexcept
Definition: Math.h:1786
size_t size_type
Definition: Defs.h:609
double Median(const T *__restrict__ i, const T *__restrict__ j, double eps=0)
Definition: Math.h:2917
#define INIT_THREAD_MONITOR()
Declares and initializes local variables used for synchronization of thread status monitoring.
#define UPDATE_THREAD_MONITOR(N)
Synchronized status monitoring of a set of image processing threads.
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
Thread synchronization data for status monitoring of parallel image processing tasks.
Auxiliary structure for data sanitization.