Version: SMASH-2.0
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_SMASH_EXPERIMENT_H_
8 #define SRC_INCLUDE_SMASH_EXPERIMENT_H_
9 
10 #include <algorithm>
11 #include <limits>
12 #include <memory>
13 #include <string>
14 #include <utility>
15 #include <vector>
16 
17 #include "actionfinderfactory.h"
18 #include "actions.h"
19 #include "bremsstrahlungaction.h"
20 #include "chrono.h"
21 #include "decayactionsfinder.h"
23 #include "energymomentumtensor.h"
24 #include "fourvector.h"
25 #include "grandcan_thermalizer.h"
26 #include "grid.h"
28 #include "outputparameters.h"
29 #include "pauliblocking.h"
30 #include "potential_globals.h"
31 #include "potentials.h"
32 #include "propagation.h"
33 #include "quantumnumbers.h"
34 #include "scatteractionphoton.h"
35 #include "scatteractionsfinder.h"
36 #include "stringprocess.h"
37 #include "thermalizationaction.h"
38 // Output
39 #include "binaryoutput.h"
40 #ifdef SMASH_USE_HEPMC
41 #include "hepmcoutput.h"
42 #endif
43 #include "icoutput.h"
44 #include "oscaroutput.h"
45 #include "thermodynamicoutput.h"
46 #ifdef SMASH_USE_ROOT
47 #include "rootoutput.h"
48 #endif
49 #include "vtkoutput.h"
50 #include "wallcrossingaction.h"
51 
52 namespace std {
63 template <typename T, typename Ratio>
64 static ostream &operator<<(ostream &out,
65  const chrono::duration<T, Ratio> &seconds) {
66  using Seconds = chrono::duration<double>;
67  using Minutes = chrono::duration<double, std::ratio<60>>;
68  using Hours = chrono::duration<double, std::ratio<60 * 60>>;
69  constexpr Minutes threshold_for_minutes{10};
70  constexpr Hours threshold_for_hours{3};
71  if (seconds < threshold_for_minutes) {
72  return out << Seconds(seconds).count() << " [s]";
73  }
74  if (seconds < threshold_for_hours) {
75  return out << Minutes(seconds).count() << " [min]";
76  }
77  return out << Hours(seconds).count() << " [h]";
78 }
79 } // namespace std
80 
81 namespace smash {
82 static constexpr int LMain = LogArea::Main::id;
83 static constexpr int LInitialConditions = LogArea::InitialConditions::id;
84 
93  public:
94  ExperimentBase() = default;
99  virtual ~ExperimentBase() = default;
100 
121  static std::unique_ptr<ExperimentBase> create(Configuration config,
122  const bf::path &output_path);
123 
130  virtual void run() = 0;
131 
137  struct InvalidModusRequest : public std::invalid_argument {
138  using std::invalid_argument::invalid_argument;
139  };
140 };
141 
142 template <typename Modus>
144 template <typename Modus>
145 std::ostream &operator<<(std::ostream &out, const Experiment<Modus> &e);
146 
171 template <typename Modus>
172 class Experiment : public ExperimentBase {
173  friend class ExperimentBase;
174 
175  public:
182  void run() override;
183 
197  explicit Experiment(Configuration config, const bf::path &output_path);
198 
204  void initialize_new_event(int event_number);
205 
212  void run_time_evolution();
213 
215  void do_final_decays();
216 
222  void final_output(const int evt_num);
223 
229 
234  Modus *modus() { return &modus_; }
235 
236  private:
248  template <typename Container>
249  bool perform_action(Action &action,
250  const Container &particles_before_actions);
259  void create_output(const std::string &format, const std::string &content,
260  const bf::path &output_path, const OutputParameters &par);
261 
268  void propagate_and_shine(double to_time);
269 
283 
285  void intermediate_output();
286 
288  void update_potentials();
289 
298  double compute_min_cell_length(double dt) const {
299  return std::sqrt(4 * dt * dt + max_transverse_distance_sqr_);
300  }
301 
303  double next_output_time() const {
304  return parameters_.outputclock->next_time();
305  }
306 
313 
316 
321  Modus modus_;
322 
325 
330  std::unique_ptr<Potentials> potentials_;
331 
336  std::unique_ptr<PauliBlocker> pauli_blocker_;
337 
342  OutputsList outputs_;
343 
345  OutputPtr dilepton_output_;
346 
348  OutputPtr photon_output_;
349 
355  std::vector<bool> nucleon_has_interacted_ = {};
360 
366  std::vector<FourVector> beam_momentum_ = {};
367 
369  std::vector<std::unique_ptr<ActionFinderInterface>> action_finders_;
370 
372  std::unique_ptr<DecayActionsFinderDilepton> dilepton_finder_;
373 
375  std::unique_ptr<ActionFinderInterface> photon_finder_;
376 
379 
381  std::unique_ptr<DensityLattice> jmu_B_lat_;
382 
384  std::unique_ptr<DensityLattice> jmu_I3_lat_;
385 
393  std::unique_ptr<DensityLattice> jmu_custom_lat_;
394 
397 
402  std::unique_ptr<RectangularLattice<FourVector>> UB_lat_ = nullptr;
403 
408  std::unique_ptr<RectangularLattice<FourVector>> UI3_lat_ = nullptr;
409 
411  std::unique_ptr<RectangularLattice<std::pair<ThreeVector, ThreeVector>>>
413 
415  std::unique_ptr<RectangularLattice<std::pair<ThreeVector, ThreeVector>>>
417 
419  std::unique_ptr<RectangularLattice<EnergyMomentumTensor>> Tmn_;
420 
422  bool printout_tmn_ = false;
423 
425  bool printout_tmn_landau_ = false;
426 
428  bool printout_v_landau_ = false;
429 
431  bool printout_lattice_td_ = false;
432 
434  std::unique_ptr<GrandCanThermalizer> thermalizer_;
435 
441 
456  const int nevents_;
457 
459  const double end_time_;
460 
466  const double delta_time_startup_;
467 
472  const bool force_decays_;
473 
475  const bool use_grid_;
476 
479 
481  const bool dileptons_switch_;
482 
484  const bool photons_switch_;
485 
488 
490  const bool IC_output_switch_;
491 
494 
496  double max_transverse_distance_sqr_ = std::numeric_limits<double>::max();
497 
506 
512 
514  SystemTimePoint time_start_ = SystemClock::now();
515 
518 
523  uint64_t interactions_total_ = 0;
524 
530 
535  uint64_t wall_actions_total_ = 0;
536 
542 
547  uint64_t total_pauli_blocked_ = 0;
548 
554 
560 
565  double total_energy_removed_ = 0.0;
566 
568  int64_t seed_ = -1;
569 
575  friend std::ostream &operator<<<>(std::ostream &out, const Experiment &e);
576 };
577 
579 template <typename Modus>
580 std::ostream &operator<<(std::ostream &out, const Experiment<Modus> &e) {
581  out << "End time: " << e.end_time_ << " fm/c\n";
582  out << e.modus_;
583  return out;
584 }
585 
586 template <typename Modus>
587 void Experiment<Modus>::create_output(const std::string &format,
588  const std::string &content,
589  const bf::path &output_path,
590  const OutputParameters &out_par) {
591  logg[LExperiment].info() << "Adding output " << content << " of format "
592  << format << std::endl;
593 
594  if (format == "VTK" && content == "Particles") {
595  outputs_.emplace_back(
596  make_unique<VtkOutput>(output_path, content, out_par));
597  } else if (format == "Root") {
598 #ifdef SMASH_USE_ROOT
599  if (content == "Initial_Conditions") {
600  outputs_.emplace_back(
601  make_unique<RootOutput>(output_path, "SMASH_IC", out_par));
602  } else {
603  outputs_.emplace_back(
604  make_unique<RootOutput>(output_path, content, out_par));
605  }
606 #else
607  logg[LExperiment].error(
608  "Root output requested, but Root support not compiled in");
609 #endif
610  } else if (format == "Binary") {
611  if (content == "Collisions" || content == "Dileptons" ||
612  content == "Photons") {
613  outputs_.emplace_back(
614  make_unique<BinaryOutputCollisions>(output_path, content, out_par));
615  } else if (content == "Particles") {
616  outputs_.emplace_back(
617  make_unique<BinaryOutputParticles>(output_path, content, out_par));
618  } else if (content == "Initial_Conditions") {
619  outputs_.emplace_back(make_unique<BinaryOutputInitialConditions>(
620  output_path, content, out_par));
621  }
622  } else if (format == "Oscar1999" || format == "Oscar2013") {
623  outputs_.emplace_back(
624  create_oscar_output(format, content, output_path, out_par));
625  } else if (content == "Thermodynamics" && format == "ASCII") {
626  outputs_.emplace_back(
627  make_unique<ThermodynamicOutput>(output_path, content, out_par));
628  } else if (content == "Thermodynamics" && format == "VTK") {
629  printout_lattice_td_ = true;
630  outputs_.emplace_back(
631  make_unique<VtkOutput>(output_path, content, out_par));
632  } else if (content == "Initial_Conditions" && format == "ASCII") {
633  outputs_.emplace_back(
634  make_unique<ICOutput>(output_path, "SMASH_IC", out_par));
635  } else if (content == "HepMC" && format == "ASCII") {
636 #ifdef SMASH_USE_HEPMC
637  outputs_.emplace_back(make_unique<HepMcOutput>(
638  output_path, "SMASH_HepMC", out_par, modus_.total_N_number(),
639  modus_.proj_N_number()));
640 #else
641  logg[LExperiment].error(
642  "HepMC output requested, but HepMC support not compiled in");
643 #endif
644  } else {
645  logg[LExperiment].error()
646  << "Unknown combination of format (" << format << ") and content ("
647  << content << "). Fix the config.";
648  }
649 }
650 
659 
821 template <typename Modus>
822 Experiment<Modus>::Experiment(Configuration config, const bf::path &output_path)
823  : parameters_(create_experiment_parameters(config)),
824  density_param_(DensityParameters(parameters_)),
825  modus_(config["Modi"], parameters_),
826  particles_(),
827  nevents_(config.take({"General", "Nevents"})),
828  end_time_(config.take({"General", "End_Time"})),
829  delta_time_startup_(parameters_.labclock->timestep_duration()),
830  force_decays_(
831  config.take({"Collision_Term", "Force_Decays_At_End"}, true)),
832  use_grid_(config.take({"General", "Use_Grid"}, true)),
833  metric_(
834  config.take({"General", "Metric_Type"}, ExpansionMode::NoExpansion),
835  config.take({"General", "Expansion_Rate"}, 0.1)),
836  dileptons_switch_(
837  config.take({"Collision_Term", "Dileptons", "Decays"}, false)),
838  photons_switch_(config.take(
839  {"Collision_Term", "Photons", "2to2_Scatterings"}, false)),
840  bremsstrahlung_switch_(
841  config.take({"Collision_Term", "Photons", "Bremsstrahlung"}, false)),
842  IC_output_switch_(config.has_value({"Output", "Initial_Conditions"})),
843  time_step_mode_(
844  config.take({"General", "Time_Step_Mode"}, TimeStepMode::Fixed)) {
845  logg[LExperiment].info() << *this;
846 
847  if (parameters_.coll_crit == CollisionCriterion::Stochastic &&
848  time_step_mode_ != TimeStepMode::Fixed) {
849  throw std::invalid_argument(
850  "The stochastic criterion can only be employed for fixed time step "
851  "mode!");
852  }
853 
854  // create finders
855  if (dileptons_switch_) {
856  dilepton_finder_ = make_unique<DecayActionsFinderDilepton>();
857  }
858  if (photons_switch_ || bremsstrahlung_switch_) {
859  n_fractional_photons_ =
860  config.take({"Collision_Term", "Photons", "Fractional_Photons"}, 100);
861  }
862  if (parameters_.two_to_one) {
863  if (parameters_.res_lifetime_factor < 0.) {
864  throw std::invalid_argument(
865  "Resonance lifetime modifier cannot be negative!");
866  }
867  if (parameters_.res_lifetime_factor < really_small) {
868  logg[LExperiment].warn(
869  "Resonance lifetime set to zero. Make sure resonances cannot "
870  "interact",
871  "inelastically (e.g. resonance chains), else SMASH is known to "
872  "hang.");
873  }
874  action_finders_.emplace_back(
875  make_unique<DecayActionsFinder>(parameters_.res_lifetime_factor));
876  }
877  bool no_coll = config.take({"Collision_Term", "No_Collisions"}, false);
878  if ((parameters_.two_to_one || parameters_.included_2to2.any() ||
879  parameters_.included_multi.any() || parameters_.strings_switch) &&
880  !no_coll) {
881  auto scat_finder = make_unique<ScatterActionsFinder>(
882  config, parameters_, nucleon_has_interacted_, modus_.total_N_number(),
883  modus_.proj_N_number());
884  max_transverse_distance_sqr_ =
885  scat_finder->max_transverse_distance_sqr(parameters_.testparticles);
886  process_string_ptr_ = scat_finder->get_process_string_ptr();
887  action_finders_.emplace_back(std::move(scat_finder));
888  } else {
889  max_transverse_distance_sqr_ =
890  parameters_.maximum_cross_section / M_PI * fm2_mb;
891  process_string_ptr_ = NULL;
892  }
893  if (modus_.is_box()) {
894  action_finders_.emplace_back(
895  make_unique<WallCrossActionsFinder>(parameters_.box_length));
896  }
897  if (IC_output_switch_) {
898  if (!modus_.is_collider()) {
899  throw std::runtime_error(
900  "Initial conditions can only be extracted in collider modus.");
901  }
902  double proper_time;
903  if (config.has_value({"Output", "Initial_Conditions", "Proper_Time"})) {
904  // Read in proper time from config
905  proper_time =
906  config.take({"Output", "Initial_Conditions", "Proper_Time"});
907  } else {
908  // Default proper time is the passing time of the two nuclei
909  double default_proper_time = modus_.nuclei_passing_time();
910  double lower_bound =
911  config.take({"Output", "Initial_Conditions", "Lower_Bound"}, 0.5);
912  if (default_proper_time >= lower_bound) {
913  proper_time = default_proper_time;
914  } else {
915  logg[LInitialConditions].warn()
916  << "Nuclei passing time is too short, hypersurface proper time set "
917  "to tau = "
918  << lower_bound << " fm.";
919  proper_time = lower_bound;
920  }
921  }
922  action_finders_.emplace_back(
923  make_unique<HyperSurfaceCrossActionsFinder>(proper_time));
924  }
925 
926  if (config.has_value({"Collision_Term", "Pauli_Blocking"})) {
927  logg[LExperiment].info() << "Pauli blocking is ON.";
928  pauli_blocker_ = make_unique<PauliBlocker>(
929  config["Collision_Term"]["Pauli_Blocking"], parameters_);
930  }
931  // In collider setup with sqrts >= 200 GeV particles don't form continuously
933  config.take({"Collision_Term", "Power_Particle_Formation"},
934  modus_.sqrt_s_NN() >= 200. ? -1. : 1.);
935 
968  // create outputs
969  logg[LExperiment].trace(source_location, " create OutputInterface objects");
970 
971  auto output_conf = config["Output"];
1196  dens_type_ = config.take({"Output", "Density_Type"}, DensityType::None);
1197  logg[LExperiment].debug()
1198  << "Density type printed to headers: " << dens_type_;
1199 
1200  const OutputParameters output_parameters(std::move(output_conf));
1201 
1202  std::vector<std::string> output_contents = output_conf.list_upmost_nodes();
1203  for (const auto &content : output_contents) {
1204  auto this_output_conf = output_conf[content.c_str()];
1205  const std::vector<std::string> formats = this_output_conf.take({"Format"});
1206  if (output_path == "") {
1207  continue;
1208  }
1209  for (const auto &format : formats) {
1210  create_output(format, content, output_path, output_parameters);
1211  }
1212  }
1213 
1214  /* We can take away the Fermi motion flag, because the collider modus is
1215  * already initialized. We only need it when potentials are enabled, but we
1216  * always have to take it, otherwise SMASH will complain about unused
1217  * options. We have to provide a default value for modi other than Collider.
1218  */
1219  const FermiMotion motion =
1220  config.take({"Modi", "Collider", "Fermi_Motion"}, FermiMotion::Off);
1221  if (config.has_value({"Potentials"})) {
1222  if (time_step_mode_ == TimeStepMode::None) {
1223  logg[LExperiment].error() << "Potentials only work with time steps!";
1224  throw std::invalid_argument("Can't use potentials without time steps!");
1225  }
1226  if (motion == FermiMotion::Frozen) {
1227  logg[LExperiment].error()
1228  << "Potentials don't work with frozen Fermi momenta! "
1229  "Use normal Fermi motion instead.";
1230  throw std::invalid_argument(
1231  "Can't use potentials "
1232  "with frozen Fermi momenta!");
1233  }
1234  logg[LExperiment].info() << "Potentials are ON. Timestep is "
1235  << parameters_.labclock->timestep_duration();
1236  // potentials need testparticles and gaussian sigma from parameters_
1237  potentials_ = make_unique<Potentials>(config["Potentials"], parameters_);
1238  }
1239 
1295  // Create lattices
1296  if (config.has_value({"Lattice"})) {
1297  // Take lattice properties from config to assign them to all lattices
1298  const std::array<double, 3> l = config.take({"Lattice", "Sizes"});
1299  const std::array<int, 3> n = config.take({"Lattice", "Cell_Number"});
1300  const std::array<double, 3> origin = config.take({"Lattice", "Origin"});
1301  const bool periodic = config.take({"Lattice", "Periodic"});
1302 
1303  if (printout_lattice_td_) {
1304  dens_type_lattice_printout_ = output_parameters.td_dens_type;
1305  printout_tmn_ = output_parameters.td_tmn;
1306  printout_tmn_landau_ = output_parameters.td_tmn_landau;
1307  printout_v_landau_ = output_parameters.td_v_landau;
1308  }
1309  if (printout_tmn_ || printout_tmn_landau_ || printout_v_landau_) {
1310  Tmn_ = make_unique<RectangularLattice<EnergyMomentumTensor>>(
1311  l, n, origin, periodic, LatticeUpdate::AtOutput);
1312  }
1313  /* Create baryon and isospin density lattices regardless of config
1314  if potentials are on. This is because they allow to compute
1315  potentials faster */
1316  if (potentials_) {
1317  if (potentials_->use_skyrme()) {
1318  jmu_B_lat_ = make_unique<DensityLattice>(l, n, origin, periodic,
1320  UB_lat_ = make_unique<RectangularLattice<FourVector>>(
1321  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1322  FB_lat_ = make_unique<
1323  RectangularLattice<std::pair<ThreeVector, ThreeVector>>>(
1324  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1325  }
1326  if (potentials_->use_symmetry()) {
1327  jmu_I3_lat_ = make_unique<DensityLattice>(l, n, origin, periodic,
1329  UI3_lat_ = make_unique<RectangularLattice<FourVector>>(
1330  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1331  FI3_lat_ = make_unique<
1332  RectangularLattice<std::pair<ThreeVector, ThreeVector>>>(
1333  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1334  }
1335  } else {
1336  if (dens_type_lattice_printout_ == DensityType::Baryon) {
1337  jmu_B_lat_ = make_unique<DensityLattice>(l, n, origin, periodic,
1339  }
1340  if (dens_type_lattice_printout_ == DensityType::BaryonicIsospin) {
1341  jmu_I3_lat_ = make_unique<DensityLattice>(l, n, origin, periodic,
1343  }
1344  }
1345  if (dens_type_lattice_printout_ != DensityType::None &&
1346  dens_type_lattice_printout_ != DensityType::BaryonicIsospin &&
1347  dens_type_lattice_printout_ != DensityType::Baryon) {
1348  jmu_custom_lat_ = make_unique<DensityLattice>(l, n, origin, periodic,
1350  }
1351  } else if (printout_lattice_td_) {
1352  logg[LExperiment].error(
1353  "If you want Thermodynamic VTK output, configure a lattice for it.");
1354  }
1355  // Warning for the mean field calculation if lattice is not on.
1356  if ((potentials_ != nullptr) && (jmu_B_lat_ == nullptr)) {
1357  logg[LExperiment].warn() << "Lattice is NOT used. Mean fields are "
1358  << "not going to be calculated.";
1359  }
1360 
1361  // Store pointers to potential and lattice accessible for Action
1362  if (parameters_.potential_affect_threshold) {
1363  UB_lat_pointer = UB_lat_.get();
1364  UI3_lat_pointer = UI3_lat_.get();
1365  pot_pointer = potentials_.get();
1366  }
1367 
1368  // Create forced thermalizer
1369  if (config.has_value({"Forced_Thermalization"})) {
1370  Configuration &&th_conf = config["Forced_Thermalization"];
1371  thermalizer_ = modus_.create_grandcan_thermalizer(th_conf);
1372  }
1373 
1374  /* Take the seed setting only after the configuration was stored to a file
1375  * in smash.cc */
1376  seed_ = config.take({"General", "Randomseed"});
1377 }
1378 
1380 const std::string hline(113, '-');
1381 
1407 std::string format_measurements(const Particles &particles,
1408  uint64_t scatterings_this_interval,
1409  const QuantumNumbers &conserved_initial,
1410  SystemTimePoint time_start, double time,
1411  double E_mean_field,
1412  double E_mean_field_initial);
1426  const Potentials &potentials,
1427  RectangularLattice<smash::DensityOnLattice> &jmu_B_lat,
1428  const ExperimentParameters &parameters);
1429 
1445 EventInfo fill_event_info(const Particles &particles, double E_mean_field,
1446  double modus_impact_parameter,
1447  const ExperimentParameters &parameters,
1448  bool projectile_target_interact);
1449 
1450 template <typename Modus>
1452  random::set_seed(seed_);
1453  logg[LExperiment].info() << "random number seed: " << seed_;
1454  /* Set seed for the next event. It has to be positive, so it can be entered
1455  * in the config.
1456  *
1457  * We have to be careful about the minimal integer, whose absolute value
1458  * cannot be represented. */
1459  int64_t r = random::advance();
1460  while (r == INT64_MIN) {
1461  r = random::advance();
1462  }
1463  seed_ = std::abs(r);
1464  /* Set the random seed used in PYTHIA hadronization
1465  * to be same with the SMASH one.
1466  * In this way we ensure that the results are reproducible
1467  * for every event if one knows SMASH random seed. */
1468  if (process_string_ptr_ != NULL) {
1469  process_string_ptr_->init_pythia_hadron_rndm();
1470  }
1471 
1472  particles_.reset();
1473 
1474  // Sample particles according to the initial conditions
1475  double start_time = modus_.initial_conditions(&particles_, parameters_);
1476  /* For box modus make sure that particles are in the box. In principle, after
1477  * a correct initialization they should be, so this is just playing it safe.
1478  */
1479  modus_.impose_boundary_conditions(&particles_, outputs_);
1480  // Reset the simulation clock
1481  double timestep = delta_time_startup_;
1482 
1483  switch (time_step_mode_) {
1484  case TimeStepMode::Fixed:
1485  break;
1486  case TimeStepMode::None:
1487  timestep = end_time_ - start_time;
1488  // Take care of the box modus + timestepless propagation
1489  const double max_dt = modus_.max_timestep(max_transverse_distance_sqr_);
1490  if (max_dt > 0. && max_dt < timestep) {
1491  timestep = max_dt;
1492  }
1493  break;
1494  }
1495  std::unique_ptr<UniformClock> clock_for_this_event;
1496  if (modus_.is_list() && (timestep < 0.0)) {
1497  throw std::runtime_error(
1498  "Timestep for the given event is negative. \n"
1499  "This might happen if the formation times of the input particles are "
1500  "larger than the specified end time of the simulation.");
1501  }
1502  clock_for_this_event = make_unique<UniformClock>(start_time, timestep);
1503  parameters_.labclock = std::move(clock_for_this_event);
1504 
1505  // Reset the output clock
1506  parameters_.outputclock->reset(start_time, true);
1507  // remove time before starting time in case of custom output times.
1508  parameters_.outputclock->remove_times_in_past(start_time);
1509 
1510  logg[LExperiment].debug(
1511  "Lab clock: t_start = ", parameters_.labclock->current_time(),
1512  ", dt = ", parameters_.labclock->timestep_duration());
1513 
1514  /* Save the initial conserved quantum numbers and total momentum in
1515  * the system for conservation checks */
1516  conserved_initial_ = QuantumNumbers(particles_);
1517  wall_actions_total_ = 0;
1518  previous_wall_actions_total_ = 0;
1519  interactions_total_ = 0;
1520  previous_interactions_total_ = 0;
1521  discarded_interactions_total_ = 0;
1522  total_pauli_blocked_ = 0;
1523  projectile_target_interact_ = false;
1524  total_hypersurface_crossing_actions_ = 0;
1525  total_energy_removed_ = 0.0;
1526  // Print output headers
1527  logg[LExperiment].info() << hline;
1528  logg[LExperiment].info() << "Time[fm] Ekin[GeV] E_MF[GeV] ETotal[GeV] "
1529  << "ETot/N[GeV] D(ETot/N)[GeV] Scatt&Decays "
1530  << "Particles Comp.Time";
1531  logg[LExperiment].info() << hline;
1532  double E_mean_field = 0.0;
1533  if (potentials_) {
1534  update_potentials();
1535  // using the lattice is necessary
1536  if ((jmu_B_lat_ != nullptr)) {
1537  E_mean_field =
1538  calculate_mean_field_energy(*potentials_, *jmu_B_lat_, parameters_);
1539  }
1540  }
1541  initial_mean_field_energy_ = E_mean_field;
1543  particles_, 0u, conserved_initial_, time_start_,
1544  parameters_.labclock->current_time(), E_mean_field,
1545  initial_mean_field_energy_);
1546 
1547  auto event_info =
1548  fill_event_info(particles_, E_mean_field, modus_.impact_parameter(),
1549  parameters_, projectile_target_interact_);
1550 
1551  // Output at event start
1552  for (const auto &output : outputs_) {
1553  output->at_eventstart(particles_, event_number, event_info);
1554  }
1555 }
1556 
1557 template <typename Modus>
1558 template <typename Container>
1560  Action &action, const Container &particles_before_actions) {
1561  // Make sure to skip invalid and Pauli-blocked actions.
1562  if (!action.is_valid(particles_)) {
1563  discarded_interactions_total_++;
1564  logg[LExperiment].debug(~einhard::DRed(), "✘ ", action,
1565  " (discarded: invalid)");
1566  return false;
1567  }
1568  action.generate_final_state();
1569  logg[LExperiment].debug("Process Type is: ", action.get_type());
1570  if (pauli_blocker_ && action.is_pauli_blocked(particles_, *pauli_blocker_)) {
1571  total_pauli_blocked_++;
1572  return false;
1573  }
1574  if (modus_.is_collider()) {
1575  /* Mark incoming nucleons as interacted - now they are permitted
1576  * to collide with nucleons from their native nucleus */
1577  bool incoming_projectile = false;
1578  bool incoming_target = false;
1579  for (const auto &incoming : action.incoming_particles()) {
1580  assert(incoming.id() >= 0);
1581  if (incoming.id() < modus_.total_N_number()) {
1582  nucleon_has_interacted_[incoming.id()] = true;
1583  }
1584  if (incoming.id() < modus_.proj_N_number()) {
1585  incoming_projectile = true;
1586  }
1587  if (incoming.id() >= modus_.proj_N_number() &&
1588  incoming.id() < modus_.total_N_number()) {
1589  incoming_target = true;
1590  }
1591  }
1592  // Check whether particles from different nuclei interacted.
1593  if (incoming_projectile & incoming_target) {
1594  projectile_target_interact_ = true;
1595  }
1596  }
1597  /* Make sure to pick a non-zero integer, because 0 is reserved for "no
1598  * interaction yet". */
1599  const auto id_process = static_cast<uint32_t>(interactions_total_ + 1);
1600  action.perform(&particles_, id_process);
1601  interactions_total_++;
1602  if (action.get_type() == ProcessType::Wall) {
1603  wall_actions_total_++;
1604  }
1605  if (action.get_type() == ProcessType::HyperSurfaceCrossing) {
1606  total_hypersurface_crossing_actions_++;
1607  total_energy_removed_ += action.incoming_particles()[0].momentum().x0();
1608  }
1609  // Calculate Eckart rest frame density at the interaction point
1610  double rho = 0.0;
1611  if (dens_type_ != DensityType::None) {
1612  const FourVector r_interaction = action.get_interaction_point();
1613  constexpr bool compute_grad = false;
1614  const bool smearing = true;
1615  rho = std::get<0>(current_eckart(r_interaction.threevec(),
1616  particles_before_actions, density_param_,
1617  dens_type_, compute_grad, smearing));
1618  }
1634  for (const auto &output : outputs_) {
1635  if (!output->is_dilepton_output() && !output->is_photon_output()) {
1636  if (output->is_IC_output() &&
1638  output->at_interaction(action, rho);
1639  } else if (!output->is_IC_output()) {
1640  output->at_interaction(action, rho);
1641  }
1642  }
1643  }
1644 
1645  // At every collision photons can be produced.
1646  // Note: We rely here on the lazy evaluation of the arguments to if.
1647  // It may happen that in a wall-crossing-action sqrt_s raises an exception.
1648  // Therefore we first have to check if the incoming particles can undergo
1649  // an em-interaction.
1650  if (photons_switch_ &&
1653  action.sqrt_s(), action.incoming_particles())) {
1654  /* Time in the action constructor is relative to
1655  * current time of incoming */
1656  constexpr double action_time = 0.;
1657  ScatterActionPhoton photon_act(action.incoming_particles(), action_time,
1658  n_fractional_photons_,
1659  action.get_total_weight());
1660 
1671  photon_act.add_dummy_hadronic_process(action.get_total_weight());
1672 
1673  // Now add the actual photon reaction channel.
1674  photon_act.add_single_process();
1675 
1676  photon_act.perform_photons(outputs_);
1677  }
1678 
1679  if (bremsstrahlung_switch_ &&
1681  action.incoming_particles())) {
1682  /* Time in the action constructor is relative to
1683  * current time of incoming */
1684  constexpr double action_time = 0.;
1685 
1686  BremsstrahlungAction brems_act(action.incoming_particles(), action_time,
1687  n_fractional_photons_,
1688  action.get_total_weight());
1689 
1701  brems_act.add_dummy_hadronic_process(action.get_total_weight());
1702 
1703  // Now add the actual bremsstrahlung reaction channel.
1704  brems_act.add_single_process();
1705 
1706  brems_act.perform_bremsstrahlung(outputs_);
1707  }
1708 
1709  logg[LExperiment].debug(~einhard::Green(), "✔ ", action);
1710  return true;
1711 }
1712 
1713 template <typename Modus>
1715  Actions actions;
1716 
1717  while (parameters_.labclock->current_time() < end_time_) {
1718  const double t = parameters_.labclock->current_time();
1719  const double dt =
1720  std::min(parameters_.labclock->timestep_duration(), end_time_ - t);
1721  logg[LExperiment].debug("Timestepless propagation for next ", dt, " fm/c.");
1722 
1723  // Perform forced thermalization if required
1724  if (thermalizer_ &&
1725  thermalizer_->is_time_to_thermalize(parameters_.labclock)) {
1726  const bool ignore_cells_under_treshold = true;
1727  thermalizer_->update_thermalizer_lattice(particles_, density_param_,
1728  ignore_cells_under_treshold);
1729  const double current_t = parameters_.labclock->current_time();
1730  thermalizer_->thermalize(particles_, current_t,
1731  parameters_.testparticles);
1732  ThermalizationAction th_act(*thermalizer_, current_t);
1733  if (th_act.any_particles_thermalized()) {
1734  perform_action(th_act, particles_);
1735  }
1736  }
1737 
1738  if (particles_.size() > 0 && action_finders_.size() > 0) {
1739  /* (1.a) Create grid. */
1740  double min_cell_length = compute_min_cell_length(dt);
1741  logg[LExperiment].debug("Creating grid with minimal cell length ",
1742  min_cell_length);
1743  const auto &grid =
1744  use_grid_ ? modus_.create_grid(particles_, min_cell_length, dt)
1745  : modus_.create_grid(particles_, min_cell_length, dt,
1747 
1748  const double gcell_vol = grid.cell_volume();
1749 
1750  /* (1.b) Iterate over cells and find actions. */
1751  grid.iterate_cells(
1752  [&](const ParticleList &search_list) {
1753  for (const auto &finder : action_finders_) {
1754  actions.insert(finder->find_actions_in_cell(
1755  search_list, dt, gcell_vol, beam_momentum_));
1756  }
1757  },
1758  [&](const ParticleList &search_list,
1759  const ParticleList &neighbors_list) {
1760  for (const auto &finder : action_finders_) {
1761  actions.insert(finder->find_actions_with_neighbors(
1762  search_list, neighbors_list, dt, beam_momentum_));
1763  }
1764  });
1765  }
1766 
1767  /* \todo (optimizations) Adapt timestep size here */
1768 
1769  /* (2) Propagation from action to action until the end of timestep */
1770  run_time_evolution_timestepless(actions);
1771 
1772  /* (3) Update potentials (if computed on the lattice) and
1773  * compute new momenta according to equations of motion */
1774  if (potentials_) {
1775  update_potentials();
1776  update_momenta(&particles_, parameters_.labclock->timestep_duration(),
1777  *potentials_, FB_lat_.get(), FI3_lat_.get());
1778  }
1779 
1780  /* (4) Expand universe if non-minkowskian metric; updates
1781  * positions and momenta according to the selected expansion */
1782  if (metric_.mode_ != ExpansionMode::NoExpansion) {
1783  expand_space_time(&particles_, parameters_, metric_);
1784  }
1785 
1786  ++(*parameters_.labclock);
1787 
1788  /* (5) Check conservation laws.
1789  *
1790  * Check conservation of conserved quantities if potentials and string
1791  * fragmentation are off. If potentials are on then momentum is conserved
1792  * only in average. If string fragmentation is on, then energy and
1793  * momentum are only very roughly conserved in high-energy collisions. */
1794  if (!potentials_ && !parameters_.strings_switch &&
1795  metric_.mode_ == ExpansionMode::NoExpansion && !IC_output_switch_) {
1796  std::string err_msg = conserved_initial_.report_deviations(particles_);
1797  if (!err_msg.empty()) {
1798  logg[LExperiment].error() << err_msg;
1799  throw std::runtime_error("Violation of conserved quantities!");
1800  }
1801  }
1802  }
1803 
1804  if (pauli_blocker_) {
1805  logg[LExperiment].info(
1806  "Interactions: Pauli-blocked/performed = ", total_pauli_blocked_, "/",
1807  interactions_total_ - wall_actions_total_);
1808  }
1809 }
1810 
1811 template <typename Modus>
1813  const double dt =
1814  propagate_straight_line(&particles_, to_time, beam_momentum_);
1815  if (dilepton_finder_ != nullptr) {
1816  for (const auto &output : outputs_) {
1817  dilepton_finder_->shine(particles_, output.get(), dt);
1818  }
1819  }
1820 }
1821 
1829 inline void check_interactions_total(uint64_t interactions_total) {
1830  constexpr uint64_t max_uint32 = std::numeric_limits<uint32_t>::max();
1831  if (interactions_total >= max_uint32) {
1832  throw std::runtime_error("Integer overflow in total interaction number!");
1833  }
1834 }
1835 
1836 template <typename Modus>
1838  const double start_time = parameters_.labclock->current_time();
1839  const double end_time =
1840  std::min(parameters_.labclock->next_time(), end_time_);
1841  double time_left = end_time - start_time;
1842  logg[LExperiment].debug(
1843  "Timestepless propagation: ", "Actions size = ", actions.size(),
1844  ", start time = ", start_time, ", end time = ", end_time);
1845 
1846  // iterate over all actions
1847  while (!actions.is_empty()) {
1848  // get next action
1849  ActionPtr act = actions.pop();
1850  if (!act->is_valid(particles_)) {
1851  discarded_interactions_total_++;
1852  logg[LExperiment].debug(~einhard::DRed(), "✘ ", act,
1853  " (discarded: invalid)");
1854  continue;
1855  }
1856  if (act->time_of_execution() > end_time) {
1857  logg[LExperiment].error(
1858  act, " scheduled later than end time: t_action[fm/c] = ",
1859  act->time_of_execution(), ", t_end[fm/c] = ", end_time);
1860  }
1861  logg[LExperiment].debug(~einhard::Green(), "✔ ", act);
1862 
1863  while (next_output_time() <= act->time_of_execution()) {
1864  logg[LExperiment].debug("Propagating until output time: ",
1865  next_output_time());
1866  propagate_and_shine(next_output_time());
1867  ++(*parameters_.outputclock);
1868  intermediate_output();
1869  }
1870 
1871  /* (1) Propagate to the next action. */
1872  logg[LExperiment].debug("Propagating until next action ", act,
1873  ", action time = ", act->time_of_execution());
1874  propagate_and_shine(act->time_of_execution());
1875 
1876  /* (2) Perform action.
1877  *
1878  * Update the positions of the incoming particles, because the information
1879  * in the action object will be outdated as the particles have been
1880  * propagated since the construction of the action. */
1881  act->update_incoming(particles_);
1882  const bool performed = perform_action(*act, particles_);
1883 
1884  /* No need to update actions for outgoing particles
1885  * if the action is not performed. */
1886  if (!performed) {
1887  continue;
1888  }
1889 
1890  /* (3) Update actions for newly-produced particles. */
1891 
1892  time_left = end_time - act->time_of_execution();
1893  const ParticleList &outgoing_particles = act->outgoing_particles();
1894  // Grid cell volume set to zero, since there is no grid
1895  const double gcell_vol = 0.0;
1896  for (const auto &finder : action_finders_) {
1897  // Outgoing particles can still decay, cross walls...
1898  actions.insert(finder->find_actions_in_cell(outgoing_particles, time_left,
1899  gcell_vol, beam_momentum_));
1900  // ... and collide with other particles.
1901  actions.insert(finder->find_actions_with_surrounding_particles(
1902  outgoing_particles, particles_, time_left, beam_momentum_));
1903  }
1904 
1905  check_interactions_total(interactions_total_);
1906  }
1907 
1908  while (next_output_time() <= end_time) {
1909  logg[LExperiment].debug("Propagating until output time: ",
1910  next_output_time());
1911  propagate_and_shine(next_output_time());
1912  ++(*parameters_.outputclock);
1913  // Avoid duplicating printout at event end time
1914  if (parameters_.outputclock->current_time() < end_time_) {
1915  intermediate_output();
1916  }
1917  }
1918  logg[LExperiment].debug("Propagating to time ", end_time);
1919  propagate_and_shine(end_time);
1920 }
1921 
1922 template <typename Modus>
1924  const uint64_t wall_actions_this_interval =
1925  wall_actions_total_ - previous_wall_actions_total_;
1926  previous_wall_actions_total_ = wall_actions_total_;
1927  const uint64_t interactions_this_interval = interactions_total_ -
1928  previous_interactions_total_ -
1929  wall_actions_this_interval;
1930  previous_interactions_total_ = interactions_total_;
1931  double E_mean_field = 0.0;
1932  if (potentials_) {
1933  // using the lattice is necessary
1934  if ((jmu_B_lat_ != nullptr)) {
1935  E_mean_field =
1936  calculate_mean_field_energy(*potentials_, *jmu_B_lat_, parameters_);
1937  /*
1938  * Mean field calculated in a box should remain approximately constant if
1939  * the system is in equilibrium, and so deviations from its original value
1940  * may signal a phase transition or other dynamical process. This
1941  * comparison only makes sense in the Box Modus, hence the condition.
1942  */
1943  if (modus_.is_box()) {
1944  double tmp = (E_mean_field - initial_mean_field_energy_) /
1945  (E_mean_field + initial_mean_field_energy_);
1946  /*
1947  * This is displayed when the system evolves away from its initial
1948  * configuration (which is when the total mean field energy in the box
1949  * deviates from its initial value).
1950  */
1951  if (std::abs(tmp) > 0.01) {
1952  logg[LExperiment].info()
1953  << "\n\n\n\t The mean field at t = "
1954  << parameters_.outputclock->current_time()
1955  << " [fm/c] differs from the mean field at t = 0:"
1956  << "\n\t\t initial_mean_field_energy_ = "
1957  << initial_mean_field_energy_ << " [GeV]"
1958  << "\n\t\t abs[(E_MF - E_MF(t=0))/(E_MF + E_MF(t=0))] = "
1959  << std::abs(tmp)
1960  << "\n\t\t E_MF/E_MF(t=0) = "
1961  << E_mean_field / initial_mean_field_energy_ << "\n\n";
1962  }
1963  }
1964  }
1965  }
1966 
1968  particles_, interactions_this_interval, conserved_initial_, time_start_,
1969  parameters_.outputclock->current_time(), E_mean_field,
1970  initial_mean_field_energy_);
1971  const LatticeUpdate lat_upd = LatticeUpdate::AtOutput;
1972 
1973  auto event_info =
1974  fill_event_info(particles_, E_mean_field, modus_.impact_parameter(),
1975  parameters_, projectile_target_interact_);
1976  // save evolution data
1977  if (!(modus_.is_box() && parameters_.outputclock->current_time() <
1978  modus_.equilibration_time())) {
1979  for (const auto &output : outputs_) {
1980  if (output->is_dilepton_output() || output->is_photon_output() ||
1981  output->is_IC_output()) {
1982  continue;
1983  }
1984 
1985  output->at_intermediate_time(particles_, parameters_.outputclock,
1986  density_param_, event_info);
1987 
1988  // Thermodynamic output on the lattice versus time
1989  switch (dens_type_lattice_printout_) {
1990  case DensityType::Baryon:
1991  update_lattice(jmu_B_lat_.get(), lat_upd, DensityType::Baryon,
1992  density_param_, particles_, false);
1993  output->thermodynamics_output(ThermodynamicQuantity::EckartDensity,
1994  DensityType::Baryon, *jmu_B_lat_);
1995  break;
1997  update_lattice(jmu_I3_lat_.get(), lat_upd,
1998  DensityType::BaryonicIsospin, density_param_,
1999  particles_, false);
2000  output->thermodynamics_output(ThermodynamicQuantity::EckartDensity,
2002  *jmu_I3_lat_);
2003  break;
2004  case DensityType::None:
2005  break;
2006  default:
2007  update_lattice(jmu_custom_lat_.get(), lat_upd,
2008  dens_type_lattice_printout_, density_param_,
2009  particles_, false);
2010  output->thermodynamics_output(ThermodynamicQuantity::EckartDensity,
2011  dens_type_lattice_printout_,
2012  *jmu_custom_lat_);
2013  }
2014  if (printout_tmn_ || printout_tmn_landau_ || printout_v_landau_) {
2015  update_lattice(Tmn_.get(), lat_upd, dens_type_lattice_printout_,
2016  density_param_, particles_);
2017  if (printout_tmn_) {
2018  output->thermodynamics_output(ThermodynamicQuantity::Tmn,
2019  dens_type_lattice_printout_, *Tmn_);
2020  }
2021  if (printout_tmn_landau_) {
2022  output->thermodynamics_output(ThermodynamicQuantity::TmnLandau,
2023  dens_type_lattice_printout_, *Tmn_);
2024  }
2025  if (printout_v_landau_) {
2026  output->thermodynamics_output(ThermodynamicQuantity::LandauVelocity,
2027  dens_type_lattice_printout_, *Tmn_);
2028  }
2029  }
2030 
2031  if (thermalizer_) {
2032  output->thermodynamics_output(*thermalizer_);
2033  }
2034  }
2035  }
2036 }
2037 
2038 template <typename Modus>
2040  if (potentials_) {
2041  if (potentials_->use_symmetry() && jmu_I3_lat_ != nullptr) {
2042  update_lattice(jmu_I3_lat_.get(), LatticeUpdate::EveryTimestep,
2043  DensityType::BaryonicIsospin, density_param_, particles_,
2044  true);
2045  }
2046  if ((potentials_->use_skyrme() || potentials_->use_symmetry()) &&
2047  jmu_B_lat_ != nullptr) {
2048  update_lattice(jmu_B_lat_.get(), LatticeUpdate::EveryTimestep,
2049  DensityType::Baryon, density_param_, particles_, true);
2050  const size_t UBlattice_size = UB_lat_->size();
2051  for (size_t i = 0; i < UBlattice_size; i++) {
2052  auto jB = (*jmu_B_lat_)[i];
2053  const FourVector flow_four_velocity_B =
2054  std::abs(jB.density()) > really_small ? jB.jmu_net() / jB.density()
2055  : FourVector();
2056  double baryon_density = jB.density();
2057  ThreeVector baryon_grad_rho = jB.grad_rho();
2058  ThreeVector baryon_dj_dt = jB.dj_dt();
2059  ThreeVector baryon_rot_j = jB.rot_j();
2060  if (potentials_->use_skyrme()) {
2061  (*UB_lat_)[i] =
2062  flow_four_velocity_B * potentials_->skyrme_pot(baryon_density);
2063  (*FB_lat_)[i] = potentials_->skyrme_force(
2064  baryon_density, baryon_grad_rho, baryon_dj_dt, baryon_rot_j);
2065  }
2066  if (potentials_->use_symmetry() && jmu_I3_lat_ != nullptr) {
2067  auto jI3 = (*jmu_I3_lat_)[i];
2068  const FourVector flow_four_velocity_I3 =
2069  std::abs(jI3.density()) > really_small
2070  ? jI3.jmu_net() / jI3.density()
2071  : FourVector();
2072  (*UI3_lat_)[i] =
2073  flow_four_velocity_I3 *
2074  potentials_->symmetry_pot(jI3.density(), baryon_density);
2075  (*FI3_lat_)[i] = potentials_->symmetry_force(
2076  jI3.density(), jI3.grad_rho(), jI3.dj_dt(), jI3.rot_j(),
2077  baryon_density, baryon_grad_rho, baryon_dj_dt, baryon_rot_j);
2078  }
2079  }
2080  }
2081  }
2082 }
2083 
2084 template <typename Modus>
2086  /* At end of time evolution: Force all resonances to decay. In order to handle
2087  * decay chains, we need to loop until no further actions occur. */
2088  uint64_t interactions_old;
2089  const auto particles_before_actions = particles_.copy_to_vector();
2090  do {
2091  Actions actions;
2092 
2093  interactions_old = interactions_total_;
2094 
2095  // Dileptons: shining of remaining resonances
2096  if (dilepton_finder_ != nullptr) {
2097  for (const auto &output : outputs_) {
2098  dilepton_finder_->shine_final(particles_, output.get(), true);
2099  }
2100  }
2101  // Find actions.
2102  for (const auto &finder : action_finders_) {
2103  actions.insert(finder->find_final_actions(particles_));
2104  }
2105  // Perform actions.
2106  while (!actions.is_empty()) {
2107  perform_action(*actions.pop(), particles_before_actions);
2108  }
2109  // loop until no more decays occur
2110  } while (interactions_total_ > interactions_old);
2111 
2112  // Dileptons: shining of stable particles at the end
2113  if (dilepton_finder_ != nullptr) {
2114  for (const auto &output : outputs_) {
2115  dilepton_finder_->shine_final(particles_, output.get(), false);
2116  }
2117  }
2118 }
2119 
2120 template <typename Modus>
2121 void Experiment<Modus>::final_output(const int evt_num) {
2122  /* make sure the experiment actually ran (note: we should compare this
2123  * to the start time, but we don't know that. Therefore, we check that
2124  * the time is positive, which should heuristically be the same). */
2125  double E_mean_field = 0.0;
2126  if (likely(parameters_.labclock > 0)) {
2127  const uint64_t wall_actions_this_interval =
2128  wall_actions_total_ - previous_wall_actions_total_;
2129  const uint64_t interactions_this_interval = interactions_total_ -
2130  previous_interactions_total_ -
2131  wall_actions_this_interval;
2132  if (potentials_) {
2133  // using the lattice is necessary
2134  if ((jmu_B_lat_ != nullptr)) {
2135  E_mean_field =
2136  calculate_mean_field_energy(*potentials_, *jmu_B_lat_, parameters_);
2137  }
2138  }
2140  particles_, interactions_this_interval, conserved_initial_, time_start_,
2141  end_time_, E_mean_field, initial_mean_field_energy_);
2142  if (IC_output_switch_ && (particles_.size() == 0)) {
2143  // Verify there is no more energy in the system if all particles were
2144  // removed when crossing the hypersurface
2145  const double remaining_energy =
2146  conserved_initial_.momentum().x0() - total_energy_removed_;
2147  if (remaining_energy > really_small) {
2148  throw std::runtime_error(
2149  "There is remaining energy in the system although all particles "
2150  "were removed.\n"
2151  "E_remain = " +
2152  std::to_string(remaining_energy) + " [GeV]");
2153  } else {
2154  logg[LExperiment].info() << hline;
2155  logg[LExperiment].info()
2156  << "Time real: " << SystemClock::now() - time_start_;
2157  logg[LExperiment].info()
2158  << "Interactions before reaching hypersurface: "
2159  << interactions_total_ - wall_actions_total_ -
2160  total_hypersurface_crossing_actions_;
2161  logg[LExperiment].info()
2162  << "Total number of particles removed on hypersurface: "
2163  << total_hypersurface_crossing_actions_;
2164  }
2165  } else {
2166  const double precent_discarded =
2167  interactions_total_ > 0
2168  ? static_cast<double>(discarded_interactions_total_) * 100.0 /
2169  interactions_total_
2170  : 0.0;
2171  std::stringstream msg_discarded;
2172  msg_discarded
2173  << "Discarded interaction number: " << discarded_interactions_total_
2174  << " (" << precent_discarded
2175  << "% of the total interaction number including wall crossings)";
2176 
2177  logg[LExperiment].info() << hline;
2178  logg[LExperiment].info()
2179  << "Time real: " << SystemClock::now() - time_start_;
2180  logg[LExperiment].debug() << msg_discarded.str();
2181 
2182  if (parameters_.coll_crit == CollisionCriterion::Stochastic &&
2183  precent_discarded > 1.0) {
2184  // The choosen threshold of 1% is a heuristical value
2185  logg[LExperiment].warn()
2186  << msg_discarded.str()
2187  << "\nThe number of discarded interactions is large, which means "
2188  "the assumption for the stochastic criterion of\n1 interaction"
2189  "per particle per timestep is probably violated. Consider "
2190  "reducing the timestep size.";
2191  }
2192 
2193  logg[LExperiment].info() << "Final interaction number: "
2194  << interactions_total_ - wall_actions_total_;
2195  }
2196 
2197  // Check if there are unformed particles
2198  int unformed_particles_count = 0;
2199  for (const auto &particle : particles_) {
2200  if (particle.formation_time() > end_time_) {
2201  unformed_particles_count++;
2202  }
2203  }
2204  if (unformed_particles_count > 0) {
2205  logg[LExperiment].warn(
2206  "End time might be too small. ", unformed_particles_count,
2207  " unformed particles were found at the end of the evolution.");
2208  }
2209  }
2210 
2211  auto event_info =
2212  fill_event_info(particles_, E_mean_field, modus_.impact_parameter(),
2213  parameters_, projectile_target_interact_);
2214 
2215  for (const auto &output : outputs_) {
2216  output->at_eventend(particles_, evt_num, event_info);
2217  }
2218 }
2219 
2220 template <typename Modus>
2222  const auto &mainlog = logg[LMain];
2223  for (int j = 0; j < nevents_; j++) {
2224  mainlog.info() << "Event " << j;
2225 
2226  // Sample initial particles, start clock, some printout and book-keeping
2227  initialize_new_event(j);
2228  /* In the ColliderModus, if the first collisions within the same nucleus are
2229  * forbidden, 'nucleon_has_interacted_', which records whether a nucleon has
2230  * collided with another nucleon, is initialized equal to false. If allowed,
2231  * 'nucleon_has_interacted' is initialized equal to true, which means these
2232  * incoming particles have experienced some fake scatterings, they can
2233  * therefore collide with each other later on since these collisions are not
2234  * "first" to them. */
2235  if (modus_.is_collider()) {
2236  if (!modus_.cll_in_nucleus()) {
2237  nucleon_has_interacted_.assign(modus_.total_N_number(), false);
2238  } else {
2239  nucleon_has_interacted_.assign(modus_.total_N_number(), true);
2240  }
2241  }
2242  /* In the ColliderModus, if Fermi motion is frozen, assign the beam momenta
2243  * to the nucleons in both the projectile and the target. */
2244  if (modus_.is_collider() && modus_.fermi_motion() == FermiMotion::Frozen) {
2245  for (int i = 0; i < modus_.total_N_number(); i++) {
2246  const auto mass_beam = particles_.copy_to_vector()[i].effective_mass();
2247  const auto v_beam = i < modus_.proj_N_number()
2248  ? modus_.velocity_projectile()
2249  : modus_.velocity_target();
2250  const auto gamma = 1.0 / std::sqrt(1.0 - v_beam * v_beam);
2251  beam_momentum_.emplace_back(FourVector(gamma * mass_beam, 0.0, 0.0,
2252  gamma * v_beam * mass_beam));
2253  }
2254  }
2255 
2256  run_time_evolution();
2257 
2258  if (force_decays_) {
2259  do_final_decays();
2260  }
2261 
2262  // Output at event end
2263  final_output(j);
2264  }
2265 }
2266 
2267 } // namespace smash
2268 
2269 #endif // SRC_INCLUDE_SMASH_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:384
smash::Experiment::particles
Particles * particles()
Provides external access to SMASH particles.
Definition: experiment.h:228
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:493
smash::Experiment::jmu_B_lat_
std::unique_ptr< DensityLattice > jmu_B_lat_
Baryon density on the lattices.
Definition: experiment.h:381
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:1714
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:422
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:419
smash::Experiment::action_finders_
std::vector< std::unique_ptr< ActionFinderInterface > > action_finders_
The Action finder objects.
Definition: experiment.h:369
smash::Experiment::initialize_new_event
void initialize_new_event(int event_number)
This is called in the beginning of each event.
Definition: experiment.h:1451
smash::LInitialConditions
static constexpr int LInitialConditions
Definition: experiment.h:83
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:484
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:154
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:565
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:541
smash::DensityParameters
A class to pre-calculate and store parameters relevant for density calculation.
Definition: density.h:108
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:134
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:366
smash::Experiment::next_output_time
double next_output_time() const
Shortcut for next output time.
Definition: experiment.h:303
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:298
TimeStepMode
TimeStepMode
The time step mode.
Definition: forwarddeclarations.h:114
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:535
smash::Experiment::final_output
void final_output(const int evt_num)
Output at the end of an event.
Definition: experiment.h:2121
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:1829
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:402
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:336
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:34
ExpansionMode::NoExpansion
smash::Experiment::modus_
Modus modus_
Instance of the Modus template parameter.
Definition: experiment.h:321
smash::Experiment::modus
Modus * modus()
Provides external access to SMASH calculation modus.
Definition: experiment.h:234
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:490
smash::Experiment::projectile_target_interact_
bool projectile_target_interact_
Whether the projectile and the target collided.
Definition: experiment.h:359
FermiMotion::Off
Don't use fermi motion.
grid.h
smash::operator<<
std::ostream & operator<<(std::ostream &out, const ActionPtr &action)
Definition: action.h:518
energymomentumtensor.h
smash::Experiment::jmu_custom_lat_
std::unique_ptr< DensityLattice > jmu_custom_lat_
Custom density on the lattices.
Definition: experiment.h:393
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:1837
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:324
propagation.h
smash::LatticeUpdate::EveryTimestep
smash::LMain
static constexpr int LMain
Definition: experiment.h:82
smash::Experiment::run
void run() override
Runs the experiment.
Definition: experiment.h:2221
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:472
smash::Experiment::total_pauli_blocked_
uint64_t total_pauli_blocked_
Total number of Pauli-blockings for current timestep.
Definition: experiment.h:547
smash::Experiment::dileptons_switch_
const bool dileptons_switch_
This indicates whether dileptons are switched on.
Definition: experiment.h:481
smash::Experiment::dens_type_
DensityType dens_type_
Type of density to be written to collision headers.
Definition: experiment.h:517
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:348
grandcan_thermalizer.h
smash::Experiment::thermalizer_
std::unique_ptr< GrandCanThermalizer > thermalizer_
Instance of class used for forced thermalization.
Definition: experiment.h:434
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: stringprocess.h:45
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:315
smash::Experiment::nevents_
const int nevents_
Number of events.
Definition: experiment.h:456
smash::Experiment::seed_
int64_t seed_
random seed for the next event.
Definition: experiment.h:568
chrono.h
smash::fill_event_info
EventInfo fill_event_info(const Particles &particles, double E_mean_field, double modus_impact_parameter, const ExperimentParameters &parameters, bool projectile_target_interact)
Generate the EventInfo object which is passed to outputs_.
Definition: experiment.cc:648
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:243
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:307
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:769
smash::Experiment::metric_
const ExpansionProperties metric_
This struct contains information on the metric to be used.
Definition: experiment.h:478
CollisionCriterion::Stochastic
Stochastic Criteiron.
actions.h
smash::Experiment::printout_tmn_landau_
bool printout_tmn_landau_
Whether to print the energy-momentum tensor in Landau frame.
Definition: experiment.h:425
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:416
smash::create_experiment_parameters
ExperimentParameters create_experiment_parameters(Configuration config)
Gathers all general Experiment parameters.
Definition: experiment.cc:335
smash::Experiment::do_final_decays
void do_final_decays()
Performs the final decays of an event.
Definition: experiment.h:2085
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:402
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:28
smash::Experiment::use_grid_
const bool use_grid_
This indicates whether to use the grid.
Definition: experiment.h:475
smash::Experiment::interactions_total_
uint64_t interactions_total_
Total number of interactions for current timestep.
Definition: experiment.h:523
smash::ExperimentBase::InvalidModusRequest
Definition: experiment.h:137
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:1923
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:580
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:1812
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:378
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:124
stringprocess.h
smash::Experiment::delta_time_startup_
const double delta_time_startup_
The clock's timestep size at start up.
Definition: experiment.h:466
smash::Experiment::IC_output_switch_
const bool IC_output_switch_
This indicates whether the IC output is enabled.
Definition: experiment.h:490
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:148
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:529
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:36
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:2039
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:412
smash::ExperimentParameters::outputclock
std::unique_ptr< Clock > outputclock
Output clock to keep track of the next output time.
Definition: experimentparameters.h:29
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.
smash::Experiment::dens_type_lattice_printout_
DensityType dens_type_lattice_printout_
Type of density for lattice printout.
Definition: experiment.h:396
smash::Experiment::dilepton_finder_
std::unique_ptr< DecayActionsFinderDilepton > dilepton_finder_
The Dilepton Action Finder.
Definition: experiment.h:372
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:431
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:587
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:536
smash::Experiment::n_fractional_photons_
int n_fractional_photons_
Number of fractional photons produced per single reaction.
Definition: experiment.h:378
smash::Experiment::dilepton_output_
OutputPtr dilepton_output_
The Dilepton output.
Definition: experiment.h:345
smash::Experiment::parameters_
ExperimentParameters parameters_
Struct of several member variables.
Definition: experiment.h:312
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:24
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:459
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:355
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:440
smash::Action::get_interaction_point
FourVector get_interaction_point() const
Get the interaction point.
Definition: action.cc:67
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:330
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:511
hepmcoutput.h
smash::Experiment
The main class, where the simulation of an experiment is executed.
Definition: experiment.h:143
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:342
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:57
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:553
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:428
smash::Experiment::discarded_interactions_total_
uint64_t discarded_interactions_total_
Total number of discarded interactions, because they were invalidated before they could be performed.
Definition: experiment.h:559
smash::Experiment::time_start_
SystemTimePoint time_start_
system starting time of the simulation
Definition: experiment.h:514
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:408
smash::Experiment::bremsstrahlung_switch_
const bool bremsstrahlung_switch_
This indicates whether bremsstrahlung is switched on.
Definition: experiment.h:487
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:26
FermiMotion
FermiMotion
Option to use Fermi Motion.
Definition: forwarddeclarations.h:94
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:375
smash::Experiment::conserved_initial_
QuantumNumbers conserved_initial_
The conserved quantities of the system.
Definition: experiment.h:505
smash::ExperimentBase
Non-template interface to Experiment<Modus>.
Definition: experiment.h:92
smash::Experiment::max_transverse_distance_sqr_
double max_transverse_distance_sqr_
Maximal distance at which particles can interact, squared.
Definition: experiment.h:496
potentials.h