Version: SMASH-1.8
experiment.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2013-2020
3  * SMASH Team
4  *
5  * GNU General Public License (GPLv3 or later)
6  */
7 #ifndef SRC_INCLUDE_EXPERIMENT_H_
8 #define SRC_INCLUDE_EXPERIMENT_H_
9 
10 #include <limits>
11 #include <memory>
12 #include <string>
13 #include <vector>
14 
15 #include "actionfinderfactory.h"
16 #include "actions.h"
17 #include "bremsstrahlungaction.h"
18 #include "chrono.h"
19 #include "decayactionsfinder.h"
21 #include "energymomentumtensor.h"
22 #include "fourvector.h"
23 #include "grandcan_thermalizer.h"
24 #include "grid.h"
26 #include "outputparameters.h"
27 #include "pauliblocking.h"
28 #include "potential_globals.h"
29 #include "potentials.h"
30 #include "processstring.h"
31 #include "propagation.h"
32 #include "quantumnumbers.h"
33 #include "scatteractionphoton.h"
34 #include "scatteractionsfinder.h"
35 #include "thermalizationaction.h"
36 // Output
37 #include "binaryoutput.h"
38 #include "icoutput.h"
39 #include "oscaroutput.h"
40 #include "thermodynamicoutput.h"
41 #ifdef SMASH_USE_ROOT
42 #include "rootoutput.h"
43 #endif
44 #include "vtkoutput.h"
45 #include "wallcrossingaction.h"
46 
47 namespace std {
58 template <typename T, typename Ratio>
59 static ostream &operator<<(ostream &out,
60  const chrono::duration<T, Ratio> &seconds) {
61  using Seconds = chrono::duration<double>;
62  using Minutes = chrono::duration<double, std::ratio<60>>;
63  using Hours = chrono::duration<double, std::ratio<60 * 60>>;
64  constexpr Minutes threshold_for_minutes{10};
65  constexpr Hours threshold_for_hours{3};
66  if (seconds < threshold_for_minutes) {
67  return out << Seconds(seconds).count() << " [s]";
68  }
69  if (seconds < threshold_for_hours) {
70  return out << Minutes(seconds).count() << " [min]";
71  }
72  return out << Hours(seconds).count() << " [h]";
73 }
74 } // namespace std
75 
76 namespace smash {
77 static constexpr int LMain = LogArea::Main::id;
78 static constexpr int LInitialConditions = LogArea::InitialConditions::id;
79 
88  public:
89  ExperimentBase() = default;
94  virtual ~ExperimentBase() = default;
95 
116  static std::unique_ptr<ExperimentBase> create(Configuration config,
117  const bf::path &output_path);
118 
125  virtual void run() = 0;
126 
132  struct InvalidModusRequest : public std::invalid_argument {
133  using std::invalid_argument::invalid_argument;
134  };
135 };
136 
137 template <typename Modus>
139 template <typename Modus>
140 std::ostream &operator<<(std::ostream &out, const Experiment<Modus> &e);
141 
166 template <typename Modus>
167 class Experiment : public ExperimentBase {
168  friend class ExperimentBase;
169 
170  public:
177  void run() override;
178 
192  explicit Experiment(Configuration config, const bf::path &output_path);
193 
199  void initialize_new_event();
200 
207  void run_time_evolution();
208 
210  void do_final_decays();
211 
217  void final_output(const int evt_num);
218 
224 
229  Modus *modus() { return &modus_; }
230 
231  private:
243  template <typename Container>
244  bool perform_action(Action &action,
245  const Container &particles_before_actions);
254  void create_output(const std::string &format, const std::string &content,
255  const bf::path &output_path, const OutputParameters &par);
256 
263  void propagate_and_shine(double to_time);
264 
278 
280  void intermediate_output();
281 
283  void update_potentials();
284 
293  double compute_min_cell_length(double dt) const {
294  return std::sqrt(4 * dt * dt + max_transverse_distance_sqr_);
295  }
296 
298  double next_output_time() const {
299  return parameters_.outputclock->next_time();
300  }
301 
308 
311 
316  Modus modus_;
317 
320 
325  std::unique_ptr<Potentials> potentials_;
326 
331  std::unique_ptr<PauliBlocker> pauli_blocker_;
332 
337  OutputsList outputs_;
338 
340  OutputPtr dilepton_output_;
341 
343  OutputPtr photon_output_;
344 
350  std::vector<bool> nucleon_has_interacted_ = {};
355 
361  std::vector<FourVector> beam_momentum_ = {};
362 
364  std::vector<std::unique_ptr<ActionFinderInterface>> action_finders_;
365 
367  std::unique_ptr<DecayActionsFinderDilepton> dilepton_finder_;
368 
370  std::unique_ptr<ActionFinderInterface> photon_finder_;
371 
374 
376  std::unique_ptr<DensityLattice> jmu_B_lat_;
377 
379  std::unique_ptr<DensityLattice> jmu_I3_lat_;
380 
388  std::unique_ptr<DensityLattice> jmu_custom_lat_;
389 
392 
397  std::unique_ptr<RectangularLattice<FourVector>> UB_lat_ = nullptr;
398 
403  std::unique_ptr<RectangularLattice<FourVector>> UI3_lat_ = nullptr;
404 
406  std::unique_ptr<RectangularLattice<std::pair<ThreeVector, ThreeVector>>>
408 
410  std::unique_ptr<RectangularLattice<std::pair<ThreeVector, ThreeVector>>>
412 
414  std::unique_ptr<RectangularLattice<EnergyMomentumTensor>> Tmn_;
415 
417  bool printout_tmn_ = false;
418 
420  bool printout_tmn_landau_ = false;
421 
423  bool printout_v_landau_ = false;
424 
426  bool printout_lattice_td_ = false;
427 
429  std::unique_ptr<GrandCanThermalizer> thermalizer_;
430 
436 
451  const int nevents_;
452 
454  const double end_time_;
455 
461  const double delta_time_startup_;
462 
467  const bool force_decays_;
468 
470  const bool use_grid_;
471 
474 
476  const bool dileptons_switch_;
477 
479  const bool photons_switch_;
480 
483 
485  const bool IC_output_switch_;
486 
489 
491  double max_transverse_distance_sqr_ = std::numeric_limits<double>::max();
492 
501 
507 
509  SystemTimePoint time_start_ = SystemClock::now();
510 
513 
518  uint64_t interactions_total_ = 0;
519 
525 
530  uint64_t wall_actions_total_ = 0;
531 
537 
542  uint64_t total_pauli_blocked_ = 0;
543 
549 
554  double total_energy_removed_ = 0.0;
555 
557  int64_t seed_ = -1;
558 
564  friend std::ostream &operator<<<>(std::ostream &out, const Experiment &e);
565 };
566 
568 template <typename Modus>
569 std::ostream &operator<<(std::ostream &out, const Experiment<Modus> &e) {
570  out << "End time: " << e.end_time_ << " fm/c\n";
571  out << e.modus_;
572  return out;
573 }
574 
575 template <typename Modus>
576 void Experiment<Modus>::create_output(const std::string &format,
577  const std::string &content,
578  const bf::path &output_path,
579  const OutputParameters &out_par) {
580  logg[LExperiment].info() << "Adding output " << content << " of format "
581  << format << std::endl;
582 
583  if (format == "VTK" && content == "Particles") {
584  outputs_.emplace_back(
585  make_unique<VtkOutput>(output_path, content, out_par));
586  } else if (format == "Root") {
587 #ifdef SMASH_USE_ROOT
588  if (content == "Initial_Conditions") {
589  outputs_.emplace_back(
590  make_unique<RootOutput>(output_path, "SMASH_IC", out_par));
591  } else {
592  outputs_.emplace_back(
593  make_unique<RootOutput>(output_path, content, out_par));
594  }
595 #else
596  logg[LExperiment].error(
597  "Root output requested, but Root support not compiled in");
598 #endif
599  } else if (format == "Binary") {
600  if (content == "Collisions" || content == "Dileptons" ||
601  content == "Photons") {
602  outputs_.emplace_back(
603  make_unique<BinaryOutputCollisions>(output_path, content, out_par));
604  } else if (content == "Particles") {
605  outputs_.emplace_back(
606  make_unique<BinaryOutputParticles>(output_path, content, out_par));
607  } else if (content == "Initial_Conditions") {
608  outputs_.emplace_back(make_unique<BinaryOutputInitialConditions>(
609  output_path, content, out_par));
610  }
611  } else if (format == "Oscar1999" || format == "Oscar2013") {
612  outputs_.emplace_back(
613  create_oscar_output(format, content, output_path, out_par));
614  } else if (content == "Thermodynamics" && format == "ASCII") {
615  outputs_.emplace_back(
616  make_unique<ThermodynamicOutput>(output_path, content, out_par));
617  } else if (content == "Thermodynamics" && format == "VTK") {
618  printout_lattice_td_ = true;
619  outputs_.emplace_back(
620  make_unique<VtkOutput>(output_path, content, out_par));
621  } else if (content == "Initial_Conditions" && format == "ASCII") {
622  outputs_.emplace_back(
623  make_unique<ICOutput>(output_path, "SMASH_IC", out_par));
624  } else {
625  logg[LExperiment].error()
626  << "Unknown combination of format (" << format << ") and content ("
627  << content << "). Fix the config.";
628  }
629 }
630 
639 
780 template <typename Modus>
781 Experiment<Modus>::Experiment(Configuration config, const bf::path &output_path)
782  : parameters_(create_experiment_parameters(config)),
783  density_param_(DensityParameters(parameters_)),
784  modus_(config["Modi"], parameters_),
785  particles_(),
786  nevents_(config.take({"General", "Nevents"})),
787  end_time_(config.take({"General", "End_Time"})),
788  delta_time_startup_(parameters_.labclock->timestep_duration()),
789  force_decays_(
790  config.take({"Collision_Term", "Force_Decays_At_End"}, true)),
791  use_grid_(config.take({"General", "Use_Grid"}, true)),
792  metric_(
793  config.take({"General", "Metric_Type"}, ExpansionMode::NoExpansion),
794  config.take({"General", "Expansion_Rate"}, 0.1)),
795  dileptons_switch_(
796  config.take({"Collision_Term", "Dileptons", "Decays"}, false)),
797  photons_switch_(config.take(
798  {"Collision_Term", "Photons", "2to2_Scatterings"}, false)),
799  bremsstrahlung_switch_(
800  config.take({"Collision_Term", "Photons", "Bremsstrahlung"}, false)),
801  IC_output_switch_(config.has_value({"Output", "Initial_Conditions"})),
802  time_step_mode_(
803  config.take({"General", "Time_Step_Mode"}, TimeStepMode::Fixed)) {
804  logg[LExperiment].info() << *this;
805 
806  // create finders
807  if (dileptons_switch_) {
808  dilepton_finder_ = make_unique<DecayActionsFinderDilepton>();
809  }
810  if (photons_switch_ || bremsstrahlung_switch_) {
811  n_fractional_photons_ =
812  config.take({"Collision_Term", "Photons", "Fractional_Photons"}, 100);
813  }
814  if (parameters_.two_to_one) {
815  if (parameters_.res_lifetime_factor < 0.) {
816  throw std::invalid_argument(
817  "Resonance lifetime modifier cannot be negative!");
818  }
819  if (parameters_.res_lifetime_factor < really_small) {
820  logg[LExperiment].warn(
821  "Resonance lifetime set to zero. Make sure resonances cannot "
822  "interact",
823  "inelastically (e.g. resonance chains), else SMASH is known to "
824  "hang.");
825  }
826  action_finders_.emplace_back(
827  make_unique<DecayActionsFinder>(parameters_.res_lifetime_factor));
828  }
829  bool no_coll = config.take({"Collision_Term", "No_Collisions"}, false);
830  if ((parameters_.two_to_one || parameters_.included_2to2.any() ||
831  parameters_.strings_switch) &&
832  !no_coll) {
833  auto scat_finder = make_unique<ScatterActionsFinder>(
834  config, parameters_, nucleon_has_interacted_, modus_.total_N_number(),
835  modus_.proj_N_number());
836  max_transverse_distance_sqr_ =
837  scat_finder->max_transverse_distance_sqr(parameters_.testparticles);
838  process_string_ptr_ = scat_finder->get_process_string_ptr();
839  action_finders_.emplace_back(std::move(scat_finder));
840  } else {
841  max_transverse_distance_sqr_ = maximum_cross_section / M_PI * fm2_mb;
842  process_string_ptr_ = NULL;
843  }
844  const double modus_l = modus_.length();
845  if (modus_l > 0.) {
846  action_finders_.emplace_back(make_unique<WallCrossActionsFinder>(modus_l));
847  }
848  if (IC_output_switch_) {
849  if (!modus_.is_collider()) {
850  throw std::runtime_error(
851  "Initial conditions can only be extracted in collider modus.");
852  }
853  double proper_time;
854  if (config.has_value({"Output", "Initial_Conditions", "Proper_Time"})) {
855  // Read in proper time from config
856  proper_time =
857  config.take({"Output", "Initial_Conditions", "Proper_Time"});
858  } else {
859  // Default proper time is the passing time of the two nuclei
860  double default_proper_time = modus_.nuclei_passing_time();
861  double lower_bound =
862  config.take({"Output", "Initial_Conditions", "Lower_Bound"}, 0.5);
863  if (default_proper_time >= lower_bound) {
864  proper_time = default_proper_time;
865  } else {
866  logg[LInitialConditions].warn()
867  << "Nuclei passing time is too short, hypersurface proper time set "
868  "to tau = "
869  << lower_bound << " fm.";
870  proper_time = lower_bound;
871  }
872  }
873  action_finders_.emplace_back(
874  make_unique<HyperSurfaceCrossActionsFinder>(proper_time));
875  }
876 
877  if (config.has_value({"Collision_Term", "Pauli_Blocking"})) {
878  logg[LExperiment].info() << "Pauli blocking is ON.";
879  pauli_blocker_ = make_unique<PauliBlocker>(
880  config["Collision_Term"]["Pauli_Blocking"], parameters_);
881  }
883  config.take({"Collision_Term", "Power_Particle_Formation"}, 1.);
884 
917  // create outputs
918  logg[LExperiment].trace(source_location, " create OutputInterface objects");
919 
920  auto output_conf = config["Output"];
1141  dens_type_ = config.take({"Output", "Density_Type"}, DensityType::None);
1142  logg[LExperiment].debug()
1143  << "Density type printed to headers: " << dens_type_;
1144 
1145  const OutputParameters output_parameters(std::move(output_conf));
1146 
1147  std::vector<std::string> output_contents = output_conf.list_upmost_nodes();
1148  for (const auto &content : output_contents) {
1149  auto this_output_conf = output_conf[content.c_str()];
1150  const std::vector<std::string> formats = this_output_conf.take({"Format"});
1151  if (output_path == "") {
1152  continue;
1153  }
1154  for (const auto &format : formats) {
1155  create_output(format, content, output_path, output_parameters);
1156  }
1157  }
1158 
1159  /* We can take away the Fermi motion flag, because the collider modus is
1160  * already initialized. We only need it when potentials are enabled, but we
1161  * always have to take it, otherwise SMASH will complain about unused
1162  * options. We have to provide a default value for modi other than Collider.
1163  */
1164  const FermiMotion motion =
1165  config.take({"Modi", "Collider", "Fermi_Motion"}, FermiMotion::Off);
1166  if (config.has_value({"Potentials"})) {
1167  if (time_step_mode_ == TimeStepMode::None) {
1168  logg[LExperiment].error() << "Potentials only work with time steps!";
1169  throw std::invalid_argument("Can't use potentials without time steps!");
1170  }
1171  if (motion == FermiMotion::Frozen) {
1172  logg[LExperiment].error()
1173  << "Potentials don't work with frozen Fermi momenta! "
1174  "Use normal Fermi motion instead.";
1175  throw std::invalid_argument(
1176  "Can't use potentials "
1177  "with frozen Fermi momenta!");
1178  }
1179  logg[LExperiment].info() << "Potentials are ON. Timestep is "
1180  << parameters_.labclock->timestep_duration();
1181  // potentials need testparticles and gaussian sigma from parameters_
1182  potentials_ = make_unique<Potentials>(config["Potentials"], parameters_);
1183  }
1184 
1239  // Create lattices
1240  if (config.has_value({"Lattice"})) {
1241  // Take lattice properties from config to assign them to all lattices
1242  const std::array<double, 3> l = config.take({"Lattice", "Sizes"});
1243  const std::array<int, 3> n = config.take({"Lattice", "Cell_Number"});
1244  const std::array<double, 3> origin = config.take({"Lattice", "Origin"});
1245  const bool periodic = config.take({"Lattice", "Periodic"});
1246 
1247  if (printout_lattice_td_) {
1248  dens_type_lattice_printout_ = output_parameters.td_dens_type;
1249  printout_tmn_ = output_parameters.td_tmn;
1250  printout_tmn_landau_ = output_parameters.td_tmn_landau;
1251  printout_v_landau_ = output_parameters.td_v_landau;
1252  }
1253  if (printout_tmn_ || printout_tmn_landau_ || printout_v_landau_) {
1254  Tmn_ = make_unique<RectangularLattice<EnergyMomentumTensor>>(
1255  l, n, origin, periodic, LatticeUpdate::AtOutput);
1256  }
1257  /* Create baryon and isospin density lattices regardless of config
1258  if potentials are on. This is because they allow to compute
1259  potentials faster */
1260  if (potentials_) {
1261  if (potentials_->use_skyrme()) {
1262  jmu_B_lat_ = make_unique<DensityLattice>(l, n, origin, periodic,
1264  UB_lat_ = make_unique<RectangularLattice<FourVector>>(
1265  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1266  FB_lat_ = make_unique<
1267  RectangularLattice<std::pair<ThreeVector, ThreeVector>>>(
1268  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1269  }
1270  if (potentials_->use_symmetry()) {
1271  jmu_I3_lat_ = make_unique<DensityLattice>(l, n, origin, periodic,
1273  UI3_lat_ = make_unique<RectangularLattice<FourVector>>(
1274  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1275  FI3_lat_ = make_unique<
1276  RectangularLattice<std::pair<ThreeVector, ThreeVector>>>(
1277  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1278  }
1279  } else {
1280  if (dens_type_lattice_printout_ == DensityType::Baryon) {
1281  jmu_B_lat_ = make_unique<DensityLattice>(l, n, origin, periodic,
1283  }
1284  if (dens_type_lattice_printout_ == DensityType::BaryonicIsospin) {
1285  jmu_I3_lat_ = make_unique<DensityLattice>(l, n, origin, periodic,
1287  }
1288  }
1289  if (dens_type_lattice_printout_ != DensityType::None &&
1290  dens_type_lattice_printout_ != DensityType::BaryonicIsospin &&
1291  dens_type_lattice_printout_ != DensityType::Baryon) {
1292  jmu_custom_lat_ = make_unique<DensityLattice>(l, n, origin, periodic,
1294  }
1295  } else if (printout_lattice_td_) {
1296  logg[LExperiment].error(
1297  "If you want Thermodynamic VTK output, configure a lattice for it.");
1298  }
1299  // Warning for the mean field calculation if lattice is not on.
1300  if ((potentials_ != nullptr) && (jmu_B_lat_ == nullptr)) {
1301  logg[LExperiment].warn() << "Lattice is NOT used. Mean fields are "
1302  << "not going to be calculated.";
1303  }
1304 
1305  // Store pointers to potential and lattice accessible for Action
1306  if (parameters_.potential_affect_threshold) {
1307  UB_lat_pointer = UB_lat_.get();
1308  UI3_lat_pointer = UI3_lat_.get();
1309  pot_pointer = potentials_.get();
1310  }
1311 
1312  // Create forced thermalizer
1313  if (config.has_value({"Forced_Thermalization"})) {
1314  Configuration &&th_conf = config["Forced_Thermalization"];
1315  thermalizer_ = modus_.create_grandcan_thermalizer(th_conf);
1316  }
1317 
1318  /* Take the seed setting only after the configuration was stored to a file
1319  * in smash.cc */
1320  seed_ = config.take({"General", "Randomseed"});
1321 }
1322 
1324 const std::string hline(113, '-');
1325 
1351 std::string format_measurements(const Particles &particles,
1352  uint64_t scatterings_this_interval,
1353  const QuantumNumbers &conserved_initial,
1354  SystemTimePoint time_start, double time,
1355  double E_mean_field,
1356  double E_mean_field_initial);
1371  const Potentials &potentials,
1372  RectangularLattice<smash::DensityOnLattice> &jmu_B_lat,
1373  const ExperimentParameters &parameters);
1374 
1375 template <typename Modus>
1377  random::set_seed(seed_);
1378  logg[LExperiment].info() << "random number seed: " << seed_;
1379  /* Set seed for the next event. It has to be positive, so it can be entered
1380  * in the config.
1381  *
1382  * We have to be careful about the minimal integer, whose absolute value
1383  * cannot be represented. */
1384  int64_t r = random::advance();
1385  while (r == INT64_MIN) {
1386  r = random::advance();
1387  }
1388  seed_ = std::abs(r);
1389  /* Set the random seed used in PYTHIA hadronization
1390  * to be same with the SMASH one.
1391  * In this way we ensure that the results are reproducible
1392  * for every event if one knows SMASH random seed. */
1393  if (process_string_ptr_ != NULL) {
1394  process_string_ptr_->init_pythia_hadron_rndm();
1395  }
1396 
1397  particles_.reset();
1398 
1399  // Sample particles according to the initial conditions
1400  double start_time = modus_.initial_conditions(&particles_, parameters_);
1401  /* For box modus make sure that particles are in the box. In principle, after
1402  * a correct initialization they should be, so this is just playing it safe.
1403  */
1404  modus_.impose_boundary_conditions(&particles_, outputs_);
1405  // Reset the simulation clock
1406  double timestep = delta_time_startup_;
1407 
1408  switch (time_step_mode_) {
1409  case TimeStepMode::Fixed:
1410  break;
1411  case TimeStepMode::None:
1412  timestep = end_time_ - start_time;
1413  // Take care of the box modus + timestepless propagation
1414  const double max_dt = modus_.max_timestep(max_transverse_distance_sqr_);
1415  if (max_dt > 0. && max_dt < timestep) {
1416  timestep = max_dt;
1417  }
1418  break;
1419  }
1420  std::unique_ptr<UniformClock> clock_for_this_event;
1421  if (modus_.is_list() && (timestep < 0.0)) {
1422  throw std::runtime_error(
1423  "Timestep for the given event is negative. \n"
1424  "This might happen if the formation times of the input particles are "
1425  "larger than the specified end time of the simulation.");
1426  }
1427  clock_for_this_event = make_unique<UniformClock>(start_time, timestep);
1428  parameters_.labclock = std::move(clock_for_this_event);
1429 
1430  // Reset the output clock
1431  parameters_.outputclock->reset(start_time, true);
1432  // remove time before starting time in case of custom output times.
1433  parameters_.outputclock->remove_times_in_past(start_time);
1434 
1435  logg[LExperiment].debug(
1436  "Lab clock: t_start = ", parameters_.labclock->current_time(),
1437  ", dt = ", parameters_.labclock->timestep_duration());
1438 
1439  /* Save the initial conserved quantum numbers and total momentum in
1440  * the system for conservation checks */
1441  conserved_initial_ = QuantumNumbers(particles_);
1442  wall_actions_total_ = 0;
1443  previous_wall_actions_total_ = 0;
1444  interactions_total_ = 0;
1445  previous_interactions_total_ = 0;
1446  total_pauli_blocked_ = 0;
1447  projectile_target_interact_ = false;
1448  total_hypersurface_crossing_actions_ = 0;
1449  total_energy_removed_ = 0.0;
1450  // Print output headers
1451  logg[LExperiment].info() << hline;
1452  logg[LExperiment].info() << "Time[fm] Ekin[GeV] E_MF[GeV] ETotal[GeV] "
1453  << "ETot/N[GeV] D(ETot/N)[GeV] Scatt&Decays "
1454  << "Particles Comp.Time";
1455  logg[LExperiment].info() << hline;
1456  double E_mean_field = 0.0;
1457  if (potentials_) {
1458  update_potentials();
1459  // using the lattice is necessary
1460  if ((jmu_B_lat_ != nullptr)) {
1461  E_mean_field =
1462  calculate_mean_field_energy(*potentials_, *jmu_B_lat_, parameters_);
1463  }
1464  }
1465  initial_mean_field_energy_ = E_mean_field;
1467  particles_, 0u, conserved_initial_, time_start_,
1468  parameters_.labclock->current_time(), E_mean_field,
1469  initial_mean_field_energy_);
1470 }
1471 
1472 template <typename Modus>
1473 template <typename Container>
1475  Action &action, const Container &particles_before_actions) {
1476  // Make sure to skip invalid and Pauli-blocked actions.
1477  if (!action.is_valid(particles_)) {
1478  logg[LExperiment].debug(~einhard::DRed(), "✘ ", action,
1479  " (discarded: invalid)");
1480  return false;
1481  }
1482  action.generate_final_state();
1483  logg[LExperiment].debug("Process Type is: ", action.get_type());
1484  if (pauli_blocker_ && action.is_pauli_blocked(particles_, *pauli_blocker_)) {
1485  total_pauli_blocked_++;
1486  return false;
1487  }
1488  if (modus_.is_collider()) {
1489  /* Mark incoming nucleons as interacted - now they are permitted
1490  * to collide with nucleons from their native nucleus */
1491  bool incoming_projectile = false;
1492  bool incoming_target = false;
1493  for (const auto &incoming : action.incoming_particles()) {
1494  assert(incoming.id() >= 0);
1495  if (incoming.id() < modus_.total_N_number()) {
1496  nucleon_has_interacted_[incoming.id()] = true;
1497  }
1498  if (incoming.id() < modus_.proj_N_number()) {
1499  incoming_projectile = true;
1500  }
1501  if (incoming.id() >= modus_.proj_N_number() &&
1502  incoming.id() < modus_.total_N_number()) {
1503  incoming_target = true;
1504  }
1505  }
1506  // Check whether particles from different nuclei interacted.
1507  if (incoming_projectile & incoming_target) {
1508  projectile_target_interact_ = true;
1509  }
1510  }
1511  /* Make sure to pick a non-zero integer, because 0 is reserved for "no
1512  * interaction yet". */
1513  const auto id_process = static_cast<uint32_t>(interactions_total_ + 1);
1514  action.perform(&particles_, id_process);
1515  interactions_total_++;
1516  if (action.get_type() == ProcessType::Wall) {
1517  wall_actions_total_++;
1518  }
1519  if (action.get_type() == ProcessType::HyperSurfaceCrossing) {
1520  total_hypersurface_crossing_actions_++;
1521  total_energy_removed_ += action.incoming_particles()[0].momentum().x0();
1522  }
1523  // Calculate Eckart rest frame density at the interaction point
1524  double rho = 0.0;
1525  if (dens_type_ != DensityType::None) {
1526  const FourVector r_interaction = action.get_interaction_point();
1527  constexpr bool compute_grad = false;
1528  const bool smearing = true;
1529  rho = std::get<0>(current_eckart(r_interaction.threevec(),
1530  particles_before_actions, density_param_,
1531  dens_type_, compute_grad, smearing));
1532  }
1548  for (const auto &output : outputs_) {
1549  if (!output->is_dilepton_output() && !output->is_photon_output()) {
1550  if (output->is_IC_output() &&
1552  output->at_interaction(action, rho);
1553  } else if (!output->is_IC_output()) {
1554  output->at_interaction(action, rho);
1555  }
1556  }
1557  }
1558 
1559  // At every collision photons can be produced.
1560  // Note: We rely here on the lazy evaluation of the arguments to if.
1561  // It may happen that in a wall-crossing-action sqrt_s raises an exception.
1562  // Therefore we first have to check if the incoming particles can undergo
1563  // an em-interaction.
1564  if (photons_switch_ &&
1567  action.sqrt_s(), action.incoming_particles())) {
1568  /* Time in the action constructor is relative to
1569  * current time of incoming */
1570  constexpr double action_time = 0.;
1571  ScatterActionPhoton photon_act(action.incoming_particles(), action_time,
1572  n_fractional_photons_,
1573  action.get_total_weight());
1574 
1585  photon_act.add_dummy_hadronic_process(action.get_total_weight());
1586 
1587  // Now add the actual photon reaction channel.
1588  photon_act.add_single_process();
1589 
1590  photon_act.perform_photons(outputs_);
1591  }
1592 
1593  if (bremsstrahlung_switch_ &&
1595  action.incoming_particles())) {
1596  /* Time in the action constructor is relative to
1597  * current time of incoming */
1598  constexpr double action_time = 0.;
1599 
1600  BremsstrahlungAction brems_act(action.incoming_particles(), action_time,
1601  n_fractional_photons_,
1602  action.get_total_weight());
1603 
1615  brems_act.add_dummy_hadronic_process(action.get_total_weight());
1616 
1617  // Now add the actual bremsstrahlung reaction channel.
1618  brems_act.add_single_process();
1619 
1620  brems_act.perform_bremsstrahlung(outputs_);
1621  }
1622 
1623  logg[LExperiment].debug(~einhard::Green(), "✔ ", action);
1624  return true;
1625 }
1626 
1627 template <typename Modus>
1629  Actions actions;
1630 
1631  while (parameters_.labclock->current_time() < end_time_) {
1632  const double t = parameters_.labclock->current_time();
1633  const double dt =
1634  std::min(parameters_.labclock->timestep_duration(), end_time_ - t);
1635  logg[LExperiment].debug("Timestepless propagation for next ", dt, " fm/c.");
1636 
1637  // Perform forced thermalization if required
1638  if (thermalizer_ &&
1639  thermalizer_->is_time_to_thermalize(parameters_.labclock)) {
1640  const bool ignore_cells_under_treshold = true;
1641  thermalizer_->update_thermalizer_lattice(particles_, density_param_,
1642  ignore_cells_under_treshold);
1643  const double current_t = parameters_.labclock->current_time();
1644  thermalizer_->thermalize(particles_, current_t,
1645  parameters_.testparticles);
1646  ThermalizationAction th_act(*thermalizer_, current_t);
1647  if (th_act.any_particles_thermalized()) {
1648  perform_action(th_act, particles_);
1649  }
1650  }
1651 
1652  if (particles_.size() > 0 && action_finders_.size() > 0) {
1653  /* (1.a) Create grid. */
1654  double min_cell_length = compute_min_cell_length(dt);
1655  logg[LExperiment].debug("Creating grid with minimal cell length ",
1656  min_cell_length);
1657  const auto &grid =
1658  use_grid_ ? modus_.create_grid(particles_, min_cell_length, dt)
1659  : modus_.create_grid(particles_, min_cell_length, dt,
1661 
1662  const double cell_vol = grid.cell_volume();
1663 
1664  /* (1.b) Iterate over cells and find actions. */
1665  grid.iterate_cells(
1666  [&](const ParticleList &search_list) {
1667  for (const auto &finder : action_finders_) {
1668  actions.insert(finder->find_actions_in_cell(
1669  search_list, dt, cell_vol, beam_momentum_));
1670  }
1671  },
1672  [&](const ParticleList &search_list,
1673  const ParticleList &neighbors_list) {
1674  for (const auto &finder : action_finders_) {
1675  actions.insert(finder->find_actions_with_neighbors(
1676  search_list, neighbors_list, dt, beam_momentum_));
1677  }
1678  });
1679  }
1680 
1681  /* \todo (optimizations) Adapt timestep size here */
1682 
1683  /* (2) Propagation from action to action until the end of timestep */
1684  run_time_evolution_timestepless(actions);
1685 
1686  /* (3) Update potentials (if computed on the lattice) and
1687  * compute new momenta according to equations of motion */
1688  if (potentials_) {
1689  update_potentials();
1690  update_momenta(&particles_, parameters_.labclock->timestep_duration(),
1691  *potentials_, FB_lat_.get(), FI3_lat_.get());
1692  }
1693 
1694  /* (4) Expand universe if non-minkowskian metric; updates
1695  * positions and momenta according to the selected expansion */
1696  if (metric_.mode_ != ExpansionMode::NoExpansion) {
1697  expand_space_time(&particles_, parameters_, metric_);
1698  }
1699 
1700  ++(*parameters_.labclock);
1701 
1702  /* (5) Check conservation laws.
1703  *
1704  * Check conservation of conserved quantities if potentials and string
1705  * fragmentation are off. If potentials are on then momentum is conserved
1706  * only in average. If string fragmentation is on, then energy and
1707  * momentum are only very roughly conserved in high-energy collisions. */
1708  if (!potentials_ && !parameters_.strings_switch &&
1709  metric_.mode_ == ExpansionMode::NoExpansion && !IC_output_switch_) {
1710  std::string err_msg = conserved_initial_.report_deviations(particles_);
1711  if (!err_msg.empty()) {
1712  logg[LExperiment].error() << err_msg;
1713  throw std::runtime_error("Violation of conserved quantities!");
1714  }
1715  }
1716  }
1717 
1718  if (pauli_blocker_) {
1719  logg[LExperiment].info(
1720  "Interactions: Pauli-blocked/performed = ", total_pauli_blocked_, "/",
1721  interactions_total_ - wall_actions_total_);
1722  }
1723 }
1724 
1725 template <typename Modus>
1727  const double dt =
1728  propagate_straight_line(&particles_, to_time, beam_momentum_);
1729  if (dilepton_finder_ != nullptr) {
1730  for (const auto &output : outputs_) {
1731  dilepton_finder_->shine(particles_, output.get(), dt);
1732  }
1733  }
1734 }
1735 
1743 inline void check_interactions_total(uint64_t interactions_total) {
1744  constexpr uint64_t max_uint32 = std::numeric_limits<uint32_t>::max();
1745  if (interactions_total >= max_uint32) {
1746  throw std::runtime_error("Integer overflow in total interaction number!");
1747  }
1748 }
1749 
1750 template <typename Modus>
1752  const double start_time = parameters_.labclock->current_time();
1753  const double end_time =
1754  std::min(parameters_.labclock->next_time(), end_time_);
1755  double time_left = end_time - start_time;
1756  logg[LExperiment].debug(
1757  "Timestepless propagation: ", "Actions size = ", actions.size(),
1758  ", start time = ", start_time, ", end time = ", end_time);
1759 
1760  // iterate over all actions
1761  while (!actions.is_empty()) {
1762  // get next action
1763  ActionPtr act = actions.pop();
1764  if (!act->is_valid(particles_)) {
1765  logg[LExperiment].debug(~einhard::DRed(), "✘ ", act,
1766  " (discarded: invalid)");
1767  continue;
1768  }
1769  if (act->time_of_execution() > end_time) {
1770  logg[LExperiment].error(
1771  act, " scheduled later than end time: t_action[fm/c] = ",
1772  act->time_of_execution(), ", t_end[fm/c] = ", end_time);
1773  }
1774  logg[LExperiment].debug(~einhard::Green(), "✔ ", act);
1775 
1776  while (next_output_time() <= act->time_of_execution()) {
1777  logg[LExperiment].debug("Propagating until output time: ",
1778  next_output_time());
1779  propagate_and_shine(next_output_time());
1780  ++(*parameters_.outputclock);
1781  intermediate_output();
1782  }
1783 
1784  /* (1) Propagate to the next action. */
1785  logg[LExperiment].debug("Propagating until next action ", act,
1786  ", action time = ", act->time_of_execution());
1787  propagate_and_shine(act->time_of_execution());
1788 
1789  /* (2) Perform action.
1790  *
1791  * Update the positions of the incoming particles, because the information
1792  * in the action object will be outdated as the particles have been
1793  * propagated since the construction of the action. */
1794  act->update_incoming(particles_);
1795  const bool performed = perform_action(*act, particles_);
1796 
1797  /* No need to update actions for outgoing particles
1798  * if the action is not performed. */
1799  if (!performed) {
1800  continue;
1801  }
1802 
1803  /* (3) Update actions for newly-produced particles. */
1804 
1805  time_left = end_time - act->time_of_execution();
1806  const ParticleList &outgoing_particles = act->outgoing_particles();
1807  // Cell volume set to zero, since there is no grid
1808  const double cell_vol = 0.0;
1809  for (const auto &finder : action_finders_) {
1810  // Outgoing particles can still decay, cross walls...
1811  actions.insert(finder->find_actions_in_cell(outgoing_particles, time_left,
1812  cell_vol, beam_momentum_));
1813  // ... and collide with other particles.
1814  actions.insert(finder->find_actions_with_surrounding_particles(
1815  outgoing_particles, particles_, time_left, beam_momentum_));
1816  }
1817 
1818  check_interactions_total(interactions_total_);
1819  }
1820 
1821  while (next_output_time() <= end_time) {
1822  logg[LExperiment].debug("Propagating until output time: ",
1823  next_output_time());
1824  propagate_and_shine(next_output_time());
1825  ++(*parameters_.outputclock);
1826  // Avoid duplicating printout at event end time
1827  if (parameters_.outputclock->current_time() < end_time_) {
1828  intermediate_output();
1829  }
1830  }
1831  logg[LExperiment].debug("Propagating to time ", end_time);
1832  propagate_and_shine(end_time);
1833 }
1834 
1835 template <typename Modus>
1837  const uint64_t wall_actions_this_interval =
1838  wall_actions_total_ - previous_wall_actions_total_;
1839  previous_wall_actions_total_ = wall_actions_total_;
1840  const uint64_t interactions_this_interval = interactions_total_ -
1841  previous_interactions_total_ -
1842  wall_actions_this_interval;
1843  previous_interactions_total_ = interactions_total_;
1844  double E_mean_field = 0.0;
1845  if (potentials_) {
1846  // using the lattice is necessary
1847  if ((jmu_B_lat_ != nullptr)) {
1848  E_mean_field =
1849  calculate_mean_field_energy(*potentials_, *jmu_B_lat_, parameters_);
1850  /*
1851  * Mean field calculated in a box should remain approximately constant if
1852  * the system is in equilibrium, and so deviations from its original value
1853  * may signal a phase transition or other dynamical process. This
1854  * comparison only makes sense in the Box Modus, hence the condition.
1855  */
1856  if (modus_.length() > 0.0) {
1857  double tmp = (E_mean_field - initial_mean_field_energy_) /
1858  (E_mean_field + initial_mean_field_energy_);
1859  /*
1860  * This is displayed when the system evolves away from its initial
1861  * configuration (which is when the total mean field energy in the box
1862  * deviates from its initial value).
1863  */
1864  if (std::abs(tmp) > 0.01) {
1865  logg[LExperiment].info()
1866  << "\n\n\n\t The mean field at t = "
1867  << parameters_.outputclock->current_time()
1868  << " [fm/c] differs from the mean field at t = 0:"
1869  << "\n\t\t initial_mean_field_energy_ = "
1870  << initial_mean_field_energy_ << " [GeV]"
1871  << "\n\t\t abs[(E_MF - E_MF(t=0))/(E_MF + E_MF(t=0))] = "
1872  << std::abs(tmp)
1873  << "\n\t\t E_MF/E_MF(t=0) = "
1874  << E_mean_field / initial_mean_field_energy_ << "\n\n";
1875  }
1876  }
1877  }
1878  }
1879 
1881  particles_, interactions_this_interval, conserved_initial_, time_start_,
1882  parameters_.outputclock->current_time(), E_mean_field,
1883  initial_mean_field_energy_);
1884  const LatticeUpdate lat_upd = LatticeUpdate::AtOutput;
1885  // save evolution data
1886  if (!(modus_.is_box() && parameters_.outputclock->current_time() <
1887  modus_.equilibration_time())) {
1888  for (const auto &output : outputs_) {
1889  if (output->is_dilepton_output() || output->is_photon_output() ||
1890  output->is_IC_output()) {
1891  continue;
1892  }
1893 
1894  output->at_intermediate_time(particles_, parameters_.outputclock,
1895  density_param_);
1896 
1897  // Thermodynamic output on the lattice versus time
1898  switch (dens_type_lattice_printout_) {
1899  case DensityType::Baryon:
1900  update_lattice(jmu_B_lat_.get(), lat_upd, DensityType::Baryon,
1901  density_param_, particles_, false);
1902  output->thermodynamics_output(ThermodynamicQuantity::EckartDensity,
1903  DensityType::Baryon, *jmu_B_lat_);
1904  break;
1906  update_lattice(jmu_I3_lat_.get(), lat_upd,
1907  DensityType::BaryonicIsospin, density_param_,
1908  particles_, false);
1909  output->thermodynamics_output(ThermodynamicQuantity::EckartDensity,
1911  *jmu_I3_lat_);
1912  break;
1913  case DensityType::None:
1914  break;
1915  default:
1916  update_lattice(jmu_custom_lat_.get(), lat_upd,
1917  dens_type_lattice_printout_, density_param_,
1918  particles_, false);
1919  output->thermodynamics_output(ThermodynamicQuantity::EckartDensity,
1920  dens_type_lattice_printout_,
1921  *jmu_custom_lat_);
1922  }
1923  if (printout_tmn_ || printout_tmn_landau_ || printout_v_landau_) {
1924  update_lattice(Tmn_.get(), lat_upd, dens_type_lattice_printout_,
1925  density_param_, particles_);
1926  if (printout_tmn_) {
1927  output->thermodynamics_output(ThermodynamicQuantity::Tmn,
1928  dens_type_lattice_printout_, *Tmn_);
1929  }
1930  if (printout_tmn_landau_) {
1931  output->thermodynamics_output(ThermodynamicQuantity::TmnLandau,
1932  dens_type_lattice_printout_, *Tmn_);
1933  }
1934  if (printout_v_landau_) {
1935  output->thermodynamics_output(ThermodynamicQuantity::LandauVelocity,
1936  dens_type_lattice_printout_, *Tmn_);
1937  }
1938  }
1939 
1940  if (thermalizer_) {
1941  output->thermodynamics_output(*thermalizer_);
1942  }
1943  }
1944  }
1945 }
1946 
1947 template <typename Modus>
1949  if (potentials_) {
1950  if (potentials_->use_symmetry() && jmu_I3_lat_ != nullptr) {
1951  update_lattice(jmu_I3_lat_.get(), LatticeUpdate::EveryTimestep,
1952  DensityType::BaryonicIsospin, density_param_, particles_,
1953  true);
1954  }
1955  if ((potentials_->use_skyrme() || potentials_->use_symmetry()) &&
1956  jmu_B_lat_ != nullptr) {
1957  update_lattice(jmu_B_lat_.get(), LatticeUpdate::EveryTimestep,
1958  DensityType::Baryon, density_param_, particles_, true);
1959  const size_t UBlattice_size = UB_lat_->size();
1960  assert(UBlattice_size == UI3_lat_->size());
1961  for (size_t i = 0; i < UBlattice_size; i++) {
1962  auto jB = (*jmu_B_lat_)[i];
1963  const FourVector flow_four_velocity_B =
1964  std::abs(jB.density()) > really_small ? jB.jmu_net() / jB.density()
1965  : FourVector();
1966  double baryon_density = jB.density();
1967  ThreeVector baryon_grad_rho = jB.grad_rho();
1968  ThreeVector baryon_dj_dt = jB.dj_dt();
1969  ThreeVector baryon_rot_j = jB.rot_j();
1970  if (potentials_->use_skyrme()) {
1971  (*UB_lat_)[i] =
1972  flow_four_velocity_B * potentials_->skyrme_pot(baryon_density);
1973  (*FB_lat_)[i] = potentials_->skyrme_force(
1974  baryon_density, baryon_grad_rho, baryon_dj_dt, baryon_rot_j);
1975  }
1976  if (potentials_->use_symmetry() && jmu_I3_lat_ != nullptr) {
1977  auto jI3 = (*jmu_I3_lat_)[i];
1978  const FourVector flow_four_velocity_I3 =
1979  std::abs(jI3.density()) > really_small
1980  ? jI3.jmu_net() / jI3.density()
1981  : FourVector();
1982  (*UI3_lat_)[i] =
1983  flow_four_velocity_I3 *
1984  potentials_->symmetry_pot(jI3.density(), baryon_density);
1985  (*FI3_lat_)[i] = potentials_->symmetry_force(
1986  jI3.density(), jI3.grad_rho(), jI3.dj_dt(), jI3.rot_j(),
1987  baryon_density, baryon_grad_rho, baryon_dj_dt, baryon_rot_j);
1988  }
1989  }
1990  }
1991  }
1992 }
1993 
1994 template <typename Modus>
1996  /* At end of time evolution: Force all resonances to decay. In order to handle
1997  * decay chains, we need to loop until no further actions occur. */
1998  uint64_t interactions_old;
1999  const auto particles_before_actions = particles_.copy_to_vector();
2000  do {
2001  Actions actions;
2002 
2003  interactions_old = interactions_total_;
2004 
2005  // Dileptons: shining of remaining resonances
2006  if (dilepton_finder_ != nullptr) {
2007  for (const auto &output : outputs_) {
2008  dilepton_finder_->shine_final(particles_, output.get(), true);
2009  }
2010  }
2011  // Find actions.
2012  for (const auto &finder : action_finders_) {
2013  actions.insert(finder->find_final_actions(particles_));
2014  }
2015  // Perform actions.
2016  while (!actions.is_empty()) {
2017  perform_action(*actions.pop(), particles_before_actions);
2018  }
2019  // loop until no more decays occur
2020  } while (interactions_total_ > interactions_old);
2021 
2022  // Dileptons: shining of stable particles at the end
2023  if (dilepton_finder_ != nullptr) {
2024  for (const auto &output : outputs_) {
2025  dilepton_finder_->shine_final(particles_, output.get(), false);
2026  }
2027  }
2028 }
2029 
2030 template <typename Modus>
2031 void Experiment<Modus>::final_output(const int evt_num) {
2032  /* make sure the experiment actually ran (note: we should compare this
2033  * to the start time, but we don't know that. Therefore, we check that
2034  * the time is positive, which should heuristically be the same). */
2035  if (likely(parameters_.labclock > 0)) {
2036  const uint64_t wall_actions_this_interval =
2037  wall_actions_total_ - previous_wall_actions_total_;
2038  const uint64_t interactions_this_interval = interactions_total_ -
2039  previous_interactions_total_ -
2040  wall_actions_this_interval;
2041  double E_mean_field = 0.0;
2042  if (potentials_) {
2043  // using the lattice is necessary
2044  if ((jmu_B_lat_ != nullptr)) {
2045  E_mean_field =
2046  calculate_mean_field_energy(*potentials_, *jmu_B_lat_, parameters_);
2047  }
2048  }
2050  particles_, interactions_this_interval, conserved_initial_, time_start_,
2051  end_time_, E_mean_field, initial_mean_field_energy_);
2052  if (IC_output_switch_ && (particles_.size() == 0)) {
2053  // Verify there is no more energy in the system if all particles were
2054  // removed when crossing the hypersurface
2055  const double remaining_energy =
2056  conserved_initial_.momentum().x0() - total_energy_removed_;
2057  if (remaining_energy > really_small) {
2058  throw std::runtime_error(
2059  "There is remaining energy in the system although all particles "
2060  "were removed.\n"
2061  "E_remain = " +
2062  std::to_string(remaining_energy) + " [GeV]");
2063  } else {
2064  logg[LExperiment].info() << hline;
2065  logg[LExperiment].info()
2066  << "Time real: " << SystemClock::now() - time_start_;
2067  logg[LExperiment].info()
2068  << "Interactions before reaching hypersurface: "
2069  << interactions_total_ - wall_actions_total_ -
2070  total_hypersurface_crossing_actions_;
2071  logg[LExperiment].info()
2072  << "Total number of particles removed on hypersurface: "
2073  << total_hypersurface_crossing_actions_;
2074  }
2075  } else {
2076  logg[LExperiment].info() << hline;
2077  logg[LExperiment].info()
2078  << "Time real: " << SystemClock::now() - time_start_;
2079  logg[LExperiment].info() << "Final interaction number: "
2080  << interactions_total_ - wall_actions_total_;
2081  }
2082 
2083  // Check if there are unformed particles
2084  int unformed_particles_count = 0;
2085  for (const auto &particle : particles_) {
2086  if (particle.formation_time() > end_time_) {
2087  unformed_particles_count++;
2088  }
2089  }
2090  if (unformed_particles_count > 0) {
2091  logg[LExperiment].warn(
2092  "End time might be too small. ", unformed_particles_count,
2093  " unformed particles were found at the end of the evolution.");
2094  }
2095  }
2096 
2097  for (const auto &output : outputs_) {
2098  output->at_eventend(particles_, evt_num, modus_.impact_parameter(),
2099  !projectile_target_interact_);
2100  }
2101 }
2102 
2103 template <typename Modus>
2105  const auto &mainlog = logg[LMain];
2106  for (int j = 0; j < nevents_; j++) {
2107  mainlog.info() << "Event " << j;
2108 
2109  // Sample initial particles, start clock, some printout and book-keeping
2110  initialize_new_event();
2111  /* In the ColliderModus, if the first collisions within the same nucleus are
2112  * forbidden, 'nucleon_has_interacted_', which records whether a nucleon has
2113  * collided with another nucleon, is initialized equal to false. If allowed,
2114  * 'nucleon_has_interacted' is initialized equal to true, which means these
2115  * incoming particles have experienced some fake scatterings, they can
2116  * therefore collide with each other later on since these collisions are not
2117  * "first" to them. */
2118  if (modus_.is_collider()) {
2119  if (!modus_.cll_in_nucleus()) {
2120  nucleon_has_interacted_.assign(modus_.total_N_number(), false);
2121  } else {
2122  nucleon_has_interacted_.assign(modus_.total_N_number(), true);
2123  }
2124  }
2125  /* In the ColliderModus, if Fermi motion is frozen, assign the beam momenta
2126  * to the nucleons in both the projectile and the target. */
2127  if (modus_.is_collider() && modus_.fermi_motion() == FermiMotion::Frozen) {
2128  for (int i = 0; i < modus_.total_N_number(); i++) {
2129  const auto mass_beam = particles_.copy_to_vector()[i].effective_mass();
2130  const auto v_beam = i < modus_.proj_N_number()
2131  ? modus_.velocity_projectile()
2132  : modus_.velocity_target();
2133  const auto gamma = 1.0 / std::sqrt(1.0 - v_beam * v_beam);
2134  beam_momentum_.emplace_back(FourVector(gamma * mass_beam, 0.0, 0.0,
2135  gamma * v_beam * mass_beam));
2136  }
2137  }
2138 
2139  // Output at event start
2140  for (const auto &output : outputs_) {
2141  output->at_eventstart(particles_, j);
2142  }
2143 
2144  run_time_evolution();
2145 
2146  if (force_decays_) {
2147  do_final_decays();
2148  }
2149 
2150  // Output at event end
2151  final_output(j);
2152  }
2153 }
2154 
2155 } // namespace smash
2156 
2157 #endif // SRC_INCLUDE_EXPERIMENT_H_
smash::LatticeUpdate
LatticeUpdate
Enumerator option for lattice updates.
Definition: lattice.h:36
smash::Experiment::jmu_I3_lat_
std::unique_ptr< DensityLattice > jmu_I3_lat_
Isospin projection density on the lattices.
Definition: experiment.h:379
smash::Experiment::particles
Particles * particles()
Provides external access to SMASH particles.
Definition: experiment.h:223
smash::Action::get_type
virtual ProcessType get_type() const
Get the process type.
Definition: action.h:131
smash::Experiment::time_step_mode_
const TimeStepMode time_step_mode_
This indicates whether to use time steps.
Definition: experiment.h:488
smash::Experiment::jmu_B_lat_
std::unique_ptr< DensityLattice > jmu_B_lat_
Baryon density on the lattices.
Definition: experiment.h:376
smash
Definition: action.h:24
smash::Experiment::run_time_evolution
void run_time_evolution()
Runs the time evolution of an event with fixed-sized time steps or without timesteps,...
Definition: experiment.h:1628
smash::SystemTimePoint
std::chrono::time_point< std::chrono::system_clock > SystemTimePoint
Type (alias) that is used to store the current time.
Definition: chrono.h:22
smash::Experiment::printout_tmn_
bool printout_tmn_
Whether to print the energy-momentum tensor.
Definition: experiment.h:417
quantumnumbers.h
outputparameters.h
ThermodynamicQuantity::TmnLandau
smash::Experiment::Tmn_
std::unique_ptr< RectangularLattice< EnergyMomentumTensor > > Tmn_
Lattices of energy-momentum tensors for printout.
Definition: experiment.h:414
smash::Experiment::action_finders_
std::vector< std::unique_ptr< ActionFinderInterface > > action_finders_
The Action finder objects.
Definition: experiment.h:364
smash::LInitialConditions
static constexpr int LInitialConditions
Definition: experiment.h:78
smash::expand_space_time
void expand_space_time(Particles *particles, const ExperimentParameters &parameters, const ExpansionProperties &metric)
Modifies positions and momentum of all particles to account for space-time deformation.
Definition: propagation.cc:86
smash::Experiment::photons_switch_
const bool photons_switch_
This indicates whether photons are switched on.
Definition: experiment.h:479
smash::maximum_cross_section
constexpr double maximum_cross_section
The maximal cross section (in mb) for which it is guaranteed that all collisions with this cross sect...
Definition: constants.h:111
smash::ScatterActionPhoton::is_kinematically_possible
static bool is_kinematically_possible(const double s_sqrt, const ParticleList &in)
Check if CM-energy is sufficient to produce hadron in final state.
Definition: scatteractionphoton.cc:157
smash::Action::generate_final_state
virtual void generate_final_state()=0
Generate the final state for this action.
smash::Experiment::total_energy_removed_
double total_energy_removed_
Total energy removed from the system in hypersurface crossing actions.
Definition: experiment.h:554
ThermodynamicQuantity::LandauVelocity
smash::Experiment::previous_wall_actions_total_
uint64_t previous_wall_actions_total_
Total number of wall-crossings for previous timestep.
Definition: experiment.h:536
smash::DensityParameters
A class to pre-calculate and store parameters relevant for density calculation.
Definition: density.h:107
smash::BremsstrahlungAction::is_bremsstrahlung_reaction
static bool is_bremsstrahlung_reaction(const ParticleList &in)
Check if particles can undergo an implemented photon process.
Definition: bremsstrahlungaction.h:132
smash::Experiment::beam_momentum_
std::vector< FourVector > beam_momentum_
The initial nucleons in the ColliderModus propagate with beam_momentum_, if Fermi motion is frozen.
Definition: experiment.h:361
smash::Experiment::next_output_time
double next_output_time() const
Shortcut for next output time.
Definition: experiment.h:298
smash::Experiment::compute_min_cell_length
double compute_min_cell_length(double dt) const
Calculate the minimal size for the grid cells such that the ScatterActionsFinder will find all collis...
Definition: experiment.h:293
TimeStepMode
TimeStepMode
The time step mode.
Definition: forwarddeclarations.h:113
smash::Actions::insert
void insert(ActionList &&new_acts)
Insert a list of actions into this object.
Definition: actions.h:76
hypersurfacecrossingaction.h
smash::LatticeUpdate::AtOutput
smash::Experiment::wall_actions_total_
uint64_t wall_actions_total_
Total number of wall-crossings for current timestep.
Definition: experiment.h:530
smash::Experiment::final_output
void final_output(const int evt_num)
Output at the end of an event.
Definition: experiment.h:2031
smash::check_interactions_total
void check_interactions_total(uint64_t interactions_total)
Make sure interactions_total can be represented as a 32-bit integer.
Definition: experiment.h:1743
smash::Experiment::UB_lat_
std::unique_ptr< RectangularLattice< FourVector > > UB_lat_
Lattices for Skyrme potentials (evaluated in the local rest frame) times the baryon flow 4-velocity.
Definition: experiment.h:397
smash::LExperiment
static constexpr int LExperiment
Definition: outputparameters.h:19
decayactionsfinder.h
smash::Experiment::pauli_blocker_
std::unique_ptr< PauliBlocker > pauli_blocker_
An instance of PauliBlocker class that stores parameters needed for Pauli blocking calculations and c...
Definition: experiment.h:331
binaryoutput.h
smash::DensityType::BaryonicIsospin
smash::Action::is_pauli_blocked
bool is_pauli_blocked(const Particles &particles, const PauliBlocker &p_bl) const
Check if the action is Pauli-blocked.
Definition: action.cc:36
ExpansionMode::NoExpansion
smash::Experiment::modus_
Modus modus_
Instance of the Modus template parameter.
Definition: experiment.h:316
smash::Experiment::modus
Modus * modus()
Provides external access to SMASH calculation modus.
Definition: experiment.h:229
smash::Experiment::initialize_new_event
void initialize_new_event()
This is called in the beginning of each event.
Definition: experiment.h:1376
smash::hline
const std::string hline(113, '-')
String representing a horizontal line.
smash::format_measurements
std::string format_measurements(const Particles &particles, uint64_t scatterings_this_interval, const QuantumNumbers &conserved_initial, SystemTimePoint time_start, double time, double E_mean_field, double E_mean_field_initial)
Generate the tabulated string which will be printed to the screen when SMASH is running.
Definition: experiment.cc:449
smash::Experiment::projectile_target_interact_
bool projectile_target_interact_
Whether the projectile and the target collided.
Definition: experiment.h:354
FermiMotion::Off
Don't use fermi motion.
grid.h
smash::operator<<
std::ostream & operator<<(std::ostream &out, const ActionPtr &action)
Definition: action.h:463
energymomentumtensor.h
smash::Experiment::jmu_custom_lat_
std::unique_ptr< DensityLattice > jmu_custom_lat_
Custom density on the lattices.
Definition: experiment.h:388
smash::Experiment::run_time_evolution_timestepless
void run_time_evolution_timestepless(Actions &actions)
Performs all the propagations and actions during a certain time interval neglecting the influence of ...
Definition: experiment.h:1751
smash::random::set_seed
void set_seed(T &&seed)
Sets the seed of the random number engine.
Definition: random.h:71
smash::Experiment::particles_
Particles particles_
Complete particle list.
Definition: experiment.h:319
propagation.h
smash::LatticeUpdate::EveryTimestep
smash::LMain
static constexpr int LMain
Definition: experiment.h:77
smash::Experiment::run
void run() override
Runs the experiment.
Definition: experiment.h:2104
fourvector.h
smash::Experiment::force_decays_
const bool force_decays_
This indicates whether we force all resonances to decay in the last timestep.
Definition: experiment.h:467
smash::Experiment::total_pauli_blocked_
uint64_t total_pauli_blocked_
Total number of Pauli-blockings for current timestep.
Definition: experiment.h:542
smash::Experiment::dileptons_switch_
const bool dileptons_switch_
This indicates whether dileptons are switched on.
Definition: experiment.h:476
smash::Experiment::dens_type_
DensityType dens_type_
Type of density to be written to collision headers.
Definition: experiment.h:512
FermiMotion::Frozen
Use fermi motion without potentials.
smash::logg
std::array< einhard::Logger<>, std::tuple_size< LogArea::AreaTuple >::value > logg
An array that stores all pre-configured Logger objects.
Definition: logging.cc:39
smash::Experiment::photon_output_
OutputPtr photon_output_
The Photon output.
Definition: experiment.h:343
grandcan_thermalizer.h
smash::Experiment::thermalizer_
std::unique_ptr< GrandCanThermalizer > thermalizer_
Instance of class used for forced thermalization.
Definition: experiment.h:429
smash::really_small
constexpr double really_small
Numerical error tolerance.
Definition: constants.h:37
pauliblocking.h
smash::Configuration
Interface to the SMASH configuration files.
Definition: configuration.h:464
smash::StringProcess
String excitation processes used in SMASH.
Definition: processstring.h:47
smash::UI3_lat_pointer
RectangularLattice< FourVector > * UI3_lat_pointer
Pointer to the symmmetry potential on the lattice.
Definition: potential_globals.cc:16
wallcrossingaction.h
smash::Experiment::density_param_
DensityParameters density_param_
Structure to precalculate and hold parameters for density computations.
Definition: experiment.h:310
smash::Experiment::nevents_
const int nevents_
Number of events.
Definition: experiment.h:451
smash::Experiment::seed_
int64_t seed_
random seed for the next event.
Definition: experiment.h:557
chrono.h
smash::ThreeVector
Definition: threevector.h:31
source_location
#define source_location
Hackery that is required to output the location in the source code where the log statement occurs.
Definition: logging.h:240
smash::ProcessType::Wall
box wall crossing
smash::ScatterActionPhoton::is_photon_reaction
static bool is_photon_reaction(const ParticleList &in)
Check if particles can undergo an implemented photon process.
Definition: scatteractionphoton.h:145
smash::OutputParameters
Helper structure for Experiment to hold output options and parameters.
Definition: outputparameters.h:25
smash::format
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:304
einhard::Color
A stream modifier that allows to colorize the log output.
Definition: einhard.hpp:142
smash::create_oscar_output
std::unique_ptr< OutputInterface > create_oscar_output(const std::string &format, const std::string &content, const bf::path &path, const OutputParameters &out_par)
Definition: oscaroutput.cc:768
smash::Experiment::metric_
const ExpansionProperties metric_
This struct contains information on the metric to be used.
Definition: experiment.h:473
actions.h
smash::Experiment::printout_tmn_landau_
bool printout_tmn_landau_
Whether to print the energy-momentum tensor in Landau frame.
Definition: experiment.h:420
smash::Experiment::FI3_lat_
std::unique_ptr< RectangularLattice< std::pair< ThreeVector, ThreeVector > > > FI3_lat_
Lattices for the electric and magnetic component of the symmetry force.
Definition: experiment.h:411
smash::create_experiment_parameters
ExperimentParameters create_experiment_parameters(Configuration config)
Gathers all general Experiment parameters.
Definition: experiment.cc:328
smash::Experiment::do_final_decays
void do_final_decays()
Performs the final decays of an event.
Definition: experiment.h:1995
smash::update_lattice
void update_lattice(RectangularLattice< T > *lat, const LatticeUpdate update, const DensityType dens_type, const DensityParameters &par, const Particles &particles, const bool compute_gradient=false)
Updates the contents on the lattice.
Definition: density.h:401
smash::UB_lat_pointer
RectangularLattice< FourVector > * UB_lat_pointer
Pointer to the skyrme potential on the lattice.
Definition: potential_globals.cc:15
smash::Actions::pop
ActionPtr pop()
Return the first action in the list and removes it from the list.
Definition: actions.h:59
smash::make_unique
std::unique_ptr< T > make_unique(Args &&... args)
Definition for make_unique Is in C++14's standard library; necessary for older compilers.
Definition: cxx14compat.h:25
rootoutput.h
ThermodynamicQuantity::Tmn
bremsstrahlungaction.h
smash::ThermalizationAction::any_particles_thermalized
bool any_particles_thermalized() const
This method checks, if there are particles in the region to be thermalized.
Definition: thermalizationaction.h:39
smash::Action::is_valid
bool is_valid(const Particles &particles) const
Check whether the action still applies.
Definition: action.cc:30
smash::Experiment::use_grid_
const bool use_grid_
This indicates whether to use the grid.
Definition: experiment.h:470
smash::Experiment::interactions_total_
uint64_t interactions_total_
Total number of interactions for current timestep.
Definition: experiment.h:518
smash::ExperimentBase::InvalidModusRequest
Definition: experiment.h:132
smash::Action::get_total_weight
virtual double get_total_weight() const =0
Return the total weight value, which is mainly used for the weight output entry.
smash::Experiment::intermediate_output
void intermediate_output()
Intermediate output during an event.
Definition: experiment.h:1836
smash::operator<<
std::ostream & operator<<(std::ostream &out, const Experiment< Modus > &e)
Creates a verbose textual description of the setup of the Experiment.
Definition: experiment.h:569
smash::Actions::is_empty
bool is_empty() const
Definition: actions.h:52
smash::Experiment::propagate_and_shine
void propagate_and_shine(double to_time)
Propagate all particles until time to_time without any interactions and shine dileptons.
Definition: experiment.h:1726
smash::ThermalizationAction
Definition: thermalizationaction.h:24
smash::ParticleData::formation_power_
static double formation_power_
Power with which the cross section scaling factor grows in time.
Definition: particledata.h:368
actionfinderfactory.h
decayactionsfinderdilepton.h
smash::Action::perform
virtual void perform(Particles *particles, uint32_t id_process)
Actually perform the action, e.g.
Definition: action.cc:95
smash::Experiment::delta_time_startup_
const double delta_time_startup_
The clock's timestep size at start up.
Definition: experiment.h:461
smash::Experiment::IC_output_switch_
const bool IC_output_switch_
This indicates whether the IC output is enabled.
Definition: experiment.h:485
smash::current_eckart
std::tuple< double, FourVector, ThreeVector, ThreeVector, ThreeVector > current_eckart(const ThreeVector &r, const ParticleList &plist, const DensityParameters &par, DensityType dens_type, bool compute_gradient, bool smearing)
Calculates Eckart rest frame density and 4-current of a given density type and optionally the gradien...
Definition: density.cc:149
smash::Experiment::Experiment
Experiment(Configuration config, const bf::path &output_path)
Create a new Experiment.
smash::CellSizeStrategy::Largest
Make cells as large as possible.
ThermodynamicQuantity::EckartDensity
vtkoutput.h
oscaroutput.h
smash::Experiment::previous_interactions_total_
uint64_t previous_interactions_total_
Total number of interactions for previous timestep.
Definition: experiment.h:524
smash::Action::sqrt_s
double sqrt_s() const
Determine the total energy in the center-of-mass frame [GeV].
Definition: action.h:266
smash::Particles
Definition: particles.h:33
smash::DensityType
DensityType
Allows to choose which kind of density to calculate.
Definition: density.h:35
scatteractionphoton.h
smash::Experiment::perform_action
bool perform_action(Action &action, const Container &particles_before_actions)
Perform the given action.
smash::Experiment::update_potentials
void update_potentials()
Recompute potentials on lattices if necessary.
Definition: experiment.h:1948
smash::random::advance
Engine::result_type advance()
Advance the engine's state and return the generated value.
Definition: random.h:78
smash::Experiment::FB_lat_
std::unique_ptr< RectangularLattice< std::pair< ThreeVector, ThreeVector > > > FB_lat_
Lattices for the electric and magnetic components of the Skyrme force.
Definition: experiment.h:407
smash::ExperimentParameters::outputclock
std::unique_ptr< Clock > outputclock
Output clock to keep track of the next output time.
Definition: experimentparameters.h:28
smash::fm2_mb
constexpr double fm2_mb
mb <-> fm^2 conversion factor.
Definition: constants.h:28
smash::ExperimentBase::~ExperimentBase
virtual ~ExperimentBase()=default
The virtual destructor avoids undefined behavior when destroying derived objects.
processstring.h
smash::Experiment::dens_type_lattice_printout_
DensityType dens_type_lattice_printout_
Type of density for lattice printout.
Definition: experiment.h:391
smash::Experiment::dilepton_finder_
std::unique_ptr< DecayActionsFinderDilepton > dilepton_finder_
The Dilepton Action Finder.
Definition: experiment.h:367
thermalizationaction.h
smash::Actions::size
ActionList::size_type size() const
Definition: actions.h:95
smash::Experiment::printout_lattice_td_
bool printout_lattice_td_
Whether to print the thermodynamics quantities evaluated on the lattices.
Definition: experiment.h:426
smash::Action
Definition: action.h:35
smash::Experiment::create_output
void create_output(const std::string &format, const std::string &content, const bf::path &output_path, const OutputParameters &par)
Create a list of output files.
Definition: experiment.h:576
smash::calculate_mean_field_energy
double calculate_mean_field_energy(const Potentials &potentials, RectangularLattice< smash::DensityOnLattice > &jmu_B_lat, const ExperimentParameters &parameters)
Calculate the total mean field energy of the system; this will be printed to the screen when SMASH is...
Definition: experiment.cc:496
smash::Experiment::n_fractional_photons_
int n_fractional_photons_
Number of fractional photons produced per single reaction.
Definition: experiment.h:373
smash::Experiment::dilepton_output_
OutputPtr dilepton_output_
The Dilepton output.
Definition: experiment.h:340
smash::Experiment::parameters_
ExperimentParameters parameters_
Struct of several member variables.
Definition: experiment.h:307
smash::ProcessType::HyperSurfaceCrossing
Hypersurface crossing Particles are removed from the evolution and printed to a separate output to se...
smash::ExperimentParameters
Helper structure for Experiment.
Definition: experimentparameters.h:23
smash::DensityType::None
smash::FourVector
Definition: fourvector.h:33
smash::Experiment::end_time_
const double end_time_
simulation time at which the evolution is stopped.
Definition: experiment.h:454
smash::Experiment::nucleon_has_interacted_
std::vector< bool > nucleon_has_interacted_
nucleon_has_interacted_ labels whether the particles in the nuclei have experienced any collisions or...
Definition: experiment.h:350
smash::Experiment::process_string_ptr_
StringProcess * process_string_ptr_
Pointer to the string process class object, which is used to set the random seed for PYTHIA objects i...
Definition: experiment.h:435
smash::Action::get_interaction_point
FourVector get_interaction_point() const
Get the interaction point.
Definition: action.cc:69
smash::propagate_straight_line
double propagate_straight_line(Particles *particles, double to_time, const std::vector< FourVector > &beam_momentum)
Propagates the positions of all particles on a straight line to a given moment.
Definition: propagation.cc:44
smash::Experiment::potentials_
std::unique_ptr< Potentials > potentials_
An instance of potentials class, that stores parameters of potentials, calculates them and their grad...
Definition: experiment.h:325
likely
#define likely(x)
Tell the branch predictor that this expression is likely true.
Definition: macros.h:14
smash::ExperimentBase::ExperimentBase
ExperimentBase()=default
smash::Experiment::initial_mean_field_energy_
double initial_mean_field_energy_
The initial total mean field energy in the system.
Definition: experiment.h:506
smash::Experiment
The main class, where the simulation of an experiment is executed.
Definition: experiment.h:138
smash::pdg::n
constexpr int n
Neutron.
Definition: pdgcode_constants.h:30
smash::Experiment::outputs_
OutputsList outputs_
A list of output formaters.
Definition: experiment.h:337
smash::pot_pointer
Potentials * pot_pointer
Pointer to a Potential class.
Definition: potential_globals.cc:17
scatteractionsfinder.h
smash::Action::incoming_particles
const ParticleList & incoming_particles() const
Get the list of particles that go into the action.
Definition: action.cc:59
smash::Experiment::total_hypersurface_crossing_actions_
uint64_t total_hypersurface_crossing_actions_
Total number of particles removed from the evolution in hypersurface crossing actions.
Definition: experiment.h:548
smash::DensityType::Baryon
potential_globals.h
thermodynamicoutput.h
smash::Experiment::printout_v_landau_
bool printout_v_landau_
Whether to print the 4-velocity in Landau fram.
Definition: experiment.h:423
smash::Experiment::time_start_
SystemTimePoint time_start_
system starting time of the simulation
Definition: experiment.h:509
smash::ExperimentBase::create
static std::unique_ptr< ExperimentBase > create(Configuration config, const bf::path &output_path)
Factory method that creates and initializes a new Experiment<Modus>.
smash::Experiment::UI3_lat_
std::unique_ptr< RectangularLattice< FourVector > > UI3_lat_
Lattices for symmetry potentials (evaluated in the local rest frame) times the isospin flow 4-velocit...
Definition: experiment.h:403
smash::Experiment::bremsstrahlung_switch_
const bool bremsstrahlung_switch_
This indicates whether bremsstrahlung is switched on.
Definition: experiment.h:482
TimeStepMode::Fixed
Use fixed time step.
smash::QuantumNumbers
Definition: quantumnumbers.h:53
icoutput.h
smash::ExpansionProperties
Struct containing the type of the metric and the expansion parameter of the metric.
Definition: propagation.h:25
FermiMotion
FermiMotion
Option to use Fermi Motion.
Definition: forwarddeclarations.h:93
smash::ExperimentBase::run
virtual void run()=0
Runs the experiment.
smash::Actions
Definition: actions.h:29
smash::update_momenta
void update_momenta(Particles *particles, double dt, const Potentials &pot, RectangularLattice< std::pair< ThreeVector, ThreeVector >> *FB_lat, RectangularLattice< std::pair< ThreeVector, ThreeVector >> *FI3_lat)
Updates the momenta of all particles at the current time step according to the equations of motion:
Definition: propagation.cc:111
smash::Experiment::photon_finder_
std::unique_ptr< ActionFinderInterface > photon_finder_
The (Scatter) Actions Finder for Direct Photons.
Definition: experiment.h:370
smash::Experiment::conserved_initial_
QuantumNumbers conserved_initial_
The conserved quantities of the system.
Definition: experiment.h:500
smash::ExperimentBase
Non-template interface to Experiment<Modus>.
Definition: experiment.h:87
smash::Experiment::max_transverse_distance_sqr_
double max_transverse_distance_sqr_
Maximal distance at which particles can interact, squared.
Definition: experiment.h:491
potentials.h