|
GInX
Geodesics Integrator tool for particles in GRMHD using Adaptive Mesh Refinement using AMReX.
|
Barycentric Lagrange interpolation related functions definitions. More...
#include "AMReX_Box.H"#include "AMReX_Config.H"#include "AMReX_MFIter.H"#include "AMReX_REAL.H"#include "AMReX_Random.H"#include "AMReX_RandomEngine.H"#include "AMReX_Scan.H"#include "cctk_Types.h"#include <array>#include <cctk.h>#include <iostream>
Include dependency graph for Interpolator.hxx:
This graph shows which files directly or indirectly include this file:Go to the source code of this file.
Namespaces | |
| namespace | GInX |
| Interpolators namespace. | |
Functions | |
| template<int N> | |
| AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE CCTK_ATTRIBUTE_ALWAYS_INLINE void | GInX::barycentric_cubic_1d (CCTK_REAL &value, const int(&points)[N], const CCTK_REAL *weights, const CCTK_REAL(&values)[N], const CCTK_REAL &x, const CCTK_REAL &plo, const CCTK_REAL &dx) |
| Do a Barycentric interpolation in one direction. More... | |
| template<int INTERPOLATION_ORDER> | |
| AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE CCTK_ATTRIBUTE_ALWAYS_INLINE CCTK_REAL | GInX::barycentric_cubic_3d (const amrex::Array4< CCTK_REAL const > &f, const long int &i0, const long int &j0, const long int &k0, const CCTK_REAL &x, const CCTK_REAL &y, const CCTK_REAL &z, const amrex::GpuArray< double, 3 > &dx, const amrex::GpuArray< double, 3 > &plo, const int &comp) |
| Do a Barycentric interpolation in three direction for a vectorial function. More... | |
| template<int INTERPOLATION_ORDER> | |
| AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE CCTK_ATTRIBUTE_ALWAYS_INLINE CCTK_REAL | GInX::barycentric_cubic_3d (const amrex::Array4< CCTK_REAL const > &f, const long int &i0, const long int &j0, const long int &k0, const CCTK_REAL &x, const CCTK_REAL &y, const CCTK_REAL &z, const amrex::GpuArray< double, 3 > &dx, const amrex::GpuArray< double, 3 > &plo) |
| Do a Barycentric interpolation in three direction for a scalar function. More... | |
| template<int N> | |
| AMREX_GPU_DEVICE AMREX_GPU_HOST AMREX_FORCE_INLINE CCTK_ATTRIBUTE_ALWAYS_INLINE void | GInX::der_barycentric_cubic_1d (CCTK_REAL &f_x, CCTK_REAL &d_f_x, const int(&points)[N], const CCTK_REAL *weights, const CCTK_REAL(&values)[N], const CCTK_REAL &x, const CCTK_REAL &plo, const CCTK_REAL &dx) |
| Do a Barycentric interpolation and its derivative in one direction. More... | |
| template<int INTERPOLATION_ORDER> | |
| AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE CCTK_ATTRIBUTE_ALWAYS_INLINE void | GInX::barycentric_derivative_and_interpolate (CCTK_REAL &f_xyz, CCTK_REAL &df_xyz_0, CCTK_REAL &df_xyz_1, CCTK_REAL &df_xyz_2, amrex::Array4< CCTK_REAL const > const &f, const long int &i0, const long int &j0, const long int &k0, const CCTK_REAL &x, const CCTK_REAL &y, const CCTK_REAL &z, const amrex::GpuArray< double, 3 > &dx, const amrex::GpuArray< double, 3 > &plo, const int &comp) |
| Do a Barycentric interpolation and its derivatives in three directions. More... | |
| template<int INTERPOLATION_ORDER> | |
| AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE CCTK_ATTRIBUTE_ALWAYS_INLINE void | GInX::interpolate_array (amrex::GpuArray< CCTK_REAL, 6 > &array6, const amrex::Array4< CCTK_REAL const > &f, const long int &i0, const long int &j0, const long int &k0, const CCTK_REAL &x, const CCTK_REAL &y, const CCTK_REAL &z, const amrex::GpuArray< double, 3 > &dx, const amrex::GpuArray< double, 3 > &plo) |
| Do the barycentric interpolation for a 3 dimensional and symmetric matrix. More... | |
| template<int INTERPOLATION_ORDER> | |
| AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE CCTK_ATTRIBUTE_ALWAYS_INLINE void | GInX::d_interpolate_array (amrex::GpuArray< CCTK_REAL, 6 > &array6, amrex::GpuArray< amrex::GpuArray< CCTK_REAL, 6 >, 3 > &d_array6, const amrex::Array4< CCTK_REAL const > &f, const long int &i0, const long int &j0, const long int &k0, const CCTK_REAL &x, const CCTK_REAL &y, const CCTK_REAL &z, const amrex::GpuArray< double, 3 > &dx, const amrex::GpuArray< double, 3 > &plo) |
| Do the barycentric interpolation for a 3 dimensional and symmetric matrix and its derivatives. More... | |
| template<int INTERPOLATION_ORDER> | |
| AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE CCTK_ATTRIBUTE_ALWAYS_INLINE void | GInX::d_interpolate_array (amrex::GpuArray< CCTK_REAL, 3 > &array3, amrex::GpuArray< amrex::GpuArray< CCTK_REAL, 3 >, 3 > &d_array3, const amrex::Array4< CCTK_REAL const > &f, const long int &i0, const long int &j0, const long int &k0, const CCTK_REAL &x, const CCTK_REAL &y, const CCTK_REAL &z, const amrex::GpuArray< double, 3 > &dx, const amrex::GpuArray< double, 3 > &plo) |
| Do the barycentric interpolation for a 3 dimensional and symmetric vector and its derivatives. More... | |
| template<int INTERPOLATION_ORDER> | |
| AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE CCTK_ATTRIBUTE_ALWAYS_INLINE void | GInX::d_interpolate_array (CCTK_REAL &array, amrex::GpuArray< CCTK_REAL, 3 > &d_array, const amrex::Array4< CCTK_REAL const > &f, const long int &i0, const long int &j0, const long int &k0, const CCTK_REAL &x, const CCTK_REAL &y, const CCTK_REAL &z, const amrex::GpuArray< double, 3 > &dx, const amrex::GpuArray< double, 3 > &plo) |
| Do the barycentric interpolation for a scalar function its derivatives. More... | |
Barycentric Lagrange interpolation related functions definitions.
This file define the Barycentric interpolation for function values and its derivatives in 1 and 3 dimensions.