Version: SMASH-3.3
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 "fluidizationaction.h"
28 #include "fourvector.h"
29 #include "grandcan_thermalizer.h"
30 #include "grid.h"
32 #include "icparameters.h"
33 #include "numeric_cast.h"
34 #include "outputparameters.h"
35 #include "pauliblocking.h"
36 #include "potential_globals.h"
37 #include "potentials.h"
38 #include "propagation.h"
39 #include "quantumnumbers.h"
40 #include "scatteractionphoton.h"
41 #include "scatteractionsfinder.h"
42 #include "stringprocess.h"
43 #include "thermalizationaction.h"
44 // Output
45 #include "binaryoutput.h"
46 #ifdef SMASH_USE_HEPMC
47 #include "hepmcoutput.h"
48 #endif
49 #ifdef SMASH_USE_RIVET
50 #include "rivetoutput.h"
51 #endif
52 #include "icoutput.h"
53 #include "oscaroutput.h"
55 #include "thermodynamicoutput.h"
56 #ifdef SMASH_USE_ROOT
57 #include "rootoutput.h"
58 #endif
59 #include "freeforallaction.h"
60 #include "vtkoutput.h"
61 #include "wallcrossingaction.h"
62 
63 namespace std {
74 template <typename T, typename Ratio>
75 static ostream &operator<<(ostream &out,
76  const chrono::duration<T, Ratio> &seconds) {
77  using Seconds = chrono::duration<double>;
78  using Minutes = chrono::duration<double, std::ratio<60>>;
79  using Hours = chrono::duration<double, std::ratio<60 * 60>>;
80  constexpr Minutes threshold_for_minutes{10};
81  constexpr Hours threshold_for_hours{3};
82  if (seconds < threshold_for_minutes) {
83  return out << Seconds(seconds).count() << " [s]";
84  }
85  if (seconds < threshold_for_hours) {
86  return out << Minutes(seconds).count() << " [min]";
87  }
88  return out << Hours(seconds).count() << " [h]";
89 }
90 } // namespace std
91 
92 namespace smash {
93 static constexpr int LMain = LogArea::Main::id;
94 static constexpr int LInitialConditions = LogArea::InitialConditions::id;
95 
104  public:
105  ExperimentBase() = default;
110  virtual ~ExperimentBase() = default;
111 
132  static std::unique_ptr<ExperimentBase> create(
133  Configuration &config, const std::filesystem::path &output_path);
134 
141  virtual void run() = 0;
142 
148  struct InvalidModusRequest : public std::invalid_argument {
149  using std::invalid_argument::invalid_argument;
150  };
151 
157  struct NonExistingOutputPathRequest : public std::invalid_argument {
158  using std::invalid_argument::invalid_argument;
159  };
160 };
161 
162 template <typename Modus>
163 class Experiment;
164 template <typename Modus>
165 std::ostream &operator<<(std::ostream &out, const Experiment<Modus> &e);
166 
191 template <typename Modus>
192 class Experiment : public ExperimentBase {
193  friend class ExperimentBase;
194 
195  public:
202  void run() override;
203 
217  explicit Experiment(Configuration &config,
218  const std::filesystem::path &output_path);
219 
225  void initialize_new_event();
226 
244  void run_time_evolution(const double t_end, ParticleList &&add_plist = {},
245  ParticleList &&remove_plist = {});
246 
252  void do_final_interactions();
253 
255  void final_output();
256 
263  std::vector<Particles> *all_ensembles() { return &ensembles_; }
264 
269  Modus *modus() { return &modus_; }
270 
275  void increase_event_number();
276 
277  private:
290  bool perform_action(Action &action, int i_ensemble,
291  bool include_pauli_blocking = true);
300  void create_output(const std::string &format, const std::string &content,
301  const std::filesystem::path &output_path,
302  const OutputParameters &par);
303 
311  void propagate_and_shine(double to_time, Particles &particles);
312 
325  void run_time_evolution_timestepless(Actions &actions, int i_ensemble,
326  const double end_time_propagation);
327 
329  void intermediate_output();
330 
332  void update_potentials();
333 
342  double compute_min_cell_length(double dt) const {
345  }
346  return std::sqrt(4 * dt * dt + max_transverse_distance_sqr_);
347  }
348 
350  double next_output_time() const {
351  return parameters_.outputclock->next_time();
352  }
353 
359 
365  bool is_finished();
366 
373 
376 
381  Modus modus_;
382 
384  std::vector<Particles> ensembles_;
385 
390  std::unique_ptr<Potentials> potentials_;
391 
396  std::unique_ptr<PauliBlocker> pauli_blocker_;
397 
402  OutputsList outputs_;
403 
405  OutputPtr dilepton_output_;
406 
408  OutputPtr photon_output_;
409 
414  std::vector<bool> projectile_target_interact_;
415 
421  std::vector<FourVector> beam_momentum_ = {};
422 
424  std::vector<std::unique_ptr<ActionFinderInterface>> action_finders_;
425 
427  std::unique_ptr<DecayActionsFinderDilepton> dilepton_finder_;
428 
430  std::unique_ptr<ActionFinderInterface> photon_finder_;
431 
434 
436  std::unique_ptr<DensityLattice> j_QBS_lat_;
437 
439  std::unique_ptr<DensityLattice> jmu_B_lat_;
440 
442  std::unique_ptr<DensityLattice> jmu_I3_lat_;
443 
445  std::unique_ptr<DensityLattice> jmu_el_lat_;
446 
448  std::unique_ptr<FieldsLattice> fields_lat_;
449 
457  std::unique_ptr<DensityLattice> jmu_custom_lat_;
458 
461 
466  std::unique_ptr<RectangularLattice<FourVector>> UB_lat_ = nullptr;
467 
472  std::unique_ptr<RectangularLattice<FourVector>> UI3_lat_ = nullptr;
473 
478  std::unique_ptr<RectangularLattice<std::pair<ThreeVector, ThreeVector>>>
480 
482  std::unique_ptr<RectangularLattice<std::pair<ThreeVector, ThreeVector>>>
484 
486  std::unique_ptr<RectangularLattice<std::pair<ThreeVector, ThreeVector>>>
488 
490  std::unique_ptr<RectangularLattice<EnergyMomentumTensor>> Tmn_;
491 
493  std::unique_ptr<RectangularLattice<FourVector>> old_jmu_auxiliary_;
495  std::unique_ptr<RectangularLattice<FourVector>> new_jmu_auxiliary_;
497  std::unique_ptr<RectangularLattice<std::array<FourVector, 4>>>
499 
501  std::unique_ptr<RectangularLattice<FourVector>> old_fields_auxiliary_;
503  std::unique_ptr<RectangularLattice<FourVector>> new_fields_auxiliary_;
505  std::unique_ptr<RectangularLattice<std::array<FourVector, 4>>>
507 
509  bool printout_rho_eckart_ = false;
510 
512  bool printout_tmn_ = false;
513 
515  bool printout_tmn_landau_ = false;
516 
518  bool printout_v_landau_ = false;
519 
521  bool printout_j_QBS_ = false;
522 
524  bool printout_lattice_td_ = false;
525 
529 
531  std::unique_ptr<GrandCanThermalizer> thermalizer_;
532 
538 
553  int nevents_ = 0;
554 
564 
572 
574  int event_ = 0;
575 
578 
583  int max_events_ = 0;
584 
586  const double end_time_;
587 
593  const double delta_time_startup_;
594 
599  const bool force_decays_;
600 
602  const bool use_grid_;
603 
606 
608  const bool dileptons_switch_;
609 
611  const bool photons_switch_;
612 
615 
620  const bool IC_switch_;
621 
623  const bool IC_dynamic_;
624 
627 
632  double max_transverse_distance_sqr_ = std::numeric_limits<double>::max();
633 
642 
648 
650  SystemTimePoint time_start_ = SystemClock::now();
651 
654 
659  uint64_t interactions_total_ = 0;
660 
666 
671  uint64_t wall_actions_total_ = 0;
672 
678 
683  uint64_t total_pauli_blocked_ = 0;
684 
690 
696 
700  double total_energy_removed_ = 0.0;
701 
706 
709 
711  int64_t seed_ = -1;
712 
718  friend std::ostream &operator<<<>(std::ostream &out, const Experiment &e);
719 };
720 
722 template <typename Modus>
723 std::ostream &operator<<(std::ostream &out, const Experiment<Modus> &e) {
724  out << "End time: " << e.end_time_ << " fm\n";
725  out << e.modus_;
726  return out;
727 }
728 
729 template <typename Modus>
730 void Experiment<Modus>::create_output(const std::string &format,
731  const std::string &content,
732  const std::filesystem::path &output_path,
733  const OutputParameters &out_par) {
734  // Disable output which do not properly work with multiple ensembles
735  if (ensembles_.size() > 1) {
736  auto abort_because_of = [](const std::string &s) {
737  throw std::invalid_argument(
738  s + " output is not available with multiple parallel ensembles.");
739  };
740  if (content == "Initial_Conditions") {
741  abort_because_of("Initial_Conditions");
742  }
743  if ((format == "HepMC") || (format == "HepMC_asciiv3") ||
744  (format == "HepMC_treeroot")) {
745  abort_because_of("HepMC");
746  }
747  if (content == "Rivet") {
748  abort_because_of("Rivet");
749  }
750  if (content == "Collisions") {
751  logg[LExperiment].warn(
752  "Information coming from different ensembles in 'Collisions' output "
753  "is not distinguishable.\nSuch an output with multiple parallel "
754  "ensembles should only be used if later in the data analysis\nit is "
755  "not necessary to trace back which data belongs to which ensemble.");
756  }
757  }
758 
759  if (format == "VTK" && content == "Particles") {
760  outputs_.emplace_back(
761  std::make_unique<VtkOutput>(output_path, content, out_par));
762  } else if (format == "Root") {
763 #ifdef SMASH_USE_ROOT
764  if (content == "Initial_Conditions") {
765  outputs_.emplace_back(
766  std::make_unique<RootOutput>(output_path, "SMASH_IC", out_par));
767  } else {
768  outputs_.emplace_back(
769  std::make_unique<RootOutput>(output_path, content, out_par));
770  }
771 #else
772  logg[LExperiment].error(
773  "Root output requested, but Root support not compiled in");
774 #endif
775  } else if ((format == "Binary" || format == "Oscar2013_bin") &&
776  (content == "Collisions" || content == "Particles" ||
777  content == "Dileptons" || content == "Photons" ||
778  content == "Initial_Conditions")) {
779  outputs_.emplace_back(
780  create_binary_output(format, content, output_path, out_par));
781  } else if (format == "Oscar1999" || format == "Oscar2013") {
782  outputs_.emplace_back(
783  create_oscar_output(format, content, output_path, out_par));
784  } else if (format == "ASCII" &&
785  (content == "Particles" || content == "Collisions" ||
786  content == "Dileptons" || content == "Photons" ||
787  content == "Initial_Conditions")) {
788  outputs_.emplace_back(
789  create_oscar_output(format, content, output_path, out_par));
790  } else if (content == "Thermodynamics" && format == "ASCII") {
791  outputs_.emplace_back(
792  std::make_unique<ThermodynamicOutput>(output_path, content, out_par));
793  } else if (content == "Thermodynamics" &&
794  (format == "Lattice_ASCII" || format == "Lattice_Binary")) {
795  printout_full_lattice_any_td_ = true;
796  outputs_.emplace_back(std::make_unique<ThermodynamicLatticeOutput>(
797  output_path, content, out_par, format == "Lattice_ASCII",
798  format == "Lattice_Binary"));
799  } else if (content == "Thermodynamics" && format == "VTK") {
800  printout_lattice_td_ = true;
801  outputs_.emplace_back(
802  std::make_unique<VtkOutput>(output_path, content, out_par));
803  } else if (content == "Initial_Conditions" && format == "For_vHLLE") {
804  if (IC_dynamic_) {
805  throw std::invalid_argument(
806  "Dynamic initial conditions are only available in Oscar2013 and "
807  "Binary formats.");
808  }
809  outputs_.emplace_back(
810  std::make_unique<ICOutput>(output_path, "SMASH_IC_For_vHLLE", out_par));
811  } else if ((format == "HepMC") || (format == "HepMC_asciiv3") ||
812  (format == "HepMC_treeroot")) {
813 #ifdef SMASH_USE_HEPMC
814  if (content == "Particles") {
815  if ((format == "HepMC") || (format == "HepMC_asciiv3")) {
816  outputs_.emplace_back(std::make_unique<HepMcOutput>(
817  output_path, "SMASH_HepMC_particles", false, "asciiv3"));
818  } else if (format == "HepMC_treeroot") {
819 #ifdef SMASH_USE_HEPMC_ROOTIO
820  outputs_.emplace_back(std::make_unique<HepMcOutput>(
821  output_path, "SMASH_HepMC_particles", false, "root"));
822 #else
823  logg[LExperiment].error(
824  "Requested HepMC_treeroot output not available, "
825  "ROOT or HepMC3_ROOTIO missing or not found by cmake.");
826 #endif
827  }
828  } else if (content == "Collisions") {
829  if ((format == "HepMC") || (format == "HepMC_asciiv3")) {
830  outputs_.emplace_back(std::make_unique<HepMcOutput>(
831  output_path, "SMASH_HepMC_collisions", true, "asciiv3"));
832  } else if (format == "HepMC_treeroot") {
833 #ifdef SMASH_USE_HEPMC_ROOTIO
834  outputs_.emplace_back(std::make_unique<HepMcOutput>(
835  output_path, "SMASH_HepMC_collisions", true, "root"));
836 #else
837  logg[LExperiment].error(
838  "Requested HepMC_treeroot output not available, "
839  "ROOT or HepMC3_ROOTIO missing or not found by cmake.");
840 #endif
841  }
842  } else {
843  logg[LExperiment].error(
844  "HepMC only available for Particles and "
845  "Collisions content. Requested for " +
846  content + ".");
847  }
848 #else
849  logg[LExperiment].error(
850  "HepMC output requested, but HepMC support not compiled in");
851 #endif
852  } else if (content == "Coulomb" && format == "VTK") {
853  outputs_.emplace_back(
854  std::make_unique<VtkOutput>(output_path, "Fields", out_par));
855  } else if (content == "Rivet") {
856 #ifdef SMASH_USE_RIVET
857  // flag to ensure that the Rivet format has not been already assigned
858  static bool rivet_format_already_selected = false;
859  // if the next check is true, then we are trying to assign the format twice
860  if (rivet_format_already_selected) {
861  logg[LExperiment].warn(
862  "Rivet output format can only be one, either YODA or YODA-full. "
863  "Only your first valid choice will be used.");
864  return;
865  }
866  if (format == "YODA") {
867  outputs_.emplace_back(std::make_unique<RivetOutput>(
868  output_path, "SMASH_Rivet", false, out_par.rivet_parameters));
869  rivet_format_already_selected = true;
870  } else if (format == "YODA-full") {
871  outputs_.emplace_back(std::make_unique<RivetOutput>(
872  output_path, "SMASH_Rivet_full", true, out_par.rivet_parameters));
873  rivet_format_already_selected = true;
874  } else {
875  logg[LExperiment].error("Rivet format " + format +
876  "not one of YODA or YODA-full");
877  }
878 #else
879  logg[LExperiment].error(
880  "Rivet output requested, but Rivet support not compiled in");
881 #endif
882  } else {
883  logg[LExperiment].error()
884  << "Unknown combination of format (" << format << ") and content ("
885  << content << "). Fix the config.";
886  }
887 
888  logg[LExperiment].info() << "Added output " << content << " of format "
889  << format << "\n";
890 }
891 
900 
901 template <typename Modus>
903  const std::filesystem::path &output_path)
904  : parameters_(create_experiment_parameters(config)),
905  density_param_(DensityParameters(parameters_)),
906  modus_(std::invoke([&]() {
907  /* This immediately invoked lambda is a work-around to cope with the
908  * fact that the "Collisions_Within_Nucleus" key belongs to the
909  * "Collider" section, but is used by the ScatterActionsFinder through
910  * the ScatterActionsFinderParameters member. Here that key is taken
911  * from the main configuration and put there back after the "Collider"
912  * section is extracted. If this were not done in this way, the
913  * sub-configuration given to ColliderModus would be deleted not empty
914  * at the end of its constructor and this would throw an exception.*/
916  const bool restore_key = config.has_value(key);
917  const bool temporary_taken_key = config.take(key);
918  auto modus_config =
920  if (restore_key) {
921  config.set_value(key, temporary_taken_key);
922  }
923  return Modus{std::move(modus_config), parameters_};
924  })),
925  ensembles_(parameters_.n_ensembles),
926  end_time_(config.take(InputKeys::gen_endTime)),
927  delta_time_startup_(parameters_.labclock->timestep_duration()),
928  force_decays_(config.take(InputKeys::collTerm_forceDecaysAtEnd)),
929  use_grid_(config.take(InputKeys::gen_useGrid)),
930  metric_(config.take(InputKeys::gen_metricType),
932  dileptons_switch_(config.take(InputKeys::collTerm_dileptons_decays)),
933  photons_switch_(
935  bremsstrahlung_switch_(
937  IC_switch_(config.has_section(InputSections::o_initialConditions) &&
938  modus_.is_IC_for_hybrid()),
939  IC_dynamic_(IC_switch_ ? (modus_.IC_parameters().type ==
941  : false),
942  time_step_mode_(config.take(InputKeys::gen_timeStepMode)) {
943  logg[LExperiment].info() << *this;
944 
945  const bool user_wants_nevents = config.has_value(InputKeys::gen_nevents);
946  const bool user_wants_min_nonempty =
948  if (user_wants_nevents == user_wants_min_nonempty) {
949  throw std::invalid_argument(
950  "Please specify either Nevents or Minimum_Nonempty_Ensembles.");
951  }
952  if (user_wants_nevents) {
953  event_counting_ = EventCounting::FixedNumber;
954  nevents_ = config.take(InputKeys::gen_nevents);
955  } else {
956  event_counting_ = EventCounting::MinimumNonEmpty;
957  minimum_nonempty_ensembles_ =
959  int max_ensembles =
961  max_events_ = numeric_cast<int>(std::ceil(
962  static_cast<double>(max_ensembles) / parameters_.n_ensembles));
963  }
964 
965  // covariant derivatives can only be done with covariant smearing
966  if (parameters_.derivatives_mode == DerivativesMode::CovariantGaussian &&
967  parameters_.smearing_mode != SmearingMode::CovariantGaussian) {
968  throw std::invalid_argument(
969  "Covariant Gaussian derivatives only make sense for Covariant Gaussian "
970  "smearing!");
971  }
972 
973  // for triangular smearing:
974  // the weight needs to be larger than 1./7. for the center cell to contribute
975  // more than the surrounding cells
976  if (parameters_.smearing_mode == SmearingMode::Discrete &&
977  parameters_.discrete_weight < (1. / 7.)) {
978  throw std::invalid_argument(
979  "The central weight for discrete smearing should be >= 1./7.");
980  }
981 
982  if (parameters_.coll_crit == CollisionCriterion::Stochastic &&
983  (time_step_mode_ != TimeStepMode::Fixed || !use_grid_)) {
984  throw std::invalid_argument(
985  "The stochastic criterion can only be employed for fixed time step "
986  "mode and with a grid!");
987  }
988 
989  if (modus_.is_box() && (time_step_mode_ != TimeStepMode::Fixed)) {
990  throw std::invalid_argument(
991  "The box modus can only be used with the fixed time step mode!");
992  }
993 
994  logg[LExperiment].info("Using ", parameters_.testparticles,
995  " testparticles per particle.");
996  logg[LExperiment].info("Using ", parameters_.n_ensembles,
997  " parallel ensembles.");
998 
999  if (modus_.is_box() && config.read(InputKeys::collTerm_totXsStrategy) !=
1001  logg[LExperiment].warn(
1002  "To preserve detailed balance in a box simulation, it is recommended\n"
1003  "to use the bottom-up strategy for evaluating total cross sections.\n"
1004  "Consider adding the following line to the 'Collision_Term' section "
1005  "in your configuration file:\n"
1006  " Total_Cross_Section_Strategy: \"BottomUp\"");
1007  }
1008  if (modus_.is_box() && config.read(InputKeys::collTerm_pseudoresonance) !=
1010  logg[LExperiment].warn(
1011  "To preserve detailed balance in a box simulation, it is recommended "
1012  "to not include the pseudoresonances,\nas they artificially increase "
1013  "the resonance production without changing the corresponding "
1014  "decay.\nConsider adding the following line to the 'Collision_Term' "
1015  "section in your configuration file:\n Pseudoresonance: \"None\"");
1016  }
1017 
1018  const bool IC_output = config.has_section(InputSections::o_initialConditions);
1019  if (IC_output != modus_.is_IC_for_hybrid()) {
1020  throw std::invalid_argument(
1021  "The 'Initial_Conditions' subsection must be present in both 'Output' "
1022  "and 'Modi: Collider' sections.");
1023  }
1024 
1025  /* In collider setup with sqrts >= 200 GeV particles don't form continuously
1026  *
1027  * NOTE: This key has to be taken before the ScatterActionsFinder is created
1028  * because there the "String_Parameters" is extracted as sub-config and
1029  * all parameters but this one are taken. If this one is still there
1030  * the configuration temporary object will be destroyed not empty, hence
1031  * throwing an exception.
1032  */
1035  modus_.sqrt_s_NN() >= 200. ? -1. : 1.);
1036 
1037  // create finders
1038  if (dileptons_switch_) {
1039  dilepton_finder_ = std::make_unique<DecayActionsFinderDilepton>();
1040  }
1041  if (photons_switch_ || bremsstrahlung_switch_) {
1042  n_fractional_photons_ =
1044  }
1045  if (parameters_.two_to_one) {
1046  if (parameters_.res_lifetime_factor < 0.) {
1047  throw std::invalid_argument(
1048  "Resonance lifetime modifier cannot be negative!");
1049  }
1050  if (parameters_.res_lifetime_factor < really_small) {
1051  logg[LExperiment].warn(
1052  "Resonance lifetime set to zero. Make sure resonances cannot "
1053  "interact",
1054  "inelastically (e.g. resonance chains), else SMASH is known to "
1055  "hang.");
1056  }
1057  action_finders_.emplace_back(std::make_unique<DecayActionsFinder>(
1058  parameters_.res_lifetime_factor, parameters_.do_non_strong_decays,
1059  force_decays_, parameters_.spin_interaction_type));
1060  }
1061  bool no_coll = config.take(InputKeys::collTerm_noCollisions);
1062  if ((parameters_.two_to_one || parameters_.included_2to2.any() ||
1063  parameters_.included_multi.any() || parameters_.strings_switch) &&
1064  !no_coll) {
1065  parameters_.use_monash_tune_default =
1066  (modus_.is_collider() && modus_.sqrt_s_NN() >= 200.);
1067  auto scat_finder =
1068  std::make_unique<ScatterActionsFinder>(config, parameters_);
1069  max_transverse_distance_sqr_ =
1070  scat_finder->max_transverse_distance_sqr(parameters_.testparticles);
1071  process_string_ptr_ = scat_finder->get_process_string_ptr();
1072  action_finders_.emplace_back(std::move(scat_finder));
1073  } else {
1074  max_transverse_distance_sqr_ =
1075  parameters_.maximum_cross_section / M_PI * fm2_mb;
1076  process_string_ptr_ = NULL;
1077  }
1078  if (modus_.is_box()) {
1079  action_finders_.emplace_back(
1080  std::make_unique<WallCrossActionsFinder>(parameters_.box_length));
1081  }
1082 
1083  if (IC_switch_) {
1084  const InitialConditionParameters &IC_parameters = modus_.IC_parameters();
1085  if (IC_dynamic_) {
1086  // Dynamic fluidization
1087  action_finders_.emplace_back(std::make_unique<DynamicFluidizationFinder>(
1088  modus_.fluid_lattice(), modus_.fluid_background(), IC_parameters));
1089  } else {
1090  // Iso-tau hypersurface
1091  double rapidity_cut = IC_parameters.rapidity_cut.value();
1092 
1093  if (modus_.calculation_frame_is_fixed_target() && rapidity_cut != 0.0) {
1094  throw std::runtime_error(
1095  "Rapidity cut for initial conditions output is not implemented "
1096  "in the fixed target calculation frame. \nPlease use "
1097  "\"center of velocity\" or \"center of mass\" as a "
1098  "\"Calculation_Frame\" instead.");
1099  }
1100 
1101  double pT_cut = IC_parameters.pT_cut.value();
1102  if (rapidity_cut > 0.0 || pT_cut > 0.0) {
1103  kinematic_cuts_for_IC_output_ = true;
1104  }
1105 
1106  const double proper_time = std::invoke([&]() {
1107  if (IC_parameters.proper_time.has_value()) {
1108  return IC_parameters.proper_time.value();
1109  } else {
1110  // Scaling factor applied to the switching time and the lower bound
1111  const double scaling = IC_parameters.proper_time_scaling.value();
1112  // Lower bound for the switching time
1113  const double lower_bound =
1114  IC_parameters.lower_bound.value() * scaling;
1115  // Default proper time is the passing time of the two nuclei
1116  const double default_proper_time =
1117  modus_.nuclei_passing_time() * scaling;
1118  if (default_proper_time >= lower_bound) {
1119  logg[LInitialConditions].info()
1120  << "Nuclei passing time is " << default_proper_time << " fm.";
1121  return default_proper_time;
1122  } else {
1123  logg[LInitialConditions].warn()
1124  << "Nuclei passing time is too short, hypersurface proper time "
1125  << "set to tau = " << lower_bound << " fm.";
1126  return lower_bound;
1127  }
1128  }
1129  });
1130 
1131  action_finders_.emplace_back(
1132  std::make_unique<HyperSurfaceCrossActionsFinder>(
1133  proper_time, rapidity_cut, pT_cut));
1134  }
1135  }
1136 
1138  logg[LExperiment].info() << "Pauli blocking is ON.";
1139  pauli_blocker_ = std::make_unique<PauliBlocker>(
1142  parameters_);
1143  }
1144 
1439  // create outputs
1441  " create OutputInterface objects");
1442  dens_type_ = config.take(InputKeys::output_densityType);
1443  logg[LExperiment].debug()
1444  << "Density type printed to headers: " << dens_type_;
1445 
1446  /* Parse configuration about output contents and formats, doing all logical
1447  * checks about specified formats, creating all needed output objects. Note
1448  * that we first extract the output sub configuration without the "Output:"
1449  * enclosing section to easily get all output contents and then we reintroduce
1450  * it for the actual parsing (remember we parse database keys which have
1451  * labels from the top-level only).
1452  */
1453  auto output_conf = config.extract_sub_configuration(
1455  if (output_path == "") {
1456  throw std::invalid_argument(
1457  "Invalid empty output path provided to Experiment constructor.");
1458  } else if (!std::filesystem::exists(output_path)) {
1459  logg[LExperiment].fatal(
1460  "Output path \"" + output_path.string() +
1461  "\" used to create an Experiment object does not exist.");
1462  throw NonExistingOutputPathRequest("Attempt to use not existing path.");
1463  } else if (!std::filesystem::is_directory(output_path)) {
1464  logg[LExperiment].fatal("Output path \"" + output_path.string() +
1465  "\" used to create an Experiment object "
1466  "exists, but it is not a directory.");
1467  throw std::logic_error("Attempt to use invalid existing path.");
1468  }
1469  const std::vector<std::string> output_contents =
1470  output_conf.list_upmost_nodes();
1471  if (output_conf.is_empty()) {
1472  logg[LExperiment].warn() << "No \"Output\" section found in the input "
1473  "file. No output file will be produced.";
1474  } else {
1475  output_conf.enclose_into_section(InputSections::output);
1476  }
1477  std::vector<std::vector<std::string>> list_of_formats(output_contents.size());
1478  std::transform(
1479  output_contents.cbegin(), output_contents.cend(), list_of_formats.begin(),
1480  [&output_conf](std::string content) -> std::vector<std::string> {
1481  /* Note that the "Format" key has an empty list as default, although it
1482  * is a required key, because then here below the error for the user is
1483  * more informative, if the key was not given in the input file. */
1484  return output_conf.take(InputKeys::get_output_format_key(content));
1485  });
1486  auto abort_because_of_invalid_input_file = []() {
1487  throw std::invalid_argument("Invalid configuration input file.");
1488  };
1489  const OutputParameters output_parameters(std::move(output_conf));
1490  for (std::size_t i = 0; i < output_contents.size(); ++i) {
1491  if (output_contents[i] == "Particles" ||
1492  output_contents[i] == "Collisions" ||
1493  output_contents[i] == "Dileptons" || output_contents[i] == "Photons" ||
1494  output_contents[i] == "Initial_Conditions") {
1495  assert(output_parameters.quantities.count(output_contents[i]) > 0);
1496  const bool quantities_given_nonempty =
1497  !output_parameters.quantities.at(output_contents[i]).empty();
1498  auto formats_contains = [&list_of_formats, &i](const std::string &label) {
1499  return std::find(list_of_formats[i].begin(), list_of_formats[i].end(),
1500  label) != list_of_formats[i].end();
1501  };
1502  const bool custom_ascii_requested = formats_contains("ASCII");
1503  const bool custom_binary_requested = formats_contains("Binary");
1504  const bool custom_requested =
1505  custom_ascii_requested || custom_binary_requested;
1506  const bool oscar2013_requested = formats_contains("Oscar2013");
1507  const bool oscar2013_bin_requested = formats_contains("Oscar2013_bin");
1508  const bool is_extended = (output_contents[i] == "Particles")
1509  ? output_parameters.part_extended
1510  : output_parameters.coll_extended;
1511  const auto &default_quantities =
1514  const bool are_given_quantities_oscar2013_ones =
1515  output_parameters.quantities.at(output_contents[i]) ==
1516  default_quantities;
1517  if (quantities_given_nonempty != custom_requested) {
1518  logg[LExperiment].fatal()
1519  << "Non-empty \"Quantities\" and \"ASCII\"/\"Binary\" format have "
1520  << "not been specified both for " << std::quoted(output_contents[i])
1521  << " in config file.";
1522  abort_because_of_invalid_input_file();
1523  }
1524  if (custom_ascii_requested && oscar2013_requested &&
1525  are_given_quantities_oscar2013_ones) {
1526  logg[LExperiment].fatal()
1527  << "The specified \"Quantities\" for the ASCII format are the same "
1528  "as those of the requested \"Oscar2013\"\nformat for "
1529  << std::quoted(output_contents[i])
1530  << " and this would produce the same output file twice.";
1531  abort_because_of_invalid_input_file();
1532  }
1533  if (custom_binary_requested && oscar2013_bin_requested &&
1534  are_given_quantities_oscar2013_ones) {
1535  logg[LExperiment].fatal()
1536  << "The specified \"Quantities\" for the binary format are the "
1537  "same as those of the requested \"Oscar2013_bin\"\nformat for "
1538  << std::quoted(output_contents[i])
1539  << " and this would produce the same output file twice.";
1540  abort_because_of_invalid_input_file();
1541  }
1542  }
1543 
1544  if (list_of_formats[i].empty()) {
1545  logg[LExperiment].fatal()
1546  << "Empty or unspecified list of formats for "
1547  << std::quoted(output_contents[i]) << " content.";
1548  abort_because_of_invalid_input_file();
1549  } else if (std::find(list_of_formats[i].begin(), list_of_formats[i].end(),
1550  "None") != list_of_formats[i].end()) {
1551  if (list_of_formats[i].size() > 1) {
1552  logg[LExperiment].fatal()
1553  << "Use of \"None\" output format together with other formats is "
1554  "not allowed.\nInvalid \"Format\" key for "
1555  << std::quoted(output_contents[i]) << " content.";
1556  abort_because_of_invalid_input_file();
1557  } else {
1558  // Clear vector so that the for below is skipped and no output created
1559  list_of_formats[i].clear();
1560  }
1561  } else if (std::set<std::string> tmp_set(list_of_formats[i].begin(),
1562  list_of_formats[i].end());
1563  list_of_formats[i].size() != tmp_set.size()) {
1564  auto join_container = [](const auto &container) {
1565  std::string result{};
1566  std::for_each(container.cbegin(), container.cend(),
1567  [&result](const std::string s) {
1568  result += (result == "") ? s : ", " + s;
1569  });
1570  return result;
1571  };
1572  const std::string old_formats = join_container(list_of_formats[i]),
1573  new_formats = join_container(tmp_set);
1574  logg[LExperiment].warn()
1575  << "Found the same output format multiple times for "
1576  << std::quoted(output_contents[i])
1577  << " content. Duplicates will be ignored:\n 'Format: [" << old_formats
1578  << "] -> [" << new_formats << "]'";
1579  list_of_formats[i].assign(tmp_set.begin(), tmp_set.end());
1580  }
1581  }
1582 
1583  /* Repeat loop over output_contents here to create all outputs after having
1584  * validated all content specifications. This is more user-friendly. */
1585  std::size_t total_number_of_requested_formats = 0;
1586  for (std::size_t i = 0; i < output_contents.size(); ++i) {
1587  for (const auto &format : list_of_formats[i]) {
1588  create_output(format, output_contents[i], output_path, output_parameters);
1589  ++total_number_of_requested_formats;
1590  }
1591  }
1592 
1593  if (outputs_.size() != total_number_of_requested_formats) {
1594  logg[LExperiment].fatal()
1595  << "At least one invalid output format has been provided.";
1596  abort_because_of_invalid_input_file();
1597  }
1598 
1599  /* We can take away the Fermi motion flag, because the collider modus is
1600  * already initialized. We only need it when potentials are enabled, but we
1601  * always have to take it, otherwise SMASH will complain about unused
1602  * options. We have to provide a default value for modi other than Collider.
1603  */
1604  if (config.has_section(InputSections::potentials)) {
1605  if (time_step_mode_ == TimeStepMode::None) {
1606  logg[LExperiment].error() << "Potentials only work with time steps!";
1607  throw std::invalid_argument("Can't use potentials without time steps!");
1608  }
1609  if (modus_.fermi_motion() == FermiMotion::Frozen) {
1610  logg[LExperiment].error()
1611  << "Potentials don't work with frozen Fermi momenta! "
1612  "Use normal Fermi motion instead.";
1613  throw std::invalid_argument(
1614  "Can't use potentials "
1615  "with frozen Fermi momenta!");
1616  }
1617  logg[LExperiment].info() << "Potentials are ON. Timestep is "
1618  << parameters_.labclock->timestep_duration();
1619  // potentials need density calculation parameters from parameters_
1620  potentials_ = std::make_unique<Potentials>(
1622  parameters_);
1623  // make sure that vdf potentials are not used together with Skyrme
1624  // or symmetry potentials
1625  if (potentials_->use_skyrme() && potentials_->use_vdf()) {
1626  throw std::runtime_error(
1627  "Can't use Skyrme and VDF potentials at the same time!");
1628  }
1629  if (potentials_->use_symmetry() && potentials_->use_vdf()) {
1630  throw std::runtime_error(
1631  "Can't use symmetry and VDF potentials at the same time!");
1632  }
1633  if (potentials_->use_skyrme()) {
1634  logg[LExperiment].info() << "Skyrme potentials are:\n";
1635  logg[LExperiment].info()
1636  << "\t\tSkyrme_A [MeV] = " << potentials_->skyrme_a() << "\n";
1637  logg[LExperiment].info()
1638  << "\t\tSkyrme_B [MeV] = " << potentials_->skyrme_b() << "\n";
1639  logg[LExperiment].info()
1640  << "\t\t Skyrme_tau = " << potentials_->skyrme_tau() << "\n";
1641  }
1642  if (potentials_->use_symmetry()) {
1643  logg[LExperiment].info()
1644  << "Symmetry potential is:"
1645  << "\n S_pot [MeV] = " << potentials_->symmetry_S_pot() << "\n";
1646  }
1647  if (potentials_->use_vdf()) {
1648  logg[LExperiment].info() << "VDF potential parameters are:\n";
1649  logg[LExperiment].info() << "\t\tsaturation density [fm^-3] = "
1650  << potentials_->saturation_density() << "\n";
1651  for (int i = 0; i < potentials_->number_of_terms(); i++) {
1652  logg[LExperiment].info()
1653  << "\t\tCoefficient_" << i + 1 << " = "
1654  << 1000.0 * (potentials_->coeffs())[i] << " [MeV] \t Power_"
1655  << i + 1 << " = " << (potentials_->powers())[i] << "\n";
1656  }
1657  }
1658  // if potentials are on, derivatives need to be calculated
1659  if (parameters_.derivatives_mode == DerivativesMode::Off &&
1660  parameters_.field_derivatives_mode == FieldDerivativesMode::ChainRule) {
1661  throw std::invalid_argument(
1662  "Derivatives are necessary for running with potentials.\n"
1663  "Derivatives_Mode: \"Off\" only makes sense for "
1664  "Field_Derivatives_Mode: \"Direct\"!\nUse \"Covariant Gaussian\" or "
1665  "\"Finite difference\".");
1666  }
1667  // for computational efficiency, we want to turn off the derivatives of jmu
1668  // and the rest frame density derivatives if direct derivatives are used
1669  if (parameters_.field_derivatives_mode == FieldDerivativesMode::Direct) {
1670  parameters_.derivatives_mode = DerivativesMode::Off;
1671  parameters_.rho_derivatives_mode = RestFrameDensityDerivativesMode::Off;
1672  }
1673  switch (parameters_.derivatives_mode) {
1675  logg[LExperiment].info() << "Covariant Gaussian derivatives are ON";
1676  break;
1678  logg[LExperiment].info() << "Finite difference derivatives are ON";
1679  break;
1680  case DerivativesMode::Off:
1681  logg[LExperiment].info() << "Gradients of baryon current are OFF";
1682  break;
1683  }
1684  switch (parameters_.rho_derivatives_mode) {
1686  logg[LExperiment].info() << "Rest frame density derivatives are ON";
1687  break;
1689  logg[LExperiment].info() << "Rest frame density derivatives are OFF";
1690  break;
1691  }
1692  // direct or chain rule derivatives only make sense for the VDF potentials
1693  if (potentials_->use_vdf()) {
1694  switch (parameters_.field_derivatives_mode) {
1696  logg[LExperiment].info() << "Chain rule field derivatives are ON";
1697  break;
1699  logg[LExperiment].info() << "Direct field derivatives are ON";
1700  break;
1701  }
1702  }
1703  /*
1704  * Necessary safety checks
1705  */
1706  // VDF potentials need derivatives of rest frame density or fields
1707  if (potentials_->use_vdf() && (parameters_.rho_derivatives_mode ==
1709  parameters_.field_derivatives_mode ==
1711  throw std::runtime_error(
1712  "Can't use VDF potentials without rest frame density derivatives or "
1713  "direct field derivatives!");
1714  }
1715  // potentials require using gradients
1716  if (parameters_.derivatives_mode == DerivativesMode::Off &&
1717  parameters_.field_derivatives_mode == FieldDerivativesMode::ChainRule) {
1718  throw std::runtime_error(
1719  "Can't use potentials without gradients of baryon current (Skyrme, "
1720  "VDF)"
1721  " or direct field derivatives (VDF)!");
1722  }
1723  // direct field derivatives only make sense for the VDF potentials
1724  if (!(potentials_->use_vdf()) &&
1725  parameters_.field_derivatives_mode == FieldDerivativesMode::Direct) {
1726  throw std::invalid_argument(
1727  "Field_Derivatives_Mode: \"Direct\" only makes sense for the VDF "
1728  "potentials!\nUse Field_Derivatives_Mode: \"Chain Rule\" or comment "
1729  "this option out (Chain Rule is default)");
1730  }
1731  }
1732 
1733  // information about the type of smearing
1734  switch (parameters_.smearing_mode) {
1736  logg[LExperiment].info() << "Smearing type: Covariant Gaussian";
1737  break;
1739  logg[LExperiment].info() << "Smearing type: Discrete with weight = "
1740  << parameters_.discrete_weight;
1741  break;
1743  logg[LExperiment].info() << "Smearing type: Triangular with range = "
1744  << parameters_.triangular_range;
1745  break;
1746  }
1747 
1748  // Create lattices
1749  if (config.has_section(InputSections::lattice)) {
1750  bool automatic = config.take(InputKeys::lattice_automatic);
1751  bool all_geometrical_properties_specified =
1755  if (!automatic && !all_geometrical_properties_specified) {
1756  throw std::invalid_argument(
1757  "The lattice was requested to be manually generated, but some\n"
1758  "lattice geometrical property was not specified. Be sure to provide\n"
1759  "both \"Cell_Number\" and \"Origin\" and \"Sizes\".");
1760  }
1761  if (automatic && all_geometrical_properties_specified) {
1762  throw std::invalid_argument(
1763  "The lattice was requested to be automatically generated, but all\n"
1764  "lattice geometrical properties were specified. In this case you\n"
1765  "need to set \"Automatic: False\".");
1766  }
1767  bool periodic = config.take(InputKeys::lattice_periodic, modus_.is_box());
1768  const auto [l, n, origin] = [&config, automatic, this]() {
1769  if (!automatic) {
1770  return std::make_tuple<std::array<double, 3>, std::array<int, 3>,
1771  std::array<double, 3>>(
1775  } else {
1776  std::array<double, 3> l_default{20., 20., 20.};
1777  std::array<int, 3> n_default{10, 10, 10};
1778  std::array<double, 3> origin_default{-20., -20., -20.};
1779  if (modus_.is_list() && !modus_.is_box()) {
1780  logg[LExperiment].fatal(
1781  "The lattice in List modus should be manually specified.");
1782  throw std::invalid_argument("Invalid Lattice setup.");
1783  } else if (modus_.is_collider()) {
1784  // Estimates on how far particles could get in x, y, z. The
1785  // default lattice is currently not contracted for afterburner runs
1786  const double gamma = modus_.sqrt_s_NN() / (2.0 * nucleon_mass);
1787  const double max_z = 5.0 / gamma + end_time_;
1788  const double estimated_max_transverse_velocity = 0.7;
1789  const double max_xy =
1790  5.0 + estimated_max_transverse_velocity * end_time_;
1791  origin_default = {-max_xy, -max_xy, -max_z};
1792  l_default = {2 * max_xy, 2 * max_xy, 2 * max_z};
1793  // For collider modus only, impose a minimum size of 30fm since the
1794  // heuristic above for determining the lattice expects the end time to
1795  // be large compared to the nucleus size
1796  const double minimum_extension = 30.;
1797  for (auto i = std::size_t{0}; i < l_default.size(); i++) {
1798  if (l_default[i] < minimum_extension) {
1799  logg[LExperiment].debug()
1800  << "Automatic lattice extension in direction " << i
1801  << " heuristically determined as " << l_default[i]
1802  << " fm is smaller than " << minimum_extension
1803  << " fm. Imposing minimum size.";
1804  l_default[i] = minimum_extension;
1805  origin_default[i] = -0.5 * minimum_extension;
1806  }
1807  }
1808  // Go for approximately 0.8 fm cell size and contract lattice in z by
1809  // gamma factor in case of smearing where smearing length is bound to
1810  // the lattice cell length
1811  const int n_xy = numeric_cast<int>(std::ceil(l_default[0] / 0.8));
1812  const bool to_be_contracted =
1813  (parameters_.smearing_mode == SmearingMode::Discrete ||
1814  parameters_.smearing_mode == SmearingMode::Triangular);
1815  const double contraction_factor = (to_be_contracted) ? gamma : 1.0;
1816  const int nz = numeric_cast<int>(
1817  std::ceil(l_default[2] / 0.8 * contraction_factor));
1818  n_default = {n_xy, n_xy, nz};
1819  } else if (modus_.is_box()) {
1820  origin_default = {0., 0., 0.};
1821  const double bl = modus_.length();
1822  l_default = {bl, bl, bl};
1823  const int n_xyz = numeric_cast<int>(std::ceil(bl / 0.5));
1824  n_default = {n_xyz, n_xyz, n_xyz};
1825  } else if (modus_.is_sphere()) {
1826  // Maximal distance from (0, 0, 0) at which a particle
1827  // may be found at the end of the simulation
1828  const double max_d = modus_.radius() + end_time_;
1829  origin_default = {-max_d, -max_d, -max_d};
1830  l_default = {2 * max_d, 2 * max_d, 2 * max_d};
1831  // Go for approximately 0.8 fm cell size
1832  const int n_xyz = numeric_cast<int>(std::ceil(2 * max_d / 0.8));
1833  n_default = {n_xyz, n_xyz, n_xyz};
1834  }
1835  // Take lattice properties from config to assign them to all lattices
1836  return std::make_tuple<std::array<double, 3>, std::array<int, 3>,
1837  std::array<double, 3>>(
1838  config.take(InputKeys::lattice_sizes, l_default),
1839  config.take(InputKeys::lattice_cellNumber, n_default),
1840  config.take(InputKeys::lattice_origin, origin_default));
1841  }
1842  }();
1843 
1844  logg[LExperiment].info()
1845  << "Lattice is ON. Origin = (" << origin[0] << "," << origin[1] << ","
1846  << origin[2] << "), sizes = (" << l[0] << "," << l[1] << "," << l[2]
1847  << "), number of cells = (" << n[0] << "," << n[1] << "," << n[2]
1848  << "), periodic = " << std::boolalpha << periodic;
1849 
1850  if (printout_lattice_td_ || printout_full_lattice_any_td_) {
1851  dens_type_lattice_printout_ = output_parameters.td_dens_type;
1852  printout_rho_eckart_ = output_parameters.td_rho_eckart;
1853  printout_tmn_ = output_parameters.td_tmn;
1854  printout_tmn_landau_ = output_parameters.td_tmn_landau;
1855  printout_v_landau_ = output_parameters.td_v_landau;
1856  printout_j_QBS_ = output_parameters.td_jQBS;
1857  }
1858  if (printout_tmn_ || printout_tmn_landau_ || printout_v_landau_) {
1859  Tmn_ = std::make_unique<RectangularLattice<EnergyMomentumTensor>>(
1860  l, n, origin, periodic, LatticeUpdate::AtOutput);
1861  }
1862  if (printout_j_QBS_) {
1863  j_QBS_lat_ = std::make_unique<DensityLattice>(l, n, origin, periodic,
1865  }
1866  /* Create baryon and isospin density lattices regardless of config
1867  if potentials are on. This is because they allow to compute
1868  potentials faster */
1869  if (potentials_) {
1870  // Create auxiliary lattices for baryon four-current calculation
1871  old_jmu_auxiliary_ = std::make_unique<RectangularLattice<FourVector>>(
1872  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1873  new_jmu_auxiliary_ = std::make_unique<RectangularLattice<FourVector>>(
1874  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1875  four_gradient_auxiliary_ =
1876  std::make_unique<RectangularLattice<std::array<FourVector, 4>>>(
1877  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1878 
1879  if (potentials_->use_skyrme()) {
1880  jmu_B_lat_ = std::make_unique<DensityLattice>(
1881  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1882  UB_lat_ = std::make_unique<RectangularLattice<FourVector>>(
1883  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1884  FB_lat_ = std::make_unique<
1885  RectangularLattice<std::pair<ThreeVector, ThreeVector>>>(
1886  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1887  }
1888  if (potentials_->use_symmetry()) {
1889  jmu_I3_lat_ = std::make_unique<DensityLattice>(
1890  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1891  UI3_lat_ = std::make_unique<RectangularLattice<FourVector>>(
1892  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1893  FI3_lat_ = std::make_unique<
1894  RectangularLattice<std::pair<ThreeVector, ThreeVector>>>(
1895  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1896  }
1897  if (potentials_->use_coulomb()) {
1898  jmu_el_lat_ = std::make_unique<DensityLattice>(
1899  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1900  EM_lat_ = std::make_unique<
1901  RectangularLattice<std::pair<ThreeVector, ThreeVector>>>(
1902  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1903  }
1904  if (potentials_->use_vdf()) {
1905  jmu_B_lat_ = std::make_unique<DensityLattice>(
1906  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1907  UB_lat_ = std::make_unique<RectangularLattice<FourVector>>(
1908  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1909  FB_lat_ = std::make_unique<
1910  RectangularLattice<std::pair<ThreeVector, ThreeVector>>>(
1911  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1912  }
1913  if (parameters_.field_derivatives_mode == FieldDerivativesMode::Direct) {
1914  // Create auxiliary lattices for field calculation
1915  old_fields_auxiliary_ =
1916  std::make_unique<RectangularLattice<FourVector>>(
1917  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1918  new_fields_auxiliary_ =
1919  std::make_unique<RectangularLattice<FourVector>>(
1920  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1921  fields_four_gradient_auxiliary_ =
1922  std::make_unique<RectangularLattice<std::array<FourVector, 4>>>(
1923  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1924 
1925  // Create the fields lattice
1926  fields_lat_ = std::make_unique<FieldsLattice>(
1927  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1928  }
1929  }
1930  if (dens_type_lattice_printout_ == DensityType::Baryon && !jmu_B_lat_) {
1931  jmu_B_lat_ = std::make_unique<DensityLattice>(l, n, origin, periodic,
1933  }
1934  if (dens_type_lattice_printout_ == DensityType::BaryonicIsospin &&
1935  !jmu_I3_lat_) {
1936  jmu_I3_lat_ = std::make_unique<DensityLattice>(l, n, origin, periodic,
1938  }
1939  if (dens_type_lattice_printout_ != DensityType::None &&
1940  dens_type_lattice_printout_ != DensityType::BaryonicIsospin &&
1941  dens_type_lattice_printout_ != DensityType::Baryon) {
1942  jmu_custom_lat_ = std::make_unique<DensityLattice>(
1943  l, n, origin, periodic, LatticeUpdate::AtOutput);
1944  }
1945  } else if (printout_lattice_td_ || printout_full_lattice_any_td_) {
1946  logg[LExperiment].error(
1947  "If you want Therm. VTK or Lattice output, configure a lattice for "
1948  "it.");
1949  } else if (potentials_ && potentials_->use_coulomb()) {
1950  logg[LExperiment].error(
1951  "Coulomb potential requires a lattice. Please add one to the "
1952  "configuration");
1953  }
1954 
1955  // Warning for the mean field calculation if lattice is not on.
1956  if ((potentials_ != nullptr) && (jmu_B_lat_ == nullptr)) {
1957  logg[LExperiment].warn() << "Lattice is NOT used. Mean-field energy is "
1958  << "not going to be calculated.";
1959  }
1960 
1961  // Store pointers to potential and lattice accessible for Action
1962  if (parameters_.potential_affect_threshold) {
1963  UB_lat_pointer = UB_lat_.get();
1964  UI3_lat_pointer = UI3_lat_.get();
1965  pot_pointer = potentials_.get();
1966  }
1967 
1968  // Throw fatal if DerivativesMode == FiniteDifference and lattice is not on.
1969  if ((parameters_.derivatives_mode == DerivativesMode::FiniteDifference) &&
1970  (jmu_B_lat_ == nullptr)) {
1971  throw std::runtime_error(
1972  "Lattice is necessary to calculate finite difference gradients.");
1973  }
1974 
1975  // Create forced thermalizer
1977  Configuration th_conf = config.extract_complete_sub_configuration(
1979  thermalizer_ = modus_.create_grandcan_thermalizer(th_conf);
1980  }
1981 
1982  /* Take the seed setting only after the configuration was stored to a file
1983  * in smash.cc */
1984  seed_ = config.take(InputKeys::gen_randomseed);
1985 }
1986 
1988 const std::string hline(113, '-');
1989 
2014 std::string format_measurements(const std::vector<Particles> &ensembles,
2015  uint64_t scatterings_this_interval,
2016  const QuantumNumbers &conserved_initial,
2017  SystemTimePoint time_start, double time,
2018  double E_mean_field,
2019  double E_mean_field_initial);
2034  const Potentials &potentials,
2036  RectangularLattice<std::pair<ThreeVector, ThreeVector>> *em_lattice,
2037  const ExperimentParameters &parameters);
2038 
2055 EventInfo fill_event_info(const std::vector<Particles> &ensembles,
2056  double E_mean_field, double modus_impact_parameter,
2057  const ExperimentParameters &parameters,
2058  bool projectile_target_interact,
2059  bool kinematic_cut_for_SMASH_IC);
2060 
2061 template <typename Modus>
2063  random::set_seed(seed_);
2064  logg[LExperiment].info() << "random number seed: " << seed_;
2065  /* Set seed for the next event. It has to be positive, so it can be entered
2066  * in the config.
2067  *
2068  * We have to be careful about the minimal integer, whose absolute value
2069  * cannot be represented. */
2070  int64_t r = random::advance();
2071  while (r == INT64_MIN) {
2072  r = random::advance();
2073  }
2074  seed_ = std::abs(r);
2075  /* Set the random seed used in PYTHIA hadronization
2076  * to be same with the SMASH one.
2077  * In this way we ensure that the results are reproducible
2078  * for every event if one knows SMASH random seed. */
2079  if (process_string_ptr_ != NULL) {
2080  process_string_ptr_->init_pythia_hadron_rndm();
2081  }
2082 
2083  for (Particles &particles : ensembles_) {
2084  particles.reset();
2085  }
2086 
2087  // Sample particles according to the initial conditions
2088  double start_time = -1.0;
2089 
2090  // Sample impact parameter only once per all ensembles
2091  // It should be the same for all ensembles
2092  if (modus_.is_collider()) {
2093  modus_.sample_impact();
2094  logg[LExperiment].info("Impact parameter = ", modus_.impact_parameter(),
2095  " fm");
2096  }
2097  for (Particles &particles : ensembles_) {
2098  start_time = modus_.initial_conditions(&particles, parameters_);
2099  }
2100  /* For box modus make sure that particles are in the box. In principle, after
2101  * a correct initialization they should be, so this is just playing it safe.
2102  */
2103  for (Particles &particles : ensembles_) {
2104  modus_.impose_boundary_conditions(&particles, outputs_);
2105  }
2106  // Reset the simulation clock
2107  double timestep = delta_time_startup_;
2108 
2109  switch (time_step_mode_) {
2110  case TimeStepMode::Fixed:
2111  break;
2112  case TimeStepMode::None:
2113  timestep = end_time_ - start_time;
2114  // Take care of the box modus + timestepless propagation
2115  const double max_dt = modus_.max_timestep(max_transverse_distance_sqr_);
2116  if (max_dt > 0. && max_dt < timestep) {
2117  timestep = max_dt;
2118  }
2119  break;
2120  }
2121  std::unique_ptr<UniformClock> clock_for_this_event;
2122  if (modus_.is_list() && (timestep < 0.0)) {
2123  throw std::runtime_error(
2124  "Timestep for the given event is negative. \n"
2125  "This might happen if the formation times of the input particles are "
2126  "larger than the specified end time of the simulation.");
2127  }
2128  clock_for_this_event =
2129  std::make_unique<UniformClock>(start_time, timestep, end_time_);
2130  parameters_.labclock = std::move(clock_for_this_event);
2131 
2132  // Reset the output clock
2133  parameters_.outputclock->reset(start_time, true);
2134  // remove time before starting time in case of custom output times.
2135  parameters_.outputclock->remove_times_in_past(start_time);
2136 
2137  logg[LExperiment].debug(
2138  "Lab clock: t_start = ", parameters_.labclock->current_time(),
2139  ", dt = ", parameters_.labclock->timestep_duration());
2140 
2141  /* Save the initial conserved quantum numbers and total momentum in
2142  * the system for conservation checks */
2143  conserved_initial_ = QuantumNumbers(ensembles_);
2144  wall_actions_total_ = 0;
2145  previous_wall_actions_total_ = 0;
2146  interactions_total_ = 0;
2147  previous_interactions_total_ = 0;
2148  discarded_interactions_total_ = 0;
2149  total_pauli_blocked_ = 0;
2150  projectile_target_interact_.assign(parameters_.n_ensembles, false);
2151  total_hypersurface_crossing_actions_ = 0;
2152  total_energy_removed_ = 0.0;
2153  total_energy_violated_by_Pythia_ = 0.0;
2154  // Print output headers
2155  logg[LExperiment].info() << hline;
2156  logg[LExperiment].info() << "Time[fm] Ekin[GeV] E_MF[GeV] ETotal[GeV] "
2157  << "ETot/N[GeV] D(ETot/N)[GeV] Scatt&Decays "
2158  << "Particles Comp.Time";
2159  logg[LExperiment].info() << hline;
2160  double E_mean_field = 0.0;
2161  if (potentials_) {
2162  // update_potentials();
2163  // if (parameters.outputclock->current_time() == 0.0 )
2164  // using the lattice is necessary
2165  if ((jmu_B_lat_ != nullptr)) {
2166  update_lattice(jmu_B_lat_.get(), old_jmu_auxiliary_.get(),
2167  new_jmu_auxiliary_.get(), four_gradient_auxiliary_.get(),
2168  LatticeUpdate::EveryTimestep, DensityType::Baryon,
2169  density_param_, ensembles_,
2170  parameters_.labclock->timestep_duration(), true);
2171  // Because there was no lattice at t=-Delta_t, the time derivatives
2172  // drho_dt and dj^mu/dt at t=0 are huge, while they shouldn't be; we
2173  // overwrite the time derivative to zero by hand.
2174  for (auto &node : *jmu_B_lat_) {
2175  node.overwrite_drho_dt_to_zero();
2176  node.overwrite_djmu_dt_to_zero();
2177  }
2178  E_mean_field = calculate_mean_field_energy(*potentials_, *jmu_B_lat_,
2179  EM_lat_.get(), parameters_);
2180  }
2181  }
2182  initial_mean_field_energy_ = E_mean_field;
2184  ensembles_, 0u, conserved_initial_, time_start_,
2185  parameters_.labclock->current_time(), E_mean_field,
2186  initial_mean_field_energy_);
2187 
2188  // Output at event start
2189  for (const auto &output : outputs_) {
2190  for (int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
2191  auto event_info = fill_event_info(
2192  ensembles_, E_mean_field, modus_.impact_parameter(), parameters_,
2193  projectile_target_interact_[i_ens], kinematic_cuts_for_IC_output_);
2194  output->at_eventstart(ensembles_[i_ens], {event_, i_ens}, event_info);
2195  }
2196  // For thermodynamic output
2197  output->at_eventstart(ensembles_, event_);
2198  // For thermodynamic lattice output
2199  if (printout_full_lattice_any_td_) {
2200  if (printout_rho_eckart_) {
2201  switch (dens_type_lattice_printout_) {
2202  case DensityType::Baryon:
2203  output->at_eventstart(event_, ThermodynamicQuantity::EckartDensity,
2204  DensityType::Baryon, *jmu_B_lat_);
2205  break;
2207  output->at_eventstart(event_, ThermodynamicQuantity::EckartDensity,
2208  DensityType::BaryonicIsospin, *jmu_I3_lat_);
2209  break;
2210  case DensityType::None:
2211  break;
2212  default:
2213  output->at_eventstart(event_, ThermodynamicQuantity::EckartDensity,
2215  *jmu_custom_lat_);
2216  }
2217  }
2218  if (printout_tmn_) {
2219  output->at_eventstart(event_, ThermodynamicQuantity::Tmn,
2220  dens_type_lattice_printout_, *Tmn_);
2221  }
2222  if (printout_tmn_landau_) {
2223  output->at_eventstart(event_, ThermodynamicQuantity::TmnLandau,
2224  dens_type_lattice_printout_, *Tmn_);
2225  }
2226  if (printout_v_landau_) {
2227  output->at_eventstart(event_, ThermodynamicQuantity::LandauVelocity,
2228  dens_type_lattice_printout_, *Tmn_);
2229  }
2230  if (printout_j_QBS_) {
2231  output->at_eventstart(event_, ThermodynamicQuantity::j_QBS,
2232  dens_type_lattice_printout_, *j_QBS_lat_);
2233  }
2234  }
2235  }
2236 
2237  /* In the ColliderModus, if Fermi motion is frozen, assign the beam momenta
2238  * to the nucleons in both the projectile and the target. Every ensemble
2239  * gets the same beam momenta, so no need to create beam_momenta_ vector
2240  * for every ensemble.
2241  */
2242  if (modus_.is_collider() && modus_.fermi_motion() == FermiMotion::Frozen) {
2243  for (ParticleData &particle : ensembles_[0]) {
2244  const double m = particle.effective_mass();
2245  double v_beam = 0.0;
2246  if (particle.belongs_to() == BelongsTo::Projectile) {
2247  v_beam = modus_.velocity_projectile();
2248  } else if (particle.belongs_to() == BelongsTo::Target) {
2249  v_beam = modus_.velocity_target();
2250  }
2251  const double gamma = 1.0 / std::sqrt(1.0 - v_beam * v_beam);
2252  beam_momentum_.emplace_back(
2253  FourVector(gamma * m, 0.0, 0.0, gamma * v_beam * m));
2254  } // loop over particles
2255  }
2256 }
2257 
2258 template <typename Modus>
2259 bool Experiment<Modus>::perform_action(Action &action, int i_ensemble,
2260  bool include_pauli_blocking) {
2261  Particles &particles = ensembles_[i_ensemble];
2262  auto &incoming = action.incoming_particles();
2263  // Make sure to skip invalid and Pauli-blocked actions.
2264  if (!action.is_valid(particles)) {
2265  discarded_interactions_total_++;
2266  logg[LExperiment].debug(~einhard::DRed(), "✘ ", action,
2267  " (discarded: invalid)");
2268  return false;
2269  }
2270  const bool core_in_incoming =
2271  std::any_of(incoming.begin(), incoming.end(),
2272  [](const ParticleData &p) { return p.is_core(); });
2273  if (core_in_incoming) {
2274  if (action.get_type() == ProcessType::FluidizationNoRemoval) {
2275  // If the incoming particle is already core, the action should not happen.
2276  logg[LExperiment].debug() << "Discarding " << incoming[0].id();
2277  return false;
2278  } else if (action.get_type() != ProcessType::Elastic) {
2279  /* Only elastic collisions can happen between core and corona particles
2280  * (1→N can still happen) */
2281  const bool all_core_in_incoming =
2282  std::all_of(incoming.begin(), incoming.end(),
2283  [](const ParticleData &p) { return p.is_core(); });
2284  if (!all_core_in_incoming) {
2285  return false;
2286  }
2287  }
2288  }
2289  try {
2290  action.generate_final_state();
2291  } catch (Action::StochasticBelowEnergyThreshold &) {
2292  return false;
2293  }
2294  logg[LExperiment].debug("Process Type is: ", action.get_type());
2295  if (include_pauli_blocking && pauli_blocker_ &&
2296  action.is_pauli_blocked(ensembles_, *pauli_blocker_)) {
2297  total_pauli_blocked_++;
2298  return false;
2299  }
2300 
2301  // Prepare projectile_target_interact_, it's used for output
2302  // to signal that there was some interaction in this event
2303  if (modus_.is_collider()) {
2304  int count_target = 0, count_projectile = 0;
2305  for (const auto &p : incoming) {
2306  if (p.belongs_to() == BelongsTo::Projectile) {
2307  count_projectile++;
2308  } else if (p.belongs_to() == BelongsTo::Target) {
2309  count_target++;
2310  }
2311  }
2312  if (count_target > 0 && count_projectile > 0) {
2313  projectile_target_interact_[i_ensemble] = true;
2314  }
2315  }
2316 
2317  /* Make sure to pick a non-zero integer, because 0 is reserved for "no
2318  * interaction yet". */
2319  const auto id_process = static_cast<uint32_t>(interactions_total_ + 1);
2320  // we perform the action and collect possible energy violations by Pythia
2321  total_energy_violated_by_Pythia_ += action.perform(&particles, id_process);
2322 
2323  interactions_total_++;
2324  if (action.get_type() == ProcessType::Wall) {
2325  wall_actions_total_++;
2326  }
2327  if (action.get_type() == ProcessType::Fluidization) {
2328  total_hypersurface_crossing_actions_++;
2329  total_energy_removed_ += action.incoming_particles()[0].momentum().x0();
2330  }
2331  // Calculate Eckart rest frame density at the interaction point
2332  double rho = 0.0;
2333  if (dens_type_ != DensityType::None) {
2334  const FourVector r_interaction = action.get_interaction_point();
2335  constexpr bool compute_grad = false;
2336  const bool smearing = true;
2337  // todo(oliiny): it's a rough density estimate from a single ensemble.
2338  // It might actually be appropriate for output. Discuss.
2339  rho = std::get<0>(current_eckart(r_interaction.threevec(), particles,
2340  density_param_, dens_type_, compute_grad,
2341  smearing));
2342  }
2358  for (const auto &output : outputs_) {
2359  if (output->is_dilepton_output() || output->is_photon_output()) {
2360  continue;
2361  }
2362  if (output->is_IC_output()) {
2363  if (action.get_type() == ProcessType::Fluidization ||
2364  action.get_type() == ProcessType::FluidizationNoRemoval) {
2365  output->at_interaction(action, rho);
2366  }
2367  } else {
2368  output->at_interaction(action, rho);
2369  }
2370  }
2371 
2372  // At every collision photons can be produced.
2373  // Note: We rely here on the lazy evaluation of the arguments to if.
2374  // It may happen that in a wall-crossing-action sqrt_s raises an exception.
2375  // Therefore we first have to check if the incoming particles can undergo
2376  // an em-interaction.
2377  if (photons_switch_ &&
2378  ScatterActionPhoton::is_photon_reaction(action.incoming_particles()) &&
2379  ScatterActionPhoton::is_kinematically_possible(
2380  action.sqrt_s(), action.incoming_particles())) {
2381  /* Time in the action constructor is relative to
2382  * current time of incoming */
2383  constexpr double action_time = 0.;
2384  ScatterActionPhoton photon_act(
2385  action.incoming_particles(), action_time, n_fractional_photons_,
2386  action.get_total_weight(), parameters_.spin_interaction_type);
2387 
2398  photon_act.add_dummy_hadronic_process(action.get_total_weight());
2399 
2400  // Now add the actual photon reaction channel.
2401  photon_act.add_single_process();
2402 
2403  photon_act.perform_photons(outputs_);
2404  }
2405 
2406  if (bremsstrahlung_switch_ &&
2407  BremsstrahlungAction::is_bremsstrahlung_reaction(
2408  action.incoming_particles())) {
2409  /* Time in the action constructor is relative to
2410  * current time of incoming */
2411  constexpr double action_time = 0.;
2412 
2413  BremsstrahlungAction brems_act(
2414  action.incoming_particles(), action_time, n_fractional_photons_,
2415  action.get_total_weight(), parameters_.spin_interaction_type);
2416 
2428  brems_act.add_dummy_hadronic_process(action.get_total_weight());
2429 
2430  // Now add the actual bremsstrahlung reaction channel.
2431  brems_act.add_single_process();
2432 
2433  brems_act.perform_bremsstrahlung(outputs_);
2434  }
2435 
2436  logg[LExperiment].debug(~einhard::Green(), "✔ ", action);
2437  return true;
2438 }
2439 
2450 void validate_and_adjust_particle_list(ParticleList &particle_list);
2451 
2452 template <typename Modus>
2454  ParticleList &&add_plist,
2455  ParticleList &&remove_plist) {
2456  if (!add_plist.empty() || !remove_plist.empty()) {
2457  if (ensembles_.size() > 1) {
2458  throw std::runtime_error(
2459  "Adding or removing particles from SMASH is only possible when one "
2460  "ensemble is used.");
2461  }
2462  const double action_time = parameters_.labclock->current_time();
2463  /* Use two if statements. The first one is to check if the particles are
2464  * valid. Since this might remove all particles, a second if statement is
2465  * needed to avoid executing the action in that case.*/
2466  if (!add_plist.empty()) {
2468  }
2469  if (!add_plist.empty()) {
2470  // Create and perform action to add particle(s)
2471  auto action_add_particles = std::make_unique<FreeforallAction>(
2472  ParticleList{}, add_plist, action_time);
2473  perform_action(*action_add_particles, 0);
2474  }
2475  // Also here 2 if statements are needed as above.
2476  if (!remove_plist.empty()) {
2477  validate_and_adjust_particle_list(remove_plist);
2478  }
2479  if (!remove_plist.empty()) {
2480  ParticleList found_particles_to_remove;
2481  for (const auto &particle_to_remove : remove_plist) {
2482  const auto iterator_to_particle_to_be_removed_in_ensemble =
2483  std::find_if(
2484  ensembles_[0].begin(), ensembles_[0].end(),
2485  [&particle_to_remove, &action_time](const ParticleData &p) {
2487  particle_to_remove, p, action_time);
2488  });
2489  if (iterator_to_particle_to_be_removed_in_ensemble !=
2490  ensembles_[0].end())
2491  found_particles_to_remove.push_back(
2492  *iterator_to_particle_to_be_removed_in_ensemble);
2493  }
2494  // Sort the particles found to be removed according to their id and look
2495  // for duplicates (sorting is needed to call std::adjacent_find).
2496  std::sort(found_particles_to_remove.begin(),
2497  found_particles_to_remove.end(),
2498  [](const ParticleData &p1, const ParticleData &p2) {
2499  return p1.id() < p2.id();
2500  });
2501  const auto iterator_to_first_duplicate = std::adjacent_find(
2502  found_particles_to_remove.begin(), found_particles_to_remove.end(),
2503  [](const ParticleData &p1, const ParticleData &p2) {
2504  return p1.id() == p2.id();
2505  });
2506  if (iterator_to_first_duplicate != found_particles_to_remove.end()) {
2507  logg[LExperiment].error() << "The same particle has been asked to be "
2508  "removed multiple times:\n"
2509  << *iterator_to_first_duplicate;
2510  throw std::logic_error("Particle cannot be removed twice!");
2511  }
2512  if (auto delta = remove_plist.size() - found_particles_to_remove.size();
2513  delta > 0) {
2514  logg[LExperiment].warn(
2515  "When trying to remove particle(s) at the beginning ",
2516  "of the system evolution,\n", delta,
2517  " particle(s) could not be found and will be ignored.");
2518  }
2519  if (!found_particles_to_remove.empty()) {
2520  [[maybe_unused]] const auto number_particles_before_removal =
2521  ensembles_[0].size();
2522  // Create and perform action to remove particles
2523  auto action_remove_particles = std::make_unique<FreeforallAction>(
2524  found_particles_to_remove, ParticleList{}, action_time);
2525  perform_action(*action_remove_particles, 0);
2526 
2527  assert(number_particles_before_removal -
2528  found_particles_to_remove.size() ==
2529  ensembles_[0].size());
2530  }
2531  }
2532  }
2533 
2534  if (t_end > end_time_) {
2535  logg[LExperiment].fatal()
2536  << "Evolution asked to be run until " << t_end << " > " << end_time_
2537  << " and this cannot be done (because of how the clock works).";
2538  throw std::logic_error(
2539  "Experiment cannot evolve the system beyond End_Time.");
2540  }
2541  while (*(parameters_.labclock) < t_end) {
2542  const double dt = parameters_.labclock->timestep_duration();
2543  logg[LExperiment].debug("Timestepless propagation for next ", dt, " fm.");
2544 
2545  // Perform forced thermalization if required
2546  if (thermalizer_ &&
2547  thermalizer_->is_time_to_thermalize(parameters_.labclock)) {
2548  const bool ignore_cells_under_treshold = true;
2549  // Thermodynamics in thermalizer is computed from all ensembles,
2550  // but thermalization actions act on each ensemble independently
2551  thermalizer_->update_thermalizer_lattice(ensembles_, density_param_,
2552  ignore_cells_under_treshold);
2553  const double current_t = parameters_.labclock->current_time();
2554  for (int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
2555  thermalizer_->thermalize(ensembles_[i_ens], current_t,
2556  parameters_.testparticles);
2557  ThermalizationAction th_act(*thermalizer_, current_t);
2558  if (th_act.any_particles_thermalized()) {
2559  perform_action(th_act, i_ens);
2560  }
2561  }
2562  }
2563 
2564  if (IC_dynamic_) {
2565  modus_.build_fluidization_lattice(parameters_.labclock->current_time(),
2566  ensembles_, density_param_);
2567  }
2568 
2569  std::vector<Actions> actions(parameters_.n_ensembles);
2570  for (int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
2571  actions[i_ens].clear();
2572  if (ensembles_[i_ens].size() > 0 && action_finders_.size() > 0) {
2573  /* (1.a) Create grid. */
2574  const double min_cell_length = compute_min_cell_length(dt);
2575  logg[LExperiment].debug("Creating grid with minimal cell length ",
2576  min_cell_length);
2577  /* For the hyper-surface-crossing actions also unformed particles are
2578  * searched and therefore needed on the grid. */
2579  const bool include_unformed_particles = IC_switch_;
2580  const auto &grid =
2581  use_grid_ ? modus_.create_grid(ensembles_[i_ens], min_cell_length,
2582  dt, parameters_.coll_crit,
2583  include_unformed_particles)
2584  : modus_.create_grid(ensembles_[i_ens], min_cell_length,
2585  dt, parameters_.coll_crit,
2586  include_unformed_particles,
2587  CellSizeStrategy::Largest);
2588 
2589  const double gcell_vol = grid.cell_volume();
2590  /* (1.b) Iterate over cells and find actions. */
2591  grid.iterate_cells(
2592  [&](const ParticleList &search_list) {
2593  for (const auto &finder : action_finders_) {
2594  actions[i_ens].insert(finder->find_actions_in_cell(
2595  search_list, dt, gcell_vol, beam_momentum_));
2596  }
2597  },
2598  [&](const ParticleList &search_list,
2599  const ParticleList &neighbors_list) {
2600  for (const auto &finder : action_finders_) {
2601  actions[i_ens].insert(finder->find_actions_with_neighbors(
2602  search_list, neighbors_list, dt, beam_momentum_));
2603  }
2604  });
2605  }
2606  }
2607 
2608  /* \todo (optimizations) Adapt timestep size here */
2609 
2610  /* (2) Propagate from action to action until next output or timestep end */
2611  const double end_timestep_time = parameters_.labclock->next_time();
2612  while (next_output_time() < end_timestep_time) {
2613  for (int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
2614  run_time_evolution_timestepless(actions[i_ens], i_ens,
2615  next_output_time());
2616  }
2617  ++(*parameters_.outputclock);
2618 
2619  intermediate_output();
2620  }
2621  for (int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
2622  run_time_evolution_timestepless(actions[i_ens], i_ens, end_timestep_time);
2623  }
2624 
2625  /* (3) Update potentials (if computed on the lattice) and
2626  * compute new momenta according to equations of motion */
2627  if (potentials_) {
2628  update_potentials();
2629  update_momenta(ensembles_, parameters_.labclock->timestep_duration(),
2630  *potentials_, FB_lat_.get(), FI3_lat_.get(), EM_lat_.get(),
2631  jmu_B_lat_.get());
2632  }
2633 
2634  /* (4) Expand universe if non-minkowskian metric; updates
2635  * positions and momenta according to the selected expansion */
2636  if (metric_.mode_ != ExpansionMode::NoExpansion) {
2637  for (Particles &particles : ensembles_) {
2638  expand_space_time(&particles, parameters_, metric_);
2639  }
2640  }
2641 
2642  ++(*parameters_.labclock);
2643 
2644  /* (5) Check conservation laws.
2645  *
2646  * Check conservation of conserved quantities if potentials and string
2647  * fragmentation are off. If potentials are on then momentum is conserved
2648  * only in average. If string fragmentation is on, then energy and
2649  * momentum are only very roughly conserved in high-energy collisions. */
2650  if (!potentials_ && !parameters_.strings_switch &&
2651  metric_.mode_ == ExpansionMode::NoExpansion && !IC_switch_) {
2652  std::string err_msg = conserved_initial_.report_deviations(ensembles_);
2653  if (!err_msg.empty()) {
2654  logg[LExperiment].error() << err_msg;
2655  throw std::runtime_error("Violation of conserved quantities!");
2656  }
2657  }
2658  }
2659 
2660  /* Increment once more the output clock in order to have it prepared for the
2661  * final_output() call. Once close to the end time, the while-loop above to
2662  * produce intermediate output is not entered as the next_output_time() is
2663  * never strictly smaller than end_timestep_time (they are usually equal).
2664  * Since in the final_output() function the current time of the output clock
2665  * is used to produce the output, this has to be incremented before producing
2666  * the final output and it makes sense to do it here.
2667  */
2668  ++(*parameters_.outputclock);
2669 
2670  if (pauli_blocker_) {
2671  logg[LExperiment].info(
2672  "Interactions: Pauli-blocked/performed = ", total_pauli_blocked_, "/",
2673  interactions_total_ - wall_actions_total_);
2674  }
2675 }
2676 
2677 template <typename Modus>
2679  Particles &particles) {
2680  const double dt =
2681  propagate_straight_line(&particles, to_time, beam_momentum_);
2682  if (dilepton_finder_ != nullptr) {
2683  for (const auto &output : outputs_) {
2684  dilepton_finder_->shine(particles, output.get(), dt);
2685  }
2686  }
2687 }
2688 
2696 inline void check_interactions_total(uint64_t interactions_total) {
2697  constexpr uint64_t max_uint32 = std::numeric_limits<uint32_t>::max();
2698  if (interactions_total >= max_uint32) {
2699  throw std::runtime_error("Integer overflow in total interaction number!");
2700  }
2701 }
2702 
2703 template <typename Modus>
2705  Actions &actions, int i_ensemble, const double end_time_propagation) {
2706  Particles &particles = ensembles_[i_ensemble];
2707  logg[LExperiment].debug(
2708  "Timestepless propagation: ", "Actions size = ", actions.size(),
2709  ", end time = ", end_time_propagation);
2710 
2711  // iterate over all actions
2712  while (!actions.is_empty()) {
2713  if (actions.earliest_time() > end_time_propagation) {
2714  break;
2715  }
2716  // get next action
2717  ActionPtr act = actions.pop();
2718  if (!act->is_valid(particles)) {
2719  discarded_interactions_total_++;
2720  logg[LExperiment].debug(~einhard::DRed(), "✘ ", act,
2721  " (discarded: invalid)");
2722  continue;
2723  }
2724  logg[LExperiment].debug(~einhard::Green(), "✔ ", act,
2725  ", action time = ", act->time_of_execution());
2726 
2727  /* (1) Propagate to the next action. */
2728  propagate_and_shine(act->time_of_execution(), particles);
2729 
2730  /* (2) Perform action.
2731  *
2732  * Update the positions of the incoming particles, because the information
2733  * in the action object will be outdated as the particles have been
2734  * propagated since the construction of the action. */
2735  act->update_incoming(particles);
2736  const bool performed = perform_action(*act, i_ensemble);
2737 
2738  /* No need to update actions for outgoing particles
2739  * if the action is not performed. */
2740  if (!performed) {
2741  continue;
2742  }
2743 
2744  /* (3) Update actions for newly-produced particles. */
2745 
2746  const double end_time_timestep = parameters_.labclock->next_time();
2747  // New actions are always search until the end of the current timestep
2748  const double time_left = end_time_timestep - act->time_of_execution();
2749  const ParticleList &outgoing_particles = act->outgoing_particles();
2750  // Grid cell volume set to zero, since there is no grid
2751  const double gcell_vol = 0.0;
2752  for (const auto &finder : action_finders_) {
2753  // Outgoing particles can still decay, cross walls...
2754  actions.insert(finder->find_actions_in_cell(outgoing_particles, time_left,
2755  gcell_vol, beam_momentum_));
2756  // ... and collide with other particles.
2757  actions.insert(finder->find_actions_with_surrounding_particles(
2758  outgoing_particles, particles, time_left, beam_momentum_));
2759  }
2760 
2761  check_interactions_total(interactions_total_);
2762  }
2763 
2764  propagate_and_shine(end_time_propagation, particles);
2765 }
2766 
2767 template <typename Modus>
2769  const uint64_t wall_actions_this_interval =
2770  wall_actions_total_ - previous_wall_actions_total_;
2771  previous_wall_actions_total_ = wall_actions_total_;
2772  const uint64_t interactions_this_interval = interactions_total_ -
2773  previous_interactions_total_ -
2774  wall_actions_this_interval;
2775  previous_interactions_total_ = interactions_total_;
2776  double E_mean_field = 0.0;
2779  double computational_frame_time = 0.0;
2780  if (potentials_) {
2781  // using the lattice is necessary
2782  if ((jmu_B_lat_ != nullptr)) {
2783  E_mean_field = calculate_mean_field_energy(*potentials_, *jmu_B_lat_,
2784  EM_lat_.get(), parameters_);
2785  /*
2786  * Mean field calculated in a box should remain approximately constant if
2787  * the system is in equilibrium, and so deviations from its original value
2788  * may signal a phase transition or other dynamical process. This
2789  * comparison only makes sense in the Box Modus, hence the condition.
2790  */
2791  if (modus_.is_box()) {
2792  double tmp = (E_mean_field - initial_mean_field_energy_) /
2793  (E_mean_field + initial_mean_field_energy_);
2794  /*
2795  * This is displayed when the system evolves away from its initial
2796  * configuration (which is when the total mean field energy in the box
2797  * deviates from its initial value).
2798  */
2799  if (std::abs(tmp) > 0.01) {
2800  logg[LExperiment].info()
2801  << "\n\n\n\t The mean field at t = "
2802  << parameters_.outputclock->current_time()
2803  << " [fm] differs from the mean field at t = 0:"
2804  << "\n\t\t initial_mean_field_energy_ = "
2805  << initial_mean_field_energy_ << " [GeV]"
2806  << "\n\t\t abs[(E_MF - E_MF(t=0))/(E_MF + E_MF(t=0))] = "
2807  << std::abs(tmp)
2808  << "\n\t\t E_MF/E_MF(t=0) = "
2809  << E_mean_field / initial_mean_field_energy_ << "\n\n";
2810  }
2811  }
2812  }
2813  }
2814 
2816  ensembles_, interactions_this_interval, conserved_initial_, time_start_,
2817  parameters_.outputclock->current_time(), E_mean_field,
2818  initial_mean_field_energy_);
2819  const LatticeUpdate lat_upd = LatticeUpdate::AtOutput;
2820 
2821  // save evolution data
2822  if (!(modus_.is_box() && parameters_.outputclock->current_time() <
2823  modus_.equilibration_time())) {
2824  for (const auto &output : outputs_) {
2825  if (output->is_dilepton_output() || output->is_photon_output() ||
2826  output->is_IC_output()) {
2827  continue;
2828  }
2829  for (int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
2830  auto event_info = fill_event_info(
2831  ensembles_, E_mean_field, modus_.impact_parameter(), parameters_,
2832  projectile_target_interact_[i_ens], kinematic_cuts_for_IC_output_);
2833 
2834  output->at_intermediate_time(ensembles_[i_ens], parameters_.outputclock,
2835  density_param_, {event_, i_ens},
2836  event_info);
2837  computational_frame_time = event_info.current_time;
2838  }
2839  // For thermodynamic output
2840  output->at_intermediate_time(ensembles_, parameters_.outputclock,
2841  density_param_);
2842 
2843  // Thermodynamic output on the lattice versus time
2844  if (printout_rho_eckart_) {
2845  switch (dens_type_lattice_printout_) {
2846  case DensityType::Baryon:
2848  jmu_B_lat_.get(), lat_upd, DensityType::Baryon, density_param_,
2849  ensembles_, false);
2850  output->thermodynamics_output(ThermodynamicQuantity::EckartDensity,
2851  DensityType::Baryon, *jmu_B_lat_);
2852  output->thermodynamics_lattice_output(*jmu_B_lat_,
2853  computational_frame_time);
2854  break;
2857  jmu_I3_lat_.get(), lat_upd, DensityType::BaryonicIsospin,
2858  density_param_, ensembles_, false);
2859  output->thermodynamics_output(ThermodynamicQuantity::EckartDensity,
2861  *jmu_I3_lat_);
2862  output->thermodynamics_lattice_output(*jmu_I3_lat_,
2863  computational_frame_time);
2864  break;
2865  case DensityType::None:
2866  break;
2867  default:
2869  jmu_custom_lat_.get(), lat_upd, dens_type_lattice_printout_,
2870  density_param_, ensembles_, false);
2871  output->thermodynamics_output(ThermodynamicQuantity::EckartDensity,
2872  dens_type_lattice_printout_,
2873  *jmu_custom_lat_);
2874  output->thermodynamics_lattice_output(*jmu_custom_lat_,
2875  computational_frame_time);
2876  }
2877  }
2878  if (printout_tmn_ || printout_tmn_landau_ || printout_v_landau_) {
2880  Tmn_.get(), lat_upd, dens_type_lattice_printout_, density_param_,
2881  ensembles_, false);
2882  if (printout_tmn_) {
2883  output->thermodynamics_output(ThermodynamicQuantity::Tmn,
2884  dens_type_lattice_printout_, *Tmn_);
2885  output->thermodynamics_lattice_output(
2886  ThermodynamicQuantity::Tmn, *Tmn_, computational_frame_time);
2887  }
2888  if (printout_tmn_landau_) {
2889  output->thermodynamics_output(ThermodynamicQuantity::TmnLandau,
2890  dens_type_lattice_printout_, *Tmn_);
2891  output->thermodynamics_lattice_output(
2893  computational_frame_time);
2894  }
2895  if (printout_v_landau_) {
2896  output->thermodynamics_output(ThermodynamicQuantity::LandauVelocity,
2897  dens_type_lattice_printout_, *Tmn_);
2898  output->thermodynamics_lattice_output(
2900  computational_frame_time);
2901  }
2902  }
2903  if (EM_lat_) {
2904  output->fields_output("Efield", "Bfield", *EM_lat_);
2905  }
2906  if (printout_j_QBS_) {
2907  output->thermodynamics_lattice_output(
2908  *j_QBS_lat_, computational_frame_time, ensembles_, density_param_);
2909  }
2910 
2911  if (thermalizer_) {
2912  output->thermodynamics_output(*thermalizer_);
2913  }
2914  }
2915  }
2916 }
2917 
2918 template <typename Modus>
2920  if (potentials_) {
2921  if (potentials_->use_symmetry() && jmu_I3_lat_ != nullptr) {
2922  update_lattice(jmu_I3_lat_.get(), old_jmu_auxiliary_.get(),
2923  new_jmu_auxiliary_.get(), four_gradient_auxiliary_.get(),
2924  LatticeUpdate::EveryTimestep, DensityType::BaryonicIsospin,
2925  density_param_, ensembles_,
2926  parameters_.labclock->timestep_duration(), true);
2927  }
2928  if ((potentials_->use_skyrme() || potentials_->use_symmetry()) &&
2929  jmu_B_lat_ != nullptr) {
2930  update_lattice(jmu_B_lat_.get(), old_jmu_auxiliary_.get(),
2931  new_jmu_auxiliary_.get(), four_gradient_auxiliary_.get(),
2932  LatticeUpdate::EveryTimestep, DensityType::Baryon,
2933  density_param_, ensembles_,
2934  parameters_.labclock->timestep_duration(), true);
2935  const size_t UBlattice_size = UB_lat_->size();
2936  for (size_t i = 0; i < UBlattice_size; i++) {
2937  auto jB = (*jmu_B_lat_)[i];
2938  const FourVector flow_four_velocity_B =
2939  std::abs(jB.rho()) > very_small_double ? jB.jmu_net() / jB.rho()
2940  : FourVector();
2941  double baryon_density = jB.rho();
2942  ThreeVector baryon_grad_j0 = jB.grad_j0();
2943  ThreeVector baryon_dvecj_dt = jB.dvecj_dt();
2944  ThreeVector baryon_curl_vecj = jB.curl_vecj();
2945  if (potentials_->use_skyrme()) {
2946  (*UB_lat_)[i] =
2947  flow_four_velocity_B * potentials_->skyrme_pot(baryon_density);
2948  (*FB_lat_)[i] =
2949  potentials_->skyrme_force(baryon_density, baryon_grad_j0,
2950  baryon_dvecj_dt, baryon_curl_vecj);
2951  }
2952  if (potentials_->use_symmetry() && jmu_I3_lat_ != nullptr) {
2953  auto jI3 = (*jmu_I3_lat_)[i];
2954  const FourVector flow_four_velocity_I3 =
2955  std::abs(jI3.rho()) > very_small_double
2956  ? jI3.jmu_net() / jI3.rho()
2957  : FourVector();
2958  (*UI3_lat_)[i] = flow_four_velocity_I3 *
2959  potentials_->symmetry_pot(jI3.rho(), baryon_density);
2960  (*FI3_lat_)[i] = potentials_->symmetry_force(
2961  jI3.rho(), jI3.grad_j0(), jI3.dvecj_dt(), jI3.curl_vecj(),
2962  baryon_density, baryon_grad_j0, baryon_dvecj_dt,
2963  baryon_curl_vecj);
2964  }
2965  }
2966  }
2967  if (potentials_->use_coulomb()) {
2969  jmu_el_lat_.get(), LatticeUpdate::EveryTimestep, DensityType::Charge,
2970  density_param_, ensembles_, true);
2971  for (size_t i = 0; i < EM_lat_->size(); i++) {
2972  ThreeVector electric_field = {0., 0., 0.};
2973  ThreeVector position = jmu_el_lat_->cell_center(i);
2974  jmu_el_lat_->integrate_volume(electric_field,
2975  Potentials::E_field_integrand,
2976  potentials_->coulomb_r_cut(), position);
2977  ThreeVector magnetic_field = {0., 0., 0.};
2978  jmu_el_lat_->integrate_volume(magnetic_field,
2979  Potentials::B_field_integrand,
2980  potentials_->coulomb_r_cut(), position);
2981  (*EM_lat_)[i] = std::make_pair(electric_field, magnetic_field);
2982  }
2983  } // if ((potentials_->use_skyrme() || ...
2984  if (potentials_->use_vdf() && jmu_B_lat_ != nullptr) {
2985  update_lattice(jmu_B_lat_.get(), old_jmu_auxiliary_.get(),
2986  new_jmu_auxiliary_.get(), four_gradient_auxiliary_.get(),
2987  LatticeUpdate::EveryTimestep, DensityType::Baryon,
2988  density_param_, ensembles_,
2989  parameters_.labclock->timestep_duration(), true);
2990  if (parameters_.field_derivatives_mode == FieldDerivativesMode::Direct) {
2992  fields_lat_.get(), old_fields_auxiliary_.get(),
2993  new_fields_auxiliary_.get(), fields_four_gradient_auxiliary_.get(),
2994  jmu_B_lat_.get(), LatticeUpdate::EveryTimestep, *potentials_,
2995  parameters_.labclock->timestep_duration());
2996  }
2997  const size_t UBlattice_size = UB_lat_->size();
2998  for (size_t i = 0; i < UBlattice_size; i++) {
2999  auto jB = (*jmu_B_lat_)[i];
3000  (*UB_lat_)[i] = potentials_->vdf_pot(jB.rho(), jB.jmu_net());
3001  switch (parameters_.field_derivatives_mode) {
3003  (*FB_lat_)[i] = potentials_->vdf_force(
3004  jB.rho(), jB.drho_dxnu().x0(), jB.drho_dxnu().threevec(),
3005  jB.grad_rho_cross_vecj(), jB.jmu_net().x0(), jB.grad_j0(),
3006  jB.jmu_net().threevec(), jB.dvecj_dt(), jB.curl_vecj());
3007  break;
3009  auto Amu = (*fields_lat_)[i];
3010  (*FB_lat_)[i] = potentials_->vdf_force(
3011  Amu.grad_A0(), Amu.dvecA_dt(), Amu.curl_vecA());
3012  break;
3013  }
3014  } // for (size_t i = 0; i < UBlattice_size; i++)
3015  } // if potentials_->use_vdf()
3016  }
3017 }
3018 
3019 template <typename Modus>
3021  /* At end of time evolution: Force all resonances to decay. In order to handle
3022  * decay chains, we need to loop until no further actions occur. */
3023  bool actions_performed, actions_found;
3024  uint64_t interactions_old;
3025  do {
3026  actions_found = false;
3027  interactions_old = interactions_total_;
3028  for (int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
3029  Actions actions;
3030  // Dileptons: shining of remaining resonances
3031  if (dilepton_finder_ != nullptr) {
3032  for (const auto &output : outputs_) {
3033  dilepton_finder_->shine_final(ensembles_[i_ens], output.get(), true);
3034  }
3035  }
3036  // Find actions.
3037  for (const auto &finder : action_finders_) {
3038  auto found_actions = finder->find_final_actions(ensembles_[i_ens]);
3039  if (!found_actions.empty()) {
3040  actions.insert(std::move(found_actions));
3041  actions_found = true;
3042  }
3043  }
3044  // Perform actions.
3045  while (!actions.is_empty()) {
3046  perform_action(*actions.pop(), i_ens, false);
3047  }
3048  }
3049  actions_performed = interactions_total_ > interactions_old;
3050  // Throw an error if actions were found but not performed
3051  if (actions_found && !actions_performed) {
3052  throw std::runtime_error("Final actions were found but not performed.");
3053  }
3054  // loop until no more decays occur
3055  } while (actions_performed);
3056 
3057  // Dileptons: shining of stable particles at the end
3058  if (dilepton_finder_ != nullptr) {
3059  for (const auto &output : outputs_) {
3060  for (Particles &particles : ensembles_) {
3061  dilepton_finder_->shine_final(particles, output.get(), false);
3062  }
3063  }
3064  }
3065 }
3066 
3067 template <typename Modus>
3069  /* make sure the experiment actually ran (note: we should compare this
3070  * to the start time, but we don't know that. Therefore, we check that
3071  * the time is positive, which should heuristically be the same). */
3072  double E_mean_field = 0.0;
3073  if (likely(parameters_.labclock > 0)) {
3074  const uint64_t wall_actions_this_interval =
3075  wall_actions_total_ - previous_wall_actions_total_;
3076  const uint64_t interactions_this_interval = interactions_total_ -
3077  previous_interactions_total_ -
3078  wall_actions_this_interval;
3079  if (potentials_) {
3080  // using the lattice is necessary
3081  if ((jmu_B_lat_ != nullptr)) {
3082  E_mean_field = calculate_mean_field_energy(*potentials_, *jmu_B_lat_,
3083  EM_lat_.get(), parameters_);
3084  }
3085  }
3086  if (std::abs(parameters_.labclock->current_time() - end_time_) >
3087  really_small) {
3088  logg[LExperiment].warn()
3089  << "SMASH not propagated until configured end time. Current time = "
3090  << parameters_.labclock->current_time()
3091  << "fm. End time = " << end_time_ << "fm.";
3092  } else {
3094  ensembles_, interactions_this_interval, conserved_initial_,
3095  time_start_, end_time_, E_mean_field, initial_mean_field_energy_);
3096  }
3097  int total_particles = 0;
3098  for (const Particles &particles : ensembles_) {
3099  total_particles += particles.size();
3100  }
3101  if (IC_switch_ && (total_particles == 0)) {
3102  const double initial_system_energy_plus_Pythia_violations =
3103  conserved_initial_.momentum().x0() + total_energy_violated_by_Pythia_;
3104  const double fraction_of_total_system_energy_removed =
3105  initial_system_energy_plus_Pythia_violations / total_energy_removed_;
3106  // Verify there is no more energy in the system if all particles were
3107  // removed when crossing the hypersurface
3108  if (std::fabs(fraction_of_total_system_energy_removed - 1.) >
3109  really_small) {
3110  throw std::runtime_error(
3111  "There is remaining energy in the system although all particles "
3112  "were removed.\n"
3113  "E_remain = " +
3114  std::to_string((initial_system_energy_plus_Pythia_violations -
3115  total_energy_removed_)) +
3116  " [GeV]");
3117  } else {
3118  logg[LExperiment].info() << hline;
3119  logg[LExperiment].info()
3120  << "Time real: " << SystemClock::now() - time_start_;
3121  logg[LExperiment].info()
3122  << "Interactions before reaching hypersurface: "
3123  << interactions_total_ - wall_actions_total_ -
3124  total_hypersurface_crossing_actions_;
3125  logg[LExperiment].info()
3126  << "Total number of particles removed on hypersurface: "
3127  << total_hypersurface_crossing_actions_;
3128  }
3129  } else {
3130  const double precent_discarded =
3131  interactions_total_ > 0
3132  ? static_cast<double>(discarded_interactions_total_) * 100.0 /
3133  interactions_total_
3134  : 0.0;
3135  std::stringstream msg_discarded;
3136  msg_discarded
3137  << "Discarded interaction number: " << discarded_interactions_total_
3138  << " (" << precent_discarded
3139  << "% of the total interaction number including wall crossings)";
3140 
3141  logg[LExperiment].info() << hline;
3142  logg[LExperiment].info()
3143  << "Time real: " << SystemClock::now() - time_start_;
3144  logg[LExperiment].debug() << msg_discarded.str();
3145 
3146  if (parameters_.coll_crit == CollisionCriterion::Stochastic &&
3147  precent_discarded > 1.0) {
3148  // The chosen threshold of 1% is a heuristical value
3149  logg[LExperiment].warn()
3150  << msg_discarded.str()
3151  << "\nThe number of discarded interactions is large, which means "
3152  "the assumption for the stochastic criterion of\n1 interaction "
3153  "per particle per timestep is probably violated. Consider "
3154  "reducing the timestep size.";
3155  }
3156 
3157  logg[LExperiment].info() << "Final interaction number: "
3158  << interactions_total_ - wall_actions_total_;
3159  }
3160 
3161  // Check if there are unformed particles
3162  int unformed_particles_count = 0;
3163  for (const Particles &particles : ensembles_) {
3164  for (const ParticleData &particle : particles) {
3165  if (particle.formation_time() > end_time_) {
3166  unformed_particles_count++;
3167  }
3168  }
3169  }
3170  if (unformed_particles_count > 0) {
3171  logg[LExperiment].warn(
3172  "End time might be too small. ", unformed_particles_count,
3173  " unformed particles were found at the end of the evolution.");
3174  }
3175  }
3176 
3177  // Keep track of how many ensembles had interactions
3178  count_nonempty_ensembles();
3179 
3180  for (const auto &output : outputs_) {
3181  for (int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
3182  auto event_info = fill_event_info(
3183  ensembles_, E_mean_field, modus_.impact_parameter(), parameters_,
3184  projectile_target_interact_[i_ens], kinematic_cuts_for_IC_output_);
3185  output->at_eventend(ensembles_[i_ens], {event_, i_ens}, event_info);
3186  }
3187  // For thermodynamic output
3188  output->at_eventend(ensembles_, event_);
3189 
3190  // For thermodynamic lattice output
3191  if (printout_rho_eckart_) {
3192  if (dens_type_lattice_printout_ != DensityType::None) {
3193  output->at_eventend(ThermodynamicQuantity::EckartDensity);
3194  }
3195  }
3196  if (printout_tmn_) {
3197  output->at_eventend(ThermodynamicQuantity::Tmn);
3198  }
3199  if (printout_tmn_landau_) {
3200  output->at_eventend(ThermodynamicQuantity::TmnLandau);
3201  }
3202  if (printout_v_landau_) {
3203  output->at_eventend(ThermodynamicQuantity::LandauVelocity);
3204  }
3205  if (printout_j_QBS_) {
3206  output->at_eventend(ThermodynamicQuantity::j_QBS);
3207  }
3208  }
3209 }
3210 
3211 template <typename Modus>
3213  for (bool has_interaction : projectile_target_interact_) {
3214  if (has_interaction) {
3215  nonempty_ensembles_++;
3216  }
3217  }
3218 }
3219 
3220 template <typename Modus>
3222  if (event_counting_ == EventCounting::FixedNumber) {
3223  return event_ >= nevents_;
3224  }
3225  if (event_counting_ == EventCounting::MinimumNonEmpty) {
3226  if (nonempty_ensembles_ >= minimum_nonempty_ensembles_) {
3227  return true;
3228  }
3229  if (event_ >= max_events_) {
3230  logg[LExperiment].warn()
3231  << "Maximum number of events (" << max_events_
3232  << ") exceeded. Stopping calculation. "
3233  << "The fraction of empty ensembles is "
3234  << (1.0 - static_cast<double>(nonempty_ensembles_) /
3235  (event_ * parameters_.n_ensembles))
3236  << ". If this fraction is expected, try increasing the "
3237  "Maximum_Ensembles_Run.";
3238  return true;
3239  }
3240  return false;
3241  }
3242  throw std::runtime_error("Event counting option is invalid");
3243  return false;
3244 }
3245 
3246 template <typename Modus>
3248  event_++;
3249 }
3250 
3251 template <typename Modus>
3253  const auto &mainlog = logg[LMain];
3254  for (event_ = 0; !is_finished(); event_++) {
3255  mainlog.info() << "Event " << event_;
3256 
3257  // Sample initial particles, start clock, some printout and book-keeping
3258  initialize_new_event();
3259 
3260  run_time_evolution(end_time_);
3261 
3262  do_final_interactions();
3263 
3264  // Output at event end
3265  final_output();
3266  }
3267 }
3268 
3269 } // namespace smash
3270 
3271 #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
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:103
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:192
const bool bremsstrahlung_switch_
This indicates whether bremsstrahlung is switched on.
Definition: experiment.h:614
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:2678
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:605
std::unique_ptr< ActionFinderInterface > photon_finder_
The (Scatter) Actions Finder for Direct Photons.
Definition: experiment.h:430
const bool force_decays_
This indicates whether we force all resonances to decay in the last timestep.
Definition: experiment.h:599
double initial_mean_field_energy_
The initial total mean field energy in the system.
Definition: experiment.h:647
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:730
std::vector< std::unique_ptr< ActionFinderInterface > > action_finders_
The Action finder objects.
Definition: experiment.h:424
bool printout_tmn_
Whether to print the energy-momentum tensor.
Definition: experiment.h:512
QuantumNumbers conserved_initial_
The conserved quantities of the system.
Definition: experiment.h:641
DensityParameters density_param_
Structure to precalculate and hold parameters for density computations.
Definition: experiment.h:375
double next_output_time() const
Shortcut for next output time.
Definition: experiment.h:350
double total_energy_removed_
Total energy removed from the system in hypersurface crossing actions.
Definition: experiment.h:700
DensityType dens_type_lattice_printout_
Type of density for lattice printout.
Definition: experiment.h:460
double max_transverse_distance_sqr_
Maximal distance at which particles can interact in case of the geometric criterion,...
Definition: experiment.h:632
const bool IC_dynamic_
This indicates if the IC is dynamic.
Definition: experiment.h:623
bool printout_j_QBS_
Whether to print the Q, B, S 4-currents.
Definition: experiment.h:521
std::unique_ptr< GrandCanThermalizer > thermalizer_
Instance of class used for forced thermalization.
Definition: experiment.h:531
void count_nonempty_ensembles()
Counts the number of ensembles in wich interactions took place at the end of an event.
Definition: experiment.h:3212
const TimeStepMode time_step_mode_
This indicates whether to use time steps.
Definition: experiment.h:626
std::unique_ptr< DensityLattice > jmu_custom_lat_
Custom density on the lattices.
Definition: experiment.h:457
bool printout_full_lattice_any_td_
Whether to print the thermodynamics quantities evaluated on the lattices, point by point,...
Definition: experiment.h:528
void increase_event_number()
Increases the event number by one.
Definition: experiment.h:3247
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:2704
Particles * first_ensemble()
Provides external access to SMASH particles.
Definition: experiment.h:261
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:2453
int minimum_nonempty_ensembles_
The number of ensembles, in which interactions take place, to be calculated.
Definition: experiment.h:563
std::unique_ptr< DensityLattice > jmu_B_lat_
Baryon density on the lattice.
Definition: experiment.h:439
DensityType dens_type_
Type of density to be written to collision headers.
Definition: experiment.h:653
void intermediate_output()
Intermediate output during an event.
Definition: experiment.h:2768
std::vector< FourVector > beam_momentum_
The initial nucleons in the ColliderModus propagate with beam_momentum_, if Fermi motion is frozen.
Definition: experiment.h:421
std::unique_ptr< RectangularLattice< std::pair< ThreeVector, ThreeVector > > > FI3_lat_
Lattices for the electric and magnetic component of the symmetry force.
Definition: experiment.h:483
std::unique_ptr< RectangularLattice< FourVector > > new_jmu_auxiliary_
Auxiliary lattice for values of jmu at a time step t0 + dt.
Definition: experiment.h:495
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:342
void initialize_new_event()
This is called in the beginning of each event.
Definition: experiment.h:2062
bool is_finished()
Checks wether the desired number events have been calculated.
Definition: experiment.h:3221
int n_fractional_photons_
Number of fractional photons produced per single reaction.
Definition: experiment.h:433
const double delta_time_startup_
The clock's timestep size at start up.
Definition: experiment.h:593
std::unique_ptr< RectangularLattice< EnergyMomentumTensor > > Tmn_
Lattices of energy-momentum tensors for printout.
Definition: experiment.h:490
std::vector< Particles > ensembles_
Complete particle list, all ensembles in one vector.
Definition: experiment.h:384
std::unique_ptr< RectangularLattice< FourVector > > new_fields_auxiliary_
Auxiliary lattice for values of Amu at a time step t0 + dt.
Definition: experiment.h:503
SystemTimePoint time_start_
system starting time of the simulation
Definition: experiment.h:650
uint64_t previous_wall_actions_total_
Total number of wall-crossings for previous timestep.
Definition: experiment.h:677
const bool IC_switch_
This indicates whether the experiment will be used as initial condition for hydrodynamics.
Definition: experiment.h:620
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:479
OutputsList outputs_
A list of output formaters.
Definition: experiment.h:402
bool printout_rho_eckart_
Whether to print the Eckart rest frame density.
Definition: experiment.h:509
int event_
Current event.
Definition: experiment.h:574
std::unique_ptr< RectangularLattice< std::array< FourVector, 4 > > > fields_four_gradient_auxiliary_
Auxiliary lattice for calculating the four-gradient of Amu.
Definition: experiment.h:506
std::unique_ptr< PauliBlocker > pauli_blocker_
An instance of PauliBlocker class that stores parameters needed for Pauli blocking calculations and c...
Definition: experiment.h:396
bool printout_v_landau_
Whether to print the 4-velocity in Landau frame.
Definition: experiment.h:518
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:472
bool printout_lattice_td_
Whether to print the thermodynamics quantities evaluated on the lattices.
Definition: experiment.h:524
void do_final_interactions()
Performs the final decays of an event.
Definition: experiment.h:3020
std::unique_ptr< RectangularLattice< std::pair< ThreeVector, ThreeVector > > > EM_lat_
Lattices for electric and magnetic field in fm^-2.
Definition: experiment.h:487
std::unique_ptr< FieldsLattice > fields_lat_
Mean-field A^mu on the lattice.
Definition: experiment.h:448
int max_events_
Maximum number of events to be calculated in order obtain the desired number of non-empty events usin...
Definition: experiment.h:583
int nonempty_ensembles_
Number of ensembles containing an interaction.
Definition: experiment.h:577
OutputPtr photon_output_
The Photon output.
Definition: experiment.h:408
EventCounting event_counting_
The way in which the number of calculated events is specified.
Definition: experiment.h:571
const bool dileptons_switch_
This indicates whether dileptons are switched on.
Definition: experiment.h:608
Modus modus_
Instance of the Modus template parameter.
Definition: experiment.h:381
std::unique_ptr< RectangularLattice< FourVector > > old_jmu_auxiliary_
Auxiliary lattice for values of jmu at a time step t0.
Definition: experiment.h:493
std::unique_ptr< DecayActionsFinderDilepton > dilepton_finder_
The Dilepton Action Finder.
Definition: experiment.h:427
std::unique_ptr< DensityLattice > j_QBS_lat_
4-current for j_QBS lattice output
Definition: experiment.h:436
bool printout_tmn_landau_
Whether to print the energy-momentum tensor in Landau frame.
Definition: experiment.h:515
std::unique_ptr< RectangularLattice< std::array< FourVector, 4 > > > four_gradient_auxiliary_
Auxiliary lattice for calculating the four-gradient of jmu.
Definition: experiment.h:498
const bool photons_switch_
This indicates whether photons are switched on.
Definition: experiment.h:611
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:537
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:466
ExperimentParameters parameters_
Struct of several member variables.
Definition: experiment.h:372
std::unique_ptr< RectangularLattice< FourVector > > old_fields_auxiliary_
Auxiliary lattice for values of Amu at a time step t0.
Definition: experiment.h:501
std::unique_ptr< DensityLattice > jmu_el_lat_
Electric charge density on the lattice.
Definition: experiment.h:445
std::unique_ptr< Potentials > potentials_
An instance of potentials class, that stores parameters of potentials, calculates them and their grad...
Definition: experiment.h:390
uint64_t wall_actions_total_
Total number of wall-crossings for current timestep.
Definition: experiment.h:671
uint64_t total_hypersurface_crossing_actions_
Total number of particles removed from the evolution in hypersurface crossing actions.
Definition: experiment.h:689
uint64_t interactions_total_
Total number of interactions for current timestep.
Definition: experiment.h:659
const bool use_grid_
This indicates whether to use the grid.
Definition: experiment.h:602
std::unique_ptr< DensityLattice > jmu_I3_lat_
Isospin projection density on the lattice.
Definition: experiment.h:442
int nevents_
Number of events.
Definition: experiment.h:553
uint64_t total_pauli_blocked_
Total number of Pauli-blockings for current timestep.
Definition: experiment.h:683
void final_output()
Output at the end of an event.
Definition: experiment.h:3068
std::vector< bool > projectile_target_interact_
Whether the projectile and the target collided.
Definition: experiment.h:414
Modus * modus()
Provides external access to SMASH calculation modus.
Definition: experiment.h:269
void update_potentials()
Recompute potentials on lattices if necessary.
Definition: experiment.h:2919
OutputPtr dilepton_output_
The Dilepton output.
Definition: experiment.h:405
int64_t seed_
random seed for the next event.
Definition: experiment.h:711
std::vector< Particles > * all_ensembles()
Getter for all ensembles.
Definition: experiment.h:263
uint64_t previous_interactions_total_
Total number of interactions for previous timestep.
Definition: experiment.h:665
uint64_t discarded_interactions_total_
Total number of discarded interactions, because they were invalidated before they could be performed.
Definition: experiment.h:695
void run() override
Runs the experiment.
Definition: experiment.h:3252
bool kinematic_cuts_for_IC_output_
This indicates whether kinematic cuts are enabled for the IC output.
Definition: experiment.h:708
const double end_time_
simulation time at which the evolution is stopped.
Definition: experiment.h:586
double total_energy_violated_by_Pythia_
Total energy violation introduced by Pythia.
Definition: experiment.h:705
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:59
static double formation_power_
Power with which the cross section scaling factor grows in time.
Definition: particledata.h:471
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
A container for storing conserved values.
A container class to hold all the arrays on the lattice and access them.
Definition: lattice.h:49
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:555
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:918
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
static constexpr int LInitialConditions
Definition: experiment.h:94
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:131
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:601
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:333
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:106
void check_interactions_total(uint64_t interactions_total)
Make sure interactions_total can be represented as a 32-bit integer.
Definition: experiment.h:2696
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:137
constexpr double very_small_double
A very small double, used to avoid division by zero.
Definition: constants.h:44
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:381
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:623
bool all_of(Container &&c, UnaryPredicate &&p)
Convenience wrapper for std::all_of that operates on a complete container.
Definition: algorithms.h:80
std::string to_string(ThermodynamicQuantity quantity)
Convert a ThermodynamicQuantity enum value to its corresponding string.
Definition: stringify.cc:26
constexpr double nucleon_mass
Nucleon mass in GeV.
Definition: constants.h:62
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:723
Potentials * pot_pointer
Pointer to a Potential class.
static constexpr int LMain
Definition: experiment.h:93
constexpr double really_small
Numerical error tolerance.
Definition: constants.h:41
RectangularLattice< FourVector > * UB_lat_pointer
Pointer to the skyrme potential on the lattice.
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:32
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:148
Exception class that is thrown if the requested output path in the Experiment factory is not existing...
Definition: experiment.h:157
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:2943
static const Key< bool > collTerm_photons_twoToTwoScatterings
See user guide description for more information.
Definition: input_keys.h:3141
static const Key< DensityType > output_densityType
See user guide description for more information.
Definition: input_keys.h:4901
static const Key< bool > modi_collider_collisionWithinNucleus
See user guide description for more information.
Definition: input_keys.h:3273
static const Key< bool > collTerm_noCollisions
See user guide description for more information.
Definition: input_keys.h:2458
static const Key< int > collTerm_photons_fractionalPhotons
See user guide description for more information.
Definition: input_keys.h:3166
static const Key< int > gen_minNonEmptyEnsembles_number
See user guide description for more information.
Definition: input_keys.h:1222
static const Key< int64_t > gen_randomseed
See user guide description for more information.
Definition: input_keys.h:1197
static const Key< bool > lattice_periodic
See user guide description for more information.
Definition: input_keys.h:5720
static const Key< int > gen_minNonEmptyEnsembles_maximumEnsembles
See user guide description for more information.
Definition: input_keys.h:1210
static const Key< double > gen_expansionRate
See user guide description for more information.
Definition: input_keys.h:1341
static const Key< bool > collTerm_dileptons_decays
See user guide description for more information.
Definition: input_keys.h:3129
static const Key< bool > collTerm_forceDecaysAtEnd
See user guide description for more information.
Definition: input_keys.h:2234
static const Key< double > gen_endTime
See user guide description for more information.
Definition: input_keys.h:1142
static const Key< bool > collTerm_photons_bremsstrahlung
See user guide description for more information.
Definition: input_keys.h:3153
static const Key< PseudoResonance > collTerm_pseudoresonance
See user guide description for more information.
Definition: input_keys.h:2508
static const Key< int > gen_nevents
See user guide description for more information.
Definition: input_keys.h:1184
static const Key< ExpansionMode > gen_metricType
See user guide description for more information.
Definition: input_keys.h:1415
static const Key< std::array< double, 3 > > lattice_origin
See user guide description for more information.
Definition: input_keys.h:5703
static const Key< std::array< double, 3 > > lattice_sizes
See user guide description for more information.
Definition: input_keys.h:5747
static const Key< TotalCrossSectionStrategy > collTerm_totXsStrategy
See user guide description for more information.
Definition: input_keys.h:2629
static const Key< std::array< int, 3 > > lattice_cellNumber
See user guide description for more information.
Definition: input_keys.h:5684
static const Key< bool > gen_useGrid
See user guide description for more information.
Definition: input_keys.h:1583
static const Key< bool > lattice_automatic
See user guide description for more information.
Definition: input_keys.h:5670
static const Key< TimeStepMode > gen_timeStepMode
See user guide description for more information.
Definition: input_keys.h:1556
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.