|
CoastalME (Coastal Modelling Environment)
Simulates the long-term behaviour of complex coastlines
|
#include <spatial_interpolation.h>
Public Member Functions | |
| DualSpatialInterpolator (std::vector< Point2D > const &points, std::vector< double > const &values_x, std::vector< double > const &values_y, int k_neighbors=12, double power=2.0) | |
| ~DualSpatialInterpolator () | |
| void | Interpolate (std::vector< Point2D > const &query_points, std::vector< double > &results_x, std::vector< double > &results_y) const |
Private Member Functions | |
| void | InterpolatePoint (double x, double y, double &result_x, double &result_y, std::vector< unsigned int > &indices, std::vector< double > &sq_dists) const |
Private Attributes | |
| PointCloud | m_cloud |
| std::vector< double > | m_values_x |
| std::vector< double > | m_values_y |
| KDTree * | m_kdtree |
| int | m_k_neighbors |
| double | m_power |
Static Private Attributes | |
| static constexpr double | EPSILON = 1e-10 |
Definition at line 129 of file spatial_interpolation.h.
| DualSpatialInterpolator::DualSpatialInterpolator | ( | std::vector< Point2D > const & | points, |
| std::vector< double > const & | values_x, | ||
| std::vector< double > const & | values_y, | ||
| int | k_neighbors = 12, | ||
| double | power = 2.0 ) |
DualSpatialInterpolator: Optimized interpolation for paired X/Y values
This class is optimized for interpolating wave properties that have both X and Y components (e.g., wave height X and wave height Y). It's more efficient than using two separate interpolators because:
DualSpatialInterpolator interp(points, values_x, values_y, 12, 2.0); interp.Interpolate(query_points, results_x, results_y); Constructor: Build dual interpolator for X and Y values
| points | Input point coordinates (x, y) |
| values_x | X-component values at those points |
| values_y | Y-component values at those points |
| k_neighbors | Number of nearest neighbors to use (default: 12) |
| power | IDW power parameter (default: 2.0) |
Definition at line 267 of file spatial_interpolation.cpp.
| DualSpatialInterpolator::~DualSpatialInterpolator | ( | ) |
Definition at line 288 of file spatial_interpolation.cpp.
| void DualSpatialInterpolator::Interpolate | ( | std::vector< Point2D > const & | query_points, |
| std::vector< double > & | results_x, | ||
| std::vector< double > & | results_y ) const |
Definition at line 350 of file spatial_interpolation.cpp.
Referenced by CSimulation::nInterpolateWavesToPolygonCells().
|
private |
Definition at line 293 of file spatial_interpolation.cpp.
Referenced by Interpolate().
|
staticconstexprprivate |
Definition at line 153 of file spatial_interpolation.h.
Referenced by InterpolatePoint().
|
private |
Definition at line 146 of file spatial_interpolation.h.
Referenced by DualSpatialInterpolator(), Interpolate(), and InterpolatePoint().
|
private |
Definition at line 150 of file spatial_interpolation.h.
Referenced by DualSpatialInterpolator(), Interpolate(), and InterpolatePoint().
|
private |
Definition at line 149 of file spatial_interpolation.h.
Referenced by DualSpatialInterpolator(), InterpolatePoint(), and ~DualSpatialInterpolator().
|
private |
Definition at line 151 of file spatial_interpolation.h.
Referenced by DualSpatialInterpolator(), and InterpolatePoint().
|
private |
Definition at line 147 of file spatial_interpolation.h.
Referenced by DualSpatialInterpolator(), and InterpolatePoint().
|
private |
Definition at line 148 of file spatial_interpolation.h.
Referenced by DualSpatialInterpolator(), and InterpolatePoint().