PCL
SurfaceSpline.h
Go to the documentation of this file.
1 // ____ ______ __
2 // / __ \ / ____// /
3 // / /_/ // / / /
4 // / ____// /___ / /___ PixInsight Class Library
5 // /_/ \____//_____/ PCL 2.8.5
6 // ----------------------------------------------------------------------------
7 // pcl/SurfaceSpline.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_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 }
201 
202 // ----------------------------------------------------------------------------
203 
208 class PCL_CLASS SurfaceSplineBase
209 {
210 public:
211 
216  using rbf_type = RadialBasisFunction::value_type;
217 
218 protected:
219 
223  SurfaceSplineBase() = default;
224 
229 
234 
239  {
240  }
241 
245  SurfaceSplineBase& operator =( const SurfaceSplineBase& ) = default;
246 
250  SurfaceSplineBase& operator =( SurfaceSplineBase&& ) = default;
251 
255  static void Generate( float* __restrict__ c, void** handle,
256  rbf_type, double e2, bool polynomial,
257  const float* __restrict__ x, const float* __restrict__ y, const float* __restrict__ z,
258  int n, int m, float r, const float* __restrict__ w );
259 
263  static void Generate( double* __restrict__ c, void** handle,
264  rbf_type, double e2, bool polynomial,
265  const double* __restrict__ x, const double* __restrict__ y, const double* __restrict__ z,
266  int n, int m, float r, const float* __restrict__ w );
267 
271  static double EvaluateHandle( const void* handle, double x, double y );
272 
277  static void EvaluateHandle( const void* handle, float* z, const float* x, const float* y, double x0, double y0, double r, size_type n );
278 
283  static void EvaluateHandle( const void* handle, double* z, const double* x, const double* y, double x0, double y0, double r, size_type n );
284 
288  static void SerializeHandle( IsoString& data, const void* handle );
289 
293  static void* DeserializeHandle( const IsoString& data );
294 
298  static void DestroyHandle( void* handle );
299 
303  static void* DuplicateHandle( const void* handle );
304 
305  friend class SplineWorldTransformation;
306 };
307 
308 // ----------------------------------------------------------------------------
309 
357 template <typename T>
358 class PCL_CLASS SurfaceSpline : private SurfaceSplineBase
359 {
360 public:
361 
366  using scalar = T;
367 
373 
378  using weight = float;
379 
384 
389  using rbf_type = SurfaceSplineBase::rbf_type;
390 
395  SurfaceSpline() = default;
396 
401  : m_rbf( S.m_rbf )
402  , m_havePolynomial( S.m_havePolynomial )
403  , m_eps( S.m_eps )
404  , m_eps2( S.m_eps2 )
405  , m_x( S.m_x )
406  , m_y( S.m_y )
407  , m_z( S.m_z )
408  , m_r0( S.m_r0 )
409  , m_x0( S.m_x0 )
410  , m_y0( S.m_y0 )
411  , m_order( S.m_order )
412  , m_smoothing( S.m_smoothing )
413  , m_weights( S.m_weights )
414  , m_spline( S.m_spline )
415  , m_serialization( S.m_serialization )
416  {
417  m_handle = DuplicateHandle( S.m_handle );
418  }
419 
424  : m_rbf( S.m_rbf )
425  , m_havePolynomial( S.m_havePolynomial )
426  , m_eps( S.m_eps )
427  , m_eps2( S.m_eps2 )
428  , m_x( std::move( S.m_x ) )
429  , m_y( std::move( S.m_y ) )
430  , m_z( std::move( S.m_z ) )
431  , m_r0( S.m_r0 )
432  , m_x0( S.m_x0 )
433  , m_y0( S.m_y0 )
434  , m_order( S.m_order )
435  , m_smoothing( S.m_smoothing )
436  , m_weights( std::move( S.m_weights ) )
437  , m_spline( std::move( S.m_spline ) )
438  , m_serialization( S.m_serialization )
439  {
440  m_handle = S.m_handle;
441  S.m_handle = 0;
442  }
443 
447  ~SurfaceSpline() override
448  {
449  DestroyHandle( m_handle );
450  m_handle = 0;
451  }
452 
456  SurfaceSpline& operator =( const SurfaceSpline& S )
457  {
458  if ( &S != this )
459  {
460  Clear();
461  m_rbf = S.m_rbf;
462  m_havePolynomial = S.m_havePolynomial;
463  m_eps = S.m_eps;
464  m_eps2 = S.m_eps2;
465  m_x = S.m_x;
466  m_y = S.m_y;
467  m_z = S.m_z;
468  m_r0 = S.m_r0;
469  m_x0 = S.m_x0;
470  m_y0 = S.m_y0;
471  m_order = S.m_order;
472  m_smoothing = S.m_smoothing;
473  m_weights = S.m_weights;
474  m_spline = S.m_spline;
475  m_serialization = S.m_serialization;
476  m_handle = DuplicateHandle( S.m_handle );
477  }
478  return *this;
479  }
480 
484  SurfaceSpline& operator =( SurfaceSpline&& S )
485  {
486  if ( &S != this )
487  {
488  Clear();
489  m_rbf = S.m_rbf;
490  m_havePolynomial = S.m_havePolynomial;
491  m_eps = S.m_eps;
492  m_eps2 = S.m_eps2;
493  m_x = std::move( S.m_x );
494  m_y = std::move( S.m_y );
495  m_z = std::move( S.m_z );
496  m_r0 = S.m_r0;
497  m_x0 = S.m_x0;
498  m_y0 = S.m_y0;
499  m_order = S.m_order;
500  m_smoothing = S.m_smoothing;
501  m_weights = std::move( S.m_weights );
502  m_spline = std::move( S.m_spline );
503  m_serialization = S.m_serialization;
504  m_handle = S.m_handle;
505  S.m_handle = 0;
506  }
507  return *this;
508  }
509 
514  bool IsValid() const
515  {
516  return m_handle != 0 || !m_spline.IsEmpty();
517  }
518 
522  int Length() const
523  {
524  return m_x.Length();
525  }
526 
532  vector X() const
533  {
534  vector x( m_x.Length() );
535  if ( IsValid() )
536  for ( int i = 0; i < m_x.Length(); ++i )
537  x[i] = m_x0 + m_x[i]/m_r0;
538  return x;
539  }
540 
546  vector Y() const
547  {
548  vector y( m_y.Length() );
549  if ( IsValid() )
550  for ( int i = 0; i < m_y.Length(); ++i )
551  y[i] = m_y0 + m_y[i]/m_r0;
552  return y;
553  }
554 
560  const vector& Z() const
561  {
562  return m_z;
563  }
564 
571  const weight_vector& Weights() const
572  {
573  return m_weights;
574  }
575 
581  rbf_type RBFType() const
582  {
583  return m_rbf;
584  }
585 
594  void SetRBFType( rbf_type rbf )
595  {
596  Clear();
597  m_rbf = rbf;
598  }
599 
608  double ShapeParameter() const
609  {
610  return Sqrt( m_eps2 );
611  }
612 
628  void SetShapeParameter( double eps )
629  {
630  Clear();
631  m_eps = Abs( eps );
632  m_eps2 = m_eps*m_eps;
633  }
634 
638  int Order() const
639  {
640  return m_order;
641  }
642 
666  void SetOrder( int order )
667  {
668  PCL_PRECONDITION( order > 1 )
669  Clear();
670  m_order = pcl::Max( 2, order );
671  }
672 
681  bool IsPolynomialEnabled() const
682  {
683  return m_havePolynomial;
684  }
685 
693  void EnablePolynomial( bool enable = true )
694  {
695  Clear();
696  m_havePolynomial = enable;
697  }
698 
706  void DisablePolynomial( bool disable = true )
707  {
708  EnablePolynomial( !disable );
709  }
710 
715  float Smoothing() const
716  {
717  return m_smoothing;
718  }
719 
737  void SetSmoothing( float s )
738  {
739  PCL_PRECONDITION( s >= 0 )
740  Clear();
741  m_smoothing = pcl::Max( 0.0F, s );
742  }
743 
794  void Initialize( const scalar* x, const scalar* y, const scalar* z, int n, const weight* w = nullptr )
795  {
796  PCL_PRECONDITION( x != nullptr && y != nullptr && z != nullptr )
797  PCL_PRECONDITION( n > 2 )
798 
799  if ( x == nullptr || y == nullptr || z == nullptr )
800  throw Error( "SurfaceSpline::Initialize(): Null vector pointer(s)." );
801 
802  if ( n < 3 )
803  throw Error( "SurfaceSpline::Initialize(): At least three input nodes must be specified." );
804 
805  Clear();
806 
807  if ( m_smoothing <= 0 )
808  w = nullptr;
809 
810  try
811  {
812  /*
813  * Find mean coordinates.
814  */
815  m_x0 = m_y0 = 0;
816  for ( int i = 0; i < n; ++i )
817  {
818  m_x0 += x[i];
819  m_y0 += y[i];
820  }
821  m_x0 /= n;
822  m_y0 /= n;
823 
824  /*
825  * Find radius of largest containing circle.
826  */
827  m_r0 = 0;
828  for ( int i = 0; i < n; ++i )
829  {
830  double dx = x[i] - m_x0;
831  double dy = y[i] - m_y0;
832  double r = Sqrt( dx*dx + dy*dy );
833  if ( r > m_r0 )
834  m_r0 = r;
835  }
836  if ( 1 + m_r0 == 1 )
837  throw Error( "SurfaceSpline::Initialize(): Empty or insignificant interpolation space." );
838  m_r0 = 1/m_r0;
839 
840  /*
841  * If requested, compute an optimal shape parameter.
842  */
843  if ( m_eps == 0 )
844  {
845  Array<double> R2;
846  for ( int i = 0; i < n; ++i )
847  for ( int j = i; ++j < n; )
848  {
849  double dx = x[i] - x[j];
850  double dy = y[i] - y[j];
851  R2 << dx*dx + dy*dy;
852  }
853  m_eps2 = 1/(m_r0*4*Sqrt( Median( R2.Begin(), R2.End() ) ));
854  }
855  else
856  m_eps2 = m_r0*m_eps;
857  m_eps2 *= m_eps2;
858 
859  /*
860  * Build point list with normalized node coordinates.
861  */
862  Array<NodeData> P;
863  for ( int i = 0; i < n; ++i )
864  P << NodeData( m_r0*(x[i] - m_x0),
865  m_r0*(y[i] - m_y0),
866  z[i],
867  (w != nullptr && w[i] > 0) ? w[i] : 1.0F );
868 
869  /*
870  * Find duplicate input nodes. Two nodes are considered equal if their
871  * (normalized) coordinates don't differ more than the machine epsilon
872  * for the floating point type T.
873  */
874  P.Sort();
875  Array<int> remove;
876  for ( int i = 0, j = 1; j < n; ++i, ++j )
877  if ( Abs( P[i].x - P[j].x ) <= std::numeric_limits<scalar>::epsilon() )
878  if ( Abs( P[i].y - P[j].y ) <= std::numeric_limits<scalar>::epsilon() )
879  remove << i;
880 
881  /*
882  * Build working vectors, excluding duplicate input nodes.
883  */
884  int N = n - int( remove.Length() );
885  if ( N < 3 )
886  throw Error( "SurfaceSpline::Initialize(): Less than three input nodes left after sanitization." );
887  m_x = vector( N );
888  m_y = vector( N );
889  m_z = vector( N );
890  if ( w != nullptr )
891  m_weights = weight_vector( N );
892  int i = 0, k = 0;
893  for ( int j : remove )
894  {
895  for ( ; i < j; ++i, ++k )
896  {
897  m_x[k] = P[i].x;
898  m_y[k] = P[i].y;
899  m_z[k] = P[i].z;
900  if ( w != nullptr )
901  m_weights[k] = P[i].w;
902  }
903  ++i;
904  }
905  for ( ; i < n; ++i, ++k )
906  {
907  m_x[k] = P[i].x;
908  m_y[k] = P[i].y;
909  m_z[k] = P[i].z;
910  if ( w != nullptr )
911  m_weights[k] = P[i].w;
912  }
913 
914  if ( !RadialBasisFunction::HasCoreImplementation( m_rbf ) )
915  m_spline = vector( scalar( 0 ), N + (m_havePolynomial ? ((m_order*(m_order + 1)) >> 1) : 0) );
916 
917  Generate( m_spline.Begin(), &m_handle, m_rbf, m_eps2, m_havePolynomial,
918  m_x.Begin(), m_y.Begin(), m_z.Begin(), N, m_order,
919  m_smoothing, m_weights.Begin() );
920  }
921  catch ( ... )
922  {
923  Clear();
924  throw;
925  }
926  }
927 
932  void Clear()
933  {
934  m_x.Clear();
935  m_y.Clear();
936  m_z.Clear();
937  m_weights.Clear();
938  m_spline.Clear();
939  DestroyHandle( m_handle );
940  m_handle = 0;
941  }
942 
957  {
958  IsoString data;
959  SerializeHandle( data, m_handle );
960  return data;
961  }
962 
972  scalar operator ()( double x, double y ) const
973  {
974  PCL_PRECONDITION( !m_x.IsEmpty() && !m_y.IsEmpty() )
975  PCL_CHECK( m_order >= 2 )
976  PCL_CHECK( !m_spline.IsEmpty() || m_handle != 0 )
977 
978  /*
979  * Normalized interpolation coordinates.
980  */
981  x = m_r0*(x - m_x0);
982  y = m_r0*(y - m_y0);
983 
984  /*
985  * Core RBF implementations
986  */
987  if ( m_handle != 0 )
988  return scalar( EvaluateHandle( m_handle, x, y ) );
989 
990  /*
991  * Add polynomial part of the surface spline.
992  */
993  int n = m_x.Length();
994  double z = 0;
995  if ( m_havePolynomial )
996  {
997  z += m_spline[n];
998  switch ( m_order )
999  {
1000  case 2:
1001  z += m_spline[n+1]*x + m_spline[n+2]*y;
1002  break;
1003  case 3:
1004  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;
1005  break;
1006  default:
1007  for ( int i = 1, j = n+1, i1 = (m_order*(m_order + 1))>>1, ix = 0, iy = 0; i < i1; ++i, ++j )
1008  if ( ix == 0 )
1009  {
1010  ix = iy + 1;
1011  iy = 0;
1012  z += m_spline[j] * PowI( x, ix );
1013  }
1014  else
1015  {
1016  --ix;
1017  ++iy;
1018  z += m_spline[j] * PowI( x, ix ) * PowI( y, iy );
1019  }
1020  break;
1021  }
1022  }
1023 
1024  /*
1025  * Add radial basis functions.
1026  */
1027  switch ( m_rbf )
1028  {
1029  case RadialBasisFunction::VariableOrder:
1030  for ( int i = 0; i < n; ++i )
1031  {
1032  double dx = m_x[i] - x;
1033  double dy = m_y[i] - y;
1034  double r2 = dx*dx + dy*dy;
1035  if ( r2 != 0 )
1036  {
1037  double r2m1 = r2;
1038  for ( int j = m_order; --j > 1; )
1039  r2m1 *= r2;
1040  z += m_spline[i] * r2m1 * pcl::Ln( r2 );
1041  }
1042  }
1043  break;
1044  case RadialBasisFunction::ThinPlateSpline:
1045  for ( int i = 0; i < n; ++i )
1046  {
1047  double dx = m_x[i] - x;
1048  double dy = m_y[i] - y;
1049  double r2 = dx*dx + dy*dy;
1050  if ( r2 != 0 )
1051  z += m_spline[i] * r2 * pcl::Ln( Sqrt( r2 ) );
1052  }
1053  break;
1054  case RadialBasisFunction::Gaussian:
1055  for ( int i = 0; i < n; ++i )
1056  {
1057  double dx = m_x[i] - x;
1058  double dy = m_y[i] - y;
1059  z += m_spline[i] * pcl::Exp( -m_eps2 * (dx*dx + dy*dy) );
1060  }
1061  break;
1062  case RadialBasisFunction::Multiquadric:
1063  for ( int i = 0; i < n; ++i )
1064  {
1065  double dx = m_x[i] - x;
1066  double dy = m_y[i] - y;
1067  z += m_spline[i] * pcl::Sqrt( 1 + m_eps2 * (dx*dx + dy*dy) );
1068  }
1069  break;
1070  case RadialBasisFunction::InverseMultiquadric:
1071  for ( int i = 0; i < n; ++i )
1072  {
1073  double dx = m_x[i] - x;
1074  double dy = m_y[i] - y;
1075  z += m_spline[i] / pcl::Sqrt( 1 + m_eps2 * (dx*dx + dy*dy) );
1076  }
1077  break;
1078  case RadialBasisFunction::InverseQuadratic:
1079  for ( int i = 0; i < n; ++i )
1080  {
1081  double dx = m_x[i] - x;
1082  double dy = m_y[i] - y;
1083  z += m_spline[i] / (1 + m_eps2 * (dx*dx + dy*dy));
1084  }
1085  break;
1086  default:
1087  break;
1088  }
1089 
1090  return scalar( z );
1091  }
1092 
1098  template <typename Tp>
1099  scalar operator ()( const GenericPoint<Tp>& p ) const
1100  {
1101  return operator ()( double( p.x ), double( p.y ) );
1102  }
1103 
1118  void Evaluate( scalar* Z, const scalar* X, const scalar* Y, size_type n ) const
1119  {
1120  PCL_PRECONDITION( !m_x.IsEmpty() && !m_y.IsEmpty() )
1121  PCL_PRECONDITION( X != nullptr && Y != nullptr && Z != nullptr )
1122  PCL_PRECONDITION( n > 0 )
1123  PCL_CHECK( m_order >= 2 )
1124  PCL_CHECK( !m_spline.IsEmpty() || m_handle != 0 )
1125 
1126  if ( m_handle != 0 )
1127  EvaluateHandle( m_handle, Z, X, Y, m_x0, m_y0, m_r0, n );
1128  else
1129  for ( size_type i = 0; i < n; ++i )
1130  Z[i] = operator()( X[i], Y[i] );
1131  }
1132 
1141  vector Evaluate( const vector& X, const vector& Y ) const
1142  {
1143  PCL_PRECONDITION( !m_x.IsEmpty() && !m_y.IsEmpty() )
1144  PCL_CHECK( m_order >= 2 )
1145  PCL_CHECK( !m_spline.IsEmpty() || m_handle != 0 )
1146 
1147  int n = Min( X.Length(), Y.Length() );
1148  vector Z( n );
1149  Evaluate( Z.Begin(), X.Begin(), Y.Begin(), size_type( n ) );
1150  return Z;
1151  }
1152 
1163  {
1164  return m_handle != 0;
1165  }
1166 
1167 protected:
1168 
1169  rbf_type m_rbf = RadialBasisFunction::Default;
1170  bool m_havePolynomial = true; // use 1st order polynomial for stabilization
1171  double m_eps = 0; // shape parameter, 0 -> find optimal value automatically
1172  double m_eps2 = 0; // squared shape parameter
1173  vector m_x; // vector of normalized X node coordinates
1174  vector m_y; // vector of normalized Y node coordinates
1175  vector m_z; // vector of function values - for reference only
1176  double m_r0 = 1; // scaling factor for normalization of node coordinates
1177  double m_x0 = 0; // zero offset for normalization of X node coordinates
1178  double m_y0 = 0; // zero offset for normalization of Y node coordinates
1179  int m_order = 2; // derivative order >= 2
1180  float m_smoothing = 0; // <= 0 -> interpolating spline, > 0 -> smoothing factor of approximating spline
1181  weight_vector m_weights; // optional node weights for approximating spline
1182  vector m_spline; // coefficients of the 2-D surface spline, polynomial coeffs. at tail
1183  IsoString m_serialization; // internal use: SplineWorldTransformation for delayed deserialization
1184  void* m_handle = 0; // handle to a core-implemented RBF interpolation/smoothing object
1185 
1191  struct NodeData
1192  {
1193  scalar x, y, z;
1194  weight w;
1195 
1196  NodeData( scalar a_x, scalar a_y, scalar a_z, weight a_w )
1197  : x( a_x ), y( a_y ), z( a_z ), w( a_w )
1198  {
1199  }
1200 
1201  bool operator ==( const NodeData& p ) const
1202  {
1203  return x == p.x && y == p.y;
1204  }
1205 
1206  bool operator <( const NodeData& p ) const
1207  {
1208  return (x != p.x) ? x < p.x : y < p.y;
1209  }
1210  };
1211 
1212  friend class DrizzleData;
1213  friend class DrizzleDataDecoder;
1214  friend class SplineWorldTransformation;
1215 };
1216 
1217 // ----------------------------------------------------------------------------
1218 
1232 class PCL_CLASS PointSurfaceSpline
1233 {
1234 public:
1235 
1240 
1245  using rbf_type = spline::rbf_type;
1246 
1251  PointSurfaceSpline() = default;
1252 
1257 
1262 
1270  template <class point_list1, class point_list2, class weight_vector = FVector>
1271  PointSurfaceSpline( const point_list1& P1, const point_list2& P2,
1272  float smoothness = 0, int order = 2,
1273  const weight_vector& W = weight_vector(),
1274  rbf_type rbf = RadialBasisFunction::Default,
1275  double eps = 0,
1276  bool polynomial = true )
1277  {
1278  Initialize( P1, P2, smoothness, W, order, rbf, eps, polynomial );
1279  }
1280 
1288  PointSurfaceSpline( const spline& Sx, const spline& Sy )
1289  {
1290  Initialize( Sx, Sy );
1291  }
1292 
1296  PointSurfaceSpline& operator =( const PointSurfaceSpline& other ) = default;
1297 
1301  PointSurfaceSpline& operator =( PointSurfaceSpline&& ) = default;
1302 
1380  template <class point_list1, class point_list2, class weight_vector = FVector>
1381  void Initialize( const point_list1& P1, const point_list2& P2,
1382  float smoothness = 0, const weight_vector& W = weight_vector(), int order = 2,
1383  rbf_type rbf = RadialBasisFunction::Default,
1384  double eps = 0,
1385  bool polynomial = true )
1386  {
1387  PCL_PRECONDITION( P1.Length() >= 3 )
1388  PCL_PRECONDITION( P1.Length() <= P2.Length() )
1389  PCL_PRECONDITION( order >= 2 )
1390  PCL_PRECONDITION( W.IsEmpty() || P1.Length() <= size_type( W.Length() ) )
1391 
1392  Clear();
1393 
1394  int N = int( P1.Length() );
1395  if ( N < 3 || P2.Length() < 3 )
1396  throw Error( "PointSurfaceSpline::Initialize(): At least three input nodes must be specified." );
1397  if ( int( P2.Length() ) < N || !W.IsEmpty() && int( W.Length() ) < N )
1398  throw Error( "PointSurfaceSpline::Initialize(): Insufficient data." );
1399 
1400  m_Sx.SetRBFType( rbf );
1401  m_Sx.SetOrder( order );
1402  m_Sx.SetShapeParameter( eps );
1403  m_Sx.EnablePolynomial( polynomial );
1404  m_Sx.SetSmoothing( smoothness );
1405 
1406  m_Sy.SetRBFType( rbf );
1407  m_Sy.SetOrder( order );
1408  m_Sy.SetShapeParameter( eps );
1409  m_Sy.EnablePolynomial( polynomial );
1410  m_Sy.SetSmoothing( smoothness );
1411 
1412  DVector X( N ), Y( N ), ZX( N ), ZY( N );
1413  double zxMin = std::numeric_limits<double>::max();
1414  double zxMax = -std::numeric_limits<double>::max();
1415  double zyMin = std::numeric_limits<double>::max();
1416  double zyMax = -std::numeric_limits<double>::max();
1417  for ( int i = 0; i < N; ++i )
1418  {
1419  X[i] = double( P1[i].x );
1420  Y[i] = double( P1[i].y );
1421 
1422  if ( m_incremental )
1423  {
1424  ZX[i] = X[i];
1425  ZY[i] = Y[i];
1426  m_H.Apply( ZX[i], ZY[i] );
1427  ZX[i] = double( P2[i].x ) - ZX[i];
1428  ZY[i] = double( P2[i].y ) - ZY[i];
1429  }
1430  else
1431  {
1432  ZX[i] = double( P2[i].x );
1433  ZY[i] = double( P2[i].y );
1434  }
1435 
1436  if ( ZX[i] < zxMin )
1437  zxMin = ZX[i];
1438  if ( ZX[i] > zxMax )
1439  zxMax = ZX[i];
1440  if ( ZY[i] < zyMin )
1441  zyMin = ZY[i];
1442  if ( ZY[i] > zyMax )
1443  zyMax = ZY[i];
1444  }
1445 
1446  if ( m_useSimplifiers )
1447  {
1448  SurfaceSimplifier SS;
1449  SS.EnableRejection();
1450  SS.SetRejectFraction( m_simplifierRejectFraction );
1452 
1453  // Simplified surface, X axis.
1454  DVector XXS, YXS, ZXS;
1455  m_epsX = (zxMax - zxMin)/100;
1456  double epsLow = 0, epsHigh = 10*m_epsX;
1457  for ( int i = 0; i < 200; ++i )
1458  {
1459  SS.SetTolerance( m_epsX );
1460  SS.Simplify( XXS, YXS, ZXS, X, Y, ZX );
1461  if ( XXS.Length() > m_maxSplinePoints )
1462  {
1463  epsLow = m_epsX;
1464  if ( epsHigh - epsLow <= 4*std::numeric_limits<double>::epsilon() )
1465  epsHigh *= 2;
1466  }
1467  else
1468  {
1469  if ( epsHigh - epsLow <= 2*std::numeric_limits<double>::epsilon() )
1470  break;
1471  epsHigh = m_epsX;
1472  }
1473  m_epsX = (epsLow + epsHigh)/2;
1474  }
1475  if ( XXS.Length() > m_maxSplinePoints )
1476  {
1477  m_truncatedX = true;
1478  XXS = DVector( XXS.Begin(), m_maxSplinePoints );
1479  YXS = DVector( YXS.Begin(), m_maxSplinePoints );
1480  ZXS = DVector( ZXS.Begin(), m_maxSplinePoints );
1481  }
1482 
1483  // Simplified surface, Y axis.
1484  DVector XYS, YYS, ZYS;
1485  m_epsY = (zyMax - zyMin)/100;
1486  epsLow = 0, epsHigh = 10*m_epsY;
1487  for ( int i = 0; i < 200; ++i )
1488  {
1489  SS.SetTolerance( m_epsY );
1490  SS.Simplify( XYS, YYS, ZYS, X, Y, ZY );
1491  if ( XYS.Length() > m_maxSplinePoints )
1492  {
1493  epsLow = m_epsY;
1494  if ( epsHigh - epsLow <= 4*std::numeric_limits<double>::epsilon() )
1495  epsHigh *= 2;
1496  }
1497  else
1498  {
1499  if ( epsHigh - epsLow <= 2*std::numeric_limits<double>::epsilon() )
1500  break;
1501  epsHigh = m_epsY;
1502  }
1503  m_epsY = (epsLow + epsHigh)/2;
1504  }
1505  if ( XYS.Length() > m_maxSplinePoints )
1506  {
1507  m_truncatedY = true;
1508  XYS = DVector( XYS.Begin(), m_maxSplinePoints );
1509  YYS = DVector( YYS.Begin(), m_maxSplinePoints );
1510  ZYS = DVector( ZYS.Begin(), m_maxSplinePoints );
1511  }
1512 
1513  SplineGenerationThread<weight_vector> Tx( m_Sx, XXS, YXS, ZXS, XXS.Length(), weight_vector() );
1514  SplineGenerationThread<weight_vector> Ty( m_Sy, XYS, YYS, ZYS, XYS.Length(), weight_vector() );
1515  Tx.Start( ThreadPriority::DefaultMax );
1516  Ty.Start( ThreadPriority::DefaultMax );
1517  Tx.Wait();
1518  Ty.Wait();
1519  Tx.ValidateOrThrow();
1520  Ty.ValidateOrThrow();
1521  }
1522  else
1523  {
1524  m_truncatedX = m_truncatedY = N > m_maxSplinePoints;
1525  SplineGenerationThread<weight_vector> Tx( m_Sx, X, Y, ZX, Min( N, m_maxSplinePoints ), W );
1526  SplineGenerationThread<weight_vector> Ty( m_Sy, X, Y, ZY, Min( N, m_maxSplinePoints ), W );
1527  Tx.Start( ThreadPriority::DefaultMax );
1528  Ty.Start( ThreadPriority::DefaultMax );
1529  Tx.Wait();
1530  Ty.Wait();
1531  Tx.ValidateOrThrow();
1532  Ty.ValidateOrThrow();
1533  }
1534  }
1535 
1552  void Initialize( const spline& Sx, const spline& Sy )
1553  {
1554  Clear();
1555  m_useSimplifiers = m_incremental = false;
1556  m_H = Homography();
1557  if ( Sx.IsValid() && Sy.IsValid() )
1558  {
1559  m_Sx = Sx;
1560  m_Sy = Sy;
1561  }
1562  }
1563 
1568  void Clear()
1569  {
1570  m_Sx.Clear();
1571  m_Sy.Clear();
1572  m_epsX = m_epsY = 0;
1573  m_truncatedX = m_truncatedY = false;
1574  }
1575 
1580  bool IsValid() const
1581  {
1582  return m_Sx.IsValid() && m_Sy.IsValid() && (!m_incremental || m_H.IsValid());
1583  }
1584 
1589  const spline& SplineX() const
1590  {
1591  return m_Sx;
1592  }
1593 
1598  const spline& SplineY() const
1599  {
1600  return m_Sy;
1601  }
1602 
1607  int MaxSplinePoints() const
1608  {
1609  return m_maxSplinePoints;
1610  }
1611 
1620  void SetMaxSplinePoints( int n )
1621  {
1622  PCL_PRECONDITION( n >= 3 )
1623  m_maxSplinePoints = Max( 3, n );
1624  }
1625 
1632  bool SimplifiersEnabled() const
1633  {
1634  return m_useSimplifiers;
1635  }
1636 
1644  void EnableSimplifiers( bool enable = true )
1645  {
1646  m_useSimplifiers = enable;
1647  }
1648 
1656  void DisableSimplifiers( bool disable = true )
1657  {
1658  EnableSimplifiers( !disable );
1659  }
1660 
1667  {
1668  return m_simplifierRejectFraction;
1669  }
1670 
1677  void SetSimplifierRejectFraction( float rejectFraction )
1678  {
1679  PCL_PRECONDITION( rejectFraction > 0 && rejectFraction < 1 )
1680  m_simplifierRejectFraction = Range( rejectFraction, 0.0F, 1.0F );
1681  }
1682 
1687  double ErrorX() const
1688  {
1689  return m_epsX;
1690  }
1691 
1696  double ErrorY() const
1697  {
1698  return m_epsY;
1699  }
1700 
1706  bool TruncatedX() const
1707  {
1708  return m_truncatedX;
1709  }
1710 
1716  bool TruncatedY() const
1717  {
1718  return m_truncatedY;
1719  }
1720 
1726  bool Truncated() const
1727  {
1728  return TruncatedX() || TruncatedY();
1729  }
1730 
1738  const Matrix& LinearFunction() const
1739  {
1740  return m_H;
1741  }
1742 
1750  void SetLinearFunction( const Matrix& H )
1751  {
1752  m_H = H;
1753  }
1754 
1769  {
1770  return m_incremental;
1771  }
1772 
1778  void EnableIncrementalFunction( bool enable = true )
1779  {
1780  m_incremental = enable;
1781  }
1782 
1788  void DisableIncrementalFunction( bool disable = true )
1789  {
1790  EnableIncrementalFunction( !disable );
1791  }
1792 
1796  DPoint operator ()( double x, double y ) const
1797  {
1798  PCL_PRECONDITION( ISValid() )
1799  DPoint p( m_Sx( x, y ), m_Sy( x, y ) );
1800  if ( m_incremental )
1801  p += m_H( x, y );
1802  return p;
1803  }
1804 
1808  template <typename T>
1809  DPoint operator ()( const GenericPoint<T>& p ) const
1810  {
1811  return operator ()( double( p.x ), double( p.y ) );
1812  }
1813 
1829  template <typename T>
1830  void Evaluate( T* ZX, T* ZY, const T* X, const T* Y, size_type n ) const
1831  {
1832  PCL_PRECONDITION( ISValid() )
1833  m_Sx.Evaluate( ZX, X, Y, n );
1834  m_Sy.Evaluate( ZY, X, Y, n );
1835  if ( m_incremental )
1836  for ( size_type i = 0; i < n; ++i )
1837  {
1838  DPoint dxy = m_H( X[i], Y[i] );
1839  ZX[i] += dxy.x;
1840  ZY[i] += dxy.y;
1841  }
1842  }
1843 
1855  template <class V>
1856  Array<DPoint> Evaluate( const V& X, const V& Y ) const
1857  {
1858  PCL_PRECONDITION( ISValid() )
1859  size_type n = Min( X.Length(), Y.Length() );
1860  Array<double> ZX( n ), ZY( n );
1861  Evaluate( ZX.Begin(), ZY.Begin(), X.Begin(), Y.Begin(), n );
1862  Array<DPoint> Z( n );
1863  for ( size_type i = 0; i < n; ++i )
1864  {
1865  Z[i].x = ZX[i];
1866  Z[i].y = ZY[i];
1867  }
1868  return Z;
1869  }
1870 
1883  template <class PV>
1884  Array<DPoint> Evaluate( const PV& P ) const
1885  {
1886  PCL_PRECONDITION( ISValid() )
1887  size_type n = P.Length();
1888  Array<double> X( n ), Y( n );
1889  for ( size_type i = 0; i < n; ++i )
1890  {
1891  X[i] = P[i].x;
1892  Y[i] = P[i].y;
1893  }
1894  return Evaluate( X, Y );
1895  }
1896 
1907  {
1908  return m_Sx.HasFastVectorEvaluation() && m_Sy.HasFastVectorEvaluation();
1909  }
1910 
1911 private:
1912 
1913  /*
1914  * The surface splines in the X and Y plane directions.
1915  */
1916  spline m_Sx, m_Sy;
1917  int m_maxSplinePoints = __PCL_PSSPLINE_DEFAULT_MAX_POINTS;
1918 
1919  /*
1920  * Incremental surface splines.
1921  */
1922  bool m_incremental = false; // true => fit differences w.r.t. a projective transformation
1923  Homography m_H; // base projective transformation when m_incremental = true
1924 
1925  /*
1926  * Surface simplification.
1927  */
1928  bool m_useSimplifiers = false;
1929  float m_simplifierRejectFraction = 0.1F;
1930  double m_epsX = 0; // simplification error on the X axis
1931  double m_epsY = 0; // simplification error on the Y axis
1932  bool m_truncatedX = false; // true => truncated vectors in the X axis
1933  bool m_truncatedY = false; // true => truncated vectors in the Y axis
1934 
1935  template <class weight_vector>
1936  class SplineGenerationThread : public Thread
1937  {
1938  public:
1939 
1940  SplineGenerationThread( spline& S, const DVector& X, const DVector& Y, const DVector& Z, int N, const weight_vector& W )
1941  : m_S( S ), m_X( X ), m_Y( Y ), m_Z( Z ), m_N( N ), m_W( W )
1942  {
1943  }
1944 
1945  PCL_HOT_FUNCTION void Run() override
1946  {
1947  try
1948  {
1949  m_S.Initialize( m_X.Begin(), m_Y.Begin(), m_Z.Begin(), m_N, m_W.Begin() );
1950  return;
1951  }
1952  catch ( const Exception& x )
1953  {
1954  m_exception = x;
1955  }
1956  catch ( ... )
1957  {
1958  m_exception = Error( "Unknown exception" );
1959  }
1960  m_S.Clear();
1961  }
1962 
1963  void ValidateOrThrow() const
1964  {
1965  if ( !m_S.IsValid() )
1966  throw m_exception;
1967  }
1968 
1969  private:
1970 
1971  spline& m_S;
1972  const DVector& m_X;
1973  const DVector& m_Y;
1974  const DVector& m_Z;
1975  int m_N;
1976  weight_vector m_W; // ### N.B. do not use a reference: the W ctor. argument can be a temporary object.
1977  Exception m_exception;
1978  };
1979 
1980  friend class DrizzleData;
1981  friend class DrizzleDataDecoder;
1982  friend class SplineWorldTransformation;
1983 };
1984 
1985 // ----------------------------------------------------------------------------
1986 
2005 {
2006 public:
2007 
2013 
2021 
2029 
2034 
2039 
2047  template <class point_list1, class point_list2, class weight_vector = FVector>
2048  RecursivePointSurfaceSpline( const point_list1& P1, const point_list2& P2,
2049  float smoothness = 0,
2050  int order = 2,
2051  const weight_vector& W = weight_vector(),
2052  bool allowExtrapolation = __PCL_RSSPLINE_DEFAULT_ALLOW_EXTRAPOLATION,
2053  int maxSplineLength = __PCL_RSSPLINE_DEFAULT_SPLINE_MAX_LENGTH,
2054  int bucketCapacity = __PCL_RSSPLINE_DEFAULT_TREE_BUCKET_CAPACITY,
2055  bool verbose = true )
2056  {
2057  Initialize( P1, P2, smoothness, W, order, allowExtrapolation, maxSplineLength, bucketCapacity, verbose );
2058  }
2059 
2064  {
2065  Clear();
2066  }
2067 
2144  template <class point_list1, class point_list2, class weight_vector = FVector>
2145  void Initialize( const point_list1& P1, const point_list2& P2,
2146  float smoothness = 0,
2147  const weight_vector& W = weight_vector(),
2148  int order = 2,
2149  bool allowExtrapolation = __PCL_RSSPLINE_DEFAULT_ALLOW_EXTRAPOLATION,
2150  int maxSplineLength = __PCL_RSSPLINE_DEFAULT_SPLINE_MAX_LENGTH,
2151  int bucketCapacity = __PCL_RSSPLINE_DEFAULT_TREE_BUCKET_CAPACITY,
2152  bool verbose = true )
2153  {
2154  PCL_PRECONDITION( P1.Length() >= 3 )
2155  PCL_PRECONDITION( P1.Length() <= P2.Length() )
2156  PCL_PRECONDITION( order >= 2 )
2157  PCL_PRECONDITION( W.IsEmpty() || P1.Length() <= W.Length() )
2158 
2159  Clear();
2160 
2161  if ( P1.Length() < 3 || P2.Length() < 3 )
2162  throw Error( "RecursivePointSurfaceSpline::Initialize(): At least three input nodes must be specified." );
2163 
2164  bool weighted = smoothness > 0 && !W.IsEmpty();
2165 
2166  if ( P1.Length() > P2.Length() || weighted && P1.Length() > size_type( W.Length() ) )
2167  throw Error( "RecursivePointSurfaceSpline::Initialize(): Insufficient data." );
2168 
2169  m_extrapolate = allowExtrapolation;
2170 
2171  if ( P1.Length() <= size_type( maxSplineLength ) )
2172  {
2173  StatusMonitor monitor;
2174 #ifndef __PCL_BUILDING_PIXINSIGHT_APPLICATION
2175  StandardStatus status;
2176  if ( verbose )
2177  {
2178  monitor.SetCallback( &status );
2179  monitor.Initialize( "Building surface subsplines", 100 );
2180  }
2181 #endif
2182  for ( const auto& p : P1 )
2183  {
2184  double x = double( p.x );
2185  double y = double( p.y );
2186  if ( x < m_rect.x0 )
2187  m_rect.x0 = x;
2188  else if ( x > m_rect.x1 )
2189  m_rect.x1 = x;
2190  if ( y < m_rect.y0 )
2191  m_rect.y0 = y;
2192  else if ( y > m_rect.y1 )
2193  m_rect.y1 = y;
2194  }
2195 
2196  m_spline.Initialize( P1, P2, smoothness, W, order,
2197  (order == 2) ? RadialBasisFunction::ThinPlateSpline
2198  : RadialBasisFunction::VariableOrder );
2199 
2200 #ifndef __PCL_BUILDING_PIXINSIGHT_APPLICATION
2201  if ( verbose )
2202  monitor.Complete();
2203 #endif
2204  }
2205  else
2206  {
2207  /*
2208  * Find the total interpolation region.
2209  */
2210  search_rectangle rect = search_coordinate( 0 );
2211  node_list data;
2212  for ( size_type i = 0; i < P1.Length(); ++i )
2213  {
2214  const auto& p1 = P1[i];
2215  const auto& p2 = P2[i];
2216  data << (weighted ? Node( p1, p2, W[i] ) : Node( p1, p2 ));
2217  double x = p1.x;
2218  double y = p1.y;
2219  if ( x < rect.x0 )
2220  rect.x0 = x;
2221  else if ( x > rect.x1 )
2222  rect.x1 = x;
2223  if ( y < rect.y0 )
2224  rect.y0 = y;
2225  else if ( y > rect.y1 )
2226  rect.y1 = y;
2227  }
2228  // Force a square interpolation region for optimal quadtree behavior.
2229  if ( rect.Width() < rect.Height() )
2230  rect.InflateBy( (rect.Height() - rect.Width())/2, search_coordinate( 0 ) );
2231  else
2232  rect.InflateBy( search_coordinate( 0 ), (rect.Width() - rect.Height())/2 );
2233 
2234  /*
2235  * Build the interpolation quadtree.
2236  */
2237  m_tree.Build( rect, data, bucketCapacity );
2238 
2239  /*
2240  * Build subspline interpolation data.
2241  */
2242  Array<SubsplineData> subsplineData;
2243  m_tree.Traverse(
2244  [&]( const search_rectangle& rect, const node_list& points, void*& D )
2245  {
2246  /*
2247  * Ensure redundancy for all subsplines by gathering a set of
2248  * interpolation points in an extended rectangular region around
2249  * each quadtree node.
2250  */
2251  search_rectangle r = rect.InflatedBy( 1.5*Max( rect.Width(), rect.Height() ) );
2252  node_list N = m_tree.Search( r );
2253 
2254  /*
2255  * Sort the cloud of interpolation points by proximity to a
2256  * circle centered at the argument of the minimum RBF value
2257  * (approximately 0.61 in normalized coordinates).
2258  */
2259  if ( N.Length() > size_type( maxSplineLength ) )
2260  {
2261  DPoint c = rect.Center();
2262  double d = 0.61 * (Max( r.Width(), r.Height() )/2 + Max( rect.Width(), rect.Height() )/4);
2263  N.Sort(
2264  [&]( const auto& a, const auto& b )
2265  {
2266  return Abs( a.position.DistanceTo( c ) - d ) < Abs( b.position.DistanceTo( c ) - d );
2267  } );
2268  }
2269 
2270  /*
2271  * Get the optimal subset of at most maxSplineLength redundant
2272  * interpolation points for this quadtree node.
2273  */
2274  Array<DPoint> P1, P2;
2275  Array<float> PW;
2276  for ( int j = 0, l = Min( int( N.Length() ), maxSplineLength ); j < l; ++j )
2277  {
2278  P1 << N[j].position;
2279  P2 << N[j].value;
2280  if ( weighted )
2281  PW << N[j].weight;
2282  }
2283 
2284  subsplineData << SubsplineData( P1, P2, PW, D );
2285  } );
2286 
2287  /*
2288  * Generate all subsplines.
2289  */
2290  StatusMonitor monitor;
2291 #ifndef __PCL_BUILDING_PIXINSIGHT_APPLICATION
2292  StandardStatus status;
2293  if ( verbose )
2294  {
2295  monitor.SetCallback( &status );
2296  monitor.Initialize( "Building recursive surface subsplines", subsplineData.Length() );
2297  }
2298 #endif
2299  Array<size_type> L = Thread::OptimalThreadLoads( subsplineData.Length(),
2300  1/*overheadLimit*/,
2301  m_parallel ? m_maxProcessors : 1 );
2302  AbstractImage::ThreadData threadData( monitor, subsplineData.Length() );
2304  for ( int i = 0, n = 0; i < int( L.Length() ); n += int( L[i++] ) )
2305  threads.Add( new SubsplineGenerationThread( threadData,
2306  subsplineData,
2307  smoothness,
2308  order,
2309  allowExtrapolation,
2310  maxSplineLength,
2311  bucketCapacity,
2312  n,
2313  n + int( L[i] ) ) );
2314  AbstractImage::RunThreads( threads, threadData );
2315  threads.Destroy();
2316  }
2317  }
2318 
2323  void Clear()
2324  {
2325  m_tree.Traverse(
2326  []( const search_rectangle&, node_list&, void*& D )
2327  {
2328  delete reinterpret_cast<RecursivePointSurfaceSpline*&>( D );
2329  } );
2330  m_tree.Clear();
2331  m_spline.Clear();
2332  m_rect = search_coordinate( 0 );
2333  }
2334 
2345  bool IsRecursive() const
2346  {
2347  return !m_tree.IsEmpty();
2348  }
2349 
2354  bool IsValid() const
2355  {
2356  return IsRecursive() || m_spline.IsValid();
2357  }
2358 
2368  DPoint operator ()( double x, double y ) const
2369  {
2370  if ( m_spline.IsValid() )
2371  {
2372  if ( m_extrapolate || m_rect.IncludesFast( x, y ) )
2373  return m_spline( x, y );
2374  }
2375  else
2376  {
2377  const typename tree::Node* node = m_tree.NodeAt( search_point( x, y ) );
2378  if ( node == nullptr )
2379  {
2380  if ( !m_extrapolate )
2381  return 0.0;
2382 
2383  search_rectangle r0 = m_tree.Root()->rect;
2384  if ( x <= r0.x0 )
2385  {
2386  if ( y <= r0.y0 )
2387  node = m_tree.NodeAt( search_point( r0.x0 + SearchDelta, r0.y0 + SearchDelta ) );
2388  else if ( y >= r0.y1 )
2389  node = m_tree.NodeAt( search_point( r0.x0 + SearchDelta, r0.y1 - SearchDelta ) );
2390  else
2391  node = m_tree.NodeAt( search_point( r0.x0 + SearchDelta, search_coordinate( y ) ) );
2392  }
2393  else if ( x >= r0.x1 )
2394  {
2395  if ( y <= r0.y0 )
2396  node = m_tree.NodeAt( search_point( r0.x1 - SearchDelta, r0.y0 + SearchDelta ) );
2397  else if ( y >= r0.y1 )
2398  node = m_tree.NodeAt( search_point( r0.x1 - SearchDelta, r0.y1 - SearchDelta ) );
2399  else
2400  node = m_tree.NodeAt( search_point( r0.x1 - SearchDelta, search_coordinate( y ) ) );
2401  }
2402  else if ( y <= r0.y0 )
2403  node = m_tree.NodeAt( search_point( search_coordinate( x ), r0.y0 + SearchDelta ) );
2404  else // y >= r0.y1
2405  node = m_tree.NodeAt( search_point( search_coordinate( x ), r0.y1 - SearchDelta ) );
2406 
2407  if ( node == nullptr ) // ?!
2408  return 0.0;
2409  }
2410 
2411  if ( node->IsLeaf() )
2412  {
2413  const typename tree::LeafNode* leaf = static_cast<const typename tree::LeafNode*>( node );
2414  if ( leaf->data == nullptr ) // ?!
2415  return 0.0;
2416  return reinterpret_cast<RecursivePointSurfaceSpline*>( leaf->data )->operator()( x, y );
2417  }
2418 
2419  DPoint p = 0.0;
2420  int n = 0;
2421  if ( node->nw != nullptr )
2422  {
2423  const typename tree::LeafNode* leaf;
2424  if ( y <= node->nw->rect.y1 )
2425  leaf = m_tree.LeafNodeAt( search_point( node->nw->rect.x1 - SearchDelta, search_coordinate( y ) ) );
2426  else if ( x <= node->nw->rect.x1 )
2427  leaf = m_tree.LeafNodeAt( search_point( search_coordinate( x ), node->nw->rect.y1 - SearchDelta ) );
2428  else
2429  leaf = m_tree.LeafNodeAt( search_point( node->nw->rect.x1 - SearchDelta, node->nw->rect.y1 - SearchDelta ) );
2430 
2431  if ( leaf != nullptr )
2432  {
2433  p += reinterpret_cast<RecursivePointSurfaceSpline*>( leaf->data )->operator()( x, y );
2434  ++n;
2435  }
2436  }
2437  if ( node->ne != nullptr )
2438  {
2439  const typename tree::LeafNode* leaf;
2440  if ( y <= node->ne->rect.y1 )
2441  leaf = m_tree.LeafNodeAt( search_point( node->ne->rect.x0 + SearchDelta, search_coordinate( y ) ) );
2442  else if ( x <= node->ne->rect.x0 )
2443  leaf = m_tree.LeafNodeAt( search_point( node->ne->rect.x0 + SearchDelta, node->ne->rect.y1 - SearchDelta ) );
2444  else
2445  leaf = m_tree.LeafNodeAt( search_point( search_coordinate( x ), node->ne->rect.y1 - SearchDelta ) );
2446 
2447  if ( leaf != nullptr )
2448  {
2449  p += reinterpret_cast<RecursivePointSurfaceSpline*>( leaf->data )->operator()( x, y );
2450  ++n;
2451  }
2452  }
2453  if ( node->sw != nullptr )
2454  {
2455  const typename tree::LeafNode* leaf;
2456  if ( y >= node->sw->rect.y0 )
2457  leaf = m_tree.LeafNodeAt( search_point( node->sw->rect.x1 - SearchDelta, search_coordinate( y ) ) );
2458  else if ( x <= node->sw->rect.x1 )
2459  leaf = m_tree.LeafNodeAt( search_point( search_coordinate( x ), node->sw->rect.y0 + SearchDelta ) );
2460  else
2461  leaf = m_tree.LeafNodeAt( search_point( node->sw->rect.x1 - SearchDelta, node->sw->rect.y0 + SearchDelta ) );
2462 
2463  if ( leaf != nullptr )
2464  {
2465  p += reinterpret_cast<RecursivePointSurfaceSpline*>( leaf->data )->operator()( x, y );
2466  ++n;
2467  }
2468  }
2469  if ( node->se != nullptr )
2470  {
2471  const typename tree::LeafNode* leaf;
2472  if ( y >= node->se->rect.y0 )
2473  leaf = m_tree.LeafNodeAt( search_point( node->se->rect.x0 + SearchDelta, search_coordinate( y ) ) );
2474  else if ( x <= node->se->rect.x0 )
2475  leaf = m_tree.LeafNodeAt( search_point( node->se->rect.x0 + SearchDelta, node->se->rect.y0 + SearchDelta ) );
2476  else
2477  leaf = m_tree.LeafNodeAt( search_point( search_coordinate( x ), node->se->rect.y0 + SearchDelta ) );
2478 
2479  if ( leaf != nullptr )
2480  {
2481  p += reinterpret_cast<RecursivePointSurfaceSpline*>( leaf->data )->operator()( x, y );
2482  ++n;
2483  }
2484  }
2485 
2486  if ( n > 0 )
2487  return p/double( n );
2488  }
2489 
2490  return 0.0;
2491  }
2492 
2502  template <typename T>
2503  DPoint operator ()( const GenericPoint<T>& p ) const
2504  {
2505  return operator ()( double( p.x ), double( p.y ) );
2506  }
2507 
2508 private:
2509 
2510  /*
2511  * Interpolation data point, QuadTree-compatible.
2512  */
2513  struct Node
2514  {
2515  DPoint position, value;
2516  float weight;
2517 
2518  using component = DPoint::component;
2519 
2520  template <class point1, class point2>
2521  Node( const point1& p, const point2& v )
2522  : position( double( p.x ), double( p.y ) )
2523  , value( double( v.x ), double( v.y ) )
2524  {
2525  }
2526 
2527  template <class point1, class point2>
2528  Node( const point1& p, const point2& v, float w )
2529  : position( double( p.x ), double( p.y ) )
2530  , value( double( v.x ), double( v.y ) )
2531  , weight( w )
2532  {
2533  }
2534 
2535  Node( const Node& ) = default;
2536 
2537  double& operator []( int i )
2538  {
2539  return position[i];
2540  }
2541 
2542  double operator []( int i ) const
2543  {
2544  return position[i];
2545  }
2546  };
2547 
2548  using tree = QuadTree<Node>;
2549  using node_list = typename tree::point_list;
2550  using search_rectangle = typename tree::rectangle;
2551  using search_coordinate = typename tree::coordinate;
2552  using search_point = GenericPoint<search_coordinate>;
2553 
2554  tree m_tree; // the tree of subsplines
2555  PointSurfaceSpline m_spline; // final point spline if there is no further recursion
2556  search_rectangle m_rect = search_coordinate( 0 ); // the interpolation region for this subspline
2557  bool m_extrapolate = __PCL_RSSPLINE_DEFAULT_ALLOW_EXTRAPOLATION;
2558 
2559  static constexpr search_coordinate SearchDelta =
2560  2 * std::numeric_limits<search_coordinate>::epsilon();
2561 
2562  /*
2563  * Parallel subspline generation.
2564  */
2565  struct SubsplineData
2566  {
2567  Array<DPoint> P1, P2;
2568  Array<float> PW;
2569  mutable void** nodeData;
2570 
2571  SubsplineData( const Array<DPoint>& p1, const Array<DPoint>& p2, const Array<float>& pw, void*& nd )
2572  : P1( p1 )
2573  , P2( p2 )
2574  , PW( pw )
2575  , nodeData( &nd )
2576  {
2577  }
2578  };
2579 
2580  class SubsplineGenerationThread : public Thread
2581  {
2582  public:
2583 
2584  SubsplineGenerationThread( const AbstractImage::ThreadData& data,
2585  const Array<SubsplineData>& subsplineData,
2586  float smoothness,
2587  int order,
2588  bool allowExtrapolation,
2589  int maxSplineLength,
2590  int bucketCapacity,
2591  int startIndex, int endIndex )
2592  : m_data( data )
2593  , m_subsplineData( subsplineData )
2594  , m_smoothness( smoothness )
2595  , m_order( order )
2596  , m_allowExtrapolation( allowExtrapolation )
2597  , m_maxSplineLength( maxSplineLength )
2598  , m_bucketCapacity( bucketCapacity )
2599  , m_startIndex( startIndex )
2600  , m_endIndex( endIndex )
2601  {
2602  }
2603 
2604  PCL_HOT_FUNCTION void Run() override
2605  {
2607 
2608  for ( int i = m_startIndex; i < m_endIndex; ++i )
2609  {
2610  const SubsplineData& d = m_subsplineData[i];
2611  AutoPointer<RecursivePointSurfaceSpline> s(
2612  new RecursivePointSurfaceSpline( d.P1, d.P2, m_smoothness, m_order, d.PW,
2613  m_allowExtrapolation,
2614  m_maxSplineLength,
2615  m_bucketCapacity,
2616  false/*verbose*/ ) );
2617  *reinterpret_cast<RecursivePointSurfaceSpline**>( d.nodeData ) = s.Release();
2618 
2620  }
2621 
2622  m_success = true;
2623  }
2624 
2625  operator bool() const
2626  {
2627  return m_success;
2628  }
2629 
2630  private:
2631 
2632  const AbstractImage::ThreadData& m_data;
2633  const Array<SubsplineData>& m_subsplineData;
2634  float m_smoothness;
2635  int m_order;
2636  bool m_allowExtrapolation;
2637  int m_maxSplineLength;
2638  int m_bucketCapacity;
2639  int m_startIndex, m_endIndex;
2640  bool m_success = false;
2641  };
2642 
2643  friend class DrizzleData;
2644  friend class DrizzleDataDecoder;
2645 };
2646 
2647 // ----------------------------------------------------------------------------
2648 
2649 } // pcl
2650 
2651 #endif // __PCL_SurfaceSpline_h
2652 
2653 // ----------------------------------------------------------------------------
2654 // EOF pcl/SurfaceSpline.h - Released 2024-12-28T16:53:48Z
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:239
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:5443
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:126
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.