  | 
  
     Version: SMASH-1.8 
   | 
           
 | 
 
 
 
 
#include <collidermodus.h>
ColliderModus: Provides a modus for colliding nuclei.
To use this modus, choose 
 in the configuration file.
Options for ColliderModus go in the "Modi"→"Collider" section of the configuration.
The following configuration options are understood: Collider 
Definition at line 42 of file collidermodus.h.
 
◆ ColliderModus()
Constructor. 
Takes all there is to take from the (truncated!) configuration object (only contains configuration for this modus).
- Parameters
 - 
  
    | [in] | modus_config | The configuration object that sets all initial conditions of the experiment.  | 
    | [in] | parameters | Unused, but necessary because of templated initialization  | 
  
   
- Exceptions
 - 
  
    | ColliderEmpty | if projectile or nucleus are empty (i.e. do not contain particles)  | 
    | InvalidEnergy | if sqrts from config is not large enough to support the colliding masses of the nuclei, or if E_kin or P_lab are negative  | 
    | domain_error | if more or less than exactly one of the input energy options is specified, or if custom impact parameter Values and Yields are improperly supplied  | 
  
   
- Todo:
 - include a check that only one method of specifying impact is used 
 
Definition at line 264 of file collidermodus.cc.
  266   Configuration modus_cfg = modus_config[
"Collider"];
 
  268   if (modus_cfg.has_value({
"Calculation_Frame"})) {
 
  269     frame_ = modus_cfg.take({
"Calculation_Frame"});
 
  273   if (modus_cfg.has_value({
"Collisions_Within_Nucleus"})) {
 
  276   Configuration proj_cfg = modus_cfg[
"Projectile"];
 
  277   Configuration targ_cfg = modus_cfg[
"Target"];
 
  280   bool same_file = 
false;
 
  282   if (proj_cfg.has_value({
"Deformed"})) {
 
  285   } 
else if (proj_cfg.has_value({
"Custom"})) {
 
  288         make_unique<CustomNucleus>(proj_cfg, params.testparticles, same_file);
 
  290     projectile_ = make_unique<Nucleus>(proj_cfg, params.testparticles);
 
  293     throw ColliderEmpty(
"Input Error: Projectile nucleus is empty.");
 
  297   if (targ_cfg.has_value({
"Deformed"})) {
 
  299   } 
else if (targ_cfg.has_value({
"Custom"})) {
 
  301         make_unique<CustomNucleus>(targ_cfg, params.testparticles, same_file);
 
  303     target_ = make_unique<Nucleus>(targ_cfg, params.testparticles);
 
  306     throw ColliderEmpty(
"Input Error: Target nucleus is empty.");
 
  310   if (modus_cfg.has_value({
"Fermi_Motion"})) {
 
  319   int energy_input = 0;
 
  321   const double mass_target = 
target_->mass();
 
  323   const double mass_a =
 
  325   const double mass_b = 
target_->mass() / 
target_->number_of_particles();
 
  327   if (modus_cfg.has_value({
"Sqrtsnn"})) {
 
  331       throw ModusDefault::InvalidEnergy(
 
  332           "Input Error: sqrt(s_NN) is not larger than masses:\n" +
 
  333           std::to_string(
sqrt_s_NN_) + 
" GeV <= " + std::to_string(mass_a) +
 
  334           " GeV + " + std::to_string(mass_b) + 
" GeV.");
 
  338                    mass_projec * mass_target / (mass_a * mass_b) +
 
  339                mass_projec * mass_projec + mass_target * mass_target;
 
  344   if (modus_cfg.has_value({
"E_Kin"})) {
 
  345     const double e_kin = modus_cfg.take({
"E_Kin"});
 
  347       throw ModusDefault::InvalidEnergy(
 
  349           "E_Kin must be nonnegative.");
 
  353                            mass_projec, mass_target);
 
  358   if (modus_cfg.has_value({
"P_Lab"})) {
 
  359     const double p_lab = modus_cfg.take({
"P_Lab"});
 
  361       throw ModusDefault::InvalidEnergy(
 
  363           "P_Lab must be nonnegative.");
 
  367                            mass_projec, mass_target);
 
  371   if (energy_input == 0) {
 
  372     throw std::domain_error(
 
  373         "Input Error: Non-existent collision energy. " 
  374         "Please provide one of Sqrtsnn/E_Kin/P_Lab.");
 
  376   if (energy_input > 1) {
 
  377     throw std::domain_error(
 
  378         "Input Error: Redundant collision energy. " 
  379         "Please provide only one of Sqrtsnn/E_Kin/P_Lab.");
 
  384   if (modus_cfg.has_value({
"Impact", 
"Value"})) {
 
  385     impact_ = modus_cfg.take({
"Impact", 
"Value"});
 
  390     if (modus_cfg.has_value({
"Impact", 
"Sample"})) {
 
  391       sampling_ = modus_cfg.take({
"Impact", 
"Sample"});
 
  393         if (!(modus_cfg.has_value({
"Impact", 
"Values"}) ||
 
  394               modus_cfg.has_value({
"Impact", 
"Yields"}))) {
 
  395           throw std::domain_error(
 
  396               "Input Error: Need impact parameter spectrum for custom " 
  398               "Please provide Values and Yields.");
 
  400         const std::vector<double> impacts =
 
  401             modus_cfg.take({
"Impact", 
"Values"});
 
  402         const std::vector<double> yields = modus_cfg.take({
"Impact", 
"Yields"});
 
  403         if (impacts.size() != yields.size()) {
 
  404           throw std::domain_error(
 
  405               "Input Error: Need as many impact parameter values as yields. " 
  406               "Please make sure that Values and Yields have the same length.");
 
  409             InterpolateDataLinear<double>(impacts, yields));
 
  411         const auto imp_minmax =
 
  412             std::minmax_element(impacts.begin(), impacts.end());
 
  415         yield_max_ = *std::max_element(yields.begin(), yields.end());
 
  418     if (modus_cfg.has_value({
"Impact", 
"Range"})) {
 
  419       const std::array<double, 2> range = modus_cfg.take({
"Impact", 
"Range"});
 
  423     if (modus_cfg.has_value({
"Impact", 
"Max"})) {
 
  425       imp_max_ = modus_cfg.take({
"Impact", 
"Max"});
 
  431       modus_cfg.take({
"Impact", 
"Random_Reaction_Plane"}, 
false);
 
  433   if (modus_cfg.has_value({
"Initial_Distance"})) {
 
 
 
 
◆ custom_file_path()
      
        
          | std::string smash::ColliderModus::custom_file_path  | 
          ( | 
          const std::string &  | 
          file_directory,  | 
        
        
           | 
           | 
          const std::string &  | 
          file_name  | 
        
        
           | 
          ) | 
           |  | 
        
      
 
Creates full path string consisting of file_directory and file_name Needed to initialize a customnucleus. 
- Parameters
 - 
  
    | [in] | file_directory | is the path to the external file  | 
    | [in] | file_name | is the name of the external file  | 
  
   
Definition at line 620 of file collidermodus.cc.
  623   if (file_directory.back() == 
'/') {
 
  624     return file_directory + file_name;
 
  626     return file_directory + 
'/' + file_name;
 
 
 
 
◆ initial_conditions()
Generates initial state of the particles in the system. 
In particular, it initializes the momenta and positions of nucleons withing the colliding nuclei.
- Parameters
 - 
  
    | [out] | particles | An empty list that gets filled up by this function  | 
    | [in] | parameters | The initialization parameters of the system  | 
  
   
- Returns
 - The starting time of the simulation (negative, so that nuclei collide exactly at t=0) 
 
- Exceptions
 - 
  
    | domain_error | if the velocities of each nucleus are >= 1, or if input for Fermi motion is invalid  | 
  
   
Definition at line 470 of file collidermodus.cc.
  487   if (v_a >= 1.0 || v_b >= 1.0) {
 
  488     throw std::domain_error(
 
  489         "Found velocity equal to or larger than 1 in " 
  490         "ColliderModus::initial_conditions.\nConsider using " 
  491         "the center of velocity reference frame.");
 
  507     target_->generate_fermi_momenta();
 
  517     throw std::domain_error(
"Invalid Fermi_Motion input.");
 
  527   const double d_a = std::max(0., 
projectile_->get_diffusiveness());
 
  528   const double d_b = std::max(0., 
target_->get_diffusiveness());
 
  529   const double r_a = 
projectile_->get_nuclear_radius();
 
  530   const double r_b = 
target_->get_nuclear_radius();
 
  533   const double simulation_time = -dz / std::abs(v_a);
 
  534   const double proj_z = -dz - std::sqrt(1.0 - v_a * v_a) * (r_a + d_a);
 
  535   const double targ_z =
 
  536       +dz * std::abs(v_b / v_a) + std::sqrt(1.0 - v_b * v_b) * (r_b + d_b);
 
  546   target_->copy_particles(particles);
 
  548   return simulation_time;
 
 
 
 
◆ total_N_number()
  
  
      
        
          | int smash::ColliderModus::total_N_number  | 
          ( | 
           | ) | 
           const | 
         
       
   | 
  
inline   | 
  
 
- Returns
 - The total number of test particles in the initial nuclei 
 
Definition at line 89 of file collidermodus.h.
 
 
◆ proj_N_number()
  
  
      
        
          | int smash::ColliderModus::proj_N_number  | 
          ( | 
           | ) | 
           const | 
         
       
   | 
  
inline   | 
  
 
- Returns
 - The number of test particles in the projectile nucleus 
 
Definition at line 91 of file collidermodus.h.
 
 
◆ nuclei_passing_time()
  
  
      
        
          | double smash::ColliderModus::nuclei_passing_time  | 
          ( | 
           | ) | 
           const | 
         
       
   | 
  
inline   | 
  
 
Time until nuclei have passed through each other. 
Definition at line 94 of file collidermodus.h.
   95     const double passing_distance =
 
   97     const double passing_time =
 
 
 
 
◆ velocity_projectile()
  
  
      
        
          | double smash::ColliderModus::velocity_projectile  | 
          ( | 
           | ) | 
           const | 
         
       
   | 
  
inline   | 
  
 
- Returns
 - the beam velocity of the projectile, which will be used to calculate the beam momenta in experiment.cc if Fermi motion is frozen. 
 
Definition at line 109 of file collidermodus.h.
 
 
◆ velocity_target()
  
  
      
        
          | double smash::ColliderModus::velocity_target  | 
          ( | 
           | ) | 
           const | 
         
       
   | 
  
inline   | 
  
 
- Returns
 - the beam velocity of the target, which will be used to calculate the beam momenta in experiment.cc if Fermi motion is frozen. 
 
Definition at line 114 of file collidermodus.h.
 
 
◆ cll_in_nucleus()
  
  
      
        
          | bool smash::ColliderModus::cll_in_nucleus  | 
          ( | 
           | ) | 
           | 
         
       
   | 
  
inline   | 
  
 
- Returns
 - A flag: whether to allow first collisions within the same nucleus. 
 
Definition at line 118 of file collidermodus.h.
 
 
◆ fermi_motion()
◆ is_collider()
  
  
      
        
          | bool smash::ColliderModus::is_collider  | 
          ( | 
           | ) | 
           const | 
         
       
   | 
  
inline   | 
  
 
- Returns
 - whether the modus is collider (which is, yes, trivially true) 
 
Definition at line 122 of file collidermodus.h.
 
 
◆ impact_parameter()
  
  
      
        
          | double smash::ColliderModus::impact_parameter  | 
          ( | 
           | ) | 
           const | 
         
       
   | 
  
inline   | 
  
 
 
◆ create_deformed_nucleus()
  
  
      
        
          | std::unique_ptr< DeformedNucleus > smash::ColliderModus::create_deformed_nucleus  | 
          ( | 
          Configuration &  | 
          nucleus_cfg,  | 
         
        
           | 
           | 
          const int  | 
          ntest,  | 
         
        
           | 
           | 
          const std::string &  | 
          nucleus_type  | 
         
        
           | 
          ) | 
           |  | 
         
       
   | 
  
staticprivate   | 
  
 
Configure Deformed Nucleus. 
Sets up a deformed nucleus object based on the input parameters in the configuration file. 
- Parameters
 - 
  
    | [in] | nucleus_cfg | Subset of configuration, projectile or target section.  | 
    | [in] | ntest | Number of test particles  | 
    | [in] | nucleus_type | String 'projectile' or 'target'. To display an appropriate error message.  | 
  
   
- Returns
 - Pointer to the created deformed nucleus object. 
 
Definition at line 451 of file collidermodus.cc.
  453   bool auto_deform = nucleus_cfg.take({
"Deformed", 
"Automatic"});
 
  454   bool is_beta2 = nucleus_cfg.has_value({
"Deformed", 
"Beta_2"}) ? 
true : 
false;
 
  455   bool is_beta4 = nucleus_cfg.has_value({
"Deformed", 
"Beta_4"}) ? 
true : 
false;
 
  456   std::unique_ptr<DeformedNucleus> nucleus;
 
  458   if ((auto_deform && (!is_beta2 && !is_beta4)) ||
 
  459       (!auto_deform && (is_beta2 && is_beta4))) {
 
  460     nucleus = make_unique<DeformedNucleus>(nucleus_cfg, ntest, auto_deform);
 
  463     throw std::domain_error(
"Deformation of " + nucleus_type +
 
  464                             " nucleus not configured " 
  465                             "properly, please check whether all necessary " 
  466                             "parameters are set.");
 
 
 
 
◆ same_inputfile()
Checks if target and projectile are read from the same external file if they are both initialized as a customnucleus. 
Function is only called if, projectile is customnucleus. /param[in] proj_config Configuration of projectile nucleus /param[in] targ_config Configuration of target nucleus 
Definition at line 630 of file collidermodus.cc.
  635   if (!targ_config.has_value({
"Custom"})) {
 
  638   std::string projectile_file_directory =
 
  639       proj_config.read({
"Custom", 
"File_Directory"});
 
  640   std::string target_file_directory =
 
  641       targ_config.read({
"Custom", 
"File_Directory"});
 
  642   std::string projectile_file_name = proj_config.read({
"Custom", 
"File_Name"});
 
  643   std::string target_file_name = targ_config.read({
"Custom", 
"File_Name"});
 
  645   std::string proj_path =
 
  647   std::string targ_path =
 
  649   if (proj_path == targ_path) {
 
 
 
 
◆ rotate_reaction_plane()
  
  
      
        
          | void smash::ColliderModus::rotate_reaction_plane  | 
          ( | 
          double  | 
          phi,  | 
         
        
           | 
           | 
          Particles *  | 
          particles  | 
         
        
           | 
          ) | 
           |  | 
         
       
   | 
  
private   | 
  
 
Rotate the reaction plane about the angle phi. 
- Parameters
 - 
  
    | [in] | phi | Angle about which to rotate  | 
    | [in] | particles | Particles, whose position is rotated  | 
  
   
Definition at line 551 of file collidermodus.cc.
  552   for (ParticleData &
p : *particles) {
 
  553     ThreeVector pos = 
p.position().threevec();
 
  554     ThreeVector mom = 
p.momentum().threevec();
 
  555     pos.rotate_around_z(phi);
 
  556     mom.rotate_around_z(phi);
 
  557     p.set_3position(pos);
 
  558     p.set_3momentum(mom);
 
 
 
 
◆ sample_impact()
  
  
      
        
          | void smash::ColliderModus::sample_impact  | 
          ( | 
           | ) | 
           | 
         
       
   | 
  
private   | 
  
 
Sample impact parameter. 
Samples the impact parameter from values between imp_min_ and imp_max_, if linear or quadratic sampling is used. By specifying impact parameters and corresponding yields, custom sampling can be used. This depends on the value of sampling_.
Note that imp_max_ less than imp_min_ also works fine. 
Definition at line 562 of file collidermodus.cc.
  575       double probability_random = 1;
 
  576       double probability = 0;
 
  578       while (probability_random > probability) {
 
  580         probability = (*impact_interpolation_)(b) / 
yield_max_;
 
  581         assert(probability < 1.);
 
 
 
 
◆ get_velocities()
  
  
      
        
          | std::pair< double, double > smash::ColliderModus::get_velocities  | 
          ( | 
          double  | 
          mandelstam_s,  | 
         
        
           | 
           | 
          double  | 
          m_a,  | 
         
        
           | 
           | 
          double  | 
          m_b  | 
         
        
           | 
          ) | 
           |  | 
         
       
   | 
  
private   | 
  
 
Get the frame dependent velocity for each nucleus, using the current reference frame. 
- See also
 - frame_
 
- Parameters
 - 
  
    | [in] | mandelstam_s | The total center-of-mass energy of the system.  | 
    | [in] | m_a | The (positive) mass of the projectile.  | 
    | [in] | m_b | The (positive) mass of the target.  | 
  
   
- Returns
 - A pair < v_a, v_b > containing the velocities of the nuclei. 
 
- Exceptions
 - 
  
    | domain_error | if the reference frame is not properly specified  | 
  
   
Definition at line 593 of file collidermodus.cc.
  606       v_a = 
pCM / std::sqrt(m_a * m_a + 
pCM * 
pCM);
 
  607       v_b = -
pCM / std::sqrt(m_b * m_b + 
pCM * 
pCM);
 
  613       throw std::domain_error(
 
  614           "Invalid reference frame in " 
  615           "ColliderModus::get_velocities.");
 
  617   return std::make_pair(v_a, v_b);
 
 
 
 
◆ projectile_
  
  
      
        
          | std::unique_ptr<Nucleus> smash::ColliderModus::projectile_ | 
         
       
   | 
  
private   | 
  
 
Projectile. 
The object that goes from negative z-values to positive z-values with positive velocity. 
Definition at line 140 of file collidermodus.h.
 
 
◆ target_
  
  
      
        
          | std::unique_ptr<Nucleus> smash::ColliderModus::target_ | 
         
       
   | 
  
private   | 
  
 
Target. 
The object that goes from positive z-values to negative z-values with negative velocity. In fixed target experiments, the target is at rest. 
Definition at line 148 of file collidermodus.h.
 
 
◆ total_s_
  
  
      
        
          | double smash::ColliderModus::total_s_ | 
         
       
   | 
  
private   | 
  
 
Center-of-mass energy squared of the nucleus-nucleus collision. 
Needs to be double to allow for calculations at LHC energies 
Definition at line 154 of file collidermodus.h.
 
 
◆ sqrt_s_NN_
  
  
      
        
          | double smash::ColliderModus::sqrt_s_NN_ | 
         
       
   | 
  
private   | 
  
 
Center-of-mass energy of a nucleon-nucleon collision. 
Needs to be double to allow for calculations at LHC energies 
Definition at line 160 of file collidermodus.h.
 
 
◆ impact_
  
  
      
        
          | double smash::ColliderModus::impact_ = 0. | 
         
       
   | 
  
private   | 
  
 
Impact parameter. 
The nuclei projectile_ and target_ will be shifted along the x-axis so that their centers move on antiparallel lines that are this distance apart from each other. 
Definition at line 190 of file collidermodus.h.
 
 
◆ random_reaction_plane_
  
  
      
        
          | bool smash::ColliderModus::random_reaction_plane_ | 
         
       
   | 
  
private   | 
  
 
Whether the reaction plane should be randomized. 
Definition at line 192 of file collidermodus.h.
 
 
◆ sampling_
Method used for sampling of impact parameter. 
Definition at line 194 of file collidermodus.h.
 
 
◆ imp_min_
  
  
      
        
          | double smash::ColliderModus::imp_min_ = 0.0 | 
         
       
   | 
  
private   | 
  
 
 
◆ imp_max_
  
  
      
        
          | double smash::ColliderModus::imp_max_ = 0.0 | 
         
       
   | 
  
private   | 
  
 
 
◆ yield_max_
  
  
      
        
          | double smash::ColliderModus::yield_max_ = 0.0 | 
         
       
   | 
  
private   | 
  
 
Maximum value of yield. Needed for custom impact parameter sampling. 
Definition at line 200 of file collidermodus.h.
 
 
◆ impact_interpolation_
Initial value:
Pointer to the impact parameter interpolation. 
Definition at line 202 of file collidermodus.h.
 
 
◆ initial_z_displacement_
  
  
      
        
          | double smash::ColliderModus::initial_z_displacement_ = 2.0 | 
         
       
   | 
  
private   | 
  
 
Initial z-displacement of nuclei. 
Projectile is shifted on -(this value) in z-direction and target on +(this value)*v_target/v_projectile. In this way projectile and target touch at t=0 in z=0. 
Definition at line 230 of file collidermodus.h.
 
 
◆ frame_
Reference frame for the system, as specified from config. 
Definition at line 234 of file collidermodus.h.
 
 
◆ fermi_motion_
An option to include Fermi motion ("off", "on", "frozen") 
Definition at line 238 of file collidermodus.h.
 
 
◆ cll_in_nucleus_
  
  
      
        
          | bool smash::ColliderModus::cll_in_nucleus_ = false | 
         
       
   | 
  
private   | 
  
 
An option to accept first collisions within the same nucleus. 
Definition at line 242 of file collidermodus.h.
 
 
◆ velocity_projectile_
  
  
      
        
          | double smash::ColliderModus::velocity_projectile_ = 0.0 | 
         
       
   | 
  
private   | 
  
 
 
◆ velocity_target_
  
  
      
        
          | double smash::ColliderModus::velocity_target_ = 0.0 | 
         
       
   | 
  
private   | 
  
 
 
The documentation for this class was generated from the following files:
 
 
std::unique_ptr< InterpolateDataLinear< double > > impact_interpolation_
Pointer to the impact parameter interpolation.
 
double initial_z_displacement_
Initial z-displacement of nuclei.
 
double imp_min_
Minimum value of impact parameter.
 
std::string custom_file_path(const std::string &file_directory, const std::string &file_name)
Creates full path string consisting of file_directory and file_name Needed to initialize a customnucl...
 
Sample from areal / quadratic distribution.
 
bool cll_in_nucleus_
An option to accept first collisions within the same nucleus.
 
FermiMotion fermi_motion_
An option to include Fermi motion ("off", "on", "frozen")
 
bool random_reaction_plane_
Whether the reaction plane should be randomized.
 
Sampling sampling_
Method used for sampling of impact parameter.
 
double sqrt_s_NN_
Center-of-mass energy of a nucleon-nucleon collision.
 
double s_from_plab(double plab, double m_P, double m_T)
Convert p_lab to Mandelstam-s for a fixed-target setup, with a projectile of mass m_P and momentum pl...
 
std::unique_ptr< Nucleus > target_
Target.
 
double fixed_target_projectile_v(double s, double ma, double mb)
 
bool same_inputfile(Configuration &proj_config, Configuration &targ_config)
Checks if target and projectile are read from the same external file if they are both initialized as ...
 
void sample_impact()
Sample impact parameter.
 
double s_from_Ekin(double e_kin, double m_P, double m_T)
Convert E_kin to Mandelstam-s for a fixed-target setup, with a projectile of mass m_P and a kinetic e...
 
constexpr double nucleon_mass
Nucleon mass in GeV.
 
Use fermi motion without potentials.
 
std::array< einhard::Logger<>, std::tuple_size< LogArea::AreaTuple >::value > logg
An array that stores all pre-configured Logger objects.
 
Sample from custom, user-defined distribution.
 
double velocity_projectile_
Beam velocity of the projectile.
 
T pCM(const T sqrts, const T mass_a, const T mass_b) noexcept
 
double center_of_velocity_v(double s, double ma, double mb)
 
double imp_max_
Maximum value of impact parameter.
 
FormattingHelper< T > format(const T &value, const char *unit, int width=-1, int precision=-1)
Acts as a stream modifier for std::ostream to output an object with an optional suffix string and wit...
 
double total_s_
Center-of-mass energy squared of the nucleus-nucleus collision.
 
CalculationFrame frame_
Reference frame for the system, as specified from config.
 
static constexpr int LCollider
 
double yield_max_
Maximum value of yield. Needed for custom impact parameter sampling.
 
double impact_
Impact parameter.
 
Sample from uniform distribution.
 
Use fermi motion in combination with potentials.
 
std::unique_ptr< Nucleus > projectile_
Projectile.
 
double velocity_target_
Beam velocity of the target.
 
static std::unique_ptr< DeformedNucleus > create_deformed_nucleus(Configuration &nucleus_cfg, const int ntest, const std::string &nucleus_type)
Configure Deformed Nucleus.
 
void rotate_reaction_plane(double phi, Particles *particles)
Rotate the reaction plane about the angle phi.
 
std::pair< double, double > get_velocities(double mandelstam_s, double m_a, double m_b)
Get the frame dependent velocity for each nucleus, using the current reference frame.
 
T pCM_from_s(const T s, const T mass_a, const T mass_b) noexcept