Advanced Multi-Physics (AMP)
On-Line Documentation
DelaunayInterpolation.h
Go to the documentation of this file.
1#ifndef included_AMP_DelaunayInterpolation
2#define included_AMP_DelaunayInterpolation
3
4#include <stdlib.h>
5#include <tuple>
6
7#include "AMP/utils/Array.h"
8#include "AMP/utils/kdtree.h"
9
10
11namespace AMP {
12
13
18template<class TYPE>
20{
21public:
24
25 // Deleted constructors
28
31
32
34
37 size_t get_N_tri() const;
38
39
42
43
46
47
49
54
55
57
64 void create_tessellation( size_t N, const TYPE *x, const TYPE *y, const TYPE *z );
65
66
68
78 void create_tessellation( const Array<TYPE> &x, const Array<int> &tri );
79
80
82
85 std::tuple<AMP::Array<TYPE>, AMP::Array<int>> copy_tessellation() const;
86
87
89
97 template<class TYPE2>
99
100
102
113 template<class TYPE2>
114 Array<int> find_tri( const Array<TYPE2> &xi, bool extrap = false ) const;
115
116
118
149 void calc_node_gradient( const double *f,
150 const int method,
151 double *grad,
152 const int n_it = 20 ) const;
153
154
156
163 template<class TYPE2>
165 const Array<TYPE2> &xi,
166 const Array<size_t> &nearest ) const;
167
168
170
184 template<class TYPE2>
185 std::tuple<AMP::Array<double>, AMP::Array<double>> interp_linear( const AMP::Array<double> &f,
186 const AMP::Array<TYPE2> &xi,
187 const AMP::Array<int> &index,
188 bool extrap = false ) const;
189
190
192
210 template<class TYPE2>
211 std::tuple<AMP::Array<double>, AMP::Array<double>> interp_cubic( const AMP::Array<double> &f,
212 const AMP::Array<double> &g,
213 const AMP::Array<TYPE2> &xi,
214 const AMP::Array<int> &index,
215 int extrap = 0 ) const;
216
217
219
227 static void compute_Barycentric( const int ndim, const double *x, const double *xi, double *L );
228
229
231 void clear();
232
233
234private:
237
240
242 void create_node_tri() const;
243
245 void create_kdtree() const;
246
247 // Subroutine to perform cubic interpolation for a single point
248 void interp_cubic_single( const double f[],
249 const double g[],
250 const double xi[],
251 const int index,
252 double &fi,
253 double *gi,
254 int extrap ) const;
255
256
257private: // Internal Data
258 Array<TYPE> d_x; // Pointer to the coordinates (ndim x N)
259 Array<int> d_tri; // Pointer to the coordinates (ndim+1 x N_tri)
260 mutable Array<int> d_tri_nab; // List of neighbor triangles ( ndim+1 x N_tri )
261 mutable size_t d_N_node_sum; // The sum of the number of node neighbors
262 mutable unsigned *d_N_node; // The number of neighbor nodes for each node (1xN)
263 mutable unsigned **d_node_list; // The list of neighbor nodes for each node (1xN)
264 // Note: The first element points to an array of size N_node_sum
265 mutable int *d_node_tri; // For each node, a triangle that contains that node
266 mutable kdtree *d_tree; // Nearest neighbor search tree
267};
268
269
270} // namespace AMP
271
272#endif
void create_node_tri() const
Subroutine to get the starting triangle for each node.
void interp_cubic_single(const double f[], const double g[], const double xi[], const int index, double &fi, double *gi, int extrap) const
void create_node_neighbors() const
Subroutine to get the list of nodes that are connected to each node.
void create_tri_neighbors() const
Subroutine to get the list of triangle that are connected to each triangle.
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()
Empty destructor.
DelaunayInterpolation()
Empty constructor.
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.
Array< int > find_tri(const Array< TYPE2 > &xi, bool extrap=false) const
Subroutine to find the triangle that contains the point.
void create_tessellation(const AMP::Array< TYPE > &x)
Function to construct the tessellation.
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.
void create_kdtree() const
Subroutine to get the list of triangle that are connected to each triangle.
AMP::Array< int > get_tri_nab() const
Function to return the triangles neignbors.
void clear()
Clear the data.
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.
DelaunayInterpolation(const DelaunayInterpolation &)=delete
DelaunayInterpolation & operator=(const DelaunayInterpolation &)=delete
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.
std::tuple< AMP::Array< TYPE >, AMP::Array< int > > copy_tessellation() const
Function to copy the tessellation.
Array< size_t > find_nearest(const Array< TYPE2 > &xi) const
Subroutine to find the nearest neighbor to a point.
static void compute_Barycentric(const int ndim, const double *x, const double *xi, double *L)
Subroutine to compute the Barycentric coordinates.
AMP::Array< int > get_tri() const
Function to return the triangles in the tessellation.
size_t get_N_tri() const
Function to return the number of triangles in the tessellation.
A class used to to perform kd-tree based operations.
Definition kdtree.h:33



Advanced Multi-Physics (AMP)
Oak Ridge National Laboratory
Idaho National Laboratory
Los Alamos National Laboratory
This page automatically produced from the
source code by doxygen
Last updated: Tue Mar 10 2026 13:06:41.
Comments on this page