52 #ifndef __PCL_SurfaceSpline_h
53 #define __PCL_SurfaceSpline_h
58 #include <pcl/Diagnostics.h>
71 #ifndef __PCL_BUILDING_PIXINSIGHT_APPLICATION
84 #define __PCL_PSSPLINE_DEFAULT_MAX_POINTS 2100
89 #define __PCL_RSSPLINE_DEFAULT_TREE_BUCKET_CAPACITY 64
95 #define __PCL_RSSPLINE_DEFAULT_SPLINE_MAX_LENGTH 1600
103 #define __PCL_RSSPLINE_DEFAULT_ALLOW_EXTRAPOLATION false
141 namespace RadialBasisFunction
152 __number_of_pcl_implemented_rbfs__,
153 DDMVariableOrder = 100,
154 DDMThinPlateSpline = 101,
155 DDMMultiquadric = 102,
156 Default = ThinPlateSpline
163 inline static bool Validate(
int rbf )
165 return rbf >= 0 && rbf < __number_of_pcl_implemented_rbfs__
166 || rbf == DDMThinPlateSpline
167 || rbf == DDMVariableOrder
168 || rbf == DDMMultiquadric;
175 inline static bool HasCoreImplementation(
int rbf )
177 return rbf > __number_of_pcl_implemented_rbfs__;
185 inline static bool HasDDMImplementation(
int rbf )
187 return rbf == DDMThinPlateSpline
188 || rbf == DDMVariableOrder
189 || rbf == DDMMultiquadric;
196 inline static bool IsVariableOrder(
int rbf )
198 return rbf == DDMVariableOrder || rbf == VariableOrder;
216 using rbf_type = RadialBasisFunction::value_type;
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 );
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 );
277 static void EvaluateHandle(
const void* handle,
float* z,
const float* x,
const float* y,
double x0,
double y0,
double r,
size_type n );
283 static void EvaluateHandle(
const void* handle,
double* z,
const double* x,
const double* y,
double x0,
double y0,
double r,
size_type n );
357 template <
typename T>
402 , m_havePolynomial( S.m_havePolynomial )
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 )
417 m_handle = DuplicateHandle( S.m_handle );
425 , m_havePolynomial( S.m_havePolynomial )
428 , m_x( std::move( S.m_x ) )
429 , m_y( std::move( S.m_y ) )
430 , m_z( std::move( S.m_z ) )
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 )
440 m_handle = S.m_handle;
449 DestroyHandle( m_handle );
462 m_havePolynomial = S.m_havePolynomial;
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 );
490 m_havePolynomial = S.m_havePolynomial;
493 m_x = std::move( S.m_x );
494 m_y = std::move( S.m_y );
495 m_z = std::move( S.m_z );
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;
516 return m_handle != 0 || !m_spline.IsEmpty();
536 for (
int i = 0; i < m_x.Length(); ++i )
537 x[i] = m_x0 + m_x[i]/m_r0;
550 for (
int i = 0; i < m_y.Length(); ++i )
551 y[i] = m_y0 + m_y[i]/m_r0;
610 return Sqrt( m_eps2 );
632 m_eps2 = m_eps*m_eps;
668 PCL_PRECONDITION( order > 1 )
683 return m_havePolynomial;
696 m_havePolynomial = enable;
708 EnablePolynomial( !disable );
739 PCL_PRECONDITION( s >= 0 )
796 PCL_PRECONDITION( x !=
nullptr && y !=
nullptr && z !=
nullptr )
797 PCL_PRECONDITION( n > 2 )
799 if ( x ==
nullptr || y ==
nullptr || z ==
nullptr )
800 throw Error(
"SurfaceSpline::Initialize(): Null vector pointer(s)." );
803 throw Error(
"SurfaceSpline::Initialize(): At least three input nodes must be specified." );
807 if ( m_smoothing <= 0 )
816 for (
int i = 0; i < n; ++i )
828 for (
int i = 0; i < n; ++i )
830 double dx = x[i] - m_x0;
831 double dy = y[i] - m_y0;
832 double r =
Sqrt( dx*dx + dy*dy );
837 throw Error(
"SurfaceSpline::Initialize(): Empty or insignificant interpolation space." );
846 for (
int i = 0; i < n; ++i )
847 for (
int j = i; ++j < n; )
849 double dx = x[i] - x[j];
850 double dy = y[i] - y[j];
863 for (
int i = 0; i < n; ++i )
867 (w !=
nullptr && w[i] > 0) ? w[i] : 1.0F );
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() )
884 int N = n - int( remove.
Length() );
886 throw Error(
"SurfaceSpline::Initialize(): Less than three input nodes left after sanitization." );
893 for (
int j : remove )
895 for ( ; i < j; ++i, ++k )
901 m_weights[k] = P[i].w;
905 for ( ; i < n; ++i, ++k )
911 m_weights[k] = P[i].w;
914 if ( !RadialBasisFunction::HasCoreImplementation( m_rbf ) )
915 m_spline =
vector(
scalar( 0 ), N + (m_havePolynomial ? ((m_order*(m_order + 1)) >> 1) : 0) );
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() );
939 DestroyHandle( m_handle );
959 SerializeHandle( data, m_handle );
972 scalar operator ()(
double x,
double y )
const
974 PCL_PRECONDITION( !m_x.IsEmpty() && !m_y.IsEmpty() )
975 PCL_CHECK( m_order >= 2 )
976 PCL_CHECK( !m_spline.IsEmpty() || m_handle != 0 )
988 return scalar( EvaluateHandle( m_handle, x, y ) );
993 int n = m_x.Length();
995 if ( m_havePolynomial )
1001 z += m_spline[n+1]*x + m_spline[n+2]*y;
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;
1007 for (
int i = 1, j = n+1, i1 = (m_order*(m_order + 1))>>1, ix = 0, iy = 0; i < i1; ++i, ++j )
1012 z += m_spline[j] *
PowI( x, ix );
1018 z += m_spline[j] *
PowI( x, ix ) *
PowI( y, iy );
1029 case RadialBasisFunction::VariableOrder:
1030 for (
int i = 0; i < n; ++i )
1032 double dx = m_x[i] - x;
1033 double dy = m_y[i] - y;
1034 double r2 = dx*dx + dy*dy;
1038 for (
int j = m_order; --j > 1; )
1040 z += m_spline[i] * r2m1 *
pcl::Ln( r2 );
1044 case RadialBasisFunction::ThinPlateSpline:
1045 for (
int i = 0; i < n; ++i )
1047 double dx = m_x[i] - x;
1048 double dy = m_y[i] - y;
1049 double r2 = dx*dx + dy*dy;
1054 case RadialBasisFunction::Gaussian:
1055 for (
int i = 0; i < n; ++i )
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) );
1062 case RadialBasisFunction::Multiquadric:
1063 for (
int i = 0; i < n; ++i )
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) );
1070 case RadialBasisFunction::InverseMultiquadric:
1071 for (
int i = 0; i < n; ++i )
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) );
1078 case RadialBasisFunction::InverseQuadratic:
1079 for (
int i = 0; i < n; ++i )
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));
1098 template <
typename Tp>
1101 return operator ()(
double( p.
x ),
double( p.
y ) );
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 )
1126 if ( m_handle != 0 )
1127 EvaluateHandle( m_handle, Z, X, Y, m_x0, m_y0, m_r0, n );
1130 Z[i] =
operator()( X[i], Y[i] );
1143 PCL_PRECONDITION( !m_x.IsEmpty() && !m_y.IsEmpty() )
1144 PCL_CHECK( m_order >= 2 )
1145 PCL_CHECK( !m_spline.IsEmpty() || m_handle != 0 )
1147 int n =
Min( X.Length(), Y.Length() );
1149 Evaluate( Z.Begin(), X.Begin(), Y.Begin(),
size_type( n ) );
1164 return m_handle != 0;
1169 rbf_type m_rbf = RadialBasisFunction::Default;
1170 bool m_havePolynomial =
true;
1180 float m_smoothing = 0;
1181 weight_vector m_weights;
1197 : x( a_x ), y( a_y ), z( a_z ), w( a_w )
1203 return x == p.x && y == p.y;
1208 return (x != p.x) ? x < p.x : y < p.y;
1213 friend class DrizzleDataDecoder;
1245 using rbf_type = spline::rbf_type;
1270 template <
class po
int_list1,
class po
int_list2,
class weight_vector = FVector>
1272 float smoothness = 0,
int order = 2,
1273 const weight_vector& W = weight_vector(),
1274 rbf_type rbf = RadialBasisFunction::Default,
1276 bool polynomial =
true )
1278 Initialize( P1, P2, smoothness, W, order, rbf, eps, polynomial );
1290 Initialize( Sx, Sy );
1380 template <
class po
int_list1,
class po
int_list2,
class weight_vector = FVector>
1382 float smoothness = 0,
const weight_vector& W = weight_vector(),
int order = 2,
1383 rbf_type rbf = RadialBasisFunction::Default,
1385 bool polynomial =
true )
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() ) )
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." );
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 );
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 );
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 )
1419 X[i] = double( P1[i].x );
1420 Y[i] = double( P1[i].y );
1422 if ( m_incremental )
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];
1432 ZX[i] = double( P2[i].x );
1433 ZY[i] = double( P2[i].y );
1436 if ( ZX[i] < zxMin )
1438 if ( ZX[i] > zxMax )
1440 if ( ZY[i] < zyMin )
1442 if ( ZY[i] > zyMax )
1446 if ( m_useSimplifiers )
1455 m_epsX = (zxMax - zxMin)/100;
1456 double epsLow = 0, epsHigh = 10*m_epsX;
1457 for (
int i = 0; i < 200; ++i )
1460 SS.
Simplify( XXS, YXS, ZXS, X, Y, ZX );
1461 if ( XXS.
Length() > m_maxSplinePoints )
1464 if ( epsHigh - epsLow <= 4*std::numeric_limits<double>::epsilon() )
1469 if ( epsHigh - epsLow <= 2*std::numeric_limits<double>::epsilon() )
1473 m_epsX = (epsLow + epsHigh)/2;
1475 if ( XXS.
Length() > m_maxSplinePoints )
1477 m_truncatedX =
true;
1485 m_epsY = (zyMax - zyMin)/100;
1486 epsLow = 0, epsHigh = 10*m_epsY;
1487 for (
int i = 0; i < 200; ++i )
1490 SS.
Simplify( XYS, YYS, ZYS, X, Y, ZY );
1491 if ( XYS.
Length() > m_maxSplinePoints )
1494 if ( epsHigh - epsLow <= 4*std::numeric_limits<double>::epsilon() )
1499 if ( epsHigh - epsLow <= 2*std::numeric_limits<double>::epsilon() )
1503 m_epsY = (epsLow + epsHigh)/2;
1505 if ( XYS.
Length() > m_maxSplinePoints )
1507 m_truncatedY =
true;
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 );
1519 Tx.ValidateOrThrow();
1520 Ty.ValidateOrThrow();
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 );
1531 Tx.ValidateOrThrow();
1532 Ty.ValidateOrThrow();
1555 m_useSimplifiers = m_incremental =
false;
1572 m_epsX = m_epsY = 0;
1573 m_truncatedX = m_truncatedY =
false;
1582 return m_Sx.IsValid() && m_Sy.IsValid() && (!m_incremental || m_H.IsValid());
1609 return m_maxSplinePoints;
1622 PCL_PRECONDITION( n >= 3 )
1623 m_maxSplinePoints =
Max( 3, n );
1634 return m_useSimplifiers;
1646 m_useSimplifiers = enable;
1658 EnableSimplifiers( !disable );
1668 return m_simplifierRejectFraction;
1679 PCL_PRECONDITION( rejectFraction > 0 && rejectFraction < 1 )
1680 m_simplifierRejectFraction =
Range( rejectFraction, 0.0F, 1.0F );
1708 return m_truncatedX;
1718 return m_truncatedY;
1728 return TruncatedX() || TruncatedY();
1770 return m_incremental;
1780 m_incremental = enable;
1790 EnableIncrementalFunction( !disable );
1796 DPoint operator ()(
double x,
double y )
const
1798 PCL_PRECONDITION( ISValid() )
1799 DPoint p( m_Sx( x, y ), m_Sy( x, y ) );
1800 if ( m_incremental )
1808 template <
typename T>
1811 return operator ()(
double( p.
x ),
double( p.
y ) );
1829 template <
typename T>
1832 PCL_PRECONDITION( ISValid() )
1833 m_Sx.Evaluate( ZX, X, Y, n );
1834 m_Sy.Evaluate( ZY, X, Y, n );
1835 if ( m_incremental )
1838 DPoint dxy = m_H( X[i], Y[i] );
1858 PCL_PRECONDITION( ISValid() )
1861 Evaluate( ZX.Begin(), ZY.
Begin(), X.Begin(), Y.Begin(), n );
1886 PCL_PRECONDITION( ISValid() )
1894 return Evaluate( X, Y );
1908 return m_Sx.HasFastVectorEvaluation() && m_Sy.HasFastVectorEvaluation();
1917 int m_maxSplinePoints = __PCL_PSSPLINE_DEFAULT_MAX_POINTS;
1922 bool m_incremental =
false;
1928 bool m_useSimplifiers =
false;
1929 float m_simplifierRejectFraction = 0.1F;
1932 bool m_truncatedX =
false;
1933 bool m_truncatedY =
false;
1935 template <
class weight_vector>
1936 class SplineGenerationThread :
public Thread
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 )
1945 PCL_HOT_FUNCTION
void Run()
override
1949 m_S.Initialize( m_X.Begin(), m_Y.Begin(), m_Z.Begin(), m_N, m_W.Begin() );
1952 catch (
const Exception& x )
1958 m_exception = Error(
"Unknown exception" );
1963 void ValidateOrThrow()
const
1965 if ( !m_S.IsValid() )
1977 Exception m_exception;
1980 friend class DrizzleData;
1981 friend class DrizzleDataDecoder;
1982 friend class SplineWorldTransformation;
2047 template <
class po
int_list1,
class po
int_list2,
class weight_vector = FVector>
2049 float smoothness = 0,
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 )
2057 Initialize( P1, P2, smoothness, W, order, allowExtrapolation, maxSplineLength, bucketCapacity, verbose );
2144 template <
class po
int_list1,
class po
int_list2,
class weight_vector = FVector>
2146 float smoothness = 0,
2147 const weight_vector& W = weight_vector(),
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 )
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() )
2161 if ( P1.Length() < 3 || P2.Length() < 3 )
2162 throw Error(
"RecursivePointSurfaceSpline::Initialize(): At least three input nodes must be specified." );
2164 bool weighted = smoothness > 0 && !W.IsEmpty();
2166 if ( P1.Length() > P2.Length() || weighted && P1.Length() >
size_type( W.Length() ) )
2167 throw Error(
"RecursivePointSurfaceSpline::Initialize(): Insufficient data." );
2169 m_extrapolate = allowExtrapolation;
2171 if ( P1.Length() <=
size_type( maxSplineLength ) )
2174 #ifndef __PCL_BUILDING_PIXINSIGHT_APPLICATION
2179 monitor.
Initialize(
"Building surface subsplines", 100 );
2182 for (
const auto& p : P1 )
2184 double x = double( p.x );
2185 double y = double( p.y );
2186 if ( x < m_rect.x0 )
2188 else if ( x > m_rect.x1 )
2190 if ( y < m_rect.y0 )
2192 else if ( y > m_rect.y1 )
2196 m_spline.Initialize( P1, P2, smoothness, W, order,
2197 (order == 2) ? RadialBasisFunction::ThinPlateSpline
2198 : RadialBasisFunction::VariableOrder );
2200 #ifndef __PCL_BUILDING_PIXINSIGHT_APPLICATION
2210 search_rectangle rect = search_coordinate( 0 );
2212 for (
size_type i = 0; i < P1.Length(); ++i )
2214 const auto& p1 = P1[i];
2215 const auto& p2 = P2[i];
2216 data << (weighted ? Node( p1, p2, W[i] ) : Node( p1, p2 ));
2221 else if ( x > rect.x1 )
2225 else if ( y > rect.y1 )
2229 if ( rect.Width() < rect.Height() )
2230 rect.InflateBy( (rect.Height() - rect.Width())/2, search_coordinate( 0 ) );
2232 rect.InflateBy( search_coordinate( 0 ), (rect.Width() - rect.Height())/2 );
2237 m_tree.Build( rect, data, bucketCapacity );
2244 [&](
const search_rectangle& rect,
const node_list& points,
void*& D )
2251 search_rectangle r = rect.InflatedBy( 1.5*
Max( rect.Width(), rect.Height() ) );
2252 node_list N = m_tree.Search( r );
2259 if ( N.Length() >
size_type( maxSplineLength ) )
2261 DPoint c = rect.Center();
2262 double d = 0.61 * (
Max( r.Width(), r.Height() )/2 +
Max( rect.Width(), rect.Height() )/4);
2264 [&](
const auto& a,
const auto& b )
2266 return Abs( a.position.DistanceTo( c ) - d ) <
Abs( b.position.DistanceTo( c ) - d );
2276 for (
int j = 0, l =
Min(
int( N.Length() ), maxSplineLength ); j < l; ++j )
2278 P1 << N[j].position;
2284 subsplineData << SubsplineData( P1, P2, PW, D );
2291 #ifndef __PCL_BUILDING_PIXINSIGHT_APPLICATION
2296 monitor.
Initialize(
"Building recursive surface subsplines", subsplineData.
Length() );
2301 m_parallel ? m_maxProcessors : 1 );
2304 for (
int i = 0, n = 0; i < int( L.Length() ); n += int( L[i++] ) )
2305 threads.
Add(
new SubsplineGenerationThread( threadData,
2313 n +
int( L[i] ) ) );
2326 [](
const search_rectangle&, node_list&,
void*& D )
2332 m_rect = search_coordinate( 0 );
2347 return !m_tree.IsEmpty();
2356 return IsRecursive() || m_spline.IsValid();
2368 DPoint operator ()(
double x,
double y )
const
2370 if ( m_spline.IsValid() )
2372 if ( m_extrapolate || m_rect.IncludesFast( x, y ) )
2373 return m_spline( x, y );
2377 const typename tree::Node* node = m_tree.NodeAt(
search_point( x, y ) );
2378 if ( node ==
nullptr )
2380 if ( !m_extrapolate )
2383 search_rectangle r0 = m_tree.Root()->rect;
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 ) );
2391 node = m_tree.NodeAt(
search_point( r0.x0 + SearchDelta, search_coordinate( y ) ) );
2393 else if ( x >= r0.x1 )
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 ) );
2400 node = m_tree.NodeAt(
search_point( r0.x1 - SearchDelta, search_coordinate( y ) ) );
2402 else if ( y <= r0.y0 )
2403 node = m_tree.NodeAt(
search_point( search_coordinate( x ), r0.y0 + SearchDelta ) );
2405 node = m_tree.NodeAt(
search_point( search_coordinate( x ), r0.y1 - SearchDelta ) );
2407 if ( node ==
nullptr )
2411 if ( node->IsLeaf() )
2413 const typename tree::LeafNode* leaf =
static_cast<const typename tree::LeafNode*
>( node );
2414 if ( leaf->data ==
nullptr )
2421 if ( node->nw !=
nullptr )
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 ) );
2429 leaf = m_tree.LeafNodeAt(
search_point( node->nw->rect.x1 - SearchDelta, node->nw->rect.y1 - SearchDelta ) );
2431 if ( leaf !=
nullptr )
2437 if ( node->ne !=
nullptr )
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 ) );
2445 leaf = m_tree.LeafNodeAt(
search_point( search_coordinate( x ), node->ne->rect.y1 - SearchDelta ) );
2447 if ( leaf !=
nullptr )
2453 if ( node->sw !=
nullptr )
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 ) );
2461 leaf = m_tree.LeafNodeAt(
search_point( node->sw->rect.x1 - SearchDelta, node->sw->rect.y0 + SearchDelta ) );
2463 if ( leaf !=
nullptr )
2469 if ( node->se !=
nullptr )
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 ) );
2477 leaf = m_tree.LeafNodeAt(
search_point( search_coordinate( x ), node->se->rect.y0 + SearchDelta ) );
2479 if ( leaf !=
nullptr )
2487 return p/double( n );
2502 template <
typename T>
2505 return operator ()(
double( p.
x ),
double( p.
y ) );
2520 template <
class po
int1,
class po
int2>
2521 Node(
const point1& p,
const point2& v )
2522 : position( double( p.x ), double( p.y ) )
2523 , value( double( v.x ), double( v.y ) )
2527 template <
class po
int1,
class po
int2>
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 ) )
2535 Node(
const Node& ) =
default;
2537 double& operator [](
int i )
2542 double operator [](
int i )
const
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>;
2555 PointSurfaceSpline m_spline;
2556 search_rectangle m_rect = search_coordinate( 0 );
2557 bool m_extrapolate = __PCL_RSSPLINE_DEFAULT_ALLOW_EXTRAPOLATION;
2559 static constexpr search_coordinate SearchDelta =
2560 2 * std::numeric_limits<search_coordinate>::epsilon();
2565 struct SubsplineData
2567 Array<DPoint> P1, P2;
2569 mutable void** nodeData;
2571 SubsplineData(
const Array<DPoint>& p1,
const Array<DPoint>& p2,
const Array<float>& pw,
void*& nd )
2580 class SubsplineGenerationThread :
public Thread
2584 SubsplineGenerationThread(
const AbstractImage::ThreadData& data,
2585 const Array<SubsplineData>& subsplineData,
2588 bool allowExtrapolation,
2589 int maxSplineLength,
2591 int startIndex,
int endIndex )
2593 , m_subsplineData( subsplineData )
2594 , m_smoothness( smoothness )
2596 , m_allowExtrapolation( allowExtrapolation )
2597 , m_maxSplineLength( maxSplineLength )
2598 , m_bucketCapacity( bucketCapacity )
2599 , m_startIndex( startIndex )
2600 , m_endIndex( endIndex )
2604 PCL_HOT_FUNCTION
void Run()
override
2608 for (
int i = m_startIndex; i < m_endIndex; ++i )
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,
2617 *
reinterpret_cast<RecursivePointSurfaceSpline**
>( d.nodeData ) = s.Release();
2625 operator bool()
const
2632 const AbstractImage::ThreadData& m_data;
2633 const Array<SubsplineData>& m_subsplineData;
2636 bool m_allowExtrapolation;
2637 int m_maxSplineLength;
2638 int m_bucketCapacity;
2639 int m_startIndex, m_endIndex;
2640 bool m_success =
false;
2643 friend class DrizzleData;
2644 friend class DrizzleDataDecoder;
static void RunThreads(ReferenceArray< thread > &threads, ThreadData &data, bool useAffinity=true)
size_type Length() const noexcept
64-bit floating point real vector.
Drizzle integration data parser and generator.
A simple exception with an associated error message.
Generic dynamic matrix of arbitrary dimensions.
A generic point in the two-dimensional space.
component x
Abscissa (horizontal, or X-axis coordinate).
component y
Ordinate (vertical, or Y-axis coordinate).
int Length() const noexcept
Homography geometric transformation.
Eight-bit string (ISO/IEC-8859-1 or UTF-8 string)
A process using multiple concurrent execution threads.
Vector surface spline interpolation/approximation in two dimensions.
int MaxSplinePoints() const
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)
PointSurfaceSpline()=default
Vector surface spline interpolation/approximation in two dimensions with recursive subspline generati...
RecursivePointSurfaceSpline(RecursivePointSurfaceSpline &&)=default
RecursivePointSurfaceSpline(const RecursivePointSurfaceSpline &)=delete
~RecursivePointSurfaceSpline() override
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()=default
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)
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)
virtual ~SurfaceSplineBase()
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)
SurfaceSplineBase()=default
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.
void DisablePolynomial(bool disable=true)
bool HasFastVectorEvaluation() const
bool IsPolynomialEnabled() const
IsoString CoreSerialization() const
~SurfaceSpline() override
void SetShapeParameter(double eps)
void EnablePolynomial(bool enable=true)
const weight_vector & Weights() const
SurfaceSpline(SurfaceSpline &&S)
void SetRBFType(rbf_type rbf)
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.
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
bool operator<(const Array< T, A > &x1, const Array< T, A > &x2) noexcept
Complex< T > Sqrt(const Complex< T > &c) noexcept
Complex< T > Exp(const Complex< T > &c) noexcept
Complex< T > Ln(const Complex< T > &c) noexcept
T Abs(const Complex< T > &c) noexcept
T PowI(T x, int n) noexcept
double Median(const T *__restrict__ i, const T *__restrict__ j, double eps=0)
#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
constexpr const T & Range(const T &x, const T &a, const T &b) noexcept
constexpr const T & Max(const T &a, const T &b) noexcept
Thread synchronization data for status monitoring of parallel image processing tasks.
Auxiliary structure for data sanitization.