GInX
Geodesics Integrator tool for particles in GRMHD using Adaptive Mesh Refinement using AMReX.
GInX Namespace Reference

Interpolators namespace. More...

Classes

class  BaseParticleContainer
 BaseParticleContainer abstract class definition. More...
 
class  ParticleIterator
 
class  ParticlesContainer
 ParticlesContainer class definition. More...
 
struct  PhotonsData
 Struct to managing null geodesic particles. More...
 

Functions

template<int N>
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE CCTK_ATTRIBUTE_ALWAYS_INLINE void 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 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 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 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 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 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 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 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 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...
 

Detailed Description

Interpolators namespace.

Function Documentation

◆ barycentric_cubic_1d()

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.

This function computes the interpolation of a function on just one direction using a barycentric Lagrange interpolation by doing:

\[ f(x) = \frac{\sum\limits_{i = 0}^N \frac{w_i}{x - x_i}f(x_i)}{\sum\limits_{l = 0}^N \frac{w_l}{x - x_l}} \]

where \(N\) is the order of interpolation.

Parameters
valueThe interpolated value reference.
pointsVector containing the coordinates of each point.
weightsThe weights of each datapoint.
xThe value where we are interpolating.
dxVector \(\Delta x\) on the particular direction.
ploLower coordinates in the box domain.

◆ barycentric_cubic_3d() [1/2]

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.

This function computes the interpolation of a function on three directions using a barycentric Lagrange interpolation by doing:

\[ f(x, y, z) = \frac{\sum\limits_{i,j,k = 0}^N \frac{u_i}{z - z_i}\frac{v_j}{y - y_j}\frac{w_k}{x - x_k}f(x_k, y_j, x_i)}{\sum\limits_{i,j,k = 0}^N \frac{u_i}{z - z_i}\frac{v_j}{y - y_j}\frac{w_k}{x - x_k}} = \frac{\sum\limits_{i=0}^N\frac{u_i}{z - z_i}\left(\frac{\sum\limits_{j=0}^N\frac{v_j}{y-y_j}\left(\frac{\sum\limits_{k=0}^N \frac{w_k}{x - x_k}f(x_k, y_j, x_i)}{\sum\limits_{k= 0}^N \frac{w_k}{x - x_k}}\right)}{\sum\limits_{j= 0}^N \frac{v_j}{y - y_j}}\right)}{\sum\limits_{i= 0}^N \frac{u_i}{z - z_i}} \]

where \(N\) is the order of interpolation.

See also
barycentric_cubic_1d
Parameters
fArray containing the function values.
i0Basis cell index accordingly to x.
j0Basis cell index accordingly to y.
k0Basis cell index accordingly to z.
xCoordinate x value to interpolate.
yCoordinate y value to interpolate.
zCoordinate z value to interpolate.
dxVector \(\Delta x\) with the space steps value.
ploLower values of the entire domain.
Returns
The interpolated value.

◆ barycentric_cubic_3d() [2/2]

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.

Parameters
fArray containing the function values.
i0Basis cell index accordingly to x.
j0Basis cell index accordingly to y.
k0Basis cell index accordingly to z.
xCoordinate x value to interpolate.
yCoordinate y value to interpolate.
zCoordinate z value to interpolate.
dxVector with the \(\Delta x\) on each direction.
ploLower coordinates in the domain.
compFunction's component to compute.
Returns
The interpolated value.

◆ barycentric_derivative_and_interpolate()

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.

This function computes the interpolation of a function on three directions using a barycentric Lagrange interpolation by doing:

\[ f(x, y, z) = \frac{\sum\limits_{i,j,k = 0}^N \frac{u_i}{z - z_i}\frac{v_j}{y - y_j}\frac{w_k}{x - x_k}f(x_k, y_j, x_i)}{\sum\limits_{i,j,k = 0}^N \frac{u_i}{z - z_i}\frac{v_j}{y - y_j}\frac{w_k}{x - x_k}} = \frac{\sum\limits_{i=0}^N\frac{u_i}{z - z_i}\left(\frac{\sum\limits_{j=0}^N\frac{v_j}{y-y_j}\left(\frac{\sum\limits_{k=0}^N \frac{w_k}{x - x_k}f(x_k, y_j, x_i)}{\sum\limits_{k= 0}^N \frac{w_k}{x - x_k}}\right)}{\sum\limits_{j= 0}^N \frac{v_j}{y - y_j}}\right)}{\sum\limits_{i= 0}^N \frac{u_i}{z - z_i}} \]

where \(N\) is the order of interpolation.

But at the same time computes the gradient of the same function by doing:

\[ \frac{\partial}{\partial x}f(x, y, z) = \frac{\sum\limits_{i,j,k = 0}^N \frac{u_i}{z - z_i}\frac{v_j}{y - y_j}\frac{d}{dx}\left(\frac{\frac{w_k}{x - x_k}}{\sum\limits_{l = 0}^N\frac{w_l}{x - x_l}}\right)f(x_k, y_j, x_i)}{\sum\limits_{i,j = 0}^N \frac{u_i}{z - z_i}\frac{v_j}{y - y_j}}, \]

\[ \frac{\partial}{\partial y}f(x, y, z) = \frac{\sum\limits_{i,j,k = 0}^N \frac{u_i}{z - z_i}\frac{d}{dy}\left(\frac{\frac{v_j}{y - y_j}}{\sum\limits_{l = 0}^N\frac{v_l}{y - y_l}}\right)\frac{w_k}{x - x_k}f(x_k, y_j, x_i)}{\sum\limits_{i,k = 0}^N \frac{u_i}{z - z_i}\frac{w_k}{x - x_k}}, \]

\[ \frac{\partial}{\partial z}f(x, y, z) = \frac{\sum\limits_{i,j,k = 0}^N \frac{d}{dz}\left(\frac{\frac{u_i}{z - z_i}}{\sum\limits_{l = 0}^N\frac{u_l}{z - z_l}}\right)\frac{v_j}{y - y_j}\frac{w_k}{x - x_k}f(x_k, y_j, x_i)}{\sum\limits_{i,k = 0}^N \frac{v_j}{y - y_j}\frac{w_k}{x - x_k}}. \]

by calling the function der_barycentric_cubic_1d().

See also
der_barycentric_cubic_1d
Parameters
f_xyzReference to the interpolated value, i.e., \(f(x,y,z)\).
df_xyz_0Reference to the interpolated derivative value on the direction 0, i.e., \(\partial_xf(x,y,z)\).
df_xyz_1Reference to the interpolated derivative value on the direction 1, i.e., \(\partial_yf(x,y,z)\).
df_xyz_2Reference to the interpolated derivative value on the direction 2, i.e., \(\partial_zf(x,y,z)\).
fArray containing the function values.
i0Basis cell index accordingly to x.
j0Basis cell index accordingly to y.
k0Basis cell index accordingly to z.
xCoordinate x value to interpolate.
yCoordinate y value to interpolate.
zCoordinate z value to interpolate.
dxVector \(\Delta x\) with the space steps value.
ploLower values of the entire domain.
compFunction's component to compute.

◆ d_interpolate_array() [1/3]

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.

This function computes the interpolation for a 3 dimensional and symmetric vector and its derivatives, such as the induced shift vector in just one call using the der_barycentric_cubic_1d function similarly to the barycentric_derivative_and_interpolate function.

See also
barycentric_derivative_and_interpolate
der_barycentric_cubic_1d
Parameters
array6Reference to the GpuArray of size 3 that is going to store the interpolated values.
d_array6Reference to GpuArray of size 3 x 3 that stores all the derivatives values.
fArray containing the function values.
i0Basis cell index accordingly to x.
j0Basis cell index accordingly to y.
k0Basis cell index accordingly to z.
xCoordinate x value to interpolate.
yCoordinate y value to interpolate.
zCoordinate z value to interpolate.
dxVector \(\Delta x\) with the space steps value.
ploLower values of the entire domain.

◆ d_interpolate_array() [2/3]

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.

This function computes the interpolation for a 3 dimensional and symmetric matrix and its derivatives, such as the induced metric and the extrinsic curvature tensor in just one call using the der_barycentric_cubic_1d function similarly to the barycentric_derivative_and_interpolate function.

See also
barycentric_derivative_and_interpolate
der_barycentric_cubic_1d
Parameters
array6Reference to the GpuArray of size 6 that is going to store the interpolated values.
d_array6Reference to the GpuArray of size 6 x 3 that stores all the derivatives values.
fArray containing the function values.
i0Basis cell index accordingly to x.
j0Basis cell index accordingly to y.
k0Basis cell index accordingly to z.
xCoordinate x value to interpolate.
yCoordinate y value to interpolate.
zCoordinate z value to interpolate.
dxVector \(\Delta x\) with the space steps value.
ploLower values of the entire domain.

◆ d_interpolate_array() [3/3]

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.

This function computes the interpolation for scalar function and its derivatives, such as the lapse function in just one call using the der_barycentric_cubic_1d function similar to the barycentric_derivative_and_interpolate function.

See also
barycentric_derivative_and_interpolate
der_barycentric_cubic_1d
Parameters
array6reference to the CCTK_REAL variable that is going to store the interpolated value.
d_array6Reference to the GpuArray of size 3 that stores all the derivatives values.
fArray containing the function values.
i0Basis cell index accordingly to x.
j0Basis cell index accordingly to y.
k0Basis cell index accordingly to z.
xCoordinate x value to interpolate.
yCoordinate y value to interpolate.
zCoordinate z value to interpolate.
dxVector \(\Delta x\) with the space steps value.
ploLower values of the entire domain.

◆ der_barycentric_cubic_1d()

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.

This function computes the interpolation of a function on just one direction using a barycentric Lagrange interpolation by doing:

\[ f(x) = \frac{\sum\limits_{i = 0}^N \frac{w_i}{x - x_i}f(x_i)}{\sum\limits_{i = 0}^N \frac{w_i}{x - x_i}} \]

where \(N\) is the order of interpolation.

And at the same time computes the derivative on the desired direction by doing:

\[ f'(x) = \sum\limits_{i = 0}^N f(x_i)\frac{d}{dx}\left(\frac{\frac{w_i}{x - x_i}}{\sum\limits_{i = 0}^N \frac{w_i}{x - x_i}}\right) = \frac{\sum\limits_{j = 0}^N f(x_j)\frac{w_j}{(x-x_j)}\sum\limits_{m = 0}^N\frac{w_m}{(x-x_m)^2} - \sum\limits_{j = 0}^N f(x_j)\frac{w_j}{(x-x_j)^2}\sum\limits_{m = 0}^N \frac{w_m}{x - x_m}}{\left(\sum\limits_{l = 0}^N \frac{w_l}{x - x_l}\right)^2} \]

In the case \(x = x_i\), one of the nodes, we have the expression

\[ f'(x) = \sum\limits_{j=0}^{N}\frac{w_j}{w_i}\frac{f(x_i)}{x_i-x_j}. \]

See also
barycentric_cubic_3d
barycentric_cubic_1d
Parameters
f_xReference to the final interpolated function value.
d_f_xReference to the final interpolated derivative function value.
pointsVector containing the coordinates of each point.
weightsThe weights of each datapoint.
xThe value where we are interpolating.
ploLower values of the entire domain.
dxVector \(\Delta x\) with the space steps value.

◆ interpolate_array()

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.

This function computes the interpolation for a 3 dimensional and symmetric matrix, such as the induced metric and the extrinsic curvature tensor in just one call using the barycentric_cubic_1d function similar to the barycentric_cubic_3d function.

See also
barycentric_cubic_3d
barycentric_cubic_1d
Parameters
array6Reference to GpuArray of size 6 that is going to store the interpolated values.
fArray containing the function values.
i0Basis cell index accordingly to x.
j0Basis cell index accordingly to y.
k0Basis cell index accordingly to z.
xCoordinate x value to interpolate.
yCoordinate y value to interpolate.
zCoordinate z value to interpolate.
dxVector \(\Delta x\) with the space steps value.
ploLower values of the entire domain.