CoastalME (Coastal Modelling Environment)
Simulates the long-term behaviour of complex coastlines
Loading...
Searching...
No Matches
SpatialInterpolator Class Reference

#include <spatial_interpolation.h>

Public Member Functions

 SpatialInterpolator (std::vector< Point2D > const &points, std::vector< double > const &values, int k_neighbors=12, double power=2.0)
 
 ~SpatialInterpolator ()
 
double Interpolate (double x, double y) const
 
void Interpolate (std::vector< Point2D > const &query_points, std::vector< double > &results) const
 
KDTree const * GetKDTree () const
 
PointCloud const & GetPointCloud () const
 

Private Member Functions

 SpatialInterpolator (PointCloud const &cloud, KDTree *kdtree, std::vector< double > const &values, int k_neighbors, double power)
 

Private Attributes

PointCloud m_cloud
 
std::vector< double > m_values
 
KDTreem_kdtree
 
int m_k_neighbors
 
double m_power
 
bool m_owns_kdtree
 

Static Private Attributes

static constexpr double EPSILON = 1e-10
 

Friends

class DualSpatialInterpolator
 

Detailed Description

Definition at line 86 of file spatial_interpolation.h.

Constructor & Destructor Documentation

◆ SpatialInterpolator() [1/2]

SpatialInterpolator::SpatialInterpolator ( std::vector< Point2D > const & points,
std::vector< double > const & values,
int k_neighbors = 12,
double power = 2.0 )

Spatial Interpolation Using k-Nearest Neighbors and Inverse Distance Weighting

This file implements fast spatial interpolation using:

  1. k-d tree for efficient nearest neighbor search (O(log n) per query)
  2. Inverse Distance Weighting (IDW) for smooth interpolation
  3. OpenMP parallelization for batch operations

ALGORITHM OVERVIEW:

For each query point (grid cell):

  1. Find k nearest input points using k-d tree
  2. Calculate weights based on inverse distance: w_i = 1 / dist_i^power
  3. Interpolated value = Σ(w_i * value_i) / Σ(w_i)

KEY PARAMETERS TO TUNE:

  • k_neighbors: Number of nearest points to use (typically 8-20)
    • Higher = smoother, more averaging
    • Lower = follows local variations more closely
  • power: Exponent for inverse distance weighting (typically 1.0-4.0)
    • Higher = nearby points have more influence (sharper transitions)
    • Lower = distant points have more influence (smoother transitions)
    • power=2.0 is optimized for performance (avoids pow() function)

PERFORMANCE OPTIMIZATIONS:

  • Special fast path for power=2.0 (uses squared distance directly, avoids sqrt/pow)
  • OpenMP parallelization with guided scheduling
  • Thread-local buffers to avoid allocation overhead
  • Shared k-d tree between X and Y interpolation (DualSpatialInterpolator)

USAGE EXAMPLE:

std::vector<Point2D> input_points = {{0,0}, {10,0}, {5,10}}; std::vector<double> input_values = {1.0, 2.0, 1.5}; SpatialInterpolator interp(input_points, input_values, 12, 2.0); double result = interp.Interpolate(5.0, 5.0); // Interpolate at (5,5) Constructor: Build interpolator from points and values

Parameters
pointsInput point coordinates (x, y)
valuesValues at those points
k_neighborsNumber of nearest neighbors to use (default: 12)
powerIDW power parameter (default: 2.0)

Definition at line 59 of file spatial_interpolation.cpp.

◆ ~SpatialInterpolator()

SpatialInterpolator::~SpatialInterpolator ( )

Definition at line 89 of file spatial_interpolation.cpp.

◆ SpatialInterpolator() [2/2]

SpatialInterpolator::SpatialInterpolator ( PointCloud const & cloud,
KDTree * kdtree,
std::vector< double > const & values,
int k_neighbors,
double power )
private

Definition at line 80 of file spatial_interpolation.cpp.

Member Function Documentation

◆ GetKDTree()

KDTree const * SpatialInterpolator::GetKDTree ( ) const
inline

Definition at line 104 of file spatial_interpolation.h.

◆ GetPointCloud()

PointCloud const & SpatialInterpolator::GetPointCloud ( ) const
inline

Definition at line 107 of file spatial_interpolation.h.

◆ Interpolate() [1/2]

double SpatialInterpolator::Interpolate ( double x,
double y ) const

Interpolate at a single query point

ALGORITHM:

  1. Find k nearest neighbors using k-d tree
  2. Calculate weight for each: w_i = 1 / distance_i^power
  3. Return weighted average: Σ(w_i * value_i) / Σ(w_i)

SPECIAL CASES:

  • If query point coincides with input point (dist < EPSILON): return exact value
  • If power=2.0: optimized path using squared distances (faster)
Parameters
xX coordinate of query point
yY coordinate of query point
Returns
Interpolated value

Definition at line 111 of file spatial_interpolation.cpp.

Referenced by Interpolate().

◆ Interpolate() [2/2]

void SpatialInterpolator::Interpolate ( std::vector< Point2D > const & query_points,
std::vector< double > & results ) const

Definition at line 167 of file spatial_interpolation.cpp.

Friends And Related Symbol Documentation

◆ DualSpatialInterpolator

friend class DualSpatialInterpolator
friend

Definition at line 120 of file spatial_interpolation.h.

Referenced by DualSpatialInterpolator.

Field Documentation

◆ EPSILON

double SpatialInterpolator::EPSILON = 1e-10
staticconstexprprivate

Definition at line 117 of file spatial_interpolation.h.

Referenced by Interpolate(), and Interpolate().

◆ m_cloud

PointCloud SpatialInterpolator::m_cloud
private

◆ m_k_neighbors

int SpatialInterpolator::m_k_neighbors
private

◆ m_kdtree

KDTree* SpatialInterpolator::m_kdtree
private

◆ m_owns_kdtree

bool SpatialInterpolator::m_owns_kdtree
private

◆ m_power

double SpatialInterpolator::m_power
private

◆ m_values

std::vector<double> SpatialInterpolator::m_values
private

The documentation for this class was generated from the following files: