 |
Version: SMASH-2.0
|
|
#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 43 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 266 of file collidermodus.cc.
268 Configuration modus_cfg = modus_config[
"Collider"];
270 if (modus_cfg.has_value({
"Calculation_Frame"})) {
271 frame_ = modus_cfg.take({
"Calculation_Frame"});
275 if (modus_cfg.has_value({
"Collisions_Within_Nucleus"})) {
278 Configuration proj_cfg = modus_cfg[
"Projectile"];
279 Configuration targ_cfg = modus_cfg[
"Target"];
282 bool same_file =
false;
284 if (proj_cfg.has_value({
"Deformed"})) {
287 }
else if (proj_cfg.has_value({
"Custom"})) {
290 make_unique<CustomNucleus>(proj_cfg, params.testparticles, same_file);
292 projectile_ = make_unique<Nucleus>(proj_cfg, params.testparticles);
295 throw ColliderEmpty(
"Input Error: Projectile nucleus is empty.");
299 if (targ_cfg.has_value({
"Deformed"})) {
301 }
else if (targ_cfg.has_value({
"Custom"})) {
303 make_unique<CustomNucleus>(targ_cfg, params.testparticles, same_file);
305 target_ = make_unique<Nucleus>(targ_cfg, params.testparticles);
308 throw ColliderEmpty(
"Input Error: Target nucleus is empty.");
312 if (modus_cfg.has_value({
"Fermi_Motion"})) {
321 int energy_input = 0;
323 const double mass_target =
target_->mass();
325 const double mass_a =
327 const double mass_b =
target_->mass() /
target_->number_of_particles();
329 if (modus_cfg.has_value({
"Sqrtsnn"})) {
333 throw ModusDefault::InvalidEnergy(
334 "Input Error: sqrt(s_NN) is not larger than masses:\n" +
335 std::to_string(
sqrt_s_NN_) +
" GeV <= " + std::to_string(mass_a) +
336 " GeV + " + std::to_string(mass_b) +
" GeV.");
340 mass_projec * mass_target / (mass_a * mass_b) +
341 mass_projec * mass_projec + mass_target * mass_target;
346 if (modus_cfg.has_value({
"E_Kin"})) {
347 const double e_kin = modus_cfg.take({
"E_Kin"});
349 throw ModusDefault::InvalidEnergy(
351 "E_Kin must be nonnegative.");
355 mass_projec, mass_target);
360 if (modus_cfg.has_value({
"P_Lab"})) {
361 const double p_lab = modus_cfg.take({
"P_Lab"});
363 throw ModusDefault::InvalidEnergy(
365 "P_Lab must be nonnegative.");
369 mass_projec, mass_target);
373 if (energy_input == 0) {
374 throw std::domain_error(
375 "Input Error: Non-existent collision energy. "
376 "Please provide one of Sqrtsnn/E_Kin/P_Lab.");
378 if (energy_input > 1) {
379 throw std::domain_error(
380 "Input Error: Redundant collision energy. "
381 "Please provide only one of Sqrtsnn/E_Kin/P_Lab.");
386 if (modus_cfg.has_value({
"Impact",
"Value"})) {
387 impact_ = modus_cfg.take({
"Impact",
"Value"});
392 if (modus_cfg.has_value({
"Impact",
"Sample"})) {
393 sampling_ = modus_cfg.take({
"Impact",
"Sample"});
395 if (!(modus_cfg.has_value({
"Impact",
"Values"}) ||
396 modus_cfg.has_value({
"Impact",
"Yields"}))) {
397 throw std::domain_error(
398 "Input Error: Need impact parameter spectrum for custom "
400 "Please provide Values and Yields.");
402 const std::vector<double> impacts =
403 modus_cfg.take({
"Impact",
"Values"});
404 const std::vector<double> yields = modus_cfg.take({
"Impact",
"Yields"});
405 if (impacts.size() != yields.size()) {
406 throw std::domain_error(
407 "Input Error: Need as many impact parameter values as yields. "
408 "Please make sure that Values and Yields have the same length.");
411 InterpolateDataLinear<double>(impacts, yields));
413 const auto imp_minmax =
414 std::minmax_element(impacts.begin(), impacts.end());
417 yield_max_ = *std::max_element(yields.begin(), yields.end());
420 if (modus_cfg.has_value({
"Impact",
"Range"})) {
421 const std::array<double, 2> range = modus_cfg.take({
"Impact",
"Range"});
425 if (modus_cfg.has_value({
"Impact",
"Max"})) {
427 imp_max_ = modus_cfg.take({
"Impact",
"Max"});
433 modus_cfg.take({
"Impact",
"Random_Reaction_Plane"},
false);
435 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 622 of file collidermodus.cc.
625 if (file_directory.back() ==
'/') {
626 return file_directory + file_name;
628 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 472 of file collidermodus.cc.
489 if (v_a >= 1.0 || v_b >= 1.0) {
490 throw std::domain_error(
491 "Found velocity equal to or larger than 1 in "
492 "ColliderModus::initial_conditions.\nConsider using "
493 "the center of velocity reference frame.");
509 target_->generate_fermi_momenta();
519 throw std::domain_error(
"Invalid Fermi_Motion input.");
529 const double d_a = std::max(0.,
projectile_->get_diffusiveness());
530 const double d_b = std::max(0.,
target_->get_diffusiveness());
531 const double r_a =
projectile_->get_nuclear_radius();
532 const double r_b =
target_->get_nuclear_radius();
535 const double simulation_time = -dz / std::abs(v_a);
536 const double proj_z = -dz - std::sqrt(1.0 - v_a * v_a) * (r_a + d_a);
537 const double targ_z =
538 +dz * std::abs(v_b / v_a) + std::sqrt(1.0 - v_b * v_b) * (r_b + d_b);
548 target_->copy_particles(particles);
550 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 90 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 92 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 95 of file collidermodus.h.
96 const double passing_distance =
98 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 110 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 115 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 119 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 123 of file collidermodus.h.
◆ sqrt_s_NN()
double smash::ColliderModus::sqrt_s_NN |
( |
| ) |
const |
|
inline |
- Returns
- center of mass energy per nucleon pair
Definition at line 125 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 453 of file collidermodus.cc.
455 bool auto_deform = nucleus_cfg.take({
"Deformed",
"Automatic"});
456 bool is_beta2 = nucleus_cfg.has_value({
"Deformed",
"Beta_2"}) ?
true :
false;
457 bool is_beta4 = nucleus_cfg.has_value({
"Deformed",
"Beta_4"}) ?
true :
false;
458 std::unique_ptr<DeformedNucleus> nucleus;
460 if ((auto_deform && (!is_beta2 && !is_beta4)) ||
461 (!auto_deform && (is_beta2 && is_beta4))) {
462 nucleus = make_unique<DeformedNucleus>(nucleus_cfg, ntest, auto_deform);
465 throw std::domain_error(
"Deformation of " + nucleus_type +
466 " nucleus not configured "
467 "properly, please check whether all necessary "
468 "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 632 of file collidermodus.cc.
637 if (!targ_config.has_value({
"Custom"})) {
640 std::string projectile_file_directory =
641 proj_config.read({
"Custom",
"File_Directory"});
642 std::string target_file_directory =
643 targ_config.read({
"Custom",
"File_Directory"});
644 std::string projectile_file_name = proj_config.read({
"Custom",
"File_Name"});
645 std::string target_file_name = targ_config.read({
"Custom",
"File_Name"});
647 std::string proj_path =
649 std::string targ_path =
651 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 553 of file collidermodus.cc.
554 for (ParticleData &
p : *particles) {
555 ThreeVector pos =
p.position().threevec();
556 ThreeVector mom =
p.momentum().threevec();
557 pos.rotate_around_z(phi);
558 mom.rotate_around_z(phi);
559 p.set_3position(pos);
560 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 564 of file collidermodus.cc.
577 double probability_random = 1;
578 double probability = 0;
580 while (probability_random > probability) {
582 probability = (*impact_interpolation_)(b) /
yield_max_;
583 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 595 of file collidermodus.cc.
608 v_a =
pCM / std::sqrt(m_a * m_a +
pCM *
pCM);
609 v_b = -
pCM / std::sqrt(m_b * m_b +
pCM *
pCM);
615 throw std::domain_error(
616 "Invalid reference frame in "
617 "ColliderModus::get_velocities.");
619 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 143 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 151 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 157 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 163 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 193 of file collidermodus.h.
◆ random_reaction_plane_
bool smash::ColliderModus::random_reaction_plane_ |
|
private |
Whether the reaction plane should be randomized.
Definition at line 195 of file collidermodus.h.
◆ sampling_
Method used for sampling of impact parameter.
Definition at line 197 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 203 of file collidermodus.h.
◆ impact_interpolation_
Initial value:
Pointer to the impact parameter interpolation.
Definition at line 205 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 233 of file collidermodus.h.
◆ frame_
Reference frame for the system, as specified from config.
Definition at line 237 of file collidermodus.h.
◆ fermi_motion_
An option to include Fermi motion ("off", "on", "frozen")
Definition at line 241 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 245 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