#include <DelaunayInterpolation.h>
Public Member Functions | |
| void | calc_node_gradient (const double *f, const int method, double *grad, const int n_it=20) const |
| Subroutine to calculate the gradient at each node. | |
| void | clear () |
| Clear the data. | |
| std::tuple< AMP::Array< TYPE >, AMP::Array< int > > | copy_tessellation () const |
| Function to copy the tessellation. | |
| void | create_tessellation (const AMP::Array< TYPE > &x) |
| Function to construct the tessellation. | |
| void | create_tessellation (const Array< TYPE > &x, const Array< int > &tri) |
| Function to construct the tessellation using a given tessellation. | |
| void | create_tessellation (size_t N, const TYPE *x, const TYPE *y, const TYPE *z) |
| Function to construct the tessellation. | |
| DelaunayInterpolation () | |
| Empty constructor. | |
| DelaunayInterpolation (const DelaunayInterpolation &)=delete | |
| template<class TYPE2 > | |
| Array< size_t > | find_nearest (const Array< TYPE2 > &xi) const |
| Subroutine to find the nearest neighbor to a point. | |
| template<class TYPE2 > | |
| Array< int > | find_tri (const Array< TYPE2 > &xi, bool extrap=false) const |
| Subroutine to find the triangle that contains the point. | |
| size_t | get_N_tri () const |
| Function to return the number of triangles in the tessellation. | |
| AMP::Array< int > | get_tri () const |
| Function to return the triangles in the tessellation. | |
| AMP::Array< int > | get_tri_nab () const |
| Function to return the triangles neignbors. | |
| template<class TYPE2 > | |
| std::tuple< AMP::Array< double >, AMP::Array< double > > | interp_cubic (const AMP::Array< double > &f, const AMP::Array< double > &g, const AMP::Array< TYPE2 > &xi, const AMP::Array< int > &index, int extrap=0) const |
| Subroutine to perform cubic interpoaltion. | |
| template<class TYPE2 > | |
| std::tuple< AMP::Array< double >, AMP::Array< double > > | interp_linear (const AMP::Array< double > &f, const AMP::Array< TYPE2 > &xi, const AMP::Array< int > &index, bool extrap=false) const |
| Subroutine to perform linear interpoaltion. | |
| template<class TYPE2 > | |
| Array< double > | interp_nearest (const Array< double > &f, const Array< TYPE2 > &xi, const Array< size_t > &nearest) const |
| Subroutine to perform nearest-neighbor interpoaltion. | |
| DelaunayInterpolation & | operator= (const DelaunayInterpolation &)=delete |
| ~DelaunayInterpolation () | |
| Empty destructor. | |
Static Public Member Functions | |
| static void | compute_Barycentric (const int ndim, const double *x, const double *xi, double *L) |
| Subroutine to compute the Barycentric coordinates. | |
Private Member Functions | |
| void | create_kdtree () const |
| Subroutine to get the list of triangle that are connected to each triangle. | |
| void | create_node_neighbors () const |
| Subroutine to get the list of nodes that are connected to each node. | |
| void | create_node_tri () const |
| Subroutine to get the starting triangle for each node. | |
| void | create_tri_neighbors () const |
| Subroutine to get the list of triangle that are connected to each triangle. | |
| void | interp_cubic_single (const double f[], const double g[], const double xi[], const int index, double &fi, double *gi, int extrap) const |
Private Attributes | |
| unsigned * | d_N_node |
| size_t | d_N_node_sum |
| unsigned ** | d_node_list |
| int * | d_node_tri |
| kdtree * | d_tree |
| Array< int > | d_tri |
| Array< int > | d_tri_nab |
| Array< TYPE > | d_x |
This class provides Delaunay based N-dimensional simplex interpolation.
Definition at line 19 of file DelaunayInterpolation.h.
| AMP::DelaunayInterpolation< TYPE >::DelaunayInterpolation | ( | ) |
Empty constructor.
|
delete |
| AMP::DelaunayInterpolation< TYPE >::~DelaunayInterpolation | ( | ) |
Empty destructor.
| void AMP::DelaunayInterpolation< TYPE >::calc_node_gradient | ( | const double * | f, |
| const int | method, | ||
| double * | grad, | ||
| const int | n_it = 20 |
||
| ) | const |
Subroutine to calculate the gradient at each node.
This function gets a list of the nodes that connect to each node
| f | Function values at the vertices( ndim ) |
| method | Gradient method to use 1 - Use a simple least squares method using only the local nodes. This method is relatively fast, but only first order in the gradient, causing the truncation error of the interpolate to be O(x^3). Note that it will still give a better cubic interpolation than MATLAB's cubic in griddata. 2 - Least squares method that solves a sparse system. This method is second order in the gradient yielding an interpolant that has truncation error O(x^4), but requires soving a sparse 2nx2n system. This can be reasonably fast for most systems, but is memory limited for large systems (it can grow as O(n^2*ndim^2)). Currently this method is NOT implemented in the C++ version. 3 - Least squares method that uses the matrix computed by method 2 and block Gauss-Seidel iteration to improve the gradient calculated by method 1. Usually only a finite number of iterations are needed to significantly reduce the error. Technically, this method has a trunction error that is still O(x^3), but for most purposes, this term is reduced so it is less than the O(x^4) term. For most systems, this method is not significantly faster than method 2, but it does not have the same memory limitations and is O(n*ndim^2) in memory. Typically 10-20 iterations are all that is necessary. 4 - This is the same as method 3, but does not store any internal data. This saves us from creating a large temporary structure ~10*ndim^2*N at a cost of ~2x in performance. |
| grad | (Output) Calculated gradient at the nodes( ndim x N ) |
| n_it | Optional argument specifying the number of Gauss-Seidel iterations. Only used if method = 3 or 4. |
| void AMP::DelaunayInterpolation< TYPE >::clear | ( | ) |
Clear the data.
|
static |
Subroutine to compute the Barycentric coordinates.
This function computes the Barycentric coordinates and the matrix T to convert from Barycentric coordinates to cartesian (x=T*L).
| ndim | Number of dimensions |
| x | Coordinates of the triangle vertices( ndim x ndim+1 ) |
| xi | Coordinates of the desired point( ndim x 1 ) |
| L | (output) The Barycentric coordinates of the point( ndim+1 x 1 ) |
| std::tuple< AMP::Array< TYPE >, AMP::Array< int > > AMP::DelaunayInterpolation< TYPE >::copy_tessellation | ( | ) | const |
Function to copy the tessellation.
This function copies the internal tessellation to a user-provided array.
|
private |
Subroutine to get the list of triangle that are connected to each triangle.
|
private |
Subroutine to get the list of nodes that are connected to each node.
|
private |
Subroutine to get the starting triangle for each node.
| void AMP::DelaunayInterpolation< TYPE >::create_tessellation | ( | const AMP::Array< TYPE > & | x | ) |
Function to construct the tessellation.
This function creates the tessellation using the given points.
| x | The coordinates of the vertices( ndim x N ) |
| void AMP::DelaunayInterpolation< TYPE >::create_tessellation | ( | const Array< TYPE > & | x, |
| const Array< int > & | tri | ||
| ) |
Function to construct the tessellation using a given tessellation.
This function sets the internal tessellation to match a provided tessellation. It does not check if the provided tessellation is valid. It allows the user to provide their own tesselation if desired. If sucessful, this routine returns 0.
| x | The coordinates of the vertices( ndim x N ) or to update the coordinate pointers before each call. See update_coordinates for more information. |
| tri | The tesselation( ndim+1 x N_tri ) |
| void AMP::DelaunayInterpolation< TYPE >::create_tessellation | ( | size_t | N, |
| const TYPE * | x, | ||
| const TYPE * | y, | ||
| const TYPE * | z | ||
| ) |
Function to construct the tessellation.
This function creates the tessellation using the given points.
| N | The number of vertices |
| x | The coordinates of the x vertices |
| y | The coordinates of the y vertices (may be NULL for 1D) |
| z | The coordinates of the z vertices (may be NULL for 1D/2D) |
|
private |
Subroutine to get the list of triangle that are connected to each triangle.
| Array< size_t > AMP::DelaunayInterpolation< TYPE >::find_nearest | ( | const Array< TYPE2 > & | xi | ) | const |
Subroutine to find the nearest neighbor to a point.
This function finds the nearest neighbor to each point. It is able to do perform the search in O(N^(1/ndim)) on average. Note: this function requires the calculate of the node lists if they are not stored (see set_storage_level)
| xi | Coordinates of the query points ( ndim x Ni ) |
| Array< int > AMP::DelaunayInterpolation< TYPE >::find_tri | ( | const Array< TYPE2 > & | xi, |
| bool | extrap = false |
||
| ) | const |
Subroutine to find the triangle that contains the point.
This function finds the triangle that contains the given points. This uses a simple search method, which may not be the most efficient. It is O(N*Ni). Note: this function requires the calculate of the node lists and triangle neighbor lists if they have not been calculated
| xi | Coordinates of the query points ( ndim x Ni ) |
| extrap | If the point is outside the convex hull, return the nearest triangle instead of -1 |
| size_t AMP::DelaunayInterpolation< TYPE >::get_N_tri | ( | ) | const |
Function to return the number of triangles in the tessellation.
This function returns the number of triangles in the tessellation.
| AMP::Array< int > AMP::DelaunayInterpolation< TYPE >::get_tri | ( | ) | const |
Function to return the triangles in the tessellation.
| AMP::Array< int > AMP::DelaunayInterpolation< TYPE >::get_tri_nab | ( | ) | const |
Function to return the triangles neignbors.
| std::tuple< AMP::Array< double >, AMP::Array< double > > AMP::DelaunayInterpolation< TYPE >::interp_cubic | ( | const AMP::Array< double > & | f, |
| const AMP::Array< double > & | g, | ||
| const AMP::Array< TYPE2 > & | xi, | ||
| const AMP::Array< int > & | index, | ||
| int | extrap = 0 |
||
| ) | const |
Subroutine to perform cubic interpoaltion.
This function performs cubic interpoaltion. Note: If the point is not contained within a triangle NaN will be returned.
| f | Function values at the triangle vertices( 1 x N ) |
| g | Gradient of f(x) at the triangle vertices( ndim x N ) (see calc_node_gradient if unknown) |
| xi | Coordinates of the query points( ndim x Ni ) |
| index | The index of the triangle containing the point (see find_tri) |
| extrap | Do we want to extrapolate from the current triangle 0: Do not extrapolate (NaNs will be used for points outside the domain) 1: Extrapolate using linear interpolation (using the nearest point and it's gradient) 2: Extrapolate using quadratic interpolation (using linear extrapolation for the gradient) |
|
private |
| std::tuple< AMP::Array< double >, AMP::Array< double > > AMP::DelaunayInterpolation< TYPE >::interp_linear | ( | const AMP::Array< double > & | f, |
| const AMP::Array< TYPE2 > & | xi, | ||
| const AMP::Array< int > & | index, | ||
| bool | extrap = false |
||
| ) | const |
Subroutine to perform linear interpoaltion.
This function performs linear interpoaltion. If a valid triangle index is not given, NaN will be returned. If extrap is false and the point is not within the triangle, NaN will be returned.
| f | Function values at the triangle vertices( 1 x N ) |
| xi | Coordinates of the query points( ndim x Ni ) |
| index | The index of the triangle containing the point (see find_tri) |
| extrap | Do we want to extrapolate from the current triangle Note: extrapolating can incure large error if sliver triangles on the boundary are present |
| Array< double > AMP::DelaunayInterpolation< TYPE >::interp_nearest | ( | const Array< double > & | f, |
| const Array< TYPE2 > & | xi, | ||
| const Array< size_t > & | nearest | ||
| ) | const |
Subroutine to perform nearest-neighbor interpoaltion.
This function performs nearest-neighbor interpoaltion.
| f | Function values at the triangle vertices( 1 x N ) |
| xi | Coordinates of the query points( ndim x Ni ) |
| nearest | The nearest-neighbor points (see find_nearest)( 1 x Ni) |
|
delete |
|
mutableprivate |
Definition at line 262 of file DelaunayInterpolation.h.
|
mutableprivate |
Definition at line 261 of file DelaunayInterpolation.h.
|
mutableprivate |
Definition at line 263 of file DelaunayInterpolation.h.
|
mutableprivate |
Definition at line 265 of file DelaunayInterpolation.h.
|
mutableprivate |
Definition at line 266 of file DelaunayInterpolation.h.
|
private |
Definition at line 259 of file DelaunayInterpolation.h.
|
mutableprivate |
Definition at line 260 of file DelaunayInterpolation.h.
|
private |
Definition at line 258 of file DelaunayInterpolation.h.
|
Advanced Multi-Physics (AMP) Oak Ridge National Laboratory Idaho National Laboratory Los Alamos National Laboratory |
This page automatically produced from the source code by Last updated: Tue Mar 10 2026 13:06:42. Comments on this page |