Version: SMASH-1.5
smash::ColliderModus Class Reference

#include <collidermodus.h>

ColliderModus: Provides a modus for colliding nuclei.

To use this modus, choose

General:
Modus: Collider

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.

Inheritance diagram for smash::ColliderModus:
[legend]
Collaboration diagram for smash::ColliderModus:
[legend]

Classes

struct  ColliderEmpty
 Thrown when either projectile_ or target_ nuclei are empty. More...
 

Public Member Functions

 ColliderModus (Configuration modus_config, const ExperimentParameters &parameters)
 Constructor. More...
 
double initial_conditions (Particles *particles, const ExperimentParameters &parameters)
 Generates initial state of the particles in the system. More...
 
int total_N_number () const
 
int proj_N_number () const
 
double velocity_projectile () const
 
double velocity_target () const
 
bool cll_in_nucleus ()
 
FermiMotion fermi_motion ()
 
bool is_collider () const
 
double impact_parameter () const
 
- Public Member Functions inherited from smash::ModusDefault
int impose_boundary_conditions (Particles *, const OutputsList &={})
 Enforces sensible positions for the particles. More...
 
int total_N_number () const
 
int proj_N_number () const
 
bool cll_in_nucleus () const
 
bool is_collider () const
 
double impact_parameter () const
 
double velocity_projectile () const
 
double velocity_target () const
 
FermiMotion fermi_motion () const
 
double max_timestep (double) const
 
double length () const
 
Grid< GridOptions::Normalcreate_grid (const Particles &particles, double min_cell_length, CellSizeStrategy strategy=CellSizeStrategy::Optimal) const
 Creates the Grid with normal boundary conditions. More...
 

Private Member Functions

void sample_impact ()
 Sample impact parameter. More...
 
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. More...
 

Static Private Member Functions

static std::unique_ptr< DeformedNucleuscreate_deformed_nucleus (Configuration &nucleus_cfg, const int ntest, const std::string &nucleus_type)
 Configure Deformed Nucleus. More...
 

Private Attributes

std::unique_ptr< Nucleusprojectile_
 Projectile. More...
 
std::unique_ptr< Nucleustarget_
 Target. More...
 
double total_s_
 Center-of-mass energy squared of the nucleus-nucleus collision. More...
 
double sqrt_s_NN_
 Center-of-mass energy of a nucleon-nucleon collision. More...
 
double impact_ = 0.
 Impact parameter. More...
 
Sampling sampling_ = Sampling::Quadratic
 Method used for sampling of impact parameter. More...
 
double imp_min_ = 0.0
 Minimum value of impact parameter. More...
 
double imp_max_ = 0.0
 Maximum value of impact parameter. More...
 
double yield_max_ = 0.0
 Maximum value of yield. Needed for custom impact parameter sampling. More...
 
std::unique_ptr< InterpolateDataLinear< double > > impact_interpolation_
 Pointer to the impact parameter interpolation. More...
 
double initial_z_displacement_ = 2.0
 Initial z-displacement of nuclei. More...
 
CalculationFrame frame_ = CalculationFrame::CenterOfVelocity
 Reference frame for the system, as specified from config. More...
 
FermiMotion fermi_motion_ = FermiMotion::Off
 An option to include Fermi motion ("off", "on", "frozen") More...
 
bool cll_in_nucleus_ = false
 An option to accept first collisions within the same nucleus. More...
 
double velocity_projectile_ = 0.0
 Beam velocity of the projectile. More...
 
double velocity_target_ = 0.0
 Beam velocity of the target. More...
 

Friends

std::ostream & operator<< (std::ostream &, const ColliderModus &)
 Writes the initial state for the ColliderModus to the output stream. More...
 

Constructor & Destructor Documentation

◆ ColliderModus()

smash::ColliderModus::ColliderModus ( Configuration  modus_config,
const ExperimentParameters parameters 
)
explicit

Constructor.

Takes all there is to take from the (truncated!) configuration object (only contains configuration for this modus).

Parameters
[in]modus_configThe configuration object that sets all initial conditions of the experiment.
[in]parametersUnused, but necessary because of templated initialization
Exceptions
ColliderEmptyif projectile or nucleus are empty (i.e. do not contain particles)
InvalidEnergyif 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_errorif more or less than exactly one of the input energy options is specified, or if custom impact parameter Values and Yields are improperly supplied

Determine whether to allow the first collisions within the same nucleus

Todo:
include a check that only one method of specifying impact is used

Definition at line 256 of file collidermodus.cc.

257  {
258  Configuration modus_cfg = modus_config["Collider"];
259 
260  // Get the reference frame for the collision calculation.
261  if (modus_cfg.has_value({"Calculation_Frame"})) {
262  frame_ = modus_cfg.take({"Calculation_Frame"});
263  }
264 
266  if (modus_cfg.has_value({"Collisions_Within_Nucleus"})) {
267  cll_in_nucleus_ = modus_cfg.take({"Collisions_Within_Nucleus"});
268  }
269 
270  // Set up the projectile nucleus
271  Configuration proj_cfg = modus_cfg["Projectile"];
272  if (proj_cfg.has_value({"Deformed"})) {
273  projectile_ =
274  create_deformed_nucleus(proj_cfg, params.testparticles, "projectile");
275  } else {
276  projectile_ = make_unique<Nucleus>(proj_cfg, params.testparticles);
277  }
278  if (projectile_->size() < 1) {
279  throw ColliderEmpty("Input Error: Projectile nucleus is empty.");
280  }
281 
282  // Set up the target nucleus
283  Configuration targ_cfg = modus_cfg["Target"];
284  if (targ_cfg.has_value({"Deformed"})) {
285  target_ = create_deformed_nucleus(targ_cfg, params.testparticles, "target");
286  } else {
287  target_ = make_unique<Nucleus>(targ_cfg, params.testparticles);
288  }
289  if (target_->size() < 1) {
290  throw ColliderEmpty("Input Error: Target nucleus is empty.");
291  }
292 
293  // Get the Fermi-Motion input (off, on, frozen)
294  if (modus_cfg.has_value({"Fermi_Motion"})) {
295  // We only read the value, because it is still required by the experiment
296  // class to make sure we don't use frozen Fermi momenta with potentials.
297  fermi_motion_ = modus_cfg.read({"Fermi_Motion"});
298  }
299 
300  // Get the total nucleus-nucleus collision energy. Since there is
301  // no meaningful choice for a default energy, we require the user to
302  // give one (and only one) energy input from the available options.
303  int energy_input = 0;
304  const double mass_projec = projectile_->mass();
305  const double mass_target = target_->mass();
306  // average mass of a particle in that nucleus
307  const double mass_a =
308  projectile_->mass() / projectile_->number_of_particles();
309  const double mass_b = target_->mass() / target_->number_of_particles();
310  // Option 1: Center of mass energy.
311  if (modus_cfg.has_value({"Sqrtsnn"})) {
312  sqrt_s_NN_ = modus_cfg.take({"Sqrtsnn"});
313  // Check that input satisfies the lower bound (everything at rest).
314  if (sqrt_s_NN_ <= mass_a + mass_b) {
315  throw ModusDefault::InvalidEnergy(
316  "Input Error: sqrt(s_NN) is not larger than masses:\n" +
317  std::to_string(sqrt_s_NN_) + " GeV <= " + std::to_string(mass_a) +
318  " GeV + " + std::to_string(mass_b) + " GeV.");
319  }
320  // Set the total nucleus-nucleus collision energy.
321  total_s_ = (sqrt_s_NN_ * sqrt_s_NN_ - mass_a * mass_a - mass_b * mass_b) *
322  mass_projec * mass_target / (mass_a * mass_b) +
323  mass_projec * mass_projec + mass_target * mass_target;
324  energy_input++;
325  }
326  /* Option 2: Kinetic energy per nucleon of the projectile nucleus
327  * (target at rest). */
328  if (modus_cfg.has_value({"E_Kin"})) {
329  const double e_kin = modus_cfg.take({"E_Kin"});
330  if (e_kin < 0) {
331  throw ModusDefault::InvalidEnergy(
332  "Input Error: "
333  "E_Kin must be nonnegative.");
334  }
335  // Set the total nucleus-nucleus collision energy.
336  total_s_ = s_from_Ekin(e_kin * projectile_->number_of_particles(),
337  mass_projec, mass_target);
338  sqrt_s_NN_ = std::sqrt(s_from_Ekin(e_kin, mass_a, mass_b));
339  energy_input++;
340  }
341  // Option 3: Momentum of the projectile nucleus (target at rest).
342  if (modus_cfg.has_value({"P_Lab"})) {
343  const double p_lab = modus_cfg.take({"P_Lab"});
344  if (p_lab < 0) {
345  throw ModusDefault::InvalidEnergy(
346  "Input Error: "
347  "P_Lab must be nonnegative.");
348  }
349  // Set the total nucleus-nucleus collision energy.
350  total_s_ = s_from_plab(p_lab * projectile_->number_of_particles(),
351  mass_projec, mass_target);
352  sqrt_s_NN_ = std::sqrt(s_from_plab(p_lab, mass_a, mass_b));
353  energy_input++;
354  }
355  if (energy_input == 0) {
356  throw std::domain_error(
357  "Input Error: Non-existent collision energy. "
358  "Please provide one of Sqrtsnn/E_Kin/P_Lab.");
359  }
360  if (energy_input > 1) {
361  throw std::domain_error(
362  "Input Error: Redundant collision energy. "
363  "Please provide only one of Sqrtsnn/E_Kin/P_Lab.");
364  }
365 
366  /* Impact parameter setting: Either "Value", "Range", "Max" or "Sample".
367  * Unspecified means 0 impact parameter.*/
368  if (modus_cfg.has_value({"Impact", "Value"})) {
369  impact_ = modus_cfg.take({"Impact", "Value"});
370  imp_min_ = impact_;
371  imp_max_ = impact_;
372  } else {
373  // If impact is not supplied by value, inspect sampling parameters:
374  if (modus_cfg.has_value({"Impact", "Sample"})) {
375  sampling_ = modus_cfg.take({"Impact", "Sample"});
376  if (sampling_ == Sampling::Custom) {
377  if (!(modus_cfg.has_value({"Impact", "Values"}) ||
378  modus_cfg.has_value({"Impact", "Yields"}))) {
379  throw std::domain_error(
380  "Input Error: Need impact parameter spectrum for custom "
381  "sampling. "
382  "Please provide Values and Yields.");
383  }
384  const std::vector<double> impacts =
385  modus_cfg.take({"Impact", "Values"});
386  const std::vector<double> yields = modus_cfg.take({"Impact", "Yields"});
387  if (impacts.size() != yields.size()) {
388  throw std::domain_error(
389  "Input Error: Need as many impact parameter values as yields. "
390  "Please make sure that Values and Yields have the same length.");
391  }
392  impact_interpolation_ = make_unique<InterpolateDataLinear<double>>(
393  InterpolateDataLinear<double>(impacts, yields));
394 
395  const auto imp_minmax =
396  std::minmax_element(impacts.begin(), impacts.end());
397  imp_min_ = *imp_minmax.first;
398  imp_max_ = *imp_minmax.second;
399  yield_max_ = *std::max_element(yields.begin(), yields.end());
400  }
401  }
402  if (modus_cfg.has_value({"Impact", "Range"})) {
403  const std::array<double, 2> range = modus_cfg.take({"Impact", "Range"});
404  imp_min_ = range[0];
405  imp_max_ = range[1];
406  }
407  if (modus_cfg.has_value({"Impact", "Max"})) {
408  imp_min_ = 0.0;
409  imp_max_ = modus_cfg.take({"Impact", "Max"});
410  }
411  }
413 
414  // Look for user-defined initial separation between nuclei.
415  if (modus_cfg.has_value({"Initial_Distance"})) {
416  initial_z_displacement_ = modus_cfg.take({"Initial_Distance"});
417  // the displacement is half the distance (both nuclei are shifted
418  // initial_z_displacement_ away from origin)
420  }
421 }
double yield_max_
Maximum value of yield. Needed for custom impact parameter sampling.
static std::unique_ptr< DeformedNucleus > create_deformed_nucleus(Configuration &nucleus_cfg, const int ntest, const std::string &nucleus_type)
Configure Deformed Nucleus.
double impact_
Impact parameter.
std::unique_ptr< InterpolateDataLinear< double > > impact_interpolation_
Pointer to the impact parameter interpolation.
double imp_min_
Minimum value of impact parameter.
std::unique_ptr< Nucleus > projectile_
Projectile.
Sampling sampling_
Method used for sampling of impact parameter.
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...
Definition: kinematics.h:224
Sample from custom, user-defined distribution.
FermiMotion fermi_motion_
An option to include Fermi motion ("off", "on", "frozen")
double initial_z_displacement_
Initial z-displacement of nuclei.
double sqrt_s_NN_
Center-of-mass energy of a nucleon-nucleon collision.
bool cll_in_nucleus_
An option to accept first collisions within the same nucleus.
std::unique_ptr< Nucleus > target_
Target.
double imp_max_
Maximum value of impact parameter.
CalculationFrame frame_
Reference frame for the system, as specified from config.
double total_s_
Center-of-mass energy squared of the nucleus-nucleus collision.
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...
Definition: kinematics.h:211
Here is the call graph for this function:

Member Function Documentation

◆ initial_conditions()

double smash::ColliderModus::initial_conditions ( Particles particles,
const ExperimentParameters parameters 
)

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]particlesAn empty list that gets filled up by this function
[in]parametersThe initialization parameters of the system
Returns
The starting time of the simulation (negative, so that nuclei collide exactly at t=0)
Exceptions
domain_errorif the velocities of each nucleus are >= 1, or if input for Fermi motion is invalid

Definition at line 452 of file collidermodus.cc.

453  {
454  const auto &log = logger<LogArea::Collider>();
455  sample_impact();
456 
457  log.info() << "Impact parameter = " << format(impact_, "fm");
458  // Populate the nuclei with appropriately distributed nucleons.
459  // If deformed, this includes rotating the nucleus.
460  projectile_->arrange_nucleons();
461  target_->arrange_nucleons();
462 
463  // Use the total mandelstam variable to get the frame-dependent velocity for
464  // each nucleus. Position a is projectile, position b is target.
465  double v_a, v_b;
466  std::tie(v_a, v_b) =
467  get_velocities(total_s_, projectile_->mass(), target_->mass());
468 
469  // If velocities are larger or equal to 1, throw an exception.
470  if (v_a >= 1.0 || v_b >= 1.0) {
471  throw std::domain_error(
472  "Found velocity equal to or larger than 1 in "
473  "ColliderModus::initial_conditions.\nConsider using "
474  "the center of velocity reference frame.");
475  }
476 
477  // Calculate the beam velocity of the projectile and the target, which will be
478  // used to calculate the beam momenta in experiment.cc
480  velocity_projectile_ = v_a;
481  velocity_target_ = v_b;
482  }
483 
484  // Generate Fermi momenta if necessary
487  // Frozen: Fermi momenta will be ignored during the propagation to
488  // avoid that the nuclei will fly apart.
489  projectile_->generate_fermi_momenta();
490  target_->generate_fermi_momenta();
492  log.info() << "Fermi motion is ON.";
493  } else {
494  log.info() << "FROZEN Fermi motion is on.";
495  }
496  } else if (fermi_motion_ == FermiMotion::Off) {
497  // No Fermi-momenta are generated in this case
498  log.info() << "Fermi motion is OFF.";
499  } else {
500  throw std::domain_error("Invalid Fermi_Motion input.");
501  }
502 
503  // Boost the nuclei to the appropriate velocity.
504  projectile_->boost(v_a);
505  target_->boost(v_b);
506 
507  // Shift the nuclei into starting positions. Contracted spheres with
508  // nuclear radii should touch exactly at t=0. Modus starts at negative
509  // time corresponding to additional initial displacement.
510  const double d_a = std::max(0., projectile_->get_diffusiveness());
511  const double d_b = std::max(0., target_->get_diffusiveness());
512  const double r_a = projectile_->get_nuclear_radius();
513  const double r_b = target_->get_nuclear_radius();
514  const double dz = initial_z_displacement_;
515 
516  const double simulation_time = -dz / std::abs(v_a);
517  const double proj_z = -dz - std::sqrt(1.0 - v_a * v_a) * (r_a + d_a);
518  const double targ_z =
519  +dz * std::abs(v_b / v_a) + std::sqrt(1.0 - v_b * v_b) * (r_b + d_b);
520  projectile_->shift(proj_z, +impact_ / 2.0, simulation_time);
521  target_->shift(targ_z, -impact_ / 2.0, simulation_time);
522 
523  // Put the particles in the nuclei into code particles.
524  projectile_->copy_particles(particles);
525  target_->copy_particles(particles);
526  return simulation_time;
527 }
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...
Definition: logging.h:310
double velocity_target_
Beam velocity of the target.
double impact_
Impact parameter.
std::unique_ptr< Nucleus > projectile_
Projectile.
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.
FermiMotion fermi_motion_
An option to include Fermi motion ("off", "on", "frozen")
double initial_z_displacement_
Initial z-displacement of nuclei.
Don&#39;t use fermi motion.
Use fermi motion without potentials.
std::unique_ptr< Nucleus > target_
Target.
double velocity_projectile_
Beam velocity of the projectile.
void sample_impact()
Sample impact parameter.
Use fermi motion in combination with potentials.
double total_s_
Center-of-mass energy squared of the nucleus-nucleus collision.
Here is the call graph for this function:

◆ 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 82 of file collidermodus.h.

82 { return target_->size() + projectile_->size(); }
std::unique_ptr< Nucleus > projectile_
Projectile.
std::unique_ptr< Nucleus > target_
Target.

◆ proj_N_number()

int smash::ColliderModus::proj_N_number ( ) const
inline
Returns
The number of test particles in the projectile nucleus

Definition at line 84 of file collidermodus.h.

84 { return projectile_->size(); }
std::unique_ptr< Nucleus > projectile_
Projectile.

◆ 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 90 of file collidermodus.h.

90 { return velocity_projectile_; }
double velocity_projectile_
Beam velocity of the projectile.

◆ 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 95 of file collidermodus.h.

95 { return velocity_target_; }
double velocity_target_
Beam velocity of the target.

◆ 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 99 of file collidermodus.h.

99 { return cll_in_nucleus_; }
bool cll_in_nucleus_
An option to accept first collisions within the same nucleus.

◆ fermi_motion()

FermiMotion smash::ColliderModus::fermi_motion ( )
inline
Returns
The Fermi motion type

Definition at line 101 of file collidermodus.h.

101 { return fermi_motion_; }
FermiMotion fermi_motion_
An option to include Fermi motion ("off", "on", "frozen")

◆ is_collider()

bool smash::ColliderModus::is_collider ( ) const
inline
Returns
whether the modus is collider (which is, yes, trivially true)

Definition at line 103 of file collidermodus.h.

103 { return true; }

◆ impact_parameter()

double smash::ColliderModus::impact_parameter ( ) const
inline
Returns
impact parameter of the collision

Definition at line 105 of file collidermodus.h.

105 { return impact_; }
double impact_
Impact parameter.

◆ 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_cfgSubset of configuration, projectile or target section.
[in]ntestNumber of test particles
[in]nucleus_typeString 'projectile' or 'target'. To display an appropriate error message.
Returns
Pointer to the created deformed nucleus object.

Definition at line 433 of file collidermodus.cc.

434  {
435  bool auto_deform = nucleus_cfg.take({"Deformed", "Automatic"});
436  bool is_beta2 = nucleus_cfg.has_value({"Deformed", "Beta_2"}) ? true : false;
437  bool is_beta4 = nucleus_cfg.has_value({"Deformed", "Beta_4"}) ? true : false;
438  std::unique_ptr<DeformedNucleus> nucleus;
439 
440  if ((auto_deform && (!is_beta2 && !is_beta4)) ||
441  (!auto_deform && (is_beta2 && is_beta4))) {
442  nucleus = make_unique<DeformedNucleus>(nucleus_cfg, ntest);
443  return nucleus;
444  } else {
445  throw std::domain_error("Deformation of " + nucleus_type +
446  " nucleus not configured "
447  "properly, please check whether all necessary "
448  "parameters are set.");
449  }
450 }
Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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 529 of file collidermodus.cc.

529  {
530  switch (sampling_) {
531  case Sampling::Quadratic: {
532  // quadratic sampling: Note that for bmin > bmax, this still yields
533  // the correct distribution (however canonical() = 0 is then the
534  // upper end, not the lower).
535  impact_ = std::sqrt(imp_min_ * imp_min_ +
538  } break;
539  case Sampling::Custom: {
540  // rejection sampling based on given distribution
541  assert(impact_interpolation_ != nullptr);
542  double probability_random = 1;
543  double probability = 0;
544  double b;
545  while (probability_random > probability) {
547  probability = (*impact_interpolation_)(b) / yield_max_;
548  assert(probability < 1.);
549  probability_random = random::uniform(0., 1.);
550  }
551  impact_ = b;
552  } break;
553  case Sampling::Uniform: {
554  // linear sampling. Still, min > max works fine.
556  }
557  }
558 }
double yield_max_
Maximum value of yield. Needed for custom impact parameter sampling.
double impact_
Impact parameter.
std::unique_ptr< InterpolateDataLinear< double > > impact_interpolation_
Pointer to the impact parameter interpolation.
double imp_min_
Minimum value of impact parameter.
T canonical()
Definition: random.h:110
Sample from uniform distribution.
Sampling sampling_
Method used for sampling of impact parameter.
Sample from custom, user-defined distribution.
Sample from areal / quadratic distribution.
T uniform(T min, T max)
Definition: random.h:85
double imp_max_
Maximum value of impact parameter.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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_sThe total center-of-mass energy of the system.
[in]m_aThe (positive) mass of the projectile.
[in]m_bThe (positive) mass of the target.
Returns
A pair < v_a, v_b > containing the velocities of the nuclei.
Exceptions
domain_errorif the reference frame is not properly specified

Definition at line 560 of file collidermodus.cc.

561  {
562  double v_a = 0.0;
563  double v_b = 0.0;
564  // Frame dependent calculations of velocities. Assume v_a >= 0, v_b <= 0.
565  switch (frame_) {
567  v_a = center_of_velocity_v(s, m_a, m_b);
568  v_b = -v_a;
569  break;
571  // Compute center of mass momentum.
572  double pCM = pCM_from_s(s, m_a, m_b);
573  v_a = pCM / std::sqrt(m_a * m_a + pCM * pCM);
574  v_b = -pCM / std::sqrt(m_b * m_b + pCM * pCM);
575  } break;
577  v_a = fixed_target_projectile_v(s, m_a, m_b);
578  break;
579  default:
580  throw std::domain_error(
581  "Invalid reference frame in "
582  "ColliderModus::get_velocities.");
583  }
584  return std::make_pair(v_a, v_b);
585 }
double fixed_target_projectile_v(double s, double ma, double mb)
Definition: kinematics.h:39
double center_of_velocity_v(double s, double ma, double mb)
Definition: kinematics.h:26
T pCM_from_s(const T s, const T mass_a, const T mass_b) noexcept
Definition: kinematics.h:66
T pCM(const T sqrts, const T mass_a, const T mass_b) noexcept
Definition: kinematics.h:79
CalculationFrame frame_
Reference frame for the system, as specified from config.
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ 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 121 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 129 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 135 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 141 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 163 of file collidermodus.h.

◆ sampling_

Sampling smash::ColliderModus::sampling_ = Sampling::Quadratic
private

Method used for sampling of impact parameter.

Definition at line 165 of file collidermodus.h.

◆ imp_min_

double smash::ColliderModus::imp_min_ = 0.0
private

Minimum value of impact parameter.

Definition at line 167 of file collidermodus.h.

◆ imp_max_

double smash::ColliderModus::imp_max_ = 0.0
private

Maximum value of impact parameter.

Definition at line 169 of file collidermodus.h.

◆ yield_max_

double smash::ColliderModus::yield_max_ = 0.0
private

Maximum value of yield. Needed for custom impact parameter sampling.

Definition at line 171 of file collidermodus.h.

◆ impact_interpolation_

std::unique_ptr<InterpolateDataLinear<double> > smash::ColliderModus::impact_interpolation_
private
Initial value:
=
nullptr

Pointer to the impact parameter interpolation.

Definition at line 173 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 192 of file collidermodus.h.

◆ frame_

CalculationFrame smash::ColliderModus::frame_ = CalculationFrame::CenterOfVelocity
private

Reference frame for the system, as specified from config.

Definition at line 196 of file collidermodus.h.

◆ fermi_motion_

FermiMotion smash::ColliderModus::fermi_motion_ = FermiMotion::Off
private

An option to include Fermi motion ("off", "on", "frozen")

Definition at line 200 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 204 of file collidermodus.h.

◆ velocity_projectile_

double smash::ColliderModus::velocity_projectile_ = 0.0
private

Beam velocity of the projectile.

Definition at line 208 of file collidermodus.h.

◆ velocity_target_

double smash::ColliderModus::velocity_target_ = 0.0
private

Beam velocity of the target.

Definition at line 212 of file collidermodus.h.


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