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;
206 inline static bool HasShapeParameter(
int rbf )
208 return rbf == Gaussian
209 || rbf == Multiquadric
210 || rbf == InverseMultiquadric
211 || rbf == InverseQuadratic;
229 using rbf_type = RadialBasisFunction::value_type;
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 );
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 );
290 static void EvaluateHandle(
const void* handle,
float* z,
const float* x,
const float* y,
double x0,
double y0,
double r,
size_type n );
296 static void EvaluateHandle(
const void* handle,
double* z,
const double* x,
const double* y,
double x0,
double y0,
double r,
size_type n );
370 template <
typename T>
415 , m_havePolynomial( S.m_havePolynomial )
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 )
430 m_handle = DuplicateHandle( S.m_handle );
438 , m_havePolynomial( S.m_havePolynomial )
441 , m_x( std::move( S.m_x ) )
442 , m_y( std::move( S.m_y ) )
443 , m_z( std::move( S.m_z ) )
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 )
453 m_handle = S.m_handle;
462 DestroyHandle( m_handle );
475 m_havePolynomial = S.m_havePolynomial;
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 );
503 m_havePolynomial = S.m_havePolynomial;
506 m_x = std::move( S.m_x );
507 m_y = std::move( S.m_y );
508 m_z = std::move( S.m_z );
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;
529 return m_handle != 0 || !m_spline.IsEmpty();
549 for (
int i = 0; i < m_x.Length(); ++i )
550 x[i] = m_x0 + m_x[i]/m_r0;
563 for (
int i = 0; i < m_y.Length(); ++i )
564 y[i] = m_y0 + m_y[i]/m_r0;
623 return Sqrt( m_eps2 );
645 m_eps2 = m_eps*m_eps;
681 PCL_PRECONDITION( order > 1 )
696 return m_havePolynomial;
709 m_havePolynomial = enable;
721 EnablePolynomial( !disable );
752 PCL_PRECONDITION( s >= 0 )
809 PCL_PRECONDITION( x !=
nullptr && y !=
nullptr && z !=
nullptr )
810 PCL_PRECONDITION( n > 2 )
812 if ( x ==
nullptr || y ==
nullptr || z ==
nullptr )
813 throw Error(
"SurfaceSpline::Initialize(): Null vector pointer(s)." );
816 throw Error(
"SurfaceSpline::Initialize(): At least three input nodes must be specified." );
820 if ( m_smoothing <= 0 )
829 for (
int i = 0; i < n; ++i )
841 for (
int i = 0; i < n; ++i )
843 double dx = x[i] - m_x0;
844 double dy = y[i] - m_y0;
845 double r =
Sqrt( dx*dx + dy*dy );
850 throw Error(
"SurfaceSpline::Initialize(): Empty or insignificant interpolation space." );
856 if ( unlikely( RadialBasisFunction::HasShapeParameter( m_rbf ) ) && m_eps == 0 )
859 for (
int i = 0; i < n; ++i )
860 for (
int j = i; ++j < n; )
862 double dx = x[i] - x[j];
863 double dy = y[i] - y[j];
876 for (
int i = 0; i < n; ++i )
880 (w !=
nullptr && w[i] > 0) ? w[i] : 1.0F );
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() )
897 int N = n - int( remove.
Length() );
899 throw Error(
"SurfaceSpline::Initialize(): Less than three input nodes left after sanitization." );
906 for (
int j : remove )
908 for ( ; i < j; ++i, ++k )
914 m_weights[k] = P[i].w;
918 for ( ; i < n; ++i, ++k )
924 m_weights[k] = P[i].w;
927 if ( !RadialBasisFunction::HasCoreImplementation( m_rbf ) )
928 m_spline =
vector(
scalar( 0 ), N + (m_havePolynomial ? ((m_order*(m_order + 1)) >> 1) : 0) );
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() );
952 DestroyHandle( m_handle );
972 SerializeHandle( data, m_handle );
985 scalar operator ()(
double x,
double y )
const
987 PCL_PRECONDITION( !m_x.IsEmpty() && !m_y.IsEmpty() )
988 PCL_CHECK( m_order >= 2 )
989 PCL_CHECK( !m_spline.IsEmpty() || m_handle != 0 )
1000 if ( m_handle != 0 )
1001 return scalar( EvaluateHandle( m_handle, x, y ) );
1006 int n = m_x.Length();
1008 if ( m_havePolynomial )
1014 z += m_spline[n+1]*x + m_spline[n+2]*y;
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;
1020 for (
int i = 1, j = n+1, i1 = (m_order*(m_order + 1))>>1, ix = 0, iy = 0; i < i1; ++i, ++j )
1025 z += m_spline[j] *
PowI( x, ix );
1031 z += m_spline[j] *
PowI( x, ix ) *
PowI( y, iy );
1042 case RadialBasisFunction::VariableOrder:
1043 for (
int i = 0; i < n; ++i )
1045 double dx = m_x[i] - x;
1046 double dy = m_y[i] - y;
1047 double r2 = dx*dx + dy*dy;
1051 for (
int j = m_order; --j > 1; )
1053 z += m_spline[i] * r2m1 *
pcl::Ln( r2 );
1057 case RadialBasisFunction::ThinPlateSpline:
1058 for (
int i = 0; i < n; ++i )
1060 double dx = m_x[i] - x;
1061 double dy = m_y[i] - y;
1062 double r2 = dx*dx + dy*dy;
1067 case RadialBasisFunction::Gaussian:
1068 for (
int i = 0; i < n; ++i )
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) );
1075 case RadialBasisFunction::Multiquadric:
1076 for (
int i = 0; i < n; ++i )
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) );
1083 case RadialBasisFunction::InverseMultiquadric:
1084 for (
int i = 0; i < n; ++i )
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) );
1091 case RadialBasisFunction::InverseQuadratic:
1092 for (
int i = 0; i < n; ++i )
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));
1111 template <
typename Tp>
1114 return operator ()(
double( p.
x ),
double( p.
y ) );
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 )
1139 if ( m_handle != 0 )
1140 EvaluateHandle( m_handle, Z, X, Y, m_x0, m_y0, m_r0, n );
1143 Z[i] =
operator()( X[i], Y[i] );
1156 PCL_PRECONDITION( !m_x.IsEmpty() && !m_y.IsEmpty() )
1157 PCL_CHECK( m_order >= 2 )
1158 PCL_CHECK( !m_spline.IsEmpty() || m_handle != 0 )
1160 int n =
Min( X.Length(), Y.Length() );
1162 Evaluate( Z.Begin(), X.Begin(), Y.Begin(),
size_type( n ) );
1177 return m_handle != 0;
1182 rbf_type m_rbf = RadialBasisFunction::Default;
1183 bool m_havePolynomial =
true;
1193 float m_smoothing = 0;
1194 weight_vector m_weights;
1210 : x( a_x ), y( a_y ), z( a_z ), w( a_w )
1216 return x == p.x && y == p.y;
1221 return (x != p.x) ? x < p.x : y < p.y;
1226 friend class DrizzleDataDecoder;
1258 using rbf_type = spline::rbf_type;
1283 template <
class po
int_list1,
class po
int_list2,
class weight_vector = FVector>
1285 float smoothness = 0,
int order = 2,
1286 const weight_vector& W = weight_vector(),
1287 rbf_type rbf = RadialBasisFunction::Default,
1289 bool polynomial =
true )
1291 Initialize( P1, P2, smoothness, W, order, rbf, eps, polynomial );
1303 Initialize( Sx, Sy );
1393 template <
class po
int_list1,
class po
int_list2,
class weight_vector = FVector>
1395 float smoothness = 0,
const weight_vector& W = weight_vector(),
int order = 2,
1396 rbf_type rbf = RadialBasisFunction::Default,
1398 bool polynomial =
true )
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() ) )
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." );
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 );
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 );
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 )
1432 X[i] = double( P1[i].x );
1433 Y[i] = double( P1[i].y );
1435 if ( m_incremental )
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];
1445 ZX[i] = double( P2[i].x );
1446 ZY[i] = double( P2[i].y );
1449 if ( ZX[i] < zxMin )
1451 if ( ZX[i] > zxMax )
1453 if ( ZY[i] < zyMin )
1455 if ( ZY[i] > zyMax )
1459 if ( m_useSimplifiers )
1468 m_epsX = (zxMax - zxMin)/100;
1469 double epsLow = 0, epsHigh = 10*m_epsX;
1470 for (
int i = 0; i < 200; ++i )
1473 SS.
Simplify( XXS, YXS, ZXS, X, Y, ZX );
1474 if ( XXS.
Length() > m_maxSplinePoints )
1477 if ( epsHigh - epsLow <= 4*std::numeric_limits<double>::epsilon() )
1482 if ( epsHigh - epsLow <= 2*std::numeric_limits<double>::epsilon() )
1486 m_epsX = (epsLow + epsHigh)/2;
1488 if ( XXS.
Length() > m_maxSplinePoints )
1490 m_truncatedX =
true;
1498 m_epsY = (zyMax - zyMin)/100;
1499 epsLow = 0, epsHigh = 10*m_epsY;
1500 for (
int i = 0; i < 200; ++i )
1503 SS.
Simplify( XYS, YYS, ZYS, X, Y, ZY );
1504 if ( XYS.
Length() > m_maxSplinePoints )
1507 if ( epsHigh - epsLow <= 4*std::numeric_limits<double>::epsilon() )
1512 if ( epsHigh - epsLow <= 2*std::numeric_limits<double>::epsilon() )
1516 m_epsY = (epsLow + epsHigh)/2;
1518 if ( XYS.
Length() > m_maxSplinePoints )
1520 m_truncatedY =
true;
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 );
1532 Tx.ValidateOrThrow();
1533 Ty.ValidateOrThrow();
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 );
1544 Tx.ValidateOrThrow();
1545 Ty.ValidateOrThrow();
1568 m_useSimplifiers = m_incremental =
false;
1585 m_epsX = m_epsY = 0;
1586 m_truncatedX = m_truncatedY =
false;
1595 return m_Sx.IsValid() && m_Sy.IsValid() && (!m_incremental || m_H.IsValid());
1622 return m_maxSplinePoints;
1635 PCL_PRECONDITION( n >= 3 )
1636 m_maxSplinePoints =
Max( 3, n );
1647 return m_useSimplifiers;
1659 m_useSimplifiers = enable;
1671 EnableSimplifiers( !disable );
1681 return m_simplifierRejectFraction;
1692 PCL_PRECONDITION( rejectFraction > 0 && rejectFraction < 1 )
1693 m_simplifierRejectFraction =
Range( rejectFraction, 0.0F, 1.0F );
1721 return m_truncatedX;
1731 return m_truncatedY;
1741 return TruncatedX() || TruncatedY();
1783 return m_incremental;
1793 m_incremental = enable;
1803 EnableIncrementalFunction( !disable );
1809 DPoint operator ()(
double x,
double y )
const
1811 PCL_PRECONDITION( ISValid() )
1812 DPoint p( m_Sx( x, y ), m_Sy( x, y ) );
1813 if ( m_incremental )
1821 template <
typename T>
1824 return operator ()(
double( p.
x ),
double( p.
y ) );
1842 template <
typename T>
1845 PCL_PRECONDITION( ISValid() )
1846 m_Sx.Evaluate( ZX, X, Y, n );
1847 m_Sy.Evaluate( ZY, X, Y, n );
1848 if ( m_incremental )
1851 DPoint dxy = m_H( X[i], Y[i] );
1871 PCL_PRECONDITION( ISValid() )
1874 Evaluate( ZX.Begin(), ZY.
Begin(), X.Begin(), Y.Begin(), n );
1899 PCL_PRECONDITION( ISValid() )
1907 return Evaluate( X, Y );
1921 return m_Sx.HasFastVectorEvaluation() && m_Sy.HasFastVectorEvaluation();
1930 int m_maxSplinePoints = __PCL_PSSPLINE_DEFAULT_MAX_POINTS;
1935 bool m_incremental =
false;
1941 bool m_useSimplifiers =
false;
1942 float m_simplifierRejectFraction = 0.1F;
1945 bool m_truncatedX =
false;
1946 bool m_truncatedY =
false;
1948 template <
class weight_vector>
1949 class SplineGenerationThread :
public Thread
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 )
1958 PCL_HOT_FUNCTION
void Run()
override
1962 m_S.Initialize( m_X.Begin(), m_Y.Begin(), m_Z.Begin(), m_N, m_W.Begin() );
1965 catch (
const Exception& x )
1971 m_exception = Error(
"Unknown exception" );
1976 void ValidateOrThrow()
const
1978 if ( !m_S.IsValid() )
1990 Exception m_exception;
1993 friend class DrizzleData;
1994 friend class DrizzleDataDecoder;
1995 friend class SplineWorldTransformation;
2060 template <
class po
int_list1,
class po
int_list2,
class weight_vector = FVector>
2062 float smoothness = 0,
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 )
2070 Initialize( P1, P2, smoothness, W, order, allowExtrapolation, maxSplineLength, bucketCapacity, verbose );
2157 template <
class po
int_list1,
class po
int_list2,
class weight_vector = FVector>
2159 float smoothness = 0,
2160 const weight_vector& W = weight_vector(),
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 )
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() )
2174 if ( P1.Length() < 3 || P2.Length() < 3 )
2175 throw Error(
"RecursivePointSurfaceSpline::Initialize(): At least three input nodes must be specified." );
2177 bool weighted = smoothness > 0 && !W.IsEmpty();
2179 if ( P1.Length() > P2.Length() || weighted && P1.Length() >
size_type( W.Length() ) )
2180 throw Error(
"RecursivePointSurfaceSpline::Initialize(): Insufficient data." );
2182 m_extrapolate = allowExtrapolation;
2184 if ( P1.Length() <=
size_type( maxSplineLength ) )
2187 #ifndef __PCL_BUILDING_PIXINSIGHT_APPLICATION
2192 monitor.
Initialize(
"Building surface subsplines", 100 );
2195 for (
const auto& p : P1 )
2197 double x = double( p.x );
2198 double y = double( p.y );
2199 if ( x < m_rect.x0 )
2201 else if ( x > m_rect.x1 )
2203 if ( y < m_rect.y0 )
2205 else if ( y > m_rect.y1 )
2209 m_spline.Initialize( P1, P2, smoothness, W, order,
2210 (order == 2) ? RadialBasisFunction::ThinPlateSpline
2211 : RadialBasisFunction::VariableOrder );
2213 #ifndef __PCL_BUILDING_PIXINSIGHT_APPLICATION
2223 search_rectangle rect = search_coordinate( 0 );
2225 for (
size_type i = 0; i < P1.Length(); ++i )
2227 const auto& p1 = P1[i];
2228 const auto& p2 = P2[i];
2229 data << (weighted ? Node( p1, p2, W[i] ) : Node( p1, p2 ));
2234 else if ( x > rect.x1 )
2238 else if ( y > rect.y1 )
2242 if ( rect.Width() < rect.Height() )
2243 rect.InflateBy( (rect.Height() - rect.Width())/2, search_coordinate( 0 ) );
2245 rect.InflateBy( search_coordinate( 0 ), (rect.Width() - rect.Height())/2 );
2250 m_tree.Build( rect, data, bucketCapacity );
2257 [&](
const search_rectangle& rect,
const node_list& points,
void*& D )
2264 search_rectangle r = rect.InflatedBy( 1.5*
Max( rect.Width(), rect.Height() ) );
2265 node_list N = m_tree.Search( r );
2272 if ( N.Length() >
size_type( maxSplineLength ) )
2274 DPoint c = rect.Center();
2275 double d = 0.61 * (
Max( r.Width(), r.Height() )/2 +
Max( rect.Width(), rect.Height() )/4);
2277 [&](
const auto& a,
const auto& b )
2279 return Abs( a.position.DistanceTo( c ) - d ) <
Abs( b.position.DistanceTo( c ) - d );
2289 for (
int j = 0, l =
Min(
int( N.Length() ), maxSplineLength ); j < l; ++j )
2291 P1 << N[j].position;
2297 subsplineData << SubsplineData( P1, P2, PW, D );
2304 #ifndef __PCL_BUILDING_PIXINSIGHT_APPLICATION
2309 monitor.
Initialize(
"Building recursive surface subsplines", subsplineData.
Length() );
2314 m_parallel ? m_maxProcessors : 1 );
2317 for (
int i = 0, n = 0; i < int( L.Length() ); n += int( L[i++] ) )
2318 threads.
Add(
new SubsplineGenerationThread( threadData,
2326 n +
int( L[i] ) ) );
2339 [](
const search_rectangle&, node_list&,
void*& D )
2345 m_rect = search_coordinate( 0 );
2360 return !m_tree.IsEmpty();
2369 return IsRecursive() || m_spline.IsValid();
2381 DPoint operator ()(
double x,
double y )
const
2383 if ( m_spline.IsValid() )
2385 if ( m_extrapolate || m_rect.IncludesFast( x, y ) )
2386 return m_spline( x, y );
2390 const typename tree::Node* node = m_tree.NodeAt(
search_point( x, y ) );
2391 if ( node ==
nullptr )
2393 if ( !m_extrapolate )
2396 search_rectangle r0 = m_tree.Root()->rect;
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 ) );
2404 node = m_tree.NodeAt(
search_point( r0.x0 + SearchDelta, search_coordinate( y ) ) );
2406 else if ( x >= r0.x1 )
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 ) );
2413 node = m_tree.NodeAt(
search_point( r0.x1 - SearchDelta, search_coordinate( y ) ) );
2415 else if ( y <= r0.y0 )
2416 node = m_tree.NodeAt(
search_point( search_coordinate( x ), r0.y0 + SearchDelta ) );
2418 node = m_tree.NodeAt(
search_point( search_coordinate( x ), r0.y1 - SearchDelta ) );
2420 if ( node ==
nullptr )
2424 if ( node->IsLeaf() )
2426 const typename tree::LeafNode* leaf =
static_cast<const typename tree::LeafNode*
>( node );
2427 if ( leaf->data ==
nullptr )
2434 if ( node->nw !=
nullptr )
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 ) );
2442 leaf = m_tree.LeafNodeAt(
search_point( node->nw->rect.x1 - SearchDelta, node->nw->rect.y1 - SearchDelta ) );
2444 if ( leaf !=
nullptr )
2450 if ( node->ne !=
nullptr )
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 ) );
2458 leaf = m_tree.LeafNodeAt(
search_point( search_coordinate( x ), node->ne->rect.y1 - SearchDelta ) );
2460 if ( leaf !=
nullptr )
2466 if ( node->sw !=
nullptr )
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 ) );
2474 leaf = m_tree.LeafNodeAt(
search_point( node->sw->rect.x1 - SearchDelta, node->sw->rect.y0 + SearchDelta ) );
2476 if ( leaf !=
nullptr )
2482 if ( node->se !=
nullptr )
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 ) );
2490 leaf = m_tree.LeafNodeAt(
search_point( search_coordinate( x ), node->se->rect.y0 + SearchDelta ) );
2492 if ( leaf !=
nullptr )
2500 return p/double( n );
2515 template <
typename T>
2518 return operator ()(
double( p.
x ),
double( p.
y ) );
2533 template <
class po
int1,
class po
int2>
2534 Node(
const point1& p,
const point2& v )
2535 : position( double( p.x ), double( p.y ) )
2536 , value( double( v.x ), double( v.y ) )
2540 template <
class po
int1,
class po
int2>
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 ) )
2548 Node(
const Node& ) =
default;
2550 double& operator [](
int i )
2555 double operator [](
int i )
const
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>;
2568 PointSurfaceSpline m_spline;
2569 search_rectangle m_rect = search_coordinate( 0 );
2570 bool m_extrapolate = __PCL_RSSPLINE_DEFAULT_ALLOW_EXTRAPOLATION;
2572 static constexpr search_coordinate SearchDelta =
2573 2 * std::numeric_limits<search_coordinate>::epsilon();
2578 struct SubsplineData
2580 Array<DPoint> P1, P2;
2582 mutable void** nodeData;
2584 SubsplineData(
const Array<DPoint>& p1,
const Array<DPoint>& p2,
const Array<float>& pw,
void*& nd )
2593 class SubsplineGenerationThread :
public Thread
2597 SubsplineGenerationThread(
const AbstractImage::ThreadData& data,
2598 const Array<SubsplineData>& subsplineData,
2601 bool allowExtrapolation,
2602 int maxSplineLength,
2604 int startIndex,
int endIndex )
2606 , m_subsplineData( subsplineData )
2607 , m_smoothness( smoothness )
2609 , m_allowExtrapolation( allowExtrapolation )
2610 , m_maxSplineLength( maxSplineLength )
2611 , m_bucketCapacity( bucketCapacity )
2612 , m_startIndex( startIndex )
2613 , m_endIndex( endIndex )
2617 PCL_HOT_FUNCTION
void Run()
override
2621 for (
int i = m_startIndex; i < m_endIndex; ++i )
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,
2630 *
reinterpret_cast<RecursivePointSurfaceSpline**
>( d.nodeData ) = s.Release();
2638 operator bool()
const
2645 const AbstractImage::ThreadData& m_data;
2646 const Array<SubsplineData>& m_subsplineData;
2649 bool m_allowExtrapolation;
2650 int m_maxSplineLength;
2651 int m_bucketCapacity;
2652 int m_startIndex, m_endIndex;
2653 bool m_success =
false;
2656 friend class DrizzleData;
2657 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.