Version: SMASH-3.2
experiment.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2013-2025
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 <set>
14 #include <string>
15 #include <utility>
16 #include <vector>
17 
18 #include "actionfinderfactory.h"
19 #include "actions.h"
20 #include "bremsstrahlungaction.h"
21 #include "chrono.h"
22 #include "decayactionsfinder.h"
24 #include "dynamicfluidfinder.h"
25 #include "energymomentumtensor.h"
26 #include "fields.h"
27 #include "fourvector.h"
28 #include "grandcan_thermalizer.h"
29 #include "grid.h"
31 #include "icparameters.h"
32 #include "numeric_cast.h"
33 #include "outputparameters.h"
34 #include "pauliblocking.h"
35 #include "potential_globals.h"
36 #include "potentials.h"
37 #include "propagation.h"
38 #include "quantumnumbers.h"
39 #include "scatteractionphoton.h"
40 #include "scatteractionsfinder.h"
41 #include "stringprocess.h"
42 #include "thermalizationaction.h"
43 // Output
44 #include "binaryoutput.h"
45 #ifdef SMASH_USE_HEPMC
46 #include "hepmcoutput.h"
47 #endif
48 #ifdef SMASH_USE_RIVET
49 #include "rivetoutput.h"
50 #endif
51 #include "icoutput.h"
52 #include "oscaroutput.h"
54 #include "thermodynamicoutput.h"
55 #ifdef SMASH_USE_ROOT
56 #include "rootoutput.h"
57 #endif
58 #include "freeforallaction.h"
59 #include "vtkoutput.h"
60 #include "wallcrossingaction.h"
61 
62 namespace std {
73 template <typename T, typename Ratio>
74 static ostream &operator<<(ostream &out,
75  const chrono::duration<T, Ratio> &seconds) {
76  using Seconds = chrono::duration<double>;
77  using Minutes = chrono::duration<double, std::ratio<60>>;
78  using Hours = chrono::duration<double, std::ratio<60 * 60>>;
79  constexpr Minutes threshold_for_minutes{10};
80  constexpr Hours threshold_for_hours{3};
81  if (seconds < threshold_for_minutes) {
82  return out << Seconds(seconds).count() << " [s]";
83  }
84  if (seconds < threshold_for_hours) {
85  return out << Minutes(seconds).count() << " [min]";
86  }
87  return out << Hours(seconds).count() << " [h]";
88 }
89 } // namespace std
90 
91 namespace smash {
92 static constexpr int LMain = LogArea::Main::id;
93 static constexpr int LInitialConditions = LogArea::InitialConditions::id;
94 
103  public:
104  ExperimentBase() = default;
109  virtual ~ExperimentBase() = default;
110 
131  static std::unique_ptr<ExperimentBase> create(
132  Configuration &config, const std::filesystem::path &output_path);
133 
140  virtual void run() = 0;
141 
147  struct InvalidModusRequest : public std::invalid_argument {
148  using std::invalid_argument::invalid_argument;
149  };
150 
156  struct NonExistingOutputPathRequest : public std::invalid_argument {
157  using std::invalid_argument::invalid_argument;
158  };
159 };
160 
161 template <typename Modus>
162 class Experiment;
163 template <typename Modus>
164 std::ostream &operator<<(std::ostream &out, const Experiment<Modus> &e);
165 
190 template <typename Modus>
191 class Experiment : public ExperimentBase {
192  friend class ExperimentBase;
193 
194  public:
201  void run() override;
202 
216  explicit Experiment(Configuration &config,
217  const std::filesystem::path &output_path);
218 
224  void initialize_new_event();
225 
243  void run_time_evolution(const double t_end, ParticleList &&add_plist = {},
244  ParticleList &&remove_plist = {});
245 
251  void do_final_decays();
252 
254  void final_output();
255 
262  std::vector<Particles> *all_ensembles() { return &ensembles_; }
263 
268  Modus *modus() { return &modus_; }
269 
274  void increase_event_number();
275 
276  private:
289  bool perform_action(Action &action, int i_ensemble,
290  bool include_pauli_blocking = true);
299  void create_output(const std::string &format, const std::string &content,
300  const std::filesystem::path &output_path,
301  const OutputParameters &par);
302 
310  void propagate_and_shine(double to_time, Particles &particles);
311 
324  void run_time_evolution_timestepless(Actions &actions, int i_ensemble,
325  const double end_time_propagation);
326 
328  void intermediate_output();
329 
331  void update_potentials();
332 
341  double compute_min_cell_length(double dt) const {
344  }
345  return std::sqrt(4 * dt * dt + max_transverse_distance_sqr_);
346  }
347 
349  double next_output_time() const {
350  return parameters_.outputclock->next_time();
351  }
352 
358 
364  bool is_finished();
365 
372 
375 
380  Modus modus_;
381 
383  std::vector<Particles> ensembles_;
384 
389  std::unique_ptr<Potentials> potentials_;
390 
395  std::unique_ptr<PauliBlocker> pauli_blocker_;
396 
401  OutputsList outputs_;
402 
404  OutputPtr dilepton_output_;
405 
407  OutputPtr photon_output_;
408 
413  std::vector<bool> projectile_target_interact_;
414 
420  std::vector<FourVector> beam_momentum_ = {};
421 
423  std::vector<std::unique_ptr<ActionFinderInterface>> action_finders_;
424 
426  std::unique_ptr<DecayActionsFinderDilepton> dilepton_finder_;
427 
429  std::unique_ptr<ActionFinderInterface> photon_finder_;
430 
433 
435  std::unique_ptr<DensityLattice> j_QBS_lat_;
436 
438  std::unique_ptr<DensityLattice> jmu_B_lat_;
439 
441  std::unique_ptr<DensityLattice> jmu_I3_lat_;
442 
444  std::unique_ptr<DensityLattice> jmu_el_lat_;
445 
447  std::unique_ptr<FieldsLattice> fields_lat_;
448 
456  std::unique_ptr<DensityLattice> jmu_custom_lat_;
457 
460 
465  std::unique_ptr<RectangularLattice<FourVector>> UB_lat_ = nullptr;
466 
471  std::unique_ptr<RectangularLattice<FourVector>> UI3_lat_ = nullptr;
472 
477  std::unique_ptr<RectangularLattice<std::pair<ThreeVector, ThreeVector>>>
479 
481  std::unique_ptr<RectangularLattice<std::pair<ThreeVector, ThreeVector>>>
483 
485  std::unique_ptr<RectangularLattice<std::pair<ThreeVector, ThreeVector>>>
487 
489  std::unique_ptr<RectangularLattice<EnergyMomentumTensor>> Tmn_;
490 
492  std::unique_ptr<RectangularLattice<FourVector>> old_jmu_auxiliary_;
494  std::unique_ptr<RectangularLattice<FourVector>> new_jmu_auxiliary_;
496  std::unique_ptr<RectangularLattice<std::array<FourVector, 4>>>
498 
500  std::unique_ptr<RectangularLattice<FourVector>> old_fields_auxiliary_;
502  std::unique_ptr<RectangularLattice<FourVector>> new_fields_auxiliary_;
504  std::unique_ptr<RectangularLattice<std::array<FourVector, 4>>>
506 
508  bool printout_rho_eckart_ = false;
509 
511  bool printout_tmn_ = false;
512 
514  bool printout_tmn_landau_ = false;
515 
517  bool printout_v_landau_ = false;
518 
520  bool printout_j_QBS_ = false;
521 
523  bool printout_lattice_td_ = false;
524 
528 
530  std::unique_ptr<GrandCanThermalizer> thermalizer_;
531 
537 
552  int nevents_ = 0;
553 
563 
571 
573  int event_ = 0;
574 
577 
582  int max_events_ = 0;
583 
585  const double end_time_;
586 
592  const double delta_time_startup_;
593 
598  const bool force_decays_;
599 
601  const bool use_grid_;
602 
605 
607  const bool dileptons_switch_;
608 
610  const bool photons_switch_;
611 
614 
619  const bool IC_switch_;
620 
622  const bool IC_dynamic_;
623 
626 
631  double max_transverse_distance_sqr_ = std::numeric_limits<double>::max();
632 
641 
647 
649  SystemTimePoint time_start_ = SystemClock::now();
650 
653 
658  uint64_t interactions_total_ = 0;
659 
665 
670  uint64_t wall_actions_total_ = 0;
671 
677 
682  uint64_t total_pauli_blocked_ = 0;
683 
689 
695 
699  double total_energy_removed_ = 0.0;
700 
705 
708 
710  int64_t seed_ = -1;
711 
717  friend std::ostream &operator<<<>(std::ostream &out, const Experiment &e);
718 };
719 
721 template <typename Modus>
722 std::ostream &operator<<(std::ostream &out, const Experiment<Modus> &e) {
723  out << "End time: " << e.end_time_ << " fm\n";
724  out << e.modus_;
725  return out;
726 }
737 void validate_duplicate_IC_config(double, std::optional<double>, std::string);
738 
739 template <typename Modus>
740 void Experiment<Modus>::create_output(const std::string &format,
741  const std::string &content,
742  const std::filesystem::path &output_path,
743  const OutputParameters &out_par) {
744  // Disable output which do not properly work with multiple ensembles
745  if (ensembles_.size() > 1) {
746  auto abort_because_of = [](const std::string &s) {
747  throw std::invalid_argument(
748  s + " output is not available with multiple parallel ensembles.");
749  };
750  if (content == "Initial_Conditions") {
751  abort_because_of("Initial_Conditions");
752  }
753  if ((format == "HepMC") || (format == "HepMC_asciiv3") ||
754  (format == "HepMC_treeroot")) {
755  abort_because_of("HepMC");
756  }
757  if (content == "Rivet") {
758  abort_because_of("Rivet");
759  }
760  if (content == "Collisions") {
761  logg[LExperiment].warn(
762  "Information coming from different ensembles in 'Collisions' output "
763  "is not distinguishable.\nSuch an output with multiple parallel "
764  "ensembles should only be used if later in the data analysis\nit is "
765  "not necessary to trace back which data belongs to which ensemble.");
766  }
767  }
768 
769  if (format == "VTK" && content == "Particles") {
770  outputs_.emplace_back(
771  std::make_unique<VtkOutput>(output_path, content, out_par));
772  } else if (format == "Root") {
773 #ifdef SMASH_USE_ROOT
774  if (content == "Initial_Conditions") {
775  outputs_.emplace_back(
776  std::make_unique<RootOutput>(output_path, "SMASH_IC", out_par));
777  } else {
778  outputs_.emplace_back(
779  std::make_unique<RootOutput>(output_path, content, out_par));
780  }
781 #else
782  logg[LExperiment].error(
783  "Root output requested, but Root support not compiled in");
784 #endif
785  } else if (format == "Binary" &&
786  (content == "Collisions" || content == "Dileptons" ||
787  content == "Photons" || content == "Particles" ||
788  content == "Initial_Conditions")) {
789  outputs_.emplace_back(
790  create_binary_output(format, content, output_path, out_par));
791  } else if (format == "Oscar2013_bin" &&
792  (content == "Collisions" || content == "Particles")) {
793  outputs_.emplace_back(
794  create_binary_output(format, content, output_path, out_par));
795  } else if (format == "Oscar1999" || format == "Oscar2013") {
796  outputs_.emplace_back(
797  create_oscar_output(format, content, output_path, out_par));
798  } else if (format == "ASCII" &&
799  (content == "Particles" || content == "Collisions")) {
800  outputs_.emplace_back(
801  create_oscar_output(format, content, output_path, out_par));
802  } else if (content == "Thermodynamics" && format == "ASCII") {
803  outputs_.emplace_back(
804  std::make_unique<ThermodynamicOutput>(output_path, content, out_par));
805  } else if (content == "Thermodynamics" &&
806  (format == "Lattice_ASCII" || format == "Lattice_Binary")) {
807  printout_full_lattice_any_td_ = true;
808  outputs_.emplace_back(std::make_unique<ThermodynamicLatticeOutput>(
809  output_path, content, out_par, format == "Lattice_ASCII",
810  format == "Lattice_Binary"));
811  } else if (content == "Thermodynamics" && format == "VTK") {
812  printout_lattice_td_ = true;
813  outputs_.emplace_back(
814  std::make_unique<VtkOutput>(output_path, content, out_par));
815  } else if (content == "Initial_Conditions" && format == "ASCII") {
816  if (IC_dynamic_) {
817  throw std::invalid_argument(
818  "Dynamic initial conditions are only available in Oscar2013 and "
819  "Binary formats.");
820  }
821  outputs_.emplace_back(
822  std::make_unique<ICOutput>(output_path, "SMASH_IC", out_par));
823  } else if ((format == "HepMC") || (format == "HepMC_asciiv3") ||
824  (format == "HepMC_treeroot")) {
825 #ifdef SMASH_USE_HEPMC
826  if (content == "Particles") {
827  if ((format == "HepMC") || (format == "HepMC_asciiv3")) {
828  outputs_.emplace_back(std::make_unique<HepMcOutput>(
829  output_path, "SMASH_HepMC_particles", false, "asciiv3"));
830  } else if (format == "HepMC_treeroot") {
831 #ifdef SMASH_USE_HEPMC_ROOTIO
832  outputs_.emplace_back(std::make_unique<HepMcOutput>(
833  output_path, "SMASH_HepMC_particles", false, "root"));
834 #else
835  logg[LExperiment].error(
836  "Requested HepMC_treeroot output not available, "
837  "ROOT or HepMC3_ROOTIO missing or not found by cmake.");
838 #endif
839  }
840  } else if (content == "Collisions") {
841  if ((format == "HepMC") || (format == "HepMC_asciiv3")) {
842  outputs_.emplace_back(std::make_unique<HepMcOutput>(
843  output_path, "SMASH_HepMC_collisions", true, "asciiv3"));
844  } else if (format == "HepMC_treeroot") {
845 #ifdef SMASH_USE_HEPMC_ROOTIO
846  outputs_.emplace_back(std::make_unique<HepMcOutput>(
847  output_path, "SMASH_HepMC_collisions", true, "root"));
848 #else
849  logg[LExperiment].error(
850  "Requested HepMC_treeroot output not available, "
851  "ROOT or HepMC3_ROOTIO missing or not found by cmake.");
852 #endif
853  }
854  } else {
855  logg[LExperiment].error(
856  "HepMC only available for Particles and "
857  "Collisions content. Requested for " +
858  content + ".");
859  }
860 #else
861  logg[LExperiment].error(
862  "HepMC output requested, but HepMC support not compiled in");
863 #endif
864  } else if (content == "Coulomb" && format == "VTK") {
865  outputs_.emplace_back(
866  std::make_unique<VtkOutput>(output_path, "Fields", out_par));
867  } else if (content == "Rivet") {
868 #ifdef SMASH_USE_RIVET
869  // flag to ensure that the Rivet format has not been already assigned
870  static bool rivet_format_already_selected = false;
871  // if the next check is true, then we are trying to assign the format twice
872  if (rivet_format_already_selected) {
873  logg[LExperiment].warn(
874  "Rivet output format can only be one, either YODA or YODA-full. "
875  "Only your first valid choice will be used.");
876  return;
877  }
878  if (format == "YODA") {
879  outputs_.emplace_back(std::make_unique<RivetOutput>(
880  output_path, "SMASH_Rivet", false, out_par.rivet_parameters));
881  rivet_format_already_selected = true;
882  } else if (format == "YODA-full") {
883  outputs_.emplace_back(std::make_unique<RivetOutput>(
884  output_path, "SMASH_Rivet_full", true, out_par.rivet_parameters));
885  rivet_format_already_selected = true;
886  } else {
887  logg[LExperiment].error("Rivet format " + format +
888  "not one of YODA or YODA-full");
889  }
890 #else
891  logg[LExperiment].error(
892  "Rivet output requested, but Rivet support not compiled in");
893 #endif
894  } else {
895  logg[LExperiment].error()
896  << "Unknown combination of format (" << format << ") and content ("
897  << content << "). Fix the config.";
898  }
899 
900  logg[LExperiment].info() << "Added output " << content << " of format "
901  << format << "\n";
902 }
903 
912 
913 template <typename Modus>
915  const std::filesystem::path &output_path)
916  : parameters_(create_experiment_parameters(config)),
917  density_param_(DensityParameters(parameters_)),
918  modus_(std::invoke([&]() {
919  /* This immediately invoked lambda is a work-around to cope with the
920  * fact that the "Collisions_Within_Nucleus" key belongs to the
921  * "Collider" section, but is used by the ScatterActionsFinder through
922  * the ScatterActionsFinderParameters member. Here that key is taken
923  * from the main configuration and put there back after the "Collider"
924  * section is extracted. If this were not done in this way, the
925  * sub-configuration given to ColliderModus would be deleted not empty
926  * at the end of its constructor and this would throw an exception.*/
928  const bool restore_key = config.has_value(key);
929  const bool temporary_taken_key = config.take(key);
930  auto modus_config =
932  if (restore_key) {
933  config.set_value(key, temporary_taken_key);
934  }
935  return Modus{std::move(modus_config), parameters_};
936  })),
937  ensembles_(parameters_.n_ensembles),
938  end_time_(config.take(InputKeys::gen_endTime)),
939  delta_time_startup_(parameters_.labclock->timestep_duration()),
940  force_decays_(config.take(InputKeys::collTerm_forceDecaysAtEnd)),
941  use_grid_(config.take(InputKeys::gen_useGrid)),
942  metric_(config.take(InputKeys::gen_metricType),
944  dileptons_switch_(config.take(InputKeys::collTerm_dileptons_decays)),
945  photons_switch_(
947  bremsstrahlung_switch_(
949  IC_switch_(config.has_section(InputSections::o_initialConditions) &&
950  modus_.is_IC_for_hybrid()),
951  IC_dynamic_(IC_switch_ ? (modus_.IC_parameters().type ==
953  : false),
954  time_step_mode_(config.take(InputKeys::gen_timeStepMode)) {
955  logg[LExperiment].info() << *this;
956 
957  const bool user_wants_nevents = config.has_value(InputKeys::gen_nevents);
958  const bool user_wants_min_nonempty =
960  if (user_wants_nevents == user_wants_min_nonempty) {
961  throw std::invalid_argument(
962  "Please specify either Nevents or Minimum_Nonempty_Ensembles.");
963  }
964  if (user_wants_nevents) {
965  event_counting_ = EventCounting::FixedNumber;
966  nevents_ = config.take(InputKeys::gen_nevents);
967  } else {
968  event_counting_ = EventCounting::MinimumNonEmpty;
969  minimum_nonempty_ensembles_ =
971  int max_ensembles =
973  max_events_ = numeric_cast<int>(std::ceil(
974  static_cast<double>(max_ensembles) / parameters_.n_ensembles));
975  }
976 
977  // covariant derivatives can only be done with covariant smearing
978  if (parameters_.derivatives_mode == DerivativesMode::CovariantGaussian &&
979  parameters_.smearing_mode != SmearingMode::CovariantGaussian) {
980  throw std::invalid_argument(
981  "Covariant Gaussian derivatives only make sense for Covariant Gaussian "
982  "smearing!");
983  }
984 
985  // for triangular smearing:
986  // the weight needs to be larger than 1./7. for the center cell to contribute
987  // more than the surrounding cells
988  if (parameters_.smearing_mode == SmearingMode::Discrete &&
989  parameters_.discrete_weight < (1. / 7.)) {
990  throw std::invalid_argument(
991  "The central weight for discrete smearing should be >= 1./7.");
992  }
993 
994  if (parameters_.coll_crit == CollisionCriterion::Stochastic &&
995  (time_step_mode_ != TimeStepMode::Fixed || !use_grid_)) {
996  throw std::invalid_argument(
997  "The stochastic criterion can only be employed for fixed time step "
998  "mode and with a grid!");
999  }
1000 
1001  if (modus_.is_box() && (time_step_mode_ != TimeStepMode::Fixed)) {
1002  throw std::invalid_argument(
1003  "The box modus can only be used with the fixed time step mode!");
1004  }
1005 
1006  logg[LExperiment].info("Using ", parameters_.testparticles,
1007  " testparticles per particle.");
1008  logg[LExperiment].info("Using ", parameters_.n_ensembles,
1009  " parallel ensembles.");
1010 
1011  if (modus_.is_box() && config.read(InputKeys::collTerm_totXsStrategy) !=
1013  logg[LExperiment].warn(
1014  "To preserve detailed balance in a box simulation, it is recommended\n"
1015  "to use the bottom-up strategy for evaluating total cross sections.\n"
1016  "Consider adding the following line to the 'Collision_Term' section "
1017  "in your configuration file:\n"
1018  " Total_Cross_Section_Strategy: \"BottomUp\"");
1019  }
1020  if (modus_.is_box() && config.read(InputKeys::collTerm_pseudoresonance) !=
1022  logg[LExperiment].warn(
1023  "To preserve detailed balance in a box simulation, it is recommended "
1024  "to not include the pseudoresonances,\nas they artificially increase "
1025  "the resonance production without changing the corresponding "
1026  "decay.\nConsider adding the following line to the 'Collision_Term' "
1027  "section in your configuration file:\n Pseudoresonance: \"None\"");
1028  }
1029 
1030  const bool IC_output = config.has_section(InputSections::o_initialConditions);
1031  if (IC_output != modus_.is_IC_for_hybrid()) {
1032  throw std::invalid_argument(
1033  "The 'Initial_Conditions' subsection must be present in both 'Output' "
1034  "and 'Modi: Collider' sections.");
1035  }
1036 
1037  /* In collider setup with sqrts >= 200 GeV particles don't form continuously
1038  *
1039  * NOTE: This key has to be taken before the ScatterActionsFinder is created
1040  * because there the "String_Parameters" is extracted as sub-config and
1041  * all parameters but this one are taken. If this one is still there
1042  * the configuration temporary object will be destroyed not empty, hence
1043  * throwing an exception.
1044  */
1047  modus_.sqrt_s_NN() >= 200. ? -1. : 1.);
1048 
1049  // create finders
1050  if (dileptons_switch_) {
1051  dilepton_finder_ = std::make_unique<DecayActionsFinderDilepton>();
1052  }
1053  if (photons_switch_ || bremsstrahlung_switch_) {
1054  n_fractional_photons_ =
1056  }
1057  if (parameters_.two_to_one) {
1058  if (parameters_.res_lifetime_factor < 0.) {
1059  throw std::invalid_argument(
1060  "Resonance lifetime modifier cannot be negative!");
1061  }
1062  if (parameters_.res_lifetime_factor < really_small) {
1063  logg[LExperiment].warn(
1064  "Resonance lifetime set to zero. Make sure resonances cannot "
1065  "interact",
1066  "inelastically (e.g. resonance chains), else SMASH is known to "
1067  "hang.");
1068  }
1069  action_finders_.emplace_back(std::make_unique<DecayActionsFinder>(
1070  parameters_.res_lifetime_factor, parameters_.do_non_strong_decays));
1071  }
1072  bool no_coll = config.take(InputKeys::collTerm_noCollisions);
1073  if ((parameters_.two_to_one || parameters_.included_2to2.any() ||
1074  parameters_.included_multi.any() || parameters_.strings_switch) &&
1075  !no_coll) {
1076  parameters_.use_monash_tune_default =
1077  (modus_.is_collider() && modus_.sqrt_s_NN() >= 200.);
1078  auto scat_finder =
1079  std::make_unique<ScatterActionsFinder>(config, parameters_);
1080  max_transverse_distance_sqr_ =
1081  scat_finder->max_transverse_distance_sqr(parameters_.testparticles);
1082  process_string_ptr_ = scat_finder->get_process_string_ptr();
1083  action_finders_.emplace_back(std::move(scat_finder));
1084  } else {
1085  max_transverse_distance_sqr_ =
1086  parameters_.maximum_cross_section / M_PI * fm2_mb;
1087  process_string_ptr_ = NULL;
1088  }
1089  if (modus_.is_box()) {
1090  action_finders_.emplace_back(
1091  std::make_unique<WallCrossActionsFinder>(parameters_.box_length));
1092  }
1093 
1094  if (IC_switch_) {
1095  const InitialConditionParameters &IC_parameters = modus_.IC_parameters();
1096  if (IC_dynamic_) {
1097  // Dynamic fluidization
1098  action_finders_.emplace_back(std::make_unique<DynamicFluidizationFinder>(
1099  modus_.fluid_lattice(), modus_.fluid_background(), IC_parameters));
1100  } else {
1101  // Iso-tau hypersurface
1102  /*
1103  * Due to an ongoing refactoring, the physics inputs for Initial
1104  * Conditions the former being deprecated. If there are two inconsistent
1105  * values, SMASH will not run. Otherwise it will follow the one present in
1106  * the configuration. If none are present, the default is used. When the
1107  * deprecated way is removed, the key taking will be handled in the
1108  * constructor of ColliderModus, and the logic here will be removed.
1109  */
1110  double rapidity_cut = IC_parameters.rapidity_cut.has_value()
1111  ? IC_parameters.rapidity_cut.value()
1112  : 0.0;
1114  rapidity_cut =
1116  validate_duplicate_IC_config(rapidity_cut, IC_parameters.rapidity_cut,
1117  "Rapidity_Cut");
1118  }
1119 
1120  if (rapidity_cut < 0.0) {
1121  logg[LInitialConditions].fatal()
1122  << "Rapidity cut for initial conditions configured as abs(y) < "
1123  << rapidity_cut << " is unreasonable. \nPlease choose a positive, "
1124  << "non-zero value or employ SMASH without rapidity cut.";
1125  throw std::runtime_error(
1126  "Kinematic cut for initial conditions malconfigured.");
1127  }
1128 
1129  if (modus_.calculation_frame_is_fixed_target() && rapidity_cut != 0.0) {
1130  throw std::runtime_error(
1131  "Rapidity cut for initial conditions output is not implemented "
1132  "in the fixed target calculation frame. \nPlease use "
1133  "\"center of velocity\" or \"center of mass\" as a "
1134  "\"Calculation_Frame\" instead.");
1135  }
1136 
1137  double pT_cut =
1138  IC_parameters.pT_cut.has_value() ? IC_parameters.pT_cut.value() : 0.0;
1141  validate_duplicate_IC_config(pT_cut, IC_parameters.pT_cut, "pT_Cut");
1142  }
1143  if (pT_cut < 0.0) {
1144  logg[LInitialConditions].fatal()
1145  << "transverse momentum cut for initial conditions configured as "
1146  << "pT < " << pT_cut << " is unreasonable. \nPlease choose a "
1147  << "positive, non-zero value or employ SMASH without pT cut.";
1148  throw std::runtime_error(
1149  "Kinematic cut for initial conditions misconfigured.");
1150  }
1151 
1152  if (rapidity_cut > 0.0 || pT_cut > 0.0) {
1153  kinematic_cuts_for_IC_output_ = true;
1154  }
1155 
1156  if (rapidity_cut > 0.0 && pT_cut > 0.0) {
1157  logg[LInitialConditions].info()
1158  << "Extracting initial conditions in kinematic range: "
1159  << -rapidity_cut << " <= y <= " << rapidity_cut
1160  << "; pT <= " << pT_cut << " GeV.";
1161  } else if (rapidity_cut > 0.0) {
1162  logg[LInitialConditions].info()
1163  << "Extracting initial conditions in kinematic range: "
1164  << -rapidity_cut << " <= y <= " << rapidity_cut << ".";
1165  } else if (pT_cut > 0.0) {
1166  logg[LInitialConditions].info()
1167  << "Extracting initial conditions in kinematic range: pT <= "
1168  << pT_cut << " GeV.";
1169  } else {
1170  logg[LInitialConditions].info()
1171  << "Extracting initial conditions without kinematic cuts.";
1172  }
1173 
1174  double proper_time = std::numeric_limits<double>::quiet_NaN();
1176  // Read in proper time from config
1177  proper_time =
1179  validate_duplicate_IC_config(proper_time, IC_parameters.proper_time,
1180  "Proper_Time");
1181  } else if (IC_parameters.proper_time.has_value()) {
1182  proper_time = IC_parameters.proper_time.value();
1183  } else {
1184  double lower_bound =
1186  if (IC_parameters.lower_bound.has_value())
1187  validate_duplicate_IC_config(lower_bound, IC_parameters.lower_bound,
1188  "Lower_Bound");
1189 
1190  // Default proper time is the passing time of the two nuclei
1191  double default_proper_time = modus_.nuclei_passing_time();
1192  if (default_proper_time >= lower_bound) {
1193  proper_time = default_proper_time;
1194  } else {
1195  logg[LInitialConditions].warn()
1196  << "Nuclei passing time is too short, hypersurface proper time "
1197  << "set to tau = " << lower_bound << " fm.";
1198  proper_time = lower_bound;
1199  }
1200  }
1201 
1202  action_finders_.emplace_back(
1203  std::make_unique<HyperSurfaceCrossActionsFinder>(
1204  proper_time, rapidity_cut, pT_cut));
1205  }
1206  }
1207 
1209  logg[LExperiment].info() << "Pauli blocking is ON.";
1210  pauli_blocker_ = std::make_unique<PauliBlocker>(
1213  parameters_);
1214  }
1215 
1467  // create outputs
1469  " create OutputInterface objects");
1470  dens_type_ = config.take(InputKeys::output_densityType);
1471  logg[LExperiment].debug()
1472  << "Density type printed to headers: " << dens_type_;
1473 
1474  /* Parse configuration about output contents and formats, doing all logical
1475  * checks about specified formats, creating all needed output objects. Note
1476  * that we first extract the output sub configuration without the "Output:"
1477  * enclosing section to easily get all output contents and then we reintroduce
1478  * it for the actual parsing (remember we parse database keys which have
1479  * labels from the top-level only).
1480  */
1481  auto output_conf = config.extract_sub_configuration(
1483  if (output_path == "") {
1484  throw std::invalid_argument(
1485  "Invalid empty output path provided to Experiment constructor.");
1486  } else if (!std::filesystem::exists(output_path)) {
1487  logg[LExperiment].fatal(
1488  "Output path \"" + output_path.string() +
1489  "\" used to create an Experiment object does not exist.");
1490  throw NonExistingOutputPathRequest("Attempt to use not existing path.");
1491  } else if (!std::filesystem::is_directory(output_path)) {
1492  logg[LExperiment].fatal("Output path \"" + output_path.string() +
1493  "\" used to create an Experiment object "
1494  "exists, but it is not a directory.");
1495  throw std::logic_error("Attempt to use invalid existing path.");
1496  }
1497  const std::vector<std::string> output_contents =
1498  output_conf.list_upmost_nodes();
1499  if (output_conf.is_empty()) {
1500  logg[LExperiment].warn() << "No \"Output\" section found in the input "
1501  "file. No output file will be produced.";
1502  } else {
1503  output_conf.enclose_into_section(InputSections::output);
1504  }
1505  std::vector<std::vector<std::string>> list_of_formats(output_contents.size());
1506  std::transform(
1507  output_contents.cbegin(), output_contents.cend(), list_of_formats.begin(),
1508  [&output_conf](std::string content) -> std::vector<std::string> {
1509  /* Note that the "Format" key has an empty list as default, although it
1510  * is a required key, because then here below the error for the user is
1511  * more informative, if the key was not given in the input file. */
1512  return output_conf.take(InputKeys::get_output_format_key(content));
1513  });
1514  auto abort_because_of_invalid_input_file = []() {
1515  throw std::invalid_argument("Invalid configuration input file.");
1516  };
1517  const OutputParameters output_parameters(std::move(output_conf));
1518  for (std::size_t i = 0; i < output_contents.size(); ++i) {
1519  if (output_contents[i] == "Particles" ||
1520  output_contents[i] == "Collisions") {
1521  assert(output_parameters.quantities.count(output_contents[i]) > 0);
1522  const bool quantities_given_nonempty =
1523  !output_parameters.quantities.at(output_contents[i]).empty();
1524  auto formats_contains = [&list_of_formats, &i](const std::string &label) {
1525  return std::find(list_of_formats[i].begin(), list_of_formats[i].end(),
1526  label) != list_of_formats[i].end();
1527  };
1528  const bool custom_ascii_requested = formats_contains("ASCII");
1529  const bool custom_binary_requested = formats_contains("Binary");
1530  const bool custom_requested =
1531  custom_ascii_requested || custom_binary_requested;
1532  const bool oscar2013_requested = formats_contains("Oscar2013");
1533  const bool oscar2013_bin_requested = formats_contains("Oscar2013_bin");
1534  const bool is_extended = (output_contents[i] == "Particles")
1535  ? output_parameters.part_extended
1536  : output_parameters.coll_extended;
1537  const auto &default_quantities =
1540  const bool are_given_quantities_oscar2013_ones =
1541  output_parameters.quantities.at(output_contents[i]) ==
1542  default_quantities;
1543  if (quantities_given_nonempty != custom_requested) {
1544  logg[LExperiment].fatal()
1545  << "Non-empty \"Quantities\" and \"ASCII\"/\"Binary\" format have "
1546  << "not been specified both for " << std::quoted(output_contents[i])
1547  << " in config file.";
1548  abort_because_of_invalid_input_file();
1549  }
1550  if (custom_ascii_requested && oscar2013_requested &&
1551  are_given_quantities_oscar2013_ones) {
1552  logg[LExperiment].fatal()
1553  << "The specified \"Quantities\" for the ASCII format are the same "
1554  "as those of the requested \"Oscar2013\"\nformat for "
1555  << std::quoted(output_contents[i])
1556  << " and this would produce the same output file twice.";
1557  abort_because_of_invalid_input_file();
1558  }
1559  if (custom_binary_requested && oscar2013_bin_requested &&
1560  are_given_quantities_oscar2013_ones) {
1561  logg[LExperiment].fatal()
1562  << "The specified \"Quantities\" for the binary format are the "
1563  "same as those of the requested \"Oscar2013_bin\"\nformat for "
1564  << std::quoted(output_contents[i])
1565  << " and this would produce the same output file twice.";
1566  abort_because_of_invalid_input_file();
1567  }
1568  }
1569 
1570  if (list_of_formats[i].empty()) {
1571  logg[LExperiment].fatal()
1572  << "Empty or unspecified list of formats for "
1573  << std::quoted(output_contents[i]) << " content.";
1574  abort_because_of_invalid_input_file();
1575  } else if (std::find(list_of_formats[i].begin(), list_of_formats[i].end(),
1576  "None") != list_of_formats[i].end()) {
1577  if (list_of_formats[i].size() > 1) {
1578  logg[LExperiment].fatal()
1579  << "Use of \"None\" output format together with other formats is "
1580  "not allowed.\nInvalid \"Format\" key for "
1581  << std::quoted(output_contents[i]) << " content.";
1582  abort_because_of_invalid_input_file();
1583  } else {
1584  // Clear vector so that the for below is skipped and no output created
1585  list_of_formats[i].clear();
1586  }
1587  } else if (std::set<std::string> tmp_set(list_of_formats[i].begin(),
1588  list_of_formats[i].end());
1589  list_of_formats[i].size() != tmp_set.size()) {
1590  auto join_container = [](const auto &container) {
1591  std::string result{};
1592  std::for_each(container.cbegin(), container.cend(),
1593  [&result](const std::string s) {
1594  result += (result == "") ? s : ", " + s;
1595  });
1596  return result;
1597  };
1598  const std::string old_formats = join_container(list_of_formats[i]),
1599  new_formats = join_container(tmp_set);
1600  logg[LExperiment].warn()
1601  << "Found the same output format multiple times for "
1602  << std::quoted(output_contents[i])
1603  << " content. Duplicates will be ignored:\n 'Format: [" << old_formats
1604  << "] -> [" << new_formats << "]'";
1605  list_of_formats[i].assign(tmp_set.begin(), tmp_set.end());
1606  }
1607  }
1608 
1609  /* Repeat loop over output_contents here to create all outputs after having
1610  * validated all content specifications. This is more user-friendly. */
1611  std::size_t total_number_of_requested_formats = 0;
1612  for (std::size_t i = 0; i < output_contents.size(); ++i) {
1613  for (const auto &format : list_of_formats[i]) {
1614  create_output(format, output_contents[i], output_path, output_parameters);
1615  ++total_number_of_requested_formats;
1616  }
1617  }
1618 
1619  if (outputs_.size() != total_number_of_requested_formats) {
1620  logg[LExperiment].fatal()
1621  << "At least one invalid output format has been provided.";
1622  abort_because_of_invalid_input_file();
1623  }
1624 
1625  /* We can take away the Fermi motion flag, because the collider modus is
1626  * already initialized. We only need it when potentials are enabled, but we
1627  * always have to take it, otherwise SMASH will complain about unused
1628  * options. We have to provide a default value for modi other than Collider.
1629  */
1630  if (config.has_section(InputSections::potentials)) {
1631  if (time_step_mode_ == TimeStepMode::None) {
1632  logg[LExperiment].error() << "Potentials only work with time steps!";
1633  throw std::invalid_argument("Can't use potentials without time steps!");
1634  }
1635  if (modus_.fermi_motion() == FermiMotion::Frozen) {
1636  logg[LExperiment].error()
1637  << "Potentials don't work with frozen Fermi momenta! "
1638  "Use normal Fermi motion instead.";
1639  throw std::invalid_argument(
1640  "Can't use potentials "
1641  "with frozen Fermi momenta!");
1642  }
1643  logg[LExperiment].info() << "Potentials are ON. Timestep is "
1644  << parameters_.labclock->timestep_duration();
1645  // potentials need density calculation parameters from parameters_
1646  potentials_ = std::make_unique<Potentials>(
1648  parameters_);
1649  // make sure that vdf potentials are not used together with Skyrme
1650  // or symmetry potentials
1651  if (potentials_->use_skyrme() && potentials_->use_vdf()) {
1652  throw std::runtime_error(
1653  "Can't use Skyrme and VDF potentials at the same time!");
1654  }
1655  if (potentials_->use_symmetry() && potentials_->use_vdf()) {
1656  throw std::runtime_error(
1657  "Can't use symmetry and VDF potentials at the same time!");
1658  }
1659  if (potentials_->use_skyrme()) {
1660  logg[LExperiment].info() << "Skyrme potentials are:\n";
1661  logg[LExperiment].info()
1662  << "\t\tSkyrme_A [MeV] = " << potentials_->skyrme_a() << "\n";
1663  logg[LExperiment].info()
1664  << "\t\tSkyrme_B [MeV] = " << potentials_->skyrme_b() << "\n";
1665  logg[LExperiment].info()
1666  << "\t\t Skyrme_tau = " << potentials_->skyrme_tau() << "\n";
1667  }
1668  if (potentials_->use_symmetry()) {
1669  logg[LExperiment].info()
1670  << "Symmetry potential is:"
1671  << "\n S_pot [MeV] = " << potentials_->symmetry_S_pot() << "\n";
1672  }
1673  if (potentials_->use_vdf()) {
1674  logg[LExperiment].info() << "VDF potential parameters are:\n";
1675  logg[LExperiment].info() << "\t\tsaturation density [fm^-3] = "
1676  << potentials_->saturation_density() << "\n";
1677  for (int i = 0; i < potentials_->number_of_terms(); i++) {
1678  logg[LExperiment].info()
1679  << "\t\tCoefficient_" << i + 1 << " = "
1680  << 1000.0 * (potentials_->coeffs())[i] << " [MeV] \t Power_"
1681  << i + 1 << " = " << (potentials_->powers())[i] << "\n";
1682  }
1683  }
1684  // if potentials are on, derivatives need to be calculated
1685  if (parameters_.derivatives_mode == DerivativesMode::Off &&
1686  parameters_.field_derivatives_mode == FieldDerivativesMode::ChainRule) {
1687  throw std::invalid_argument(
1688  "Derivatives are necessary for running with potentials.\n"
1689  "Derivatives_Mode: \"Off\" only makes sense for "
1690  "Field_Derivatives_Mode: \"Direct\"!\nUse \"Covariant Gaussian\" or "
1691  "\"Finite difference\".");
1692  }
1693  // for computational efficiency, we want to turn off the derivatives of jmu
1694  // and the rest frame density derivatives if direct derivatives are used
1695  if (parameters_.field_derivatives_mode == FieldDerivativesMode::Direct) {
1696  parameters_.derivatives_mode = DerivativesMode::Off;
1697  parameters_.rho_derivatives_mode = RestFrameDensityDerivativesMode::Off;
1698  }
1699  switch (parameters_.derivatives_mode) {
1701  logg[LExperiment].info() << "Covariant Gaussian derivatives are ON";
1702  break;
1704  logg[LExperiment].info() << "Finite difference derivatives are ON";
1705  break;
1706  case DerivativesMode::Off:
1707  logg[LExperiment].info() << "Gradients of baryon current are OFF";
1708  break;
1709  }
1710  switch (parameters_.rho_derivatives_mode) {
1712  logg[LExperiment].info() << "Rest frame density derivatives are ON";
1713  break;
1715  logg[LExperiment].info() << "Rest frame density derivatives are OFF";
1716  break;
1717  }
1718  // direct or chain rule derivatives only make sense for the VDF potentials
1719  if (potentials_->use_vdf()) {
1720  switch (parameters_.field_derivatives_mode) {
1722  logg[LExperiment].info() << "Chain rule field derivatives are ON";
1723  break;
1725  logg[LExperiment].info() << "Direct field derivatives are ON";
1726  break;
1727  }
1728  }
1729  /*
1730  * Necessary safety checks
1731  */
1732  // VDF potentials need derivatives of rest frame density or fields
1733  if (potentials_->use_vdf() && (parameters_.rho_derivatives_mode ==
1735  parameters_.field_derivatives_mode ==
1737  throw std::runtime_error(
1738  "Can't use VDF potentials without rest frame density derivatives or "
1739  "direct field derivatives!");
1740  }
1741  // potentials require using gradients
1742  if (parameters_.derivatives_mode == DerivativesMode::Off &&
1743  parameters_.field_derivatives_mode == FieldDerivativesMode::ChainRule) {
1744  throw std::runtime_error(
1745  "Can't use potentials without gradients of baryon current (Skyrme, "
1746  "VDF)"
1747  " or direct field derivatives (VDF)!");
1748  }
1749  // direct field derivatives only make sense for the VDF potentials
1750  if (!(potentials_->use_vdf()) &&
1751  parameters_.field_derivatives_mode == FieldDerivativesMode::Direct) {
1752  throw std::invalid_argument(
1753  "Field_Derivatives_Mode: \"Direct\" only makes sense for the VDF "
1754  "potentials!\nUse Field_Derivatives_Mode: \"Chain Rule\" or comment "
1755  "this option out (Chain Rule is default)");
1756  }
1757  }
1758 
1759  // information about the type of smearing
1760  switch (parameters_.smearing_mode) {
1762  logg[LExperiment].info() << "Smearing type: Covariant Gaussian";
1763  break;
1765  logg[LExperiment].info() << "Smearing type: Discrete with weight = "
1766  << parameters_.discrete_weight;
1767  break;
1769  logg[LExperiment].info() << "Smearing type: Triangular with range = "
1770  << parameters_.triangular_range;
1771  break;
1772  }
1773 
1774  // Create lattices
1775  if (config.has_section(InputSections::lattice)) {
1776  bool automatic = config.take(InputKeys::lattice_automatic);
1777  bool all_geometrical_properties_specified =
1781  if (!automatic && !all_geometrical_properties_specified) {
1782  throw std::invalid_argument(
1783  "The lattice was requested to be manually generated, but some\n"
1784  "lattice geometrical property was not specified. Be sure to provide\n"
1785  "both \"Cell_Number\" and \"Origin\" and \"Sizes\".");
1786  }
1787  if (automatic && all_geometrical_properties_specified) {
1788  throw std::invalid_argument(
1789  "The lattice was requested to be automatically generated, but all\n"
1790  "lattice geometrical properties were specified. In this case you\n"
1791  "need to set \"Automatic: False\".");
1792  }
1793  bool periodic = config.take(InputKeys::lattice_periodic, modus_.is_box());
1794  const auto [l, n, origin] = [&config, automatic, this]() {
1795  if (!automatic) {
1796  return std::make_tuple<std::array<double, 3>, std::array<int, 3>,
1797  std::array<double, 3>>(
1801  } else {
1802  std::array<double, 3> l_default{20., 20., 20.};
1803  std::array<int, 3> n_default{10, 10, 10};
1804  std::array<double, 3> origin_default{-20., -20., -20.};
1805  if (modus_.is_collider() || (modus_.is_list() && !modus_.is_box())) {
1806  // Estimates on how far particles could get in x, y, z. The
1807  // default lattice is currently not contracted for afterburner runs
1808  const double gam = modus_.is_collider()
1809  ? modus_.sqrt_s_NN() / (2.0 * nucleon_mass)
1810  : 1.0;
1811  const double max_z = 5.0 / gam + end_time_;
1812  const double estimated_max_transverse_velocity = 0.7;
1813  const double max_xy =
1814  5.0 + estimated_max_transverse_velocity * end_time_;
1815  origin_default = {-max_xy, -max_xy, -max_z};
1816  l_default = {2 * max_xy, 2 * max_xy, 2 * max_z};
1817  // Go for approximately 0.8 fm cell size and contract
1818  // lattice in z by gamma factor
1819  const int n_xy = numeric_cast<int>(std::ceil(2 * max_xy / 0.8));
1820  int nz = numeric_cast<int>(std::ceil(2 * max_z / 0.8));
1821  // Contract lattice by gamma factor in case of smearing where
1822  // smearing length is bound to the lattice cell length
1823  if (parameters_.smearing_mode == SmearingMode::Discrete ||
1824  parameters_.smearing_mode == SmearingMode::Triangular) {
1825  nz = numeric_cast<int>(std::ceil(2 * max_z / 0.8 * gam));
1826  }
1827  n_default = {n_xy, n_xy, nz};
1828  } else if (modus_.is_box()) {
1829  origin_default = {0., 0., 0.};
1830  const double bl = modus_.length();
1831  l_default = {bl, bl, bl};
1832  const int n_xyz = numeric_cast<int>(std::ceil(bl / 0.5));
1833  n_default = {n_xyz, n_xyz, n_xyz};
1834  } else if (modus_.is_sphere()) {
1835  // Maximal distance from (0, 0, 0) at which a particle
1836  // may be found at the end of the simulation
1837  const double max_d = modus_.radius() + end_time_;
1838  origin_default = {-max_d, -max_d, -max_d};
1839  l_default = {2 * max_d, 2 * max_d, 2 * max_d};
1840  // Go for approximately 0.8 fm cell size
1841  const int n_xyz = numeric_cast<int>(std::ceil(2 * max_d / 0.8));
1842  n_default = {n_xyz, n_xyz, n_xyz};
1843  }
1844  // Take lattice properties from config to assign them to all lattices
1845  return std::make_tuple<std::array<double, 3>, std::array<int, 3>,
1846  std::array<double, 3>>(
1847  config.take(InputKeys::lattice_sizes, l_default),
1848  config.take(InputKeys::lattice_cellNumber, n_default),
1849  config.take(InputKeys::lattice_origin, origin_default));
1850  }
1851  }();
1852 
1853  logg[LExperiment].info()
1854  << "Lattice is ON. Origin = (" << origin[0] << "," << origin[1] << ","
1855  << origin[2] << "), sizes = (" << l[0] << "," << l[1] << "," << l[2]
1856  << "), number of cells = (" << n[0] << "," << n[1] << "," << n[2]
1857  << "), periodic = " << std::boolalpha << periodic;
1858 
1859  if (printout_lattice_td_ || printout_full_lattice_any_td_) {
1860  dens_type_lattice_printout_ = output_parameters.td_dens_type;
1861  printout_rho_eckart_ = output_parameters.td_rho_eckart;
1862  printout_tmn_ = output_parameters.td_tmn;
1863  printout_tmn_landau_ = output_parameters.td_tmn_landau;
1864  printout_v_landau_ = output_parameters.td_v_landau;
1865  printout_j_QBS_ = output_parameters.td_jQBS;
1866  }
1867  if (printout_tmn_ || printout_tmn_landau_ || printout_v_landau_) {
1868  Tmn_ = std::make_unique<RectangularLattice<EnergyMomentumTensor>>(
1869  l, n, origin, periodic, LatticeUpdate::AtOutput);
1870  }
1871  if (printout_j_QBS_) {
1872  j_QBS_lat_ = std::make_unique<DensityLattice>(l, n, origin, periodic,
1874  }
1875  /* Create baryon and isospin density lattices regardless of config
1876  if potentials are on. This is because they allow to compute
1877  potentials faster */
1878  if (potentials_) {
1879  // Create auxiliary lattices for baryon four-current calculation
1880  old_jmu_auxiliary_ = std::make_unique<RectangularLattice<FourVector>>(
1881  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1882  new_jmu_auxiliary_ = std::make_unique<RectangularLattice<FourVector>>(
1883  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1884  four_gradient_auxiliary_ =
1885  std::make_unique<RectangularLattice<std::array<FourVector, 4>>>(
1886  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1887 
1888  if (potentials_->use_skyrme()) {
1889  jmu_B_lat_ = std::make_unique<DensityLattice>(
1890  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1891  UB_lat_ = std::make_unique<RectangularLattice<FourVector>>(
1892  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1893  FB_lat_ = std::make_unique<
1894  RectangularLattice<std::pair<ThreeVector, ThreeVector>>>(
1895  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1896  }
1897  if (potentials_->use_symmetry()) {
1898  jmu_I3_lat_ = std::make_unique<DensityLattice>(
1899  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1900  UI3_lat_ = std::make_unique<RectangularLattice<FourVector>>(
1901  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1902  FI3_lat_ = std::make_unique<
1903  RectangularLattice<std::pair<ThreeVector, ThreeVector>>>(
1904  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1905  }
1906  if (potentials_->use_coulomb()) {
1907  jmu_el_lat_ = std::make_unique<DensityLattice>(
1908  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1909  EM_lat_ = std::make_unique<
1910  RectangularLattice<std::pair<ThreeVector, ThreeVector>>>(
1911  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1912  }
1913  if (potentials_->use_vdf()) {
1914  jmu_B_lat_ = std::make_unique<DensityLattice>(
1915  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1916  UB_lat_ = std::make_unique<RectangularLattice<FourVector>>(
1917  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1918  FB_lat_ = std::make_unique<
1919  RectangularLattice<std::pair<ThreeVector, ThreeVector>>>(
1920  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1921  }
1922  if (parameters_.field_derivatives_mode == FieldDerivativesMode::Direct) {
1923  // Create auxiliary lattices for field calculation
1924  old_fields_auxiliary_ =
1925  std::make_unique<RectangularLattice<FourVector>>(
1926  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1927  new_fields_auxiliary_ =
1928  std::make_unique<RectangularLattice<FourVector>>(
1929  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1930  fields_four_gradient_auxiliary_ =
1931  std::make_unique<RectangularLattice<std::array<FourVector, 4>>>(
1932  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1933 
1934  // Create the fields lattice
1935  fields_lat_ = std::make_unique<FieldsLattice>(
1936  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1937  }
1938  }
1939  if (dens_type_lattice_printout_ == DensityType::Baryon && !jmu_B_lat_) {
1940  jmu_B_lat_ = std::make_unique<DensityLattice>(l, n, origin, periodic,
1942  }
1943  if (dens_type_lattice_printout_ == DensityType::BaryonicIsospin &&
1944  !jmu_I3_lat_) {
1945  jmu_I3_lat_ = std::make_unique<DensityLattice>(l, n, origin, periodic,
1947  }
1948  if (dens_type_lattice_printout_ != DensityType::None &&
1949  dens_type_lattice_printout_ != DensityType::BaryonicIsospin &&
1950  dens_type_lattice_printout_ != DensityType::Baryon) {
1951  jmu_custom_lat_ = std::make_unique<DensityLattice>(
1952  l, n, origin, periodic, LatticeUpdate::AtOutput);
1953  }
1954  } else if (printout_lattice_td_ || printout_full_lattice_any_td_) {
1955  logg[LExperiment].error(
1956  "If you want Therm. VTK or Lattice output, configure a lattice for "
1957  "it.");
1958  } else if (potentials_ && potentials_->use_coulomb()) {
1959  logg[LExperiment].error(
1960  "Coulomb potential requires a lattice. Please add one to the "
1961  "configuration");
1962  }
1963 
1964  // Warning for the mean field calculation if lattice is not on.
1965  if ((potentials_ != nullptr) && (jmu_B_lat_ == nullptr)) {
1966  logg[LExperiment].warn() << "Lattice is NOT used. Mean-field energy is "
1967  << "not going to be calculated.";
1968  }
1969 
1970  // Store pointers to potential and lattice accessible for Action
1971  if (parameters_.potential_affect_threshold) {
1972  UB_lat_pointer = UB_lat_.get();
1973  UI3_lat_pointer = UI3_lat_.get();
1974  pot_pointer = potentials_.get();
1975  }
1976 
1977  // Throw fatal if DerivativesMode == FiniteDifference and lattice is not on.
1978  if ((parameters_.derivatives_mode == DerivativesMode::FiniteDifference) &&
1979  (jmu_B_lat_ == nullptr)) {
1980  throw std::runtime_error(
1981  "Lattice is necessary to calculate finite difference gradients.");
1982  }
1983 
1984  // Create forced thermalizer
1986  Configuration th_conf = config.extract_complete_sub_configuration(
1988  thermalizer_ = modus_.create_grandcan_thermalizer(th_conf);
1989  }
1990 
1991  /* Take the seed setting only after the configuration was stored to a file
1992  * in smash.cc */
1993  seed_ = config.take(InputKeys::gen_randomseed);
1994 }
1995 
1997 const std::string hline(113, '-');
1998 
2023 std::string format_measurements(const std::vector<Particles> &ensembles,
2024  uint64_t scatterings_this_interval,
2025  const QuantumNumbers &conserved_initial,
2026  SystemTimePoint time_start, double time,
2027  double E_mean_field,
2028  double E_mean_field_initial);
2043  const Potentials &potentials,
2045  RectangularLattice<std::pair<ThreeVector, ThreeVector>> *em_lattice,
2046  const ExperimentParameters &parameters);
2047 
2064 EventInfo fill_event_info(const std::vector<Particles> &ensembles,
2065  double E_mean_field, double modus_impact_parameter,
2066  const ExperimentParameters &parameters,
2067  bool projectile_target_interact,
2068  bool kinematic_cut_for_SMASH_IC);
2069 
2070 template <typename Modus>
2072  random::set_seed(seed_);
2073  logg[LExperiment].info() << "random number seed: " << seed_;
2074  /* Set seed for the next event. It has to be positive, so it can be entered
2075  * in the config.
2076  *
2077  * We have to be careful about the minimal integer, whose absolute value
2078  * cannot be represented. */
2079  int64_t r = random::advance();
2080  while (r == INT64_MIN) {
2081  r = random::advance();
2082  }
2083  seed_ = std::abs(r);
2084  /* Set the random seed used in PYTHIA hadronization
2085  * to be same with the SMASH one.
2086  * In this way we ensure that the results are reproducible
2087  * for every event if one knows SMASH random seed. */
2088  if (process_string_ptr_ != NULL) {
2089  process_string_ptr_->init_pythia_hadron_rndm();
2090  }
2091 
2092  for (Particles &particles : ensembles_) {
2093  particles.reset();
2094  }
2095 
2096  // Sample particles according to the initial conditions
2097  double start_time = -1.0;
2098 
2099  // Sample impact parameter only once per all ensembles
2100  // It should be the same for all ensembles
2101  if (modus_.is_collider()) {
2102  modus_.sample_impact();
2103  logg[LExperiment].info("Impact parameter = ", modus_.impact_parameter(),
2104  " fm");
2105  }
2106  for (Particles &particles : ensembles_) {
2107  start_time = modus_.initial_conditions(&particles, parameters_);
2108  }
2109  /* For box modus make sure that particles are in the box. In principle, after
2110  * a correct initialization they should be, so this is just playing it safe.
2111  */
2112  for (Particles &particles : ensembles_) {
2113  modus_.impose_boundary_conditions(&particles, outputs_);
2114  }
2115  // Reset the simulation clock
2116  double timestep = delta_time_startup_;
2117 
2118  switch (time_step_mode_) {
2119  case TimeStepMode::Fixed:
2120  break;
2121  case TimeStepMode::None:
2122  timestep = end_time_ - start_time;
2123  // Take care of the box modus + timestepless propagation
2124  const double max_dt = modus_.max_timestep(max_transverse_distance_sqr_);
2125  if (max_dt > 0. && max_dt < timestep) {
2126  timestep = max_dt;
2127  }
2128  break;
2129  }
2130  std::unique_ptr<UniformClock> clock_for_this_event;
2131  if (modus_.is_list() && (timestep < 0.0)) {
2132  throw std::runtime_error(
2133  "Timestep for the given event is negative. \n"
2134  "This might happen if the formation times of the input particles are "
2135  "larger than the specified end time of the simulation.");
2136  }
2137  clock_for_this_event =
2138  std::make_unique<UniformClock>(start_time, timestep, end_time_);
2139  parameters_.labclock = std::move(clock_for_this_event);
2140 
2141  // Reset the output clock
2142  parameters_.outputclock->reset(start_time, true);
2143  // remove time before starting time in case of custom output times.
2144  parameters_.outputclock->remove_times_in_past(start_time);
2145 
2146  logg[LExperiment].debug(
2147  "Lab clock: t_start = ", parameters_.labclock->current_time(),
2148  ", dt = ", parameters_.labclock->timestep_duration());
2149 
2150  /* Save the initial conserved quantum numbers and total momentum in
2151  * the system for conservation checks */
2152  conserved_initial_ = QuantumNumbers(ensembles_);
2153  wall_actions_total_ = 0;
2154  previous_wall_actions_total_ = 0;
2155  interactions_total_ = 0;
2156  previous_interactions_total_ = 0;
2157  discarded_interactions_total_ = 0;
2158  total_pauli_blocked_ = 0;
2159  projectile_target_interact_.assign(parameters_.n_ensembles, false);
2160  total_hypersurface_crossing_actions_ = 0;
2161  total_energy_removed_ = 0.0;
2162  total_energy_violated_by_Pythia_ = 0.0;
2163  // Print output headers
2164  logg[LExperiment].info() << hline;
2165  logg[LExperiment].info() << "Time[fm] Ekin[GeV] E_MF[GeV] ETotal[GeV] "
2166  << "ETot/N[GeV] D(ETot/N)[GeV] Scatt&Decays "
2167  << "Particles Comp.Time";
2168  logg[LExperiment].info() << hline;
2169  double E_mean_field = 0.0;
2170  if (potentials_) {
2171  // update_potentials();
2172  // if (parameters.outputclock->current_time() == 0.0 )
2173  // using the lattice is necessary
2174  if ((jmu_B_lat_ != nullptr)) {
2175  update_lattice(jmu_B_lat_.get(), old_jmu_auxiliary_.get(),
2176  new_jmu_auxiliary_.get(), four_gradient_auxiliary_.get(),
2178  density_param_, ensembles_,
2179  parameters_.labclock->timestep_duration(), true);
2180  // Because there was no lattice at t=-Delta_t, the time derivatives
2181  // drho_dt and dj^mu/dt at t=0 are huge, while they shouldn't be; we
2182  // overwrite the time derivative to zero by hand.
2183  for (auto &node : *jmu_B_lat_) {
2184  node.overwrite_drho_dt_to_zero();
2185  node.overwrite_djmu_dt_to_zero();
2186  }
2187  E_mean_field = calculate_mean_field_energy(*potentials_, *jmu_B_lat_,
2188  EM_lat_.get(), parameters_);
2189  }
2190  }
2191  initial_mean_field_energy_ = E_mean_field;
2193  ensembles_, 0u, conserved_initial_, time_start_,
2194  parameters_.labclock->current_time(), E_mean_field,
2195  initial_mean_field_energy_);
2196 
2197  // Output at event start
2198  for (const auto &output : outputs_) {
2199  for (int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
2200  auto event_info = fill_event_info(
2201  ensembles_, E_mean_field, modus_.impact_parameter(), parameters_,
2202  projectile_target_interact_[i_ens], kinematic_cuts_for_IC_output_);
2203  output->at_eventstart(ensembles_[i_ens], {event_, i_ens}, event_info);
2204  }
2205  // For thermodynamic output
2206  output->at_eventstart(ensembles_, event_);
2207  // For thermodynamic lattice output
2208  if (printout_full_lattice_any_td_) {
2209  if (printout_rho_eckart_) {
2210  switch (dens_type_lattice_printout_) {
2211  case DensityType::Baryon:
2212  output->at_eventstart(event_, ThermodynamicQuantity::EckartDensity,
2213  DensityType::Baryon, *jmu_B_lat_);
2214  break;
2216  output->at_eventstart(event_, ThermodynamicQuantity::EckartDensity,
2217  DensityType::BaryonicIsospin, *jmu_I3_lat_);
2218  break;
2219  case DensityType::None:
2220  break;
2221  default:
2222  output->at_eventstart(event_, ThermodynamicQuantity::EckartDensity,
2224  *jmu_custom_lat_);
2225  }
2226  }
2227  if (printout_tmn_) {
2228  output->at_eventstart(event_, ThermodynamicQuantity::Tmn,
2229  dens_type_lattice_printout_, *Tmn_);
2230  }
2231  if (printout_tmn_landau_) {
2232  output->at_eventstart(event_, ThermodynamicQuantity::TmnLandau,
2233  dens_type_lattice_printout_, *Tmn_);
2234  }
2235  if (printout_v_landau_) {
2236  output->at_eventstart(event_, ThermodynamicQuantity::LandauVelocity,
2237  dens_type_lattice_printout_, *Tmn_);
2238  }
2239  if (printout_j_QBS_) {
2240  output->at_eventstart(event_, ThermodynamicQuantity::j_QBS,
2241  dens_type_lattice_printout_, *j_QBS_lat_);
2242  }
2243  }
2244  }
2245 
2246  /* In the ColliderModus, if Fermi motion is frozen, assign the beam momenta
2247  * to the nucleons in both the projectile and the target. Every ensemble
2248  * gets the same beam momenta, so no need to create beam_momenta_ vector
2249  * for every ensemble.
2250  */
2251  if (modus_.is_collider() && modus_.fermi_motion() == FermiMotion::Frozen) {
2252  for (ParticleData &particle : ensembles_[0]) {
2253  const double m = particle.effective_mass();
2254  double v_beam = 0.0;
2255  if (particle.belongs_to() == BelongsTo::Projectile) {
2256  v_beam = modus_.velocity_projectile();
2257  } else if (particle.belongs_to() == BelongsTo::Target) {
2258  v_beam = modus_.velocity_target();
2259  }
2260  const double gamma = 1.0 / std::sqrt(1.0 - v_beam * v_beam);
2261  beam_momentum_.emplace_back(
2262  FourVector(gamma * m, 0.0, 0.0, gamma * v_beam * m));
2263  } // loop over particles
2264  }
2265 }
2266 
2267 template <typename Modus>
2268 bool Experiment<Modus>::perform_action(Action &action, int i_ensemble,
2269  bool include_pauli_blocking) {
2270  Particles &particles = ensembles_[i_ensemble];
2271  // Make sure to skip invalid and Pauli-blocked actions.
2272  if (!action.is_valid(particles)) {
2273  discarded_interactions_total_++;
2274  logg[LExperiment].debug(~einhard::DRed(), "✘ ", action,
2275  " (discarded: invalid)");
2276  return false;
2277  }
2278  try {
2279  action.generate_final_state();
2280  } catch (Action::StochasticBelowEnergyThreshold &) {
2281  return false;
2282  }
2283  logg[LExperiment].debug("Process Type is: ", action.get_type());
2284  if (include_pauli_blocking && pauli_blocker_ &&
2285  action.is_pauli_blocked(ensembles_, *pauli_blocker_)) {
2286  total_pauli_blocked_++;
2287  return false;
2288  }
2289 
2290  // Prepare projectile_target_interact_, it's used for output
2291  // to signal that there was some interaction in this event
2292  if (modus_.is_collider()) {
2293  int count_target = 0, count_projectile = 0;
2294  for (const auto &incoming : action.incoming_particles()) {
2295  if (incoming.belongs_to() == BelongsTo::Projectile) {
2296  count_projectile++;
2297  } else if (incoming.belongs_to() == BelongsTo::Target) {
2298  count_target++;
2299  }
2300  }
2301  if (count_target > 0 && count_projectile > 0) {
2302  projectile_target_interact_[i_ensemble] = true;
2303  }
2304  }
2305 
2306  /* Make sure to pick a non-zero integer, because 0 is reserved for "no
2307  * interaction yet". */
2308  const auto id_process = static_cast<uint32_t>(interactions_total_ + 1);
2309  // we perform the action and collect possible energy violations by Pythia
2310  total_energy_violated_by_Pythia_ += action.perform(&particles, id_process);
2311 
2312  interactions_total_++;
2313  if (action.get_type() == ProcessType::Wall) {
2314  wall_actions_total_++;
2315  }
2316  if (action.get_type() == ProcessType::HyperSurfaceCrossing) {
2317  total_hypersurface_crossing_actions_++;
2318  total_energy_removed_ += action.incoming_particles()[0].momentum().x0();
2319  }
2320  // Calculate Eckart rest frame density at the interaction point
2321  double rho = 0.0;
2322  if (dens_type_ != DensityType::None) {
2323  const FourVector r_interaction = action.get_interaction_point();
2324  constexpr bool compute_grad = false;
2325  const bool smearing = true;
2326  // todo(oliiny): it's a rough density estimate from a single ensemble.
2327  // It might actually be appropriate for output. Discuss.
2328  rho = std::get<0>(current_eckart(r_interaction.threevec(), particles,
2329  density_param_, dens_type_, compute_grad,
2330  smearing));
2331  }
2347  for (const auto &output : outputs_) {
2348  if (!output->is_dilepton_output() && !output->is_photon_output()) {
2349  if (output->is_IC_output() &&
2351  output->at_interaction(action, rho);
2352  } else if (!output->is_IC_output()) {
2353  output->at_interaction(action, rho);
2354  }
2355  }
2356  }
2357 
2358  // At every collision photons can be produced.
2359  // Note: We rely here on the lazy evaluation of the arguments to if.
2360  // It may happen that in a wall-crossing-action sqrt_s raises an exception.
2361  // Therefore we first have to check if the incoming particles can undergo
2362  // an em-interaction.
2363  if (photons_switch_ &&
2366  action.sqrt_s(), action.incoming_particles())) {
2367  /* Time in the action constructor is relative to
2368  * current time of incoming */
2369  constexpr double action_time = 0.;
2370  ScatterActionPhoton photon_act(action.incoming_particles(), action_time,
2371  n_fractional_photons_,
2372  action.get_total_weight());
2373 
2384  photon_act.add_dummy_hadronic_process(action.get_total_weight());
2385 
2386  // Now add the actual photon reaction channel.
2387  photon_act.add_single_process();
2388 
2389  photon_act.perform_photons(outputs_);
2390  }
2391 
2392  if (bremsstrahlung_switch_ &&
2394  action.incoming_particles())) {
2395  /* Time in the action constructor is relative to
2396  * current time of incoming */
2397  constexpr double action_time = 0.;
2398 
2399  BremsstrahlungAction brems_act(action.incoming_particles(), action_time,
2400  n_fractional_photons_,
2401  action.get_total_weight());
2402 
2414  brems_act.add_dummy_hadronic_process(action.get_total_weight());
2415 
2416  // Now add the actual bremsstrahlung reaction channel.
2417  brems_act.add_single_process();
2418 
2419  brems_act.perform_bremsstrahlung(outputs_);
2420  }
2421 
2422  logg[LExperiment].debug(~einhard::Green(), "✔ ", action);
2423  return true;
2424 }
2425 
2436 void validate_and_adjust_particle_list(ParticleList &particle_list);
2437 
2438 template <typename Modus>
2440  ParticleList &&add_plist,
2441  ParticleList &&remove_plist) {
2442  if (!add_plist.empty() || !remove_plist.empty()) {
2443  if (ensembles_.size() > 1) {
2444  throw std::runtime_error(
2445  "Adding or removing particles from SMASH is only possible when one "
2446  "ensemble is used.");
2447  }
2448  const double action_time = parameters_.labclock->current_time();
2449  /* Use two if statements. The first one is to check if the particles are
2450  * valid. Since this might remove all particles, a second if statement is
2451  * needed to avoid executing the action in that case.*/
2452  if (!add_plist.empty()) {
2454  }
2455  if (!add_plist.empty()) {
2456  // Create and perform action to add particle(s)
2457  auto action_add_particles = std::make_unique<FreeforallAction>(
2458  ParticleList{}, add_plist, action_time);
2459  perform_action(*action_add_particles, 0);
2460  }
2461  // Also here 2 if statements are needed as above.
2462  if (!remove_plist.empty()) {
2463  validate_and_adjust_particle_list(remove_plist);
2464  }
2465  if (!remove_plist.empty()) {
2466  ParticleList found_particles_to_remove;
2467  for (const auto &particle_to_remove : remove_plist) {
2468  const auto iterator_to_particle_to_be_removed_in_ensemble =
2469  std::find_if(
2470  ensembles_[0].begin(), ensembles_[0].end(),
2471  [&particle_to_remove, &action_time](const ParticleData &p) {
2473  particle_to_remove, p, action_time);
2474  });
2475  if (iterator_to_particle_to_be_removed_in_ensemble !=
2476  ensembles_[0].end())
2477  found_particles_to_remove.push_back(
2478  *iterator_to_particle_to_be_removed_in_ensemble);
2479  }
2480  // Sort the particles found to be removed according to their id and look
2481  // for duplicates (sorting is needed to call std::adjacent_find).
2482  std::sort(found_particles_to_remove.begin(),
2483  found_particles_to_remove.end(),
2484  [](const ParticleData &p1, const ParticleData &p2) {
2485  return p1.id() < p2.id();
2486  });
2487  const auto iterator_to_first_duplicate = std::adjacent_find(
2488  found_particles_to_remove.begin(), found_particles_to_remove.end(),
2489  [](const ParticleData &p1, const ParticleData &p2) {
2490  return p1.id() == p2.id();
2491  });
2492  if (iterator_to_first_duplicate != found_particles_to_remove.end()) {
2493  logg[LExperiment].error() << "The same particle has been asked to be "
2494  "removed multiple times:\n"
2495  << *iterator_to_first_duplicate;
2496  throw std::logic_error("Particle cannot be removed twice!");
2497  }
2498  if (auto delta = remove_plist.size() - found_particles_to_remove.size();
2499  delta > 0) {
2500  logg[LExperiment].warn(
2501  "When trying to remove particle(s) at the beginning ",
2502  "of the system evolution,\n", delta,
2503  " particle(s) could not be found and will be ignored.");
2504  }
2505  if (!found_particles_to_remove.empty()) {
2506  [[maybe_unused]] const auto number_particles_before_removal =
2507  ensembles_[0].size();
2508  // Create and perform action to remove particles
2509  auto action_remove_particles = std::make_unique<FreeforallAction>(
2510  found_particles_to_remove, ParticleList{}, action_time);
2511  perform_action(*action_remove_particles, 0);
2512 
2513  assert(number_particles_before_removal -
2514  found_particles_to_remove.size() ==
2515  ensembles_[0].size());
2516  }
2517  }
2518  }
2519 
2520  if (t_end > end_time_) {
2521  logg[LExperiment].fatal()
2522  << "Evolution asked to be run until " << t_end << " > " << end_time_
2523  << " and this cannot be done (because of how the clock works).";
2524  throw std::logic_error(
2525  "Experiment cannot evolve the system beyond End_Time.");
2526  }
2527  while (*(parameters_.labclock) < t_end) {
2528  const double dt = parameters_.labclock->timestep_duration();
2529  logg[LExperiment].debug("Timestepless propagation for next ", dt, " fm.");
2530 
2531  // Perform forced thermalization if required
2532  if (thermalizer_ &&
2533  thermalizer_->is_time_to_thermalize(parameters_.labclock)) {
2534  const bool ignore_cells_under_treshold = true;
2535  // Thermodynamics in thermalizer is computed from all ensembles,
2536  // but thermalization actions act on each ensemble independently
2537  thermalizer_->update_thermalizer_lattice(ensembles_, density_param_,
2538  ignore_cells_under_treshold);
2539  const double current_t = parameters_.labclock->current_time();
2540  for (int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
2541  thermalizer_->thermalize(ensembles_[i_ens], current_t,
2542  parameters_.testparticles);
2543  ThermalizationAction th_act(*thermalizer_, current_t);
2544  if (th_act.any_particles_thermalized()) {
2545  perform_action(th_act, i_ens);
2546  }
2547  }
2548  }
2549 
2550  if (IC_dynamic_) {
2551  modus_.build_fluidization_lattice(parameters_.labclock->current_time(),
2552  ensembles_, density_param_);
2553  }
2554 
2555  std::vector<Actions> actions(parameters_.n_ensembles);
2556  for (int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
2557  actions[i_ens].clear();
2558  if (ensembles_[i_ens].size() > 0 && action_finders_.size() > 0) {
2559  /* (1.a) Create grid. */
2560  const double min_cell_length = compute_min_cell_length(dt);
2561  logg[LExperiment].debug("Creating grid with minimal cell length ",
2562  min_cell_length);
2563  /* For the hyper-surface-crossing actions also unformed particles are
2564  * searched and therefore needed on the grid. */
2565  const bool include_unformed_particles = IC_switch_;
2566  const auto &grid =
2567  use_grid_ ? modus_.create_grid(ensembles_[i_ens], min_cell_length,
2568  dt, parameters_.coll_crit,
2569  include_unformed_particles)
2570  : modus_.create_grid(ensembles_[i_ens], min_cell_length,
2571  dt, parameters_.coll_crit,
2572  include_unformed_particles,
2574 
2575  const double gcell_vol = grid.cell_volume();
2576  /* (1.b) Iterate over cells and find actions. */
2577  grid.iterate_cells(
2578  [&](const ParticleList &search_list) {
2579  for (const auto &finder : action_finders_) {
2580  actions[i_ens].insert(finder->find_actions_in_cell(
2581  search_list, dt, gcell_vol, beam_momentum_));
2582  }
2583  },
2584  [&](const ParticleList &search_list,
2585  const ParticleList &neighbors_list) {
2586  for (const auto &finder : action_finders_) {
2587  actions[i_ens].insert(finder->find_actions_with_neighbors(
2588  search_list, neighbors_list, dt, beam_momentum_));
2589  }
2590  });
2591  }
2592  }
2593 
2594  /* \todo (optimizations) Adapt timestep size here */
2595 
2596  /* (2) Propagate from action to action until next output or timestep end */
2597  const double end_timestep_time = parameters_.labclock->next_time();
2598  while (next_output_time() < end_timestep_time) {
2599  for (int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
2600  run_time_evolution_timestepless(actions[i_ens], i_ens,
2601  next_output_time());
2602  }
2603  ++(*parameters_.outputclock);
2604 
2605  intermediate_output();
2606  }
2607  for (int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
2608  run_time_evolution_timestepless(actions[i_ens], i_ens, end_timestep_time);
2609  }
2610 
2611  /* (3) Update potentials (if computed on the lattice) and
2612  * compute new momenta according to equations of motion */
2613  if (potentials_) {
2614  update_potentials();
2615  update_momenta(ensembles_, parameters_.labclock->timestep_duration(),
2616  *potentials_, FB_lat_.get(), FI3_lat_.get(), EM_lat_.get(),
2617  jmu_B_lat_.get());
2618  }
2619 
2620  /* (4) Expand universe if non-minkowskian metric; updates
2621  * positions and momenta according to the selected expansion */
2622  if (metric_.mode_ != ExpansionMode::NoExpansion) {
2623  for (Particles &particles : ensembles_) {
2624  expand_space_time(&particles, parameters_, metric_);
2625  }
2626  }
2627 
2628  ++(*parameters_.labclock);
2629 
2630  /* (5) Check conservation laws.
2631  *
2632  * Check conservation of conserved quantities if potentials and string
2633  * fragmentation are off. If potentials are on then momentum is conserved
2634  * only in average. If string fragmentation is on, then energy and
2635  * momentum are only very roughly conserved in high-energy collisions. */
2636  if (!potentials_ && !parameters_.strings_switch &&
2637  metric_.mode_ == ExpansionMode::NoExpansion && !IC_switch_) {
2638  std::string err_msg = conserved_initial_.report_deviations(ensembles_);
2639  if (!err_msg.empty()) {
2640  logg[LExperiment].error() << err_msg;
2641  throw std::runtime_error("Violation of conserved quantities!");
2642  }
2643  }
2644  }
2645 
2646  /* Increment once more the output clock in order to have it prepared for the
2647  * final_output() call. Once close to the end time, the while-loop above to
2648  * produce intermediate output is not entered as the next_output_time() is
2649  * never strictly smaller than end_timestep_time (they are usually equal).
2650  * Since in the final_output() function the current time of the output clock
2651  * is used to produce the output, this has to be incremented before producing
2652  * the final output and it makes sense to do it here.
2653  */
2654  ++(*parameters_.outputclock);
2655 
2656  if (pauli_blocker_) {
2657  logg[LExperiment].info(
2658  "Interactions: Pauli-blocked/performed = ", total_pauli_blocked_, "/",
2659  interactions_total_ - wall_actions_total_);
2660  }
2661 }
2662 
2663 template <typename Modus>
2665  Particles &particles) {
2666  const double dt =
2667  propagate_straight_line(&particles, to_time, beam_momentum_);
2668  if (dilepton_finder_ != nullptr) {
2669  for (const auto &output : outputs_) {
2670  dilepton_finder_->shine(particles, output.get(), dt);
2671  }
2672  }
2673 }
2674 
2682 inline void check_interactions_total(uint64_t interactions_total) {
2683  constexpr uint64_t max_uint32 = std::numeric_limits<uint32_t>::max();
2684  if (interactions_total >= max_uint32) {
2685  throw std::runtime_error("Integer overflow in total interaction number!");
2686  }
2687 }
2688 
2689 template <typename Modus>
2691  Actions &actions, int i_ensemble, const double end_time_propagation) {
2692  Particles &particles = ensembles_[i_ensemble];
2693  logg[LExperiment].debug(
2694  "Timestepless propagation: ", "Actions size = ", actions.size(),
2695  ", end time = ", end_time_propagation);
2696 
2697  // iterate over all actions
2698  while (!actions.is_empty()) {
2699  if (actions.earliest_time() > end_time_propagation) {
2700  break;
2701  }
2702  // get next action
2703  ActionPtr act = actions.pop();
2704  if (!act->is_valid(particles)) {
2705  discarded_interactions_total_++;
2706  logg[LExperiment].debug(~einhard::DRed(), "✘ ", act,
2707  " (discarded: invalid)");
2708  continue;
2709  }
2710  logg[LExperiment].debug(~einhard::Green(), "✔ ", act,
2711  ", action time = ", act->time_of_execution());
2712 
2713  /* (1) Propagate to the next action. */
2714  propagate_and_shine(act->time_of_execution(), particles);
2715 
2716  /* (2) Perform action.
2717  *
2718  * Update the positions of the incoming particles, because the information
2719  * in the action object will be outdated as the particles have been
2720  * propagated since the construction of the action. */
2721  act->update_incoming(particles);
2722  const bool performed = perform_action(*act, i_ensemble);
2723 
2724  /* No need to update actions for outgoing particles
2725  * if the action is not performed. */
2726  if (!performed) {
2727  continue;
2728  }
2729 
2730  /* (3) Update actions for newly-produced particles. */
2731 
2732  const double end_time_timestep = parameters_.labclock->next_time();
2733  // New actions are always search until the end of the current timestep
2734  const double time_left = end_time_timestep - act->time_of_execution();
2735  const ParticleList &outgoing_particles = act->outgoing_particles();
2736  // Grid cell volume set to zero, since there is no grid
2737  const double gcell_vol = 0.0;
2738  for (const auto &finder : action_finders_) {
2739  // Outgoing particles can still decay, cross walls...
2740  actions.insert(finder->find_actions_in_cell(outgoing_particles, time_left,
2741  gcell_vol, beam_momentum_));
2742  // ... and collide with other particles.
2743  actions.insert(finder->find_actions_with_surrounding_particles(
2744  outgoing_particles, particles, time_left, beam_momentum_));
2745  }
2746 
2747  check_interactions_total(interactions_total_);
2748  }
2749 
2750  propagate_and_shine(end_time_propagation, particles);
2751 }
2752 
2753 template <typename Modus>
2755  const uint64_t wall_actions_this_interval =
2756  wall_actions_total_ - previous_wall_actions_total_;
2757  previous_wall_actions_total_ = wall_actions_total_;
2758  const uint64_t interactions_this_interval = interactions_total_ -
2759  previous_interactions_total_ -
2760  wall_actions_this_interval;
2761  previous_interactions_total_ = interactions_total_;
2762  double E_mean_field = 0.0;
2765  double computational_frame_time = 0.0;
2766  if (potentials_) {
2767  // using the lattice is necessary
2768  if ((jmu_B_lat_ != nullptr)) {
2769  E_mean_field = calculate_mean_field_energy(*potentials_, *jmu_B_lat_,
2770  EM_lat_.get(), parameters_);
2771  /*
2772  * Mean field calculated in a box should remain approximately constant if
2773  * the system is in equilibrium, and so deviations from its original value
2774  * may signal a phase transition or other dynamical process. This
2775  * comparison only makes sense in the Box Modus, hence the condition.
2776  */
2777  if (modus_.is_box()) {
2778  double tmp = (E_mean_field - initial_mean_field_energy_) /
2779  (E_mean_field + initial_mean_field_energy_);
2780  /*
2781  * This is displayed when the system evolves away from its initial
2782  * configuration (which is when the total mean field energy in the box
2783  * deviates from its initial value).
2784  */
2785  if (std::abs(tmp) > 0.01) {
2786  logg[LExperiment].info()
2787  << "\n\n\n\t The mean field at t = "
2788  << parameters_.outputclock->current_time()
2789  << " [fm] differs from the mean field at t = 0:"
2790  << "\n\t\t initial_mean_field_energy_ = "
2791  << initial_mean_field_energy_ << " [GeV]"
2792  << "\n\t\t abs[(E_MF - E_MF(t=0))/(E_MF + E_MF(t=0))] = "
2793  << std::abs(tmp)
2794  << "\n\t\t E_MF/E_MF(t=0) = "
2795  << E_mean_field / initial_mean_field_energy_ << "\n\n";
2796  }
2797  }
2798  }
2799  }
2800 
2802  ensembles_, interactions_this_interval, conserved_initial_, time_start_,
2803  parameters_.outputclock->current_time(), E_mean_field,
2804  initial_mean_field_energy_);
2805  const LatticeUpdate lat_upd = LatticeUpdate::AtOutput;
2806 
2807  // save evolution data
2808  if (!(modus_.is_box() && parameters_.outputclock->current_time() <
2809  modus_.equilibration_time())) {
2810  for (const auto &output : outputs_) {
2811  if (output->is_dilepton_output() || output->is_photon_output() ||
2812  output->is_IC_output()) {
2813  continue;
2814  }
2815  for (int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
2816  auto event_info = fill_event_info(
2817  ensembles_, E_mean_field, modus_.impact_parameter(), parameters_,
2818  projectile_target_interact_[i_ens], kinematic_cuts_for_IC_output_);
2819 
2820  output->at_intermediate_time(ensembles_[i_ens], parameters_.outputclock,
2821  density_param_, {event_, i_ens},
2822  event_info);
2823  computational_frame_time = event_info.current_time;
2824  }
2825  // For thermodynamic output
2826  output->at_intermediate_time(ensembles_, parameters_.outputclock,
2827  density_param_);
2828 
2829  // Thermodynamic output on the lattice versus time
2830  if (printout_rho_eckart_) {
2831  switch (dens_type_lattice_printout_) {
2832  case DensityType::Baryon:
2834  jmu_B_lat_.get(), lat_upd, DensityType::Baryon, density_param_,
2835  ensembles_, false);
2836  output->thermodynamics_output(ThermodynamicQuantity::EckartDensity,
2837  DensityType::Baryon, *jmu_B_lat_);
2838  output->thermodynamics_lattice_output(*jmu_B_lat_,
2839  computational_frame_time);
2840  break;
2843  jmu_I3_lat_.get(), lat_upd, DensityType::BaryonicIsospin,
2844  density_param_, ensembles_, false);
2845  output->thermodynamics_output(ThermodynamicQuantity::EckartDensity,
2847  *jmu_I3_lat_);
2848  output->thermodynamics_lattice_output(*jmu_I3_lat_,
2849  computational_frame_time);
2850  break;
2851  case DensityType::None:
2852  break;
2853  default:
2855  jmu_custom_lat_.get(), lat_upd, dens_type_lattice_printout_,
2856  density_param_, ensembles_, false);
2857  output->thermodynamics_output(ThermodynamicQuantity::EckartDensity,
2858  dens_type_lattice_printout_,
2859  *jmu_custom_lat_);
2860  output->thermodynamics_lattice_output(*jmu_custom_lat_,
2861  computational_frame_time);
2862  }
2863  }
2864  if (printout_tmn_ || printout_tmn_landau_ || printout_v_landau_) {
2866  Tmn_.get(), lat_upd, dens_type_lattice_printout_, density_param_,
2867  ensembles_, false);
2868  if (printout_tmn_) {
2869  output->thermodynamics_output(ThermodynamicQuantity::Tmn,
2870  dens_type_lattice_printout_, *Tmn_);
2871  output->thermodynamics_lattice_output(
2872  ThermodynamicQuantity::Tmn, *Tmn_, computational_frame_time);
2873  }
2874  if (printout_tmn_landau_) {
2875  output->thermodynamics_output(ThermodynamicQuantity::TmnLandau,
2876  dens_type_lattice_printout_, *Tmn_);
2877  output->thermodynamics_lattice_output(
2879  computational_frame_time);
2880  }
2881  if (printout_v_landau_) {
2882  output->thermodynamics_output(ThermodynamicQuantity::LandauVelocity,
2883  dens_type_lattice_printout_, *Tmn_);
2884  output->thermodynamics_lattice_output(
2886  computational_frame_time);
2887  }
2888  }
2889  if (EM_lat_) {
2890  output->fields_output("Efield", "Bfield", *EM_lat_);
2891  }
2892  if (printout_j_QBS_) {
2893  output->thermodynamics_lattice_output(
2894  *j_QBS_lat_, computational_frame_time, ensembles_, density_param_);
2895  }
2896 
2897  if (thermalizer_) {
2898  output->thermodynamics_output(*thermalizer_);
2899  }
2900  }
2901  }
2902 }
2903 
2904 template <typename Modus>
2906  if (potentials_) {
2907  if (potentials_->use_symmetry() && jmu_I3_lat_ != nullptr) {
2908  update_lattice(jmu_I3_lat_.get(), old_jmu_auxiliary_.get(),
2909  new_jmu_auxiliary_.get(), four_gradient_auxiliary_.get(),
2911  density_param_, ensembles_,
2912  parameters_.labclock->timestep_duration(), true);
2913  }
2914  if ((potentials_->use_skyrme() || potentials_->use_symmetry()) &&
2915  jmu_B_lat_ != nullptr) {
2916  update_lattice(jmu_B_lat_.get(), old_jmu_auxiliary_.get(),
2917  new_jmu_auxiliary_.get(), four_gradient_auxiliary_.get(),
2919  density_param_, ensembles_,
2920  parameters_.labclock->timestep_duration(), true);
2921  const size_t UBlattice_size = UB_lat_->size();
2922  for (size_t i = 0; i < UBlattice_size; i++) {
2923  auto jB = (*jmu_B_lat_)[i];
2924  const FourVector flow_four_velocity_B =
2925  std::abs(jB.rho()) > very_small_double ? jB.jmu_net() / jB.rho()
2926  : FourVector();
2927  double baryon_density = jB.rho();
2928  ThreeVector baryon_grad_j0 = jB.grad_j0();
2929  ThreeVector baryon_dvecj_dt = jB.dvecj_dt();
2930  ThreeVector baryon_curl_vecj = jB.curl_vecj();
2931  if (potentials_->use_skyrme()) {
2932  (*UB_lat_)[i] =
2933  flow_four_velocity_B * potentials_->skyrme_pot(baryon_density);
2934  (*FB_lat_)[i] =
2935  potentials_->skyrme_force(baryon_density, baryon_grad_j0,
2936  baryon_dvecj_dt, baryon_curl_vecj);
2937  }
2938  if (potentials_->use_symmetry() && jmu_I3_lat_ != nullptr) {
2939  auto jI3 = (*jmu_I3_lat_)[i];
2940  const FourVector flow_four_velocity_I3 =
2941  std::abs(jI3.rho()) > very_small_double
2942  ? jI3.jmu_net() / jI3.rho()
2943  : FourVector();
2944  (*UI3_lat_)[i] = flow_four_velocity_I3 *
2945  potentials_->symmetry_pot(jI3.rho(), baryon_density);
2946  (*FI3_lat_)[i] = potentials_->symmetry_force(
2947  jI3.rho(), jI3.grad_j0(), jI3.dvecj_dt(), jI3.curl_vecj(),
2948  baryon_density, baryon_grad_j0, baryon_dvecj_dt,
2949  baryon_curl_vecj);
2950  }
2951  }
2952  }
2953  if (potentials_->use_coulomb()) {
2956  density_param_, ensembles_, true);
2957  for (size_t i = 0; i < EM_lat_->size(); i++) {
2958  ThreeVector electric_field = {0., 0., 0.};
2959  ThreeVector position = jmu_el_lat_->cell_center(i);
2960  jmu_el_lat_->integrate_volume(electric_field,
2962  potentials_->coulomb_r_cut(), position);
2963  ThreeVector magnetic_field = {0., 0., 0.};
2964  jmu_el_lat_->integrate_volume(magnetic_field,
2966  potentials_->coulomb_r_cut(), position);
2967  (*EM_lat_)[i] = std::make_pair(electric_field, magnetic_field);
2968  }
2969  } // if ((potentials_->use_skyrme() || ...
2970  if (potentials_->use_vdf() && jmu_B_lat_ != nullptr) {
2971  update_lattice(jmu_B_lat_.get(), old_jmu_auxiliary_.get(),
2972  new_jmu_auxiliary_.get(), four_gradient_auxiliary_.get(),
2974  density_param_, ensembles_,
2975  parameters_.labclock->timestep_duration(), true);
2976  if (parameters_.field_derivatives_mode == FieldDerivativesMode::Direct) {
2978  fields_lat_.get(), old_fields_auxiliary_.get(),
2979  new_fields_auxiliary_.get(), fields_four_gradient_auxiliary_.get(),
2980  jmu_B_lat_.get(), LatticeUpdate::EveryTimestep, *potentials_,
2981  parameters_.labclock->timestep_duration());
2982  }
2983  const size_t UBlattice_size = UB_lat_->size();
2984  for (size_t i = 0; i < UBlattice_size; i++) {
2985  auto jB = (*jmu_B_lat_)[i];
2986  (*UB_lat_)[i] = potentials_->vdf_pot(jB.rho(), jB.jmu_net());
2987  switch (parameters_.field_derivatives_mode) {
2989  (*FB_lat_)[i] = potentials_->vdf_force(
2990  jB.rho(), jB.drho_dxnu().x0(), jB.drho_dxnu().threevec(),
2991  jB.grad_rho_cross_vecj(), jB.jmu_net().x0(), jB.grad_j0(),
2992  jB.jmu_net().threevec(), jB.dvecj_dt(), jB.curl_vecj());
2993  break;
2995  auto Amu = (*fields_lat_)[i];
2996  (*FB_lat_)[i] = potentials_->vdf_force(
2997  Amu.grad_A0(), Amu.dvecA_dt(), Amu.curl_vecA());
2998  break;
2999  }
3000  } // for (size_t i = 0; i < UBlattice_size; i++)
3001  } // if potentials_->use_vdf()
3002  }
3003 }
3004 
3005 template <typename Modus>
3007  /* At end of time evolution: Force all resonances to decay. In order to handle
3008  * decay chains, we need to loop until no further actions occur. */
3009  bool actions_performed, decays_found;
3010  uint64_t interactions_old;
3011  do {
3012  decays_found = false;
3013  interactions_old = interactions_total_;
3014  for (int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
3015  Actions actions;
3016 
3017  // Dileptons: shining of remaining resonances
3018  if (dilepton_finder_ != nullptr) {
3019  for (const auto &output : outputs_) {
3020  dilepton_finder_->shine_final(ensembles_[i_ens], output.get(), true);
3021  }
3022  }
3023  // Find actions.
3024  for (const auto &finder : action_finders_) {
3025  auto found_actions = finder->find_final_actions(ensembles_[i_ens]);
3026  if (!found_actions.empty()) {
3027  actions.insert(std::move(found_actions));
3028  decays_found = true;
3029  }
3030  }
3031  // Perform actions.
3032  while (!actions.is_empty()) {
3033  perform_action(*actions.pop(), i_ens, false);
3034  }
3035  }
3036  actions_performed = interactions_total_ > interactions_old;
3037  // Throw an error if actions were found but not performed
3038  if (decays_found && !actions_performed) {
3039  throw std::runtime_error("Final decays were found but not performed.");
3040  }
3041  // loop until no more decays occur
3042  } while (actions_performed);
3043 
3044  // Dileptons: shining of stable particles at the end
3045  if (dilepton_finder_ != nullptr) {
3046  for (const auto &output : outputs_) {
3047  for (Particles &particles : ensembles_) {
3048  dilepton_finder_->shine_final(particles, output.get(), false);
3049  }
3050  }
3051  }
3052 }
3053 
3054 template <typename Modus>
3056  /* make sure the experiment actually ran (note: we should compare this
3057  * to the start time, but we don't know that. Therefore, we check that
3058  * the time is positive, which should heuristically be the same). */
3059  double E_mean_field = 0.0;
3060  if (likely(parameters_.labclock > 0)) {
3061  const uint64_t wall_actions_this_interval =
3062  wall_actions_total_ - previous_wall_actions_total_;
3063  const uint64_t interactions_this_interval = interactions_total_ -
3064  previous_interactions_total_ -
3065  wall_actions_this_interval;
3066  if (potentials_) {
3067  // using the lattice is necessary
3068  if ((jmu_B_lat_ != nullptr)) {
3069  E_mean_field = calculate_mean_field_energy(*potentials_, *jmu_B_lat_,
3070  EM_lat_.get(), parameters_);
3071  }
3072  }
3073  if (std::abs(parameters_.labclock->current_time() - end_time_) >
3074  really_small) {
3075  logg[LExperiment].warn()
3076  << "SMASH not propagated until configured end time. Current time = "
3077  << parameters_.labclock->current_time()
3078  << "fm. End time = " << end_time_ << "fm.";
3079  } else {
3081  ensembles_, interactions_this_interval, conserved_initial_,
3082  time_start_, end_time_, E_mean_field, initial_mean_field_energy_);
3083  }
3084  int total_particles = 0;
3085  for (const Particles &particles : ensembles_) {
3086  total_particles += particles.size();
3087  }
3088  if (IC_switch_ && (total_particles == 0)) {
3089  const double initial_system_energy_plus_Pythia_violations =
3090  conserved_initial_.momentum().x0() + total_energy_violated_by_Pythia_;
3091  const double fraction_of_total_system_energy_removed =
3092  initial_system_energy_plus_Pythia_violations / total_energy_removed_;
3093  // Verify there is no more energy in the system if all particles were
3094  // removed when crossing the hypersurface
3095  if (std::fabs(fraction_of_total_system_energy_removed - 1.) >
3096  really_small) {
3097  throw std::runtime_error(
3098  "There is remaining energy in the system although all particles "
3099  "were removed.\n"
3100  "E_remain = " +
3101  std::to_string((initial_system_energy_plus_Pythia_violations -
3102  total_energy_removed_)) +
3103  " [GeV]");
3104  } else {
3105  logg[LExperiment].info() << hline;
3106  logg[LExperiment].info()
3107  << "Time real: " << SystemClock::now() - time_start_;
3108  logg[LExperiment].info()
3109  << "Interactions before reaching hypersurface: "
3110  << interactions_total_ - wall_actions_total_ -
3111  total_hypersurface_crossing_actions_;
3112  logg[LExperiment].info()
3113  << "Total number of particles removed on hypersurface: "
3114  << total_hypersurface_crossing_actions_;
3115  }
3116  } else {
3117  const double precent_discarded =
3118  interactions_total_ > 0
3119  ? static_cast<double>(discarded_interactions_total_) * 100.0 /
3120  interactions_total_
3121  : 0.0;
3122  std::stringstream msg_discarded;
3123  msg_discarded
3124  << "Discarded interaction number: " << discarded_interactions_total_
3125  << " (" << precent_discarded
3126  << "% of the total interaction number including wall crossings)";
3127 
3128  logg[LExperiment].info() << hline;
3129  logg[LExperiment].info()
3130  << "Time real: " << SystemClock::now() - time_start_;
3131  logg[LExperiment].debug() << msg_discarded.str();
3132 
3133  if (parameters_.coll_crit == CollisionCriterion::Stochastic &&
3134  precent_discarded > 1.0) {
3135  // The choosen threshold of 1% is a heuristical value
3136  logg[LExperiment].warn()
3137  << msg_discarded.str()
3138  << "\nThe number of discarded interactions is large, which means "
3139  "the assumption for the stochastic criterion of\n1 interaction "
3140  "per particle per timestep is probably violated. Consider "
3141  "reducing the timestep size.";
3142  }
3143 
3144  logg[LExperiment].info() << "Final interaction number: "
3145  << interactions_total_ - wall_actions_total_;
3146  }
3147 
3148  // Check if there are unformed particles
3149  int unformed_particles_count = 0;
3150  for (const Particles &particles : ensembles_) {
3151  for (const ParticleData &particle : particles) {
3152  if (particle.formation_time() > end_time_) {
3153  unformed_particles_count++;
3154  }
3155  }
3156  }
3157  if (unformed_particles_count > 0) {
3158  logg[LExperiment].warn(
3159  "End time might be too small. ", unformed_particles_count,
3160  " unformed particles were found at the end of the evolution.");
3161  }
3162  }
3163 
3164  // Keep track of how many ensembles had interactions
3165  count_nonempty_ensembles();
3166 
3167  for (const auto &output : outputs_) {
3168  for (int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
3169  auto event_info = fill_event_info(
3170  ensembles_, E_mean_field, modus_.impact_parameter(), parameters_,
3171  projectile_target_interact_[i_ens], kinematic_cuts_for_IC_output_);
3172  output->at_eventend(ensembles_[i_ens], {event_, i_ens}, event_info);
3173  }
3174  // For thermodynamic output
3175  output->at_eventend(ensembles_, event_);
3176 
3177  // For thermodynamic lattice output
3178  if (printout_rho_eckart_) {
3179  if (dens_type_lattice_printout_ != DensityType::None) {
3180  output->at_eventend(ThermodynamicQuantity::EckartDensity);
3181  }
3182  }
3183  if (printout_tmn_) {
3184  output->at_eventend(ThermodynamicQuantity::Tmn);
3185  }
3186  if (printout_tmn_landau_) {
3187  output->at_eventend(ThermodynamicQuantity::TmnLandau);
3188  }
3189  if (printout_v_landau_) {
3190  output->at_eventend(ThermodynamicQuantity::LandauVelocity);
3191  }
3192  if (printout_j_QBS_) {
3193  output->at_eventend(ThermodynamicQuantity::j_QBS);
3194  }
3195  }
3196 }
3197 
3198 template <typename Modus>
3200  for (bool has_interaction : projectile_target_interact_) {
3201  if (has_interaction) {
3202  nonempty_ensembles_++;
3203  }
3204  }
3205 }
3206 
3207 template <typename Modus>
3209  if (event_counting_ == EventCounting::FixedNumber) {
3210  return event_ >= nevents_;
3211  }
3212  if (event_counting_ == EventCounting::MinimumNonEmpty) {
3213  if (nonempty_ensembles_ >= minimum_nonempty_ensembles_) {
3214  return true;
3215  }
3216  if (event_ >= max_events_) {
3217  logg[LExperiment].warn()
3218  << "Maximum number of events (" << max_events_
3219  << ") exceeded. Stopping calculation. "
3220  << "The fraction of empty ensembles is "
3221  << (1.0 - static_cast<double>(nonempty_ensembles_) /
3222  (event_ * parameters_.n_ensembles))
3223  << ". If this fraction is expected, try increasing the "
3224  "Maximum_Ensembles_Run.";
3225  return true;
3226  }
3227  return false;
3228  }
3229  throw std::runtime_error("Event counting option is invalid");
3230  return false;
3231 }
3232 
3233 template <typename Modus>
3235  event_++;
3236 }
3237 
3238 template <typename Modus>
3240  const auto &mainlog = logg[LMain];
3241  for (event_ = 0; !is_finished(); event_++) {
3242  mainlog.info() << "Event " << event_;
3243 
3244  // Sample initial particles, start clock, some printout and book-keeping
3245  initialize_new_event();
3246 
3247  run_time_evolution(end_time_);
3248 
3249  if (force_decays_) {
3250  do_final_decays();
3251  }
3252 
3253  // Output at event end
3254  final_output();
3255  }
3256 }
3257 
3258 } // namespace smash
3259 
3260 #endif // SRC_INCLUDE_SMASH_EXPERIMENT_H_
Collection of useful type aliases to measure and output the (real) runtime.
A stream modifier that allows to colorize the log output.
Definition: einhard.hpp:147
Action is the base class for a generic process that takes a number of incoming particles and transfor...
Definition: action.h:35
virtual ProcessType get_type() const
Get the process type.
Definition: action.h:131
virtual double get_total_weight() const =0
Return the total weight value, which is mainly used for the weight output entry.
virtual double perform(Particles *particles, uint32_t id_process)
Actually perform the action, e.g.
Definition: action.cc:128
const ParticleList & incoming_particles() const
Get the list of particles that go into the action.
Definition: action.cc:58
virtual void generate_final_state()=0
Generate the final state for this action.
double sqrt_s() const
Determine the total energy in the center-of-mass frame [GeV].
Definition: action.h:271
FourVector get_interaction_point() const
Get the interaction point.
Definition: action.cc:68
bool is_valid(const Particles &particles) const
Check whether the action still applies.
Definition: action.cc:29
bool is_pauli_blocked(const std::vector< Particles > &ensembles, const PauliBlocker &p_bl) const
Check if the action is Pauli-blocked.
Definition: action.cc:35
The Actions class abstracts the storage and manipulation of actions.
Definition: actions.h:29
ActionPtr pop()
Return the first action in the list and removes it from the list.
Definition: actions.h:59
double earliest_time() const
Return time of execution of earliest action.
Definition: actions.h:70
ActionList::size_type size() const
Definition: actions.h:98
void insert(ActionList &&new_acts)
Insert a list of actions into this object.
Definition: actions.h:79
bool is_empty() const
Definition: actions.h:52
static bool is_bremsstrahlung_reaction(const ParticleList &in)
Check if particles can undergo an implemented photon process.
Interface to the SMASH configuration files.
void set_value(Key< U > key, T &&value)
Overwrite the value of the YAML node corresponding to the specified key.
Configuration extract_sub_configuration(KeyLabels section, Configuration::GetEmpty empty_if_not_existing=Configuration::GetEmpty::No)
Create a new configuration from a then-removed section of the present object.
T read(const Key< T > &key) const
Additional interface for SMASH to read configuration values without removing them.
bool has_value(const Key< T > &key) const
Return whether there is a non-empty value behind the requested key (which is supposed not to refer to...
bool has_section(const KeyLabels &labels) const
Return whether there is a (possibly empty) section with the given labels.
Configuration extract_complete_sub_configuration(KeyLabels section, Configuration::GetEmpty empty_if_not_existing=Configuration::GetEmpty::No)
Alternative method to extract a sub-configuration, which retains the labels from the top-level in the...
T take(const Key< T > &key)
The default interface for SMASH to read configuration values.
A class to pre-calculate and store parameters relevant for density calculation.
Definition: density.h:92
Non-template interface to Experiment<Modus>.
Definition: experiment.h:102
static std::unique_ptr< ExperimentBase > create(Configuration &config, const std::filesystem::path &output_path)
Factory method that creates and initializes a new Experiment<Modus>.
Definition: experiment.cc:22
virtual ~ExperimentBase()=default
The virtual destructor avoids undefined behavior when destroying derived objects.
virtual void run()=0
Runs the experiment.
The main class, where the simulation of an experiment is executed.
Definition: experiment.h:191
const bool bremsstrahlung_switch_
This indicates whether bremsstrahlung is switched on.
Definition: experiment.h:613
void propagate_and_shine(double to_time, Particles &particles)
Propagate all particles until time to_time without any interactions and shine dileptons.
Definition: experiment.h:2664
Experiment(Configuration &config, const std::filesystem::path &output_path)
Create a new Experiment.
const ExpansionProperties metric_
This struct contains information on the metric to be used.
Definition: experiment.h:604
std::unique_ptr< ActionFinderInterface > photon_finder_
The (Scatter) Actions Finder for Direct Photons.
Definition: experiment.h:429
const bool force_decays_
This indicates whether we force all resonances to decay in the last timestep.
Definition: experiment.h:598
double initial_mean_field_energy_
The initial total mean field energy in the system.
Definition: experiment.h:646
void create_output(const std::string &format, const std::string &content, const std::filesystem::path &output_path, const OutputParameters &par)
Create a list of output files.
Definition: experiment.h:740
std::vector< std::unique_ptr< ActionFinderInterface > > action_finders_
The Action finder objects.
Definition: experiment.h:423
bool printout_tmn_
Whether to print the energy-momentum tensor.
Definition: experiment.h:511
QuantumNumbers conserved_initial_
The conserved quantities of the system.
Definition: experiment.h:640
DensityParameters density_param_
Structure to precalculate and hold parameters for density computations.
Definition: experiment.h:374
double next_output_time() const
Shortcut for next output time.
Definition: experiment.h:349
double total_energy_removed_
Total energy removed from the system in hypersurface crossing actions.
Definition: experiment.h:699
DensityType dens_type_lattice_printout_
Type of density for lattice printout.
Definition: experiment.h:459
double max_transverse_distance_sqr_
Maximal distance at which particles can interact in case of the geometric criterion,...
Definition: experiment.h:631
const bool IC_dynamic_
This indicates if the IC is dynamic.
Definition: experiment.h:622
bool printout_j_QBS_
Whether to print the Q, B, S 4-currents.
Definition: experiment.h:520
std::unique_ptr< GrandCanThermalizer > thermalizer_
Instance of class used for forced thermalization.
Definition: experiment.h:530
void count_nonempty_ensembles()
Counts the number of ensembles in wich interactions took place at the end of an event.
Definition: experiment.h:3199
const TimeStepMode time_step_mode_
This indicates whether to use time steps.
Definition: experiment.h:625
std::unique_ptr< DensityLattice > jmu_custom_lat_
Custom density on the lattices.
Definition: experiment.h:456
bool printout_full_lattice_any_td_
Whether to print the thermodynamics quantities evaluated on the lattices, point by point,...
Definition: experiment.h:527
void increase_event_number()
Increases the event number by one.
Definition: experiment.h:3234
void run_time_evolution_timestepless(Actions &actions, int i_ensemble, const double end_time_propagation)
Performs all the propagations and actions during a certain time interval neglecting the influence of ...
Definition: experiment.h:2690
Particles * first_ensemble()
Provides external access to SMASH particles.
Definition: experiment.h:260
void run_time_evolution(const double t_end, ParticleList &&add_plist={}, ParticleList &&remove_plist={})
Runs the time evolution of an event with fixed-size time steps or without timesteps,...
Definition: experiment.h:2439
int minimum_nonempty_ensembles_
The number of ensembles, in which interactions take place, to be calculated.
Definition: experiment.h:562
std::unique_ptr< DensityLattice > jmu_B_lat_
Baryon density on the lattice.
Definition: experiment.h:438
DensityType dens_type_
Type of density to be written to collision headers.
Definition: experiment.h:652
void intermediate_output()
Intermediate output during an event.
Definition: experiment.h:2754
std::vector< FourVector > beam_momentum_
The initial nucleons in the ColliderModus propagate with beam_momentum_, if Fermi motion is frozen.
Definition: experiment.h:420
std::unique_ptr< RectangularLattice< std::pair< ThreeVector, ThreeVector > > > FI3_lat_
Lattices for the electric and magnetic component of the symmetry force.
Definition: experiment.h:482
std::unique_ptr< RectangularLattice< FourVector > > new_jmu_auxiliary_
Auxiliary lattice for values of jmu at a time step t0 + dt.
Definition: experiment.h:494
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:341
void initialize_new_event()
This is called in the beginning of each event.
Definition: experiment.h:2071
bool is_finished()
Checks wether the desired number events have been calculated.
Definition: experiment.h:3208
int n_fractional_photons_
Number of fractional photons produced per single reaction.
Definition: experiment.h:432
const double delta_time_startup_
The clock's timestep size at start up.
Definition: experiment.h:592
std::unique_ptr< RectangularLattice< EnergyMomentumTensor > > Tmn_
Lattices of energy-momentum tensors for printout.
Definition: experiment.h:489
std::vector< Particles > ensembles_
Complete particle list, all ensembles in one vector.
Definition: experiment.h:383
std::unique_ptr< RectangularLattice< FourVector > > new_fields_auxiliary_
Auxiliary lattice for values of Amu at a time step t0 + dt.
Definition: experiment.h:502
SystemTimePoint time_start_
system starting time of the simulation
Definition: experiment.h:649
uint64_t previous_wall_actions_total_
Total number of wall-crossings for previous timestep.
Definition: experiment.h:676
const bool IC_switch_
This indicates whether the experiment will be used as initial condition for hydrodynamics.
Definition: experiment.h:619
std::unique_ptr< RectangularLattice< std::pair< ThreeVector, ThreeVector > > > FB_lat_
Lattices for the electric and magnetic components of the Skyrme or VDF force.
Definition: experiment.h:478
OutputsList outputs_
A list of output formaters.
Definition: experiment.h:401
bool printout_rho_eckart_
Whether to print the Eckart rest frame density.
Definition: experiment.h:508
int event_
Current event.
Definition: experiment.h:573
std::unique_ptr< RectangularLattice< std::array< FourVector, 4 > > > fields_four_gradient_auxiliary_
Auxiliary lattice for calculating the four-gradient of Amu.
Definition: experiment.h:505
void do_final_decays()
Performs the final decays of an event.
Definition: experiment.h:3006
std::unique_ptr< PauliBlocker > pauli_blocker_
An instance of PauliBlocker class that stores parameters needed for Pauli blocking calculations and c...
Definition: experiment.h:395
bool printout_v_landau_
Whether to print the 4-velocity in Landau frame.
Definition: experiment.h:517
bool perform_action(Action &action, int i_ensemble, bool include_pauli_blocking=true)
Perform the given action.
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:471
bool printout_lattice_td_
Whether to print the thermodynamics quantities evaluated on the lattices.
Definition: experiment.h:523
std::unique_ptr< RectangularLattice< std::pair< ThreeVector, ThreeVector > > > EM_lat_
Lattices for electric and magnetic field in fm^-2.
Definition: experiment.h:486
std::unique_ptr< FieldsLattice > fields_lat_
Mean-field A^mu on the lattice.
Definition: experiment.h:447
int max_events_
Maximum number of events to be calculated in order obtain the desired number of non-empty events usin...
Definition: experiment.h:582
int nonempty_ensembles_
Number of ensembles containing an interaction.
Definition: experiment.h:576
OutputPtr photon_output_
The Photon output.
Definition: experiment.h:407
EventCounting event_counting_
The way in which the number of calculated events is specified.
Definition: experiment.h:570
const bool dileptons_switch_
This indicates whether dileptons are switched on.
Definition: experiment.h:607
Modus modus_
Instance of the Modus template parameter.
Definition: experiment.h:380
std::unique_ptr< RectangularLattice< FourVector > > old_jmu_auxiliary_
Auxiliary lattice for values of jmu at a time step t0.
Definition: experiment.h:492
std::unique_ptr< DecayActionsFinderDilepton > dilepton_finder_
The Dilepton Action Finder.
Definition: experiment.h:426
std::unique_ptr< DensityLattice > j_QBS_lat_
4-current for j_QBS lattice output
Definition: experiment.h:435
bool printout_tmn_landau_
Whether to print the energy-momentum tensor in Landau frame.
Definition: experiment.h:514
std::unique_ptr< RectangularLattice< std::array< FourVector, 4 > > > four_gradient_auxiliary_
Auxiliary lattice for calculating the four-gradient of jmu.
Definition: experiment.h:497
const bool photons_switch_
This indicates whether photons are switched on.
Definition: experiment.h:610
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:536
std::unique_ptr< RectangularLattice< FourVector > > UB_lat_
Lattices for Skyrme or VDF potentials (evaluated in the local rest frame) times the baryon flow 4-vel...
Definition: experiment.h:465
ExperimentParameters parameters_
Struct of several member variables.
Definition: experiment.h:371
std::unique_ptr< RectangularLattice< FourVector > > old_fields_auxiliary_
Auxiliary lattice for values of Amu at a time step t0.
Definition: experiment.h:500
std::unique_ptr< DensityLattice > jmu_el_lat_
Electric charge density on the lattice.
Definition: experiment.h:444
std::unique_ptr< Potentials > potentials_
An instance of potentials class, that stores parameters of potentials, calculates them and their grad...
Definition: experiment.h:389
uint64_t wall_actions_total_
Total number of wall-crossings for current timestep.
Definition: experiment.h:670
uint64_t total_hypersurface_crossing_actions_
Total number of particles removed from the evolution in hypersurface crossing actions.
Definition: experiment.h:688
uint64_t interactions_total_
Total number of interactions for current timestep.
Definition: experiment.h:658
const bool use_grid_
This indicates whether to use the grid.
Definition: experiment.h:601
std::unique_ptr< DensityLattice > jmu_I3_lat_
Isospin projection density on the lattice.
Definition: experiment.h:441
int nevents_
Number of events.
Definition: experiment.h:552
uint64_t total_pauli_blocked_
Total number of Pauli-blockings for current timestep.
Definition: experiment.h:682
void final_output()
Output at the end of an event.
Definition: experiment.h:3055
std::vector< bool > projectile_target_interact_
Whether the projectile and the target collided.
Definition: experiment.h:413
Modus * modus()
Provides external access to SMASH calculation modus.
Definition: experiment.h:268
void update_potentials()
Recompute potentials on lattices if necessary.
Definition: experiment.h:2905
OutputPtr dilepton_output_
The Dilepton output.
Definition: experiment.h:404
int64_t seed_
random seed for the next event.
Definition: experiment.h:710
std::vector< Particles > * all_ensembles()
Getter for all ensembles.
Definition: experiment.h:262
uint64_t previous_interactions_total_
Total number of interactions for previous timestep.
Definition: experiment.h:664
uint64_t discarded_interactions_total_
Total number of discarded interactions, because they were invalidated before they could be performed.
Definition: experiment.h:694
void run() override
Runs the experiment.
Definition: experiment.h:3239
bool kinematic_cuts_for_IC_output_
This indicates whether kinematic cuts are enabled for the IC output.
Definition: experiment.h:707
const double end_time_
simulation time at which the evolution is stopped.
Definition: experiment.h:585
double total_energy_violated_by_Pythia_
Total energy violation introduced by Pythia.
Definition: experiment.h:704
The FourVector class holds relevant values in Minkowski spacetime with (+, −, −, −) metric signature.
Definition: fourvector.h:33
ParticleData contains the dynamic information of a certain particle.
Definition: particledata.h:58
static double formation_power_
Power with which the cross section scaling factor grows in time.
Definition: particledata.h:395
The Particles class abstracts the storage and manipulation of particles.
Definition: particles.h:33
size_t size() const
Definition: particles.h:87
A class that stores parameters of potentials, calculates potentials and their gradients.
Definition: potentials.h:36
static ThreeVector B_field_integrand(ThreeVector pos, DensityOnLattice &charge_density, ThreeVector point)
Integrand for calculating the magnetic field using the Biot-Savart formula.
Definition: potentials.h:338
static ThreeVector E_field_integrand(ThreeVector pos, DensityOnLattice &charge_density, ThreeVector point)
Integrand for calculating the electric field.
Definition: potentials.h:320
A container for storing conserved values.
A container class to hold all the arrays on the lattice and access them.
Definition: lattice.h:49
static bool is_photon_reaction(const ParticleList &in)
Check if particles can undergo an implemented photon process.
static bool is_kinematically_possible(const double s_sqrt, const ParticleList &in)
Check if CM-energy is sufficient to produce hadron in final state.
String excitation processes used in SMASH.
Definition: stringprocess.h:45
ThermalizationAction implements forced thermalization as an Action class.
bool any_particles_thermalized() const
This method checks, if there are particles in the region to be thermalized.
The ThreeVector class represents a physical three-vector with the components .
Definition: threevector.h:31
@ Frozen
Use fermi motion without potentials.
@ Dynamic
Dynamic fluidization based on local densities.
TimeStepMode
The time step mode.
@ Fixed
Use fixed time step.
@ None
Don't use time steps; propagate from action to action.
@ EckartDensity
Density in the Eckart frame.
@ Tmn
Energy-momentum tensor in lab frame.
@ LandauVelocity
Velocity of the Landau rest frame.
@ j_QBS
Electric (Q), baryonic (B) and strange (S) currents.
@ TmnLandau
Energy-momentum tensor in Landau rest frame.
@ BottomUp
Sum the existing partial contributions.
@ Stochastic
Stochastic Criteiron.
@ None
No pseudo-resonance is created.
EventCounting
Defines how the number of events is determined.
@ Invalid
Unused, only in the code for internal logic.
@ FixedNumber
The desired number of events is simulated disregarding of whether an interaction took place.
@ MinimumNonEmpty
Events are simulated until there are at least a given number of ensembles in which an interaction too...
DensityType
Allows to choose which kind of density to calculate.
#define SMASH_SOURCE_LOCATION
Hackery that is required to output the location in the source code where the log statement occurs.
Definition: logging.h:153
std::ostream & operator<<(std::ostream &out, const ActionPtr &action)
Convenience: dereferences the ActionPtr to Action.
Definition: action.h:546
std::array< einhard::Logger<>, std::tuple_size< LogArea::AreaTuple >::value > logg
An array that stores all pre-configured Logger objects.
Definition: logging.cc:40
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:217
std::unique_ptr< OutputInterface > create_oscar_output(const std::string &format, const std::string &content, const std::filesystem::path &path, const OutputParameters &out_par)
Definition: oscaroutput.cc:889
std::unique_ptr< OutputInterface > create_binary_output(const std::string &format, const std::string &content, const std::filesystem::path &path, const OutputParameters &out_par)
Create a binary output object.
#define likely(x)
Tell the branch predictor that this expression is likely true.
Definition: macros.h:14
constexpr int p
Proton.
constexpr int n
Neutron.
Engine::result_type advance()
Advance the engine's state and return the generated value.
Definition: random.h:78
void set_seed(T &&seed)
Sets the seed of the random number engine.
Definition: random.h:71
Definition: action.h:24
void validate_duplicate_IC_config(double, std::optional< double >, std::string)
The physics inputs for Initial Conditions are currently duplicated in both Output and Collider sectio...
Definition: experiment.cc:654
static constexpr int LInitialConditions
Definition: experiment.h:93
void update_momenta(std::vector< Particles > &particles, double dt, const Potentials &pot, RectangularLattice< std::pair< ThreeVector, ThreeVector >> *FB_lat, RectangularLattice< std::pair< ThreeVector, ThreeVector >> *FI3_lat, RectangularLattice< std::pair< ThreeVector, ThreeVector >> *EM_lat, DensityLattice *jB_lat)
Updates the momenta of all particles at the current time step according to the equations of motion:
Definition: propagation.cc:111
UnaryFunction for_each(Container &&c, UnaryFunction &&f)
Convenience wrapper for std::for_each that operates on a complete container.
Definition: algorithms.h:96
EventInfo fill_event_info(const std::vector< Particles > &ensembles, double E_mean_field, double modus_impact_parameter, const ExperimentParameters &parameters, bool projectile_target_interact, bool kinematic_cut_for_SMASH_IC)
Generate the EventInfo object which is passed to outputs_.
Definition: experiment.cc:595
std::string format_measurements(const std::vector< Particles > &ensembles, uint64_t scatterings_this_interval, const QuantumNumbers &conserved_initial, SystemTimePoint time_start, double time, double E_mean_field, double E_mean_field_initial)
Generate a string which will be printed to the screen when SMASH is running.
Definition: experiment.cc:327
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
void check_interactions_total(uint64_t interactions_total)
Make sure interactions_total can be represented as a 32-bit integer.
Definition: experiment.h:2682
@ Wall
See here for a short description.
@ HyperSurfaceCrossing
See here for a short description.
std::tuple< double, FourVector, ThreeVector, ThreeVector, FourVector, FourVector, FourVector, FourVector > 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:171
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
ExperimentParameters create_experiment_parameters(Configuration &config)
Gathers all general Experiment parameters.
Definition: experiment.cc:132
constexpr double very_small_double
A very small double, used to avoid division by zero.
Definition: constants.h:40
double calculate_mean_field_energy(const Potentials &potentials, RectangularLattice< smash::DensityOnLattice > &jmu_B_lat, RectangularLattice< std::pair< ThreeVector, ThreeVector >> *em_lattice, 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:375
static constexpr int LExperiment
void update_fields_lattice(RectangularLattice< FieldsOnLattice > *fields_lat, RectangularLattice< FourVector > *old_fields, RectangularLattice< FourVector > *new_fields, RectangularLattice< std::array< FourVector, 4 >> *fields_four_grad_lattice, DensityLattice *jmu_B_lat, const LatticeUpdate fields_lat_update, const Potentials &potentials, const double time_step)
Updates the contents on the lattice of FieldsOnLattice type.
Definition: fields.cc:14
bool are_particles_identical_at_given_time(const ParticleData &p1, const ParticleData &p2, double time)
Utility function to compare two ParticleData instances with respect to their PDG code,...
void update_lattice_accumulating_ensembles(RectangularLattice< T > *lat, const LatticeUpdate update, const DensityType dens_type, const DensityParameters &par, const std::vector< Particles > &ensembles, const bool compute_gradient)
Updates the contents on the lattice when ensembles are used.
Definition: density.h:644
void validate_and_adjust_particle_list(ParticleList &particle_list)
Validate a particle list adjusting each particle to be a valid SMASH particle.
Definition: experiment.cc:617
constexpr double nucleon_mass
Nucleon mass in GeV.
Definition: constants.h:58
void update_lattice(RectangularLattice< DensityOnLattice > *lat, RectangularLattice< FourVector > *old_jmu, RectangularLattice< FourVector > *new_jmu, RectangularLattice< std::array< FourVector, 4 >> *four_grad_lattice, const LatticeUpdate update, const DensityType dens_type, const DensityParameters &par, const std::vector< Particles > &ensembles, const double time_step, const bool compute_gradient)
Updates the contents on the lattice of DensityOnLattice type.
Definition: density.cc:186
LatticeUpdate
Enumerator option for lattice updates.
Definition: lattice.h:38
std::ostream & operator<<(std::ostream &out, const Experiment< Modus > &e)
Creates a verbose textual description of the setup of the Experiment.
Definition: experiment.h:722
Potentials * pot_pointer
Pointer to a Potential class.
static constexpr int LMain
Definition: experiment.h:92
constexpr double really_small
Numerical error tolerance.
Definition: constants.h:37
RectangularLattice< FourVector > * UB_lat_pointer
Pointer to the skyrme potential on the lattice.
@ Largest
Make cells as large as possible.
std::chrono::time_point< std::chrono::system_clock > SystemTimePoint
Type (alias) that is used to store the current time.
Definition: chrono.h:22
const std::string hline(113, '-')
String representing a horizontal line.
RectangularLattice< FourVector > * UI3_lat_pointer
Pointer to the symmmetry potential on the lattice.
constexpr double fm2_mb
mb <-> fm^2 conversion factor.
Definition: constants.h:28
Structure to contain custom data for output.
Struct containing the type of the metric and the expansion parameter of the metric.
Definition: propagation.h:26
Exception class that is thrown if an invalid modus is requested from the Experiment factory.
Definition: experiment.h:147
Exception class that is thrown if the requested output path in the Experiment factory is not existing...
Definition: experiment.h:156
Helper structure for Experiment.
double fixed_min_cell_length
Fixed minimal grid cell length (in fm).
const CollisionCriterion coll_crit
Employed collision criterion.
std::unique_ptr< Clock > outputclock
Output clock to keep track of the next output time.
static const Key< double > collTerm_stringParam_powerParticleFormation
See user guide description for more information.
Definition: input_keys.h:2833
static const Key< bool > collTerm_photons_twoToTwoScatterings
See user guide description for more information.
Definition: input_keys.h:3031
static const Key< DensityType > output_densityType
See user guide description for more information.
Definition: input_keys.h:4572
static const Key< bool > modi_collider_collisionWithinNucleus
See user guide description for more information.
Definition: input_keys.h:3163
static const Key< bool > collTerm_noCollisions
See user guide description for more information.
Definition: input_keys.h:2363
static const Key< int > collTerm_photons_fractionalPhotons
See user guide description for more information.
Definition: input_keys.h:3056
static const Key< double > output_initialConditions_rapidityCut
See user guide description for more information.
Definition: input_keys.h:4943
static const Key< int > gen_minNonEmptyEnsembles_number
See user guide description for more information.
Definition: input_keys.h:1173
static const Key< int64_t > gen_randomseed
See user guide description for more information.
Definition: input_keys.h:1148
static const Key< bool > lattice_periodic
See user guide description for more information.
Definition: input_keys.h:5355
static const Key< int > gen_minNonEmptyEnsembles_maximumEnsembles
See user guide description for more information.
Definition: input_keys.h:1161
static const Key< double > gen_expansionRate
See user guide description for more information.
Definition: input_keys.h:1292
static const Key< bool > collTerm_dileptons_decays
See user guide description for more information.
Definition: input_keys.h:3019
static const Key< bool > collTerm_forceDecaysAtEnd
See user guide description for more information.
Definition: input_keys.h:2144
static const Key< double > gen_endTime
See user guide description for more information.
Definition: input_keys.h:1093
static const Key< bool > collTerm_photons_bremsstrahlung
See user guide description for more information.
Definition: input_keys.h:3043
static const Key< PseudoResonance > collTerm_pseudoresonance
See user guide description for more information.
Definition: input_keys.h:2413
static const Key< int > gen_nevents
See user guide description for more information.
Definition: input_keys.h:1135
static const Key< ExpansionMode > gen_metricType
See user guide description for more information.
Definition: input_keys.h:1366
static const Key< std::array< double, 3 > > lattice_origin
See user guide description for more information.
Definition: input_keys.h:5338
static const Key< double > output_initialConditions_properTime
See user guide description for more information.
Definition: input_keys.h:4906
static const Key< std::array< double, 3 > > lattice_sizes
See user guide description for more information.
Definition: input_keys.h:5382
static const Key< TotalCrossSectionStrategy > collTerm_totXsStrategy
See user guide description for more information.
Definition: input_keys.h:2519
static const Key< std::array< int, 3 > > lattice_cellNumber
See user guide description for more information.
Definition: input_keys.h:5319
static const Key< double > output_initialConditions_lowerBound
See user guide description for more information.
Definition: input_keys.h:4886
static const Key< bool > gen_useGrid
See user guide description for more information.
Definition: input_keys.h:1525
static const Key< double > output_initialConditions_pTCut
See user guide description for more information.
Definition: input_keys.h:4924
static const Key< bool > lattice_automatic
See user guide description for more information.
Definition: input_keys.h:5305
static const Key< TimeStepMode > gen_timeStepMode
See user guide description for more information.
Definition: input_keys.h:1498
static const Section output
Section for the output information.
Definition: input_keys.h:140
static const Section potentials
Section for the potentials information.
Definition: input_keys.h:164
static const Section c_pauliBlocking
Subsection for the Pauli blocking mechanism.
Definition: input_keys.h:58
static const Section forcedThermalization
Section for the forced thermalization.
Definition: input_keys.h:74
static const Section modi
Section for the modus specific information.
Definition: input_keys.h:89
static const Section o_initialConditions
Subsection for the output initial conditions content.
Definition: input_keys.h:149
static const Section lattice
Section for the lattice.
Definition: input_keys.h:83
static const Section g_minEnsembles
Subsection for the minimum-nonempty-ensembles mechanism.
Definition: input_keys.h:79
static const std::vector< std::string > oscar2013
Quantities output in OSCAR2013 format.
static const std::vector< std::string > oscar2013extended
Quantities output in Extended OSCAR2013 format.
Helper structure for Experiment to hold output options and parameters.
RivetOutputParameters rivet_parameters
Rivet specfic parameters.