ImpactX
Namespaces | Classes | Functions
impactx::initialization Namespace Reference

Namespaces

 details
 

Classes

class  AmrCoreData
 
struct  InitSingleParticleData
 

Functions

AmrCoreData init_amr_core ()
 
AmrCoreData amrex_amrcore_gridding ()
 
AmrCoreData one_box_per_rank ()
 
constexpr amrex::RealVect undefined_geometry_prob_lo (-1.0)
 the fake domain size (lower end), used to indicate that the user made no choice for it More...
 
constexpr amrex::RealVect undefined_geometry_prob_hi (1.0)
 the fake domain size (higher end), used to indicate that the user made no choice for it More...
 
void default_init_AMReX (int argc, char *argv[])
 
void default_init_AMReX ()
 
void set_distribution_parameters_from_phase_space_inputs (amrex::ParmParse const &pp_dist, amrex::ParticleReal &lambdax, amrex::ParticleReal &lambday, amrex::ParticleReal &lambdat, amrex::ParticleReal &lambdapx, amrex::ParticleReal &lambdapy, amrex::ParticleReal &lambdapt, amrex::ParticleReal &muxpx, amrex::ParticleReal &muypy, amrex::ParticleReal &mutpt)
 
void set_distribution_parameters_from_twiss_inputs (amrex::ParmParse const &pp_dist, amrex::ParticleReal &lambdax, amrex::ParticleReal &lambday, amrex::ParticleReal &lambdat, amrex::ParticleReal &lambdapx, amrex::ParticleReal &lambdapy, amrex::ParticleReal &lambdapt, amrex::ParticleReal &muxpx, amrex::ParticleReal &muypy, amrex::ParticleReal &mutpt)
 
amrex::Vector< amrex::Real > read_mr_prob_relative ()
 
void overwrite_amrex_parser_defaults ()
 

Function Documentation

◆ amrex_amrcore_gridding()

AmrCoreData impactx::initialization::amrex_amrcore_gridding ( )

This lets AMReX AmrCore/AmrMesh build the boxes for each MPI-rank

This uses the regular logic in AmrCore/AmrMesh to grid boxes and refinement levels.

Returns
simulation_geometry the geometry (topology) of the simulation; amr_info contains information on mesh-refinement and box/grid blocks

◆ default_init_AMReX() [1/2]

void impactx::initialization::default_init_AMReX ( )

Initialize AMReX

Initializes AMReX if not already done.

◆ default_init_AMReX() [2/2]

void impactx::initialization::default_init_AMReX ( int  argc,
char *  argv[] 
)

Initialize AMReX

Initializes AMReX if not already done.

◆ init_amr_core()

AmrCoreData impactx::initialization::init_amr_core ( )

This creates AMReX boxes for each MPI-rank

Returns
simulation_geometry the geometry (topology) of the simulation; amr_info contains information on mesh-refinement and box/grid blocks

◆ one_box_per_rank()

AmrCoreData impactx::initialization::one_box_per_rank ( )

This builds one AMReX box per MPI-rank

This is a simple decomposition for particles that we default to if we do not need space charge or load balancing.

Returns
simulation_geometry the geometry (topology) of the simulation; amr_info contains information on mesh-refinement and box/grid blocks

◆ overwrite_amrex_parser_defaults()

void impactx::initialization::overwrite_amrex_parser_defaults ( )

Overwrite defaults in AMReX Inputs

This overwrites defaults in amrex::ParamParse for inputs.

◆ read_mr_prob_relative()

amrex::Vector<amrex::Real> impactx::initialization::read_mr_prob_relative ( )

◆ set_distribution_parameters_from_phase_space_inputs()

void impactx::initialization::set_distribution_parameters_from_phase_space_inputs ( amrex::ParmParse const &  pp_dist,
amrex::ParticleReal &  lambdax,
amrex::ParticleReal &  lambday,
amrex::ParticleReal &  lambdat,
amrex::ParticleReal &  lambdapx,
amrex::ParticleReal &  lambdapy,
amrex::ParticleReal &  lambdapt,
amrex::ParticleReal &  muxpx,
amrex::ParticleReal &  muypy,
amrex::ParticleReal &  mutpt 
)

Initialize the input parameters for all distributions that read phase space ellipse parameters from the input

This function sets the distribution parameters based on the provided phase space inputs. The parameters include the phase space ellipse intersections for position (lambdaX, lambdaY), time (lambdaT), momentum (lambdaPx, lambdaPy) and energy (lambdaPt), as well as the correlation terms (muxpx, muypy, mutpt).

Parameters
pp_distThe parameter parser object.
lambdaxReference to the variable storing the axis intersection of the phase space ellipse for x position.
lambdayReference to the variable storing the axis intersection of the phase space ellipse for y position.
lambdatReference to the variable storing the axis intersection of the phase space ellipse for time.
lambdapxReference to the variable storing the axis intersection of the phase space ellipse for momentum in x direction.
lambdapyReference to the variable storing the axis intersection of the phase space ellipse for momentum in y direction.
lambdaptReference to the variable storing the axis intersection of the phase space ellipse for energy.
muxpxReference to the variable storing the mean for momentum in x direction.
muypyReference to the variable storing the mean for momentum in y direction.
mutptReference to the variable storing the mean for energy.

◆ set_distribution_parameters_from_twiss_inputs()

void impactx::initialization::set_distribution_parameters_from_twiss_inputs ( amrex::ParmParse const &  pp_dist,
amrex::ParticleReal &  lambdax,
amrex::ParticleReal &  lambday,
amrex::ParticleReal &  lambdat,
amrex::ParticleReal &  lambdapx,
amrex::ParticleReal &  lambdapy,
amrex::ParticleReal &  lambdapt,
amrex::ParticleReal &  muxpx,
amrex::ParticleReal &  muypy,
amrex::ParticleReal &  mutpt 
)

Set the distribution parameters from Twiss inputs

This function reads Courant-Snyder / Twiss parameters from the provided ParmParse object and calculates the distribution parameters for a particle beam. It sets the values for the axis intercepts of the phase space ellipse (lambdaX, lambdaY, lambdaT, and lambdaPx, lambdaPy, lambdaPt) and correlation terms (muxpx, muypy, mutpt) normally accepted as input.

Parameters
pp_distThe parameter parser object.
lambdaxReference to the variable storing the axis intersection of the phase space ellipse for x position.
lambdayReference to the variable storing the axis intersection of the phase space ellipse for y position.
lambdatReference to the variable storing the axis intersection of the phase space ellipse for time.
lambdapxReference to the variable storing the axis intersection of the phase space ellipse for momentum in x direction.
lambdapyReference to the variable storing the axis intersection of the phase space ellipse for momentum in y direction.
lambdaptReference to the variable storing the axis intersection of the phase space ellipse for energy.
muxpxReference to the variable storing the mean for momentum in x direction.
muypyReference to the variable storing the mean for momentum in y direction.
mutptReference to the variable storing the mean for energy.

◆ undefined_geometry_prob_hi()

constexpr amrex::RealVect impactx::initialization::undefined_geometry_prob_hi ( 1.  0)
constexpr

the fake domain size (higher end), used to indicate that the user made no choice for it

◆ undefined_geometry_prob_lo()

constexpr amrex::RealVect impactx::initialization::undefined_geometry_prob_lo ( -1.  0)
constexpr

the fake domain size (lower end), used to indicate that the user made no choice for it