GInX
Geodesics Integrator tool for particles in GRMHD using Adaptive Mesh Refinement using AMReX.
GInX::ParticlesContainer< StructType > Class Template Reference

ParticlesContainer class definition. More...

#include <ParticlesContainer.hxx>

+ Inheritance diagram for GInX::ParticlesContainer< StructType >:
+ Collaboration diagram for GInX::ParticlesContainer< StructType >:

Public Types

using Base = BaseParticleContainer< ParticlesContainer< StructType >, StructType >
 Using BaseParticlesContainer constructor. More...
 

Public Member Functions

 ParticlesContainer (amrex::AmrCore *amr_core, const CCTK_REAL m)
 
 ~ParticlesContainer ()=default
 
void evolve (const amrex::MultiFab &lapse, const amrex::MultiFab &shift, const amrex::MultiFab &metric, const amrex::MultiFab &curv, const CCTK_REAL &dt, const int &lev) override
 Evolving using Runge-Kutta 4. More...
 
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE CCTK_ATTRIBUTE_ALWAYS_INLINE amrex::GpuArray< CCTK_REAL, 7 > compute_rhs (const amrex::GpuArray< CCTK_REAL, 7 > &u, const CCTK_REAL &t, amrex::Array4< CCTK_REAL const > const &lapse, amrex::Array4< CCTK_REAL const > const &shift, amrex::Array4< CCTK_REAL const > const &metric, amrex::Array4< CCTK_REAL const > const &K, const CCTK_REAL dt, const amrex::GpuArray< double, 3 > &dx, const int lev, const amrex::GpuArray< double, 3 > &plo)
 Computes the right hand side of the geodesic differential equation. More...
 
void normalize_velocity (const amrex::MultiFab &metric, const int level)
 Normalize the velocity accordingly to the metric on each particle position. More...
 
CCTK_INT get_mass ()
 
void redistribute_particles ()
 
- Public Member Functions inherited from GInX::BaseParticleContainer< ParticlesContainer< StructType >, StructType >
 BaseParticleContainer (amrex::AmrCore *amr_core)
 
virtual ~BaseParticleContainer ()=default
 
void initialize (Function initializer_function, const CCTK_REAL *real_params, const CCTK_INT *int_params)
 Initialize the particles given a custom initialization function. More...
 
void initialize (Function initializer_function, const amrex::MultiFab &metric, const int &level, const CCTK_REAL *real_params, const CCTK_INT *int_params)
 Initialize the particles given a custom initialization function. More...
 
virtual void evolve (const amrex::MultiFab &lapse, const amrex::MultiFab &shift, const amrex::MultiFab &metric, const amrex::MultiFab &curv, const CCTK_REAL &dt, const int &lev)=0
 
void check_banned_zones (const int &level, const CCTK_INT4 &zones, const CCTK_REAL(&x)[10], const CCTK_REAL(&y)[10], const CCTK_REAL(&z)[10], const CCTK_REAL(&radius)[10])
 
void outputParticlesAscii (const int &it, const int &plot_every, const std::string &out_dir)
 
void outputParticlesPlot (const int &it, const int &plot_every, const std::string &out_dir)
 

Protected Attributes

CCTK_REAL mass = 0.
 

Additional Inherited Members

- Public Attributes inherited from GInX::BaseParticleContainer< ParticlesContainer< StructType >, StructType >
const std::string name
 
- Static Public Attributes inherited from GInX::BaseParticleContainer< ParticlesContainer< StructType >, StructType >
static constexpr int n_attributes
 

Detailed Description

template<typename StructType>
class GInX::ParticlesContainer< StructType >

ParticlesContainer class definition.

The following class define the needed functions to evolve the position and velocity of the photons in the simulation.

Member Typedef Documentation

◆ Base

template<typename StructType >
using GInX::ParticlesContainer< StructType >::Base = BaseParticleContainer<ParticlesContainer<StructType>, StructType>

Using BaseParticlesContainer constructor.

Constructor & Destructor Documentation

◆ ParticlesContainer()

template<typename StructType >
GInX::ParticlesContainer< StructType >::ParticlesContainer ( amrex::AmrCore *  amr_core,
const CCTK_REAL  m 
)
inline

◆ ~ParticlesContainer()

template<typename StructType >
GInX::ParticlesContainer< StructType >::~ParticlesContainer ( )
default

Member Function Documentation

◆ compute_rhs()

template<typename StructType >
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE CCTK_ATTRIBUTE_ALWAYS_INLINE amrex::GpuArray< CCTK_REAL, 7 > GInX::ParticlesContainer< StructType >::compute_rhs ( const amrex::GpuArray< CCTK_REAL, 7 > &  u,
const CCTK_REAL &  t,
amrex::Array4< CCTK_REAL const > const &  lapse,
amrex::Array4< CCTK_REAL const > const &  shift,
amrex::Array4< CCTK_REAL const > const &  metric,
amrex::Array4< CCTK_REAL const > const &  K,
const CCTK_REAL  dt,
const amrex::GpuArray< double, 3 > &  dx,
const int  lev,
const amrex::GpuArray< double, 3 > &  plo 
)

Computes the right hand side of the geodesic differential equation.

Given differential equation

\[\frac{d}{dt}U = f\left(U, \frac{dU}{dx}; t\right)\]

computes

\[f\left(U, \frac{dU}{dx}; t\right)\]

where \(U\) is a vector which contains \((x, y, z, v_x, v_y, v_z, \ln E)\). The differential equation for the particles' position is

\[\frac{d}{dt} U[i] = \alpha \gamma^{ij} U[3 + j] - \beta^i\]

Where \(i,j = 0, 1, 2\), \(\gamma\) is the induced metric, \(\alpha\) is the lapse function and \(\beta\) is the shift vector.

For the Velocity_d the differential equation is:

\begin{eqnarray*} \frac{d}{dt}U[3 + i] &= -\partial_i\alpha + \left(\gamma^{kj} U[3 + k] \partial_j\alpha - \alpha K_{jk}\gamma^{jl}\gamma^{km}U[3+l]U[3+m]\right) U[3 + i]\\ & + \frac{1}{2}\alpha\gamma^{jl}\gamma^{km}U[3+l]U[3+m]\partial_i\gamma{jk} + U[3 + j] \partial_i\beta^j \end{eqnarray*}

and finally, for the \( \ln E \) the differential equation is:

\[ \dfrac{d}{dt} U[6] = \alpha K_{jk}U[3 + l]U[3 + m]\gamma^{lj}\gamma^{mk} - U[3+l]\gamma^{lj}\partial_j\alpha\]

Where \(i, j, k, l, m = 0, 1, 2\) and \(K_{ij}\) is the extrinsic curvature. We have been using Einstein notation.

Parameters
uA GpuArray of size n_attributes + the coordinates that contains the varaibles needed to evolve.
tCurrent time t.
lapseADM lapse function.
shiftShift vector \beta^i
metric3 dimensional ADM metric.
curvExtrinsic curvature.
dtTimestep.
dxSpacestep
levAMR Level of discretization.
ploPhysical lower bounds of the whole domain.
Returns
The right hind side of the differential equation.

◆ evolve()

template<typename StructType >
void GInX::ParticlesContainer< StructType >::evolve ( const amrex::MultiFab &  lapse,
const amrex::MultiFab &  shift,
const amrex::MultiFab &  metric,
const amrex::MultiFab &  curv,
const CCTK_REAL &  dt,
const int &  lev 
)
overridevirtual

Evolving using Runge-Kutta 4.

We are solving the differential equation \(\frac{dU}{dt} = f\left(U, \frac{dU}{dx}, t\right)\) using:

\[ U_{n+1} = U_n + \frac{1}{6}\Delta t \left(f_1 + 2f_2 + 2f_3 + f_4\right) \]

where:

  • \(f_1 = f(U_n, t),\)
  • \(f_2 = f\left(U_n + \frac{\Delta t}{2} f_1, t + \frac{\Delta t}{2}\right),\)
  • \(f_3 = f\left(U_n + \frac{\Delta t}{2} f_2, t + \frac{\Delta t}{2}\right),\)
  • \(f_4 = f(U_n + \Delta t f_3, t + \Delta t),\)

While computing we are checking if the particles still in the physical domain.

See also
compute_rhs()
Parameters
lapseADM lapse function.
shiftADM shift vector.
metricADM induced metric.
curvExtrinsic curvature.
dtTimestep.
levRefinement level.

Implements GInX::BaseParticleContainer< ParticlesContainer< StructType >, StructType >.

◆ get_mass()

template<typename StructType >
CCTK_INT GInX::ParticlesContainer< StructType >::get_mass ( )
inline

◆ normalize_velocity()

template<typename StructType >
void GInX::ParticlesContainer< StructType >::normalize_velocity ( const amrex::MultiFab &  metric,
const int  level 
)

Normalize the velocity accordingly to the metric on each particle position.

This function normalizes the velocity given a random initial data. The invariant mass is described by

\[ P^\mu P_\mu = 0, \]

for null geodesics and

\[ P^\mu P_\mu = -m^2, \]

for timelike particles.

Implying

\[ V^\alpha V_\alpha = V_\alpha V_\beta \gamma^{\alpha\beta} = 1, \]

for null geodesic and

\[ V^\alpha V_\alpha = V_\alpha V_\beta \gamma^{\alpha\beta} = 1 -\frac{m^2}{E^2} \]

for timelike particles.

Parameters
metricADM 3 dimension metric.
Currentrefinement level.

◆ redistribute_particles()

template<typename StructType >
void GInX::ParticlesContainer< StructType >::redistribute_particles

Member Data Documentation

◆ mass

template<typename StructType >
CCTK_REAL GInX::ParticlesContainer< StructType >::mass = 0.
protected

The documentation for this class was generated from the following file: