Version: SMASH-3.0
experiment.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2013-2023
3  * SMASH Team
4  *
5  * GNU General Public License (GPLv3 or later)
6  */
7 #ifndef SRC_INCLUDE_SMASH_EXPERIMENT_H_
8 #define SRC_INCLUDE_SMASH_EXPERIMENT_H_
9 
10 #include <algorithm>
11 #include <limits>
12 #include <memory>
13 #include <string>
14 #include <utility>
15 #include <vector>
16 
17 #include "actionfinderfactory.h"
18 #include "actions.h"
19 #include "bremsstrahlungaction.h"
20 #include "chrono.h"
21 #include "decayactionsfinder.h"
23 #include "energymomentumtensor.h"
24 #include "fields.h"
25 #include "fourvector.h"
26 #include "grandcan_thermalizer.h"
27 #include "grid.h"
29 #include "outputparameters.h"
30 #include "pauliblocking.h"
31 #include "potential_globals.h"
32 #include "potentials.h"
33 #include "propagation.h"
34 #include "quantumnumbers.h"
35 #include "scatteractionphoton.h"
36 #include "scatteractionsfinder.h"
37 #include "stringprocess.h"
38 #include "thermalizationaction.h"
39 // Output
40 #include "binaryoutput.h"
41 #ifdef SMASH_USE_HEPMC
42 #include "hepmcoutput.h"
43 #endif
44 #ifdef SMASH_USE_RIVET
45 #include "rivetoutput.h"
46 #endif
47 #include "icoutput.h"
48 #include "oscaroutput.h"
50 #include "thermodynamicoutput.h"
51 #ifdef SMASH_USE_ROOT
52 #include "rootoutput.h"
53 #endif
54 #include "freeforallaction.h"
55 #include "vtkoutput.h"
56 #include "wallcrossingaction.h"
57 
58 namespace std {
69 template <typename T, typename Ratio>
70 static ostream &operator<<(ostream &out,
71  const chrono::duration<T, Ratio> &seconds) {
72  using Seconds = chrono::duration<double>;
73  using Minutes = chrono::duration<double, std::ratio<60>>;
74  using Hours = chrono::duration<double, std::ratio<60 * 60>>;
75  constexpr Minutes threshold_for_minutes{10};
76  constexpr Hours threshold_for_hours{3};
77  if (seconds < threshold_for_minutes) {
78  return out << Seconds(seconds).count() << " [s]";
79  }
80  if (seconds < threshold_for_hours) {
81  return out << Minutes(seconds).count() << " [min]";
82  }
83  return out << Hours(seconds).count() << " [h]";
84 }
85 } // namespace std
86 
87 namespace smash {
88 static constexpr int LMain = LogArea::Main::id;
89 static constexpr int LInitialConditions = LogArea::InitialConditions::id;
90 
99  public:
100  ExperimentBase() = default;
105  virtual ~ExperimentBase() = default;
106 
127  static std::unique_ptr<ExperimentBase> create(
128  Configuration &config, const std::filesystem::path &output_path);
129 
136  virtual void run() = 0;
137 
143  struct InvalidModusRequest : public std::invalid_argument {
144  using std::invalid_argument::invalid_argument;
145  };
146 
152  struct NonExistingOutputPathRequest : public std::invalid_argument {
153  using std::invalid_argument::invalid_argument;
154  };
155 };
156 
157 template <typename Modus>
158 class Experiment;
159 template <typename Modus>
160 std::ostream &operator<<(std::ostream &out, const Experiment<Modus> &e);
161 
186 template <typename Modus>
187 class Experiment : public ExperimentBase {
188  friend class ExperimentBase;
189 
190  public:
197  void run() override;
198 
212  explicit Experiment(Configuration &config,
213  const std::filesystem::path &output_path);
214 
220  void initialize_new_event();
221 
239  void run_time_evolution(const double t_end, ParticleList &&add_plist = {},
240  ParticleList &&remove_plist = {});
241 
247  void do_final_decays();
248 
250  void final_output();
251 
258  std::vector<Particles> *all_ensembles() { return &ensembles_; }
259 
264  Modus *modus() { return &modus_; }
265 
270  void increase_event_number();
271 
272  private:
285  bool perform_action(Action &action, int i_ensemble,
286  bool include_pauli_blocking = true);
295  void create_output(const std::string &format, const std::string &content,
296  const std::filesystem::path &output_path,
297  const OutputParameters &par);
298 
306  void propagate_and_shine(double to_time, Particles &particles);
307 
320  void run_time_evolution_timestepless(Actions &actions, int i_ensemble,
321  const double end_time_propagation);
322 
324  void intermediate_output();
325 
327  void update_potentials();
328 
337  double compute_min_cell_length(double dt) const {
340  }
341  return std::sqrt(4 * dt * dt + max_transverse_distance_sqr_);
342  }
343 
345  double next_output_time() const {
346  return parameters_.outputclock->next_time();
347  }
348 
354 
360  bool is_finished();
361 
368 
371 
376  Modus modus_;
377 
379  std::vector<Particles> ensembles_;
380 
385  std::unique_ptr<Potentials> potentials_;
386 
391  std::unique_ptr<PauliBlocker> pauli_blocker_;
392 
397  OutputsList outputs_;
398 
400  OutputPtr dilepton_output_;
401 
403  OutputPtr photon_output_;
404 
409  std::vector<bool> projectile_target_interact_;
410 
416  std::vector<FourVector> beam_momentum_ = {};
417 
419  std::vector<std::unique_ptr<ActionFinderInterface>> action_finders_;
420 
422  std::unique_ptr<DecayActionsFinderDilepton> dilepton_finder_;
423 
425  std::unique_ptr<ActionFinderInterface> photon_finder_;
426 
429 
431  std::unique_ptr<DensityLattice> j_QBS_lat_;
432 
434  std::unique_ptr<DensityLattice> jmu_B_lat_;
435 
437  std::unique_ptr<DensityLattice> jmu_I3_lat_;
438 
440  std::unique_ptr<DensityLattice> jmu_el_lat_;
441 
443  std::unique_ptr<FieldsLattice> fields_lat_;
444 
452  std::unique_ptr<DensityLattice> jmu_custom_lat_;
453 
456 
461  std::unique_ptr<RectangularLattice<FourVector>> UB_lat_ = nullptr;
462 
467  std::unique_ptr<RectangularLattice<FourVector>> UI3_lat_ = nullptr;
468 
473  std::unique_ptr<RectangularLattice<std::pair<ThreeVector, ThreeVector>>>
475 
477  std::unique_ptr<RectangularLattice<std::pair<ThreeVector, ThreeVector>>>
479 
481  std::unique_ptr<RectangularLattice<std::pair<ThreeVector, ThreeVector>>>
483 
485  std::unique_ptr<RectangularLattice<EnergyMomentumTensor>> Tmn_;
486 
488  std::unique_ptr<RectangularLattice<FourVector>> old_jmu_auxiliary_;
490  std::unique_ptr<RectangularLattice<FourVector>> new_jmu_auxiliary_;
492  std::unique_ptr<RectangularLattice<std::array<FourVector, 4>>>
494 
496  std::unique_ptr<RectangularLattice<FourVector>> old_fields_auxiliary_;
498  std::unique_ptr<RectangularLattice<FourVector>> new_fields_auxiliary_;
500  std::unique_ptr<RectangularLattice<std::array<FourVector, 4>>>
502 
504  bool printout_tmn_ = false;
505 
507  bool printout_tmn_landau_ = false;
508 
510  bool printout_v_landau_ = false;
511 
513  bool printout_j_QBS_ = false;
514 
516  bool printout_lattice_td_ = false;
517 
521 
525 
529 
531  std::unique_ptr<GrandCanThermalizer> thermalizer_;
532 
538 
553  const 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 
617  const bool IC_output_switch_;
618 
621 
626  double max_transverse_distance_sqr_ = std::numeric_limits<double>::max();
627 
636 
642 
644  SystemTimePoint time_start_ = SystemClock::now();
645 
648 
653  uint64_t interactions_total_ = 0;
654 
660 
665  uint64_t wall_actions_total_ = 0;
666 
672 
677  uint64_t total_pauli_blocked_ = 0;
678 
684 
690 
694  double total_energy_removed_ = 0.0;
695 
700 
703 
705  int64_t seed_ = -1;
706 
712  friend std::ostream &operator<<<>(std::ostream &out, const Experiment &e);
713 };
714 
716 template <typename Modus>
717 std::ostream &operator<<(std::ostream &out, const Experiment<Modus> &e) {
718  out << "End time: " << e.end_time_ << " fm\n";
719  out << e.modus_;
720  return out;
721 }
722 
723 template <typename Modus>
724 void Experiment<Modus>::create_output(const std::string &format,
725  const std::string &content,
726  const std::filesystem::path &output_path,
727  const OutputParameters &out_par) {
728  logg[LExperiment].info() << "Adding output " << content << " of format "
729  << format << std::endl;
730 
731  if (format == "VTK" && content == "Particles") {
732  outputs_.emplace_back(
733  std::make_unique<VtkOutput>(output_path, content, out_par));
734  } else if (format == "Root") {
735 #ifdef SMASH_USE_ROOT
736  if (content == "Initial_Conditions") {
737  outputs_.emplace_back(
738  std::make_unique<RootOutput>(output_path, "SMASH_IC", out_par));
739  } else {
740  outputs_.emplace_back(
741  std::make_unique<RootOutput>(output_path, content, out_par));
742  }
743 #else
744  logg[LExperiment].error(
745  "Root output requested, but Root support not compiled in");
746 #endif
747  } else if (format == "Binary") {
748  if (content == "Collisions" || content == "Dileptons" ||
749  content == "Photons") {
750  outputs_.emplace_back(std::make_unique<BinaryOutputCollisions>(
751  output_path, content, out_par));
752  } else if (content == "Particles") {
753  outputs_.emplace_back(std::make_unique<BinaryOutputParticles>(
754  output_path, content, out_par));
755  } else if (content == "Initial_Conditions") {
756  outputs_.emplace_back(std::make_unique<BinaryOutputInitialConditions>(
757  output_path, content, out_par));
758  }
759  } else if (format == "Oscar1999" || format == "Oscar2013") {
760  outputs_.emplace_back(
761  create_oscar_output(format, content, output_path, out_par));
762  } else if (content == "Thermodynamics" && format == "ASCII") {
763  outputs_.emplace_back(
764  std::make_unique<ThermodynamicOutput>(output_path, content, out_par));
765  } else if (content == "Thermodynamics" &&
766  (format == "Lattice_ASCII" || format == "Lattice_Binary")) {
767  printout_full_lattice_any_td_ = true;
768  if (format == "Lattice_ASCII") {
769  printout_full_lattice_ascii_td_ = true;
770  }
771  if (format == "Lattice_Binary") {
772  printout_full_lattice_binary_td_ = true;
773  }
774  outputs_.emplace_back(std::make_unique<ThermodynamicLatticeOutput>(
775  output_path, content, out_par, printout_full_lattice_ascii_td_,
776  printout_full_lattice_binary_td_));
777  } else if (content == "Thermodynamics" && format == "VTK") {
778  printout_lattice_td_ = true;
779  outputs_.emplace_back(
780  std::make_unique<VtkOutput>(output_path, content, out_par));
781  } else if (content == "Initial_Conditions" && format == "ASCII") {
782  outputs_.emplace_back(
783  std::make_unique<ICOutput>(output_path, "SMASH_IC", out_par));
784  } else if ((format == "HepMC") || (format == "HepMC_asciiv3") ||
785  (format == "HepMC_treeroot")) {
786 #ifdef SMASH_USE_HEPMC
787  if (content == "Particles") {
788  if ((format == "HepMC") || (format == "HepMC_asciiv3")) {
789  outputs_.emplace_back(std::make_unique<HepMcOutput>(
790  output_path, "SMASH_HepMC_particles", false, "asciiv3"));
791  } else if (format == "HepMC_treeroot") {
792 #ifdef SMASH_USE_HEPMC_ROOTIO
793  outputs_.emplace_back(std::make_unique<HepMcOutput>(
794  output_path, "SMASH_HepMC_particles", false, "root"));
795 #else
796  logg[LExperiment].error(
797  "Requested HepMC_treeroot output not available, "
798  "ROOT or HepMC3 ROOTIO missing or not found by cmake.");
799 #endif
800  }
801  } else if (content == "Collisions") {
802  if ((format == "HepMC") || (format == "HepMC_asciiv3")) {
803  outputs_.emplace_back(std::make_unique<HepMcOutput>(
804  output_path, "SMASH_HepMC_collisions", true, "asciiv3"));
805  } else if (format == "HepMC_treeroot") {
806 #ifdef SMASH_USE_HEPMC_ROOTIO
807  outputs_.emplace_back(std::make_unique<HepMcOutput>(
808  output_path, "SMASH_HepMC_collisions", true, "root"));
809 #else
810  logg[LExperiment].error(
811  "Requested HepMC_treeroot output not available, "
812  "ROOT or HepMC3 ROOTIO missing or not found by cmake.");
813 #endif
814  }
815  } else {
816  logg[LExperiment].error(
817  "HepMC only available for Particles and "
818  "Collisions content. Requested for " +
819  content + ".");
820  }
821 #else
822  logg[LExperiment].error(
823  "HepMC output requested, but HepMC support not compiled in");
824 #endif
825  } else if (content == "Coulomb" && format == "VTK") {
826  outputs_.emplace_back(
827  std::make_unique<VtkOutput>(output_path, "Fields", out_par));
828  } else if (content == "Rivet") {
829 #ifdef SMASH_USE_RIVET
830  // flag to ensure that the Rivet format has not been already assigned
831  static bool rivet_format_already_selected = false;
832  // if the next check is true, then we are trying to assign the format twice
833  if (rivet_format_already_selected) {
834  logg[LExperiment].warn(
835  "Rivet output format can only be one, either YODA or YODA-full. "
836  "Only your first valid choice will be used.");
837  return;
838  }
839  if (format == "YODA") {
840  outputs_.emplace_back(std::make_unique<RivetOutput>(
841  output_path, "SMASH_Rivet", false, out_par.rivet_parameters));
842  rivet_format_already_selected = true;
843  } else if (format == "YODA-full") {
844  outputs_.emplace_back(std::make_unique<RivetOutput>(
845  output_path, "SMASH_Rivet_full", true, out_par.rivet_parameters));
846  rivet_format_already_selected = true;
847  } else {
848  logg[LExperiment].error("Rivet format " + format +
849  "not one of YODA or YODA-full");
850  }
851 #else
852  logg[LExperiment].error(
853  "Rivet output requested, but Rivet support not compiled in");
854 #endif
855  } else {
856  logg[LExperiment].error()
857  << "Unknown combination of format (" << format << ") and content ("
858  << content << "). Fix the config.";
859  }
860 }
861 
870 
871 template <typename Modus>
873  const std::filesystem::path &output_path)
874  : parameters_(create_experiment_parameters(config)),
875  density_param_(DensityParameters(parameters_)),
876  modus_(std::invoke([&]() {
877  /* This immediately invoked lambda is a work-around to cope with the
878  * fact that the "Collisions_Within_Nucleus" key belongs to the
879  * "Collider" section, but is used by the ScatterActionsFinder through
880  * the ScatterActionsFinderParameters member. Here that key is taken
881  * from the main configuration and put there back after the "Collider"
882  * section is extracted. If this were not done in this way, the
883  * sub-configuration given to ColliderModus would be deleted at the end
884  * of its constructor not empty and this would throw an exception.*/
885  const std::initializer_list<const char *> key_labels = {
886  "Modi", "Collider", "Collisions_Within_Nucleus"};
887  const bool restore_key = config.has_value(key_labels);
888  const bool temporary_taken_key = config.take(key_labels, false);
889  auto modus_config = config.extract_sub_configuration({"Modi"});
890  if (restore_key) {
891  config.set_value(key_labels, temporary_taken_key);
892  }
893  return Modus{std::move(modus_config), parameters_};
894  })),
895  ensembles_(parameters_.n_ensembles),
896  nevents_(config.take({"General", "Nevents"}, 0)),
897  end_time_(config.take({"General", "End_Time"})),
898  delta_time_startup_(parameters_.labclock->timestep_duration()),
899  force_decays_(
900  config.take({"Collision_Term", "Force_Decays_At_End"}, true)),
901  use_grid_(config.take({"General", "Use_Grid"}, true)),
902  metric_(
903  config.take({"General", "Metric_Type"}, ExpansionMode::NoExpansion),
904  config.take({"General", "Expansion_Rate"}, 0.1)),
905  dileptons_switch_(
906  config.take({"Collision_Term", "Dileptons", "Decays"}, false)),
907  photons_switch_(config.take(
908  {"Collision_Term", "Photons", "2to2_Scatterings"}, false)),
909  bremsstrahlung_switch_(
910  config.take({"Collision_Term", "Photons", "Bremsstrahlung"}, false)),
911  IC_output_switch_(config.has_value({"Output", "Initial_Conditions"})),
912  time_step_mode_(
913  config.take({"General", "Time_Step_Mode"}, TimeStepMode::Fixed)) {
914  logg[LExperiment].info() << *this;
915 
916  if (config.has_value({"General", "Minimum_Nonempty_Ensembles"})) {
917  if (config.has_value({"General", "Nevents"})) {
918  throw std::invalid_argument(
919  "Please specify either Nevents or Minimum_Nonempty_Ensembles.");
920  }
921  event_counting_ = EventCounting::MinimumNonEmpty;
922  minimum_nonempty_ensembles_ =
923  config.take({"General", "Minimum_Nonempty_Ensembles", "Number"});
924  int max_ensembles = config.take(
925  {"General", "Minimum_Nonempty_Ensembles", "Maximum_Ensembles_Run"});
926  max_events_ =
927  std::ceil(static_cast<double>(max_ensembles) / parameters_.n_ensembles);
928  } else {
929  event_counting_ = EventCounting::FixedNumber;
930  }
931 
932  // covariant derivatives can only be done with covariant smearing
933  if (parameters_.derivatives_mode == DerivativesMode::CovariantGaussian &&
934  parameters_.smearing_mode != SmearingMode::CovariantGaussian) {
935  throw std::invalid_argument(
936  "Covariant Gaussian derivatives only make sense for Covariant Gaussian "
937  "smearing!");
938  }
939 
940  // for triangular smearing:
941  // the weight needs to be larger than 1./7. for the center cell to contribute
942  // more than the surrounding cells
943  if (parameters_.smearing_mode == SmearingMode::Discrete &&
944  parameters_.discrete_weight < (1. / 7.)) {
945  throw std::invalid_argument(
946  "The central weight for discrete smearing should be >= 1./7.");
947  }
948 
949  if (parameters_.coll_crit == CollisionCriterion::Stochastic &&
950  (time_step_mode_ != TimeStepMode::Fixed || !use_grid_)) {
951  throw std::invalid_argument(
952  "The stochastic criterion can only be employed for fixed time step "
953  "mode and with a grid!");
954  }
955 
956  if (modus_.is_box() && (time_step_mode_ != TimeStepMode::Fixed)) {
957  throw std::invalid_argument(
958  "The box modus can only be used with the fixed time step mode!");
959  }
960 
961  logg[LExperiment].info("Using ", parameters_.testparticles,
962  " testparticles per particle.");
963  logg[LExperiment].info("Using ", parameters_.n_ensembles,
964  " parallel ensembles.");
965 
966  // create finders
967  if (dileptons_switch_) {
968  dilepton_finder_ = std::make_unique<DecayActionsFinderDilepton>();
969  }
970  if (photons_switch_ || bremsstrahlung_switch_) {
971  n_fractional_photons_ =
972  config.take({"Collision_Term", "Photons", "Fractional_Photons"}, 100);
973  }
974  if (parameters_.two_to_one) {
975  if (parameters_.res_lifetime_factor < 0.) {
976  throw std::invalid_argument(
977  "Resonance lifetime modifier cannot be negative!");
978  }
979  if (parameters_.res_lifetime_factor < really_small) {
980  logg[LExperiment].warn(
981  "Resonance lifetime set to zero. Make sure resonances cannot "
982  "interact",
983  "inelastically (e.g. resonance chains), else SMASH is known to "
984  "hang.");
985  }
986  action_finders_.emplace_back(std::make_unique<DecayActionsFinder>(
987  parameters_.res_lifetime_factor, parameters_.do_weak_decays));
988  }
989  bool no_coll = config.take({"Collision_Term", "No_Collisions"}, false);
990  if ((parameters_.two_to_one || parameters_.included_2to2.any() ||
991  parameters_.included_multi.any() || parameters_.strings_switch) &&
992  !no_coll) {
993  parameters_.use_monash_tune_default =
994  (modus_.is_collider() && modus_.sqrt_s_NN() >= 200.);
995  auto scat_finder =
996  std::make_unique<ScatterActionsFinder>(config, parameters_);
997  max_transverse_distance_sqr_ =
998  scat_finder->max_transverse_distance_sqr(parameters_.testparticles);
999  process_string_ptr_ = scat_finder->get_process_string_ptr();
1000  action_finders_.emplace_back(std::move(scat_finder));
1001  } else {
1002  max_transverse_distance_sqr_ =
1003  parameters_.maximum_cross_section / M_PI * fm2_mb;
1004  process_string_ptr_ = NULL;
1005  }
1006  if (modus_.is_box()) {
1007  action_finders_.emplace_back(
1008  std::make_unique<WallCrossActionsFinder>(parameters_.box_length));
1009  }
1010  if (IC_output_switch_) {
1011  if (!modus_.is_collider()) {
1012  throw std::runtime_error(
1013  "Initial conditions can only be extracted in collider modus.");
1014  }
1015  double proper_time;
1016  if (config.has_value({"Output", "Initial_Conditions", "Proper_Time"})) {
1017  // Read in proper time from config
1018  proper_time =
1019  config.take({"Output", "Initial_Conditions", "Proper_Time"});
1020  } else {
1021  // Default proper time is the passing time of the two nuclei
1022  double default_proper_time = modus_.nuclei_passing_time();
1023  double lower_bound =
1024  config.take({"Output", "Initial_Conditions", "Lower_Bound"}, 0.5);
1025  if (default_proper_time >= lower_bound) {
1026  proper_time = default_proper_time;
1027  } else {
1028  logg[LInitialConditions].warn()
1029  << "Nuclei passing time is too short, hypersurface proper time set "
1030  "to tau = "
1031  << lower_bound << " fm.";
1032  proper_time = lower_bound;
1033  }
1034  }
1035 
1036  double rapidity_cut = 0.0;
1037  if (config.has_value({"Output", "Initial_Conditions", "Rapidity_Cut"})) {
1038  rapidity_cut =
1039  config.take({"Output", "Initial_Conditions", "Rapidity_Cut"});
1040  if (rapidity_cut <= 0.0) {
1041  logg[LInitialConditions].fatal()
1042  << "Rapidity cut for initial conditions configured as abs(y) < "
1043  << rapidity_cut << " is unreasonable. \nPlease choose a positive, "
1044  << "non-zero value or employ SMASH without rapidity cut.";
1045  throw std::runtime_error(
1046  "Kinematic cut for initial conditions malconfigured.");
1047  }
1048  }
1049 
1050  if (modus_.calculation_frame_is_fixed_target() && rapidity_cut != 0.0) {
1051  throw std::runtime_error(
1052  "Rapidity cut for initial conditions output is not implemented "
1053  "in the fixed target calculation frame. \nPlease use "
1054  "\"center of velocity\" or \"center of mass\" as a "
1055  "\"Calculation_Frame\" instead.");
1056  }
1057 
1058  double transverse_momentum_cut = 0.0;
1059  if (config.has_value({"Output", "Initial_Conditions", "pT_Cut"})) {
1060  transverse_momentum_cut =
1061  config.take({"Output", "Initial_Conditions", "pT_Cut"});
1062  if (transverse_momentum_cut <= 0.0) {
1063  logg[LInitialConditions].fatal()
1064  << "transverse momentum cut for initial conditions configured as "
1065  "pT < "
1066  << rapidity_cut << " is unreasonable. \nPlease choose a positive, "
1067  << "non-zero value or employ SMASH without pT cut.";
1068  throw std::runtime_error(
1069  "Kinematic cut for initial conditions misconfigured.");
1070  }
1071  }
1072 
1073  if (rapidity_cut > 0.0 || transverse_momentum_cut > 0.0) {
1074  kinematic_cuts_for_IC_output_ = true;
1075  }
1076 
1077  if (rapidity_cut > 0.0 && transverse_momentum_cut > 0.0) {
1078  logg[LInitialConditions].info()
1079  << "Extracting initial conditions in kinematic range: "
1080  << -rapidity_cut << " <= y <= " << rapidity_cut
1081  << "; pT <= " << transverse_momentum_cut << " GeV.";
1082  } else if (rapidity_cut > 0.0) {
1083  logg[LInitialConditions].info()
1084  << "Extracting initial conditions in kinematic range: "
1085  << -rapidity_cut << " <= y <= " << rapidity_cut << ".";
1086  } else if (transverse_momentum_cut > 0.0) {
1087  logg[LInitialConditions].info()
1088  << "Extracting initial conditions in kinematic range: pT <= "
1089  << transverse_momentum_cut << " GeV.";
1090  } else {
1091  logg[LInitialConditions].info()
1092  << "Extracting initial conditions without kinematic cuts.";
1093  }
1094 
1095  action_finders_.emplace_back(
1096  std::make_unique<HyperSurfaceCrossActionsFinder>(
1097  proper_time, rapidity_cut, transverse_momentum_cut));
1098  }
1099 
1100  if (config.has_value({"Collision_Term", "Pauli_Blocking"})) {
1101  logg[LExperiment].info() << "Pauli blocking is ON.";
1102  pauli_blocker_ = std::make_unique<PauliBlocker>(
1103  config.extract_sub_configuration({"Collision_Term", "Pauli_Blocking"}),
1104  parameters_);
1105  }
1106  // In collider setup with sqrts >= 200 GeV particles don't form continuously
1108  {"Collision_Term", "String_Parameters", "Power_Particle_Formation"},
1109  modus_.sqrt_s_NN() >= 200. ? -1. : 1.);
1110 
1350  // create outputs
1352  " create OutputInterface objects");
1353  dens_type_ = config.take({"Output", "Density_Type"}, DensityType::None);
1354  logg[LExperiment].debug()
1355  << "Density type printed to headers: " << dens_type_;
1356 
1357  /* Parse configuration about output contents and formats, doing all logical
1358  * checks about specified formats, creating all needed output objects. */
1359  auto output_conf = config.extract_sub_configuration(
1360  {"Output"}, Configuration::GetEmpty::Yes);
1361  if (output_path == "") {
1362  throw std::invalid_argument(
1363  "Invalid empty output path provided to Experiment constructor.");
1364  }
1365  if (output_conf.is_empty()) {
1366  logg[LExperiment].warn() << "No \"Output\" section found in the input "
1367  "file. No output file will be produced.";
1368  }
1369  const std::vector<std::string> output_contents =
1370  output_conf.list_upmost_nodes();
1371  std::vector<std::vector<std::string>> list_of_formats(output_contents.size());
1372  std::transform(
1373  output_contents.cbegin(), output_contents.cend(), list_of_formats.begin(),
1374  [&output_conf](std::string content) -> std::vector<std::string> {
1375  /* Use here a default value for "Format" even though it is a required
1376  * key, just because then here below the error for the user is more
1377  * informative, if the key was not given in the input file. */
1378  return output_conf.take({content.c_str(), "Format"},
1379  std::vector<std::string>{});
1380  });
1381  const OutputParameters output_parameters(std::move(output_conf));
1382  std::size_t total_number_of_requested_formats = 0;
1383  auto abort_because_of_invalid_input_file = []() {
1384  throw std::invalid_argument("Invalid configuration input file.");
1385  };
1386  for (std::size_t i = 0; i < output_contents.size(); ++i) {
1387  if (list_of_formats[i].empty()) {
1388  logg[LExperiment].fatal()
1389  << "Empty or unspecified list of formats for "
1390  << std::quoted(output_contents[i]) << " content.";
1391  abort_because_of_invalid_input_file();
1392  } else if (std::find(list_of_formats[i].begin(), list_of_formats[i].end(),
1393  "None") != list_of_formats[i].end()) {
1394  if (list_of_formats[i].size() > 1) {
1395  logg[LExperiment].fatal()
1396  << "Use of \"None\" output format together with other formats is "
1397  "not allowed.\nInvalid \"Format\" key for "
1398  << std::quoted(output_contents[i]) << " content.";
1399  abort_because_of_invalid_input_file();
1400  } else {
1401  // Clear vector so that the for below is skipped and no output created
1402  list_of_formats[i].clear();
1403  }
1404  }
1405  for (const auto &format : list_of_formats[i]) {
1406  create_output(format, output_contents[i], output_path, output_parameters);
1407  ++total_number_of_requested_formats;
1408  }
1409  }
1410  if (outputs_.size() != total_number_of_requested_formats) {
1411  logg[LExperiment].fatal()
1412  << "At least one invalid output format has been provided.";
1413  abort_because_of_invalid_input_file();
1414  }
1415 
1416  /* We can take away the Fermi motion flag, because the collider modus is
1417  * already initialized. We only need it when potentials are enabled, but we
1418  * always have to take it, otherwise SMASH will complain about unused
1419  * options. We have to provide a default value for modi other than Collider.
1420  */
1421  if (config.has_value({"Potentials"})) {
1422  if (time_step_mode_ == TimeStepMode::None) {
1423  logg[LExperiment].error() << "Potentials only work with time steps!";
1424  throw std::invalid_argument("Can't use potentials without time steps!");
1425  }
1426  if (modus_.fermi_motion() == FermiMotion::Frozen) {
1427  logg[LExperiment].error()
1428  << "Potentials don't work with frozen Fermi momenta! "
1429  "Use normal Fermi motion instead.";
1430  throw std::invalid_argument(
1431  "Can't use potentials "
1432  "with frozen Fermi momenta!");
1433  }
1434  logg[LExperiment].info() << "Potentials are ON. Timestep is "
1435  << parameters_.labclock->timestep_duration();
1436  // potentials need density calculation parameters from parameters_
1437  potentials_ = std::make_unique<Potentials>(
1438  config.extract_sub_configuration({"Potentials"}), parameters_);
1439  // make sure that vdf potentials are not used together with Skyrme
1440  // or symmetry potentials
1441  if (potentials_->use_skyrme() && potentials_->use_vdf()) {
1442  throw std::runtime_error(
1443  "Can't use Skyrme and VDF potentials at the same time!");
1444  }
1445  if (potentials_->use_symmetry() && potentials_->use_vdf()) {
1446  throw std::runtime_error(
1447  "Can't use symmetry and VDF potentials at the same time!");
1448  }
1449  if (potentials_->use_skyrme()) {
1450  logg[LExperiment].info() << "Skyrme potentials are:\n";
1451  logg[LExperiment].info()
1452  << "\t\tSkyrme_A [MeV] = " << potentials_->skyrme_a() << "\n";
1453  logg[LExperiment].info()
1454  << "\t\tSkyrme_B [MeV] = " << potentials_->skyrme_b() << "\n";
1455  logg[LExperiment].info()
1456  << "\t\t Skyrme_tau = " << potentials_->skyrme_tau() << "\n";
1457  }
1458  if (potentials_->use_symmetry()) {
1459  logg[LExperiment].info()
1460  << "Symmetry potential is:"
1461  << "\n S_pot [MeV] = " << potentials_->symmetry_S_pot() << "\n";
1462  }
1463  if (potentials_->use_vdf()) {
1464  logg[LExperiment].info() << "VDF potential parameters are:\n";
1465  logg[LExperiment].info() << "\t\tsaturation density [fm^-3] = "
1466  << potentials_->saturation_density() << "\n";
1467  for (int i = 0; i < potentials_->number_of_terms(); i++) {
1468  logg[LExperiment].info()
1469  << "\t\tCoefficient_" << i + 1 << " = "
1470  << 1000.0 * (potentials_->coeffs())[i] << " [MeV] \t Power_"
1471  << i + 1 << " = " << (potentials_->powers())[i] << "\n";
1472  }
1473  }
1474  // if potentials are on, derivatives need to be calculated
1475  if (parameters_.derivatives_mode == DerivativesMode::Off &&
1476  parameters_.field_derivatives_mode == FieldDerivativesMode::ChainRule) {
1477  throw std::invalid_argument(
1478  "Derivatives are necessary for running with potentials.\n"
1479  "Derivatives_Mode: \"Off\" only makes sense for "
1480  "Field_Derivatives_Mode: \"Direct\"!\nUse \"Covariant Gaussian\" or "
1481  "\"Finite difference\".");
1482  }
1483  // for computational efficiency, we want to turn off the derivatives of jmu
1484  // and the rest frame density derivatives if direct derivatives are used
1485  if (parameters_.field_derivatives_mode == FieldDerivativesMode::Direct) {
1486  parameters_.derivatives_mode = DerivativesMode::Off;
1487  parameters_.rho_derivatives_mode = RestFrameDensityDerivativesMode::Off;
1488  }
1489  switch (parameters_.derivatives_mode) {
1491  logg[LExperiment].info() << "Covariant Gaussian derivatives are ON";
1492  break;
1494  logg[LExperiment].info() << "Finite difference derivatives are ON";
1495  break;
1496  case DerivativesMode::Off:
1497  logg[LExperiment].info() << "Gradients of baryon current are OFF";
1498  break;
1499  }
1500  switch (parameters_.rho_derivatives_mode) {
1502  logg[LExperiment].info() << "Rest frame density derivatives are ON";
1503  break;
1505  logg[LExperiment].info() << "Rest frame density derivatives are OFF";
1506  break;
1507  }
1508  // direct or chain rule derivatives only make sense for the VDF potentials
1509  if (potentials_->use_vdf()) {
1510  switch (parameters_.field_derivatives_mode) {
1512  logg[LExperiment].info() << "Chain rule field derivatives are ON";
1513  break;
1515  logg[LExperiment].info() << "Direct field derivatives are ON";
1516  break;
1517  }
1518  }
1519  /*
1520  * Necessary safety checks
1521  */
1522  // VDF potentials need derivatives of rest frame density or fields
1523  if (potentials_->use_vdf() && (parameters_.rho_derivatives_mode ==
1525  parameters_.field_derivatives_mode ==
1527  throw std::runtime_error(
1528  "Can't use VDF potentials without rest frame density derivatives or "
1529  "direct field derivatives!");
1530  }
1531  // potentials require using gradients
1532  if (parameters_.derivatives_mode == DerivativesMode::Off &&
1533  parameters_.field_derivatives_mode == FieldDerivativesMode::ChainRule) {
1534  throw std::runtime_error(
1535  "Can't use potentials without gradients of baryon current (Skyrme, "
1536  "VDF)"
1537  " or direct field derivatives (VDF)!");
1538  }
1539  // direct field derivatives only make sense for the VDF potentials
1540  if (!(potentials_->use_vdf()) &&
1541  parameters_.field_derivatives_mode == FieldDerivativesMode::Direct) {
1542  throw std::invalid_argument(
1543  "Field_Derivatives_Mode: \"Direct\" only makes sense for the VDF "
1544  "potentials!\nUse Field_Derivatives_Mode: \"Chain Rule\" or comment "
1545  "this option out (Chain Rule is default)");
1546  }
1547  }
1548 
1549  // information about the type of smearing
1550  switch (parameters_.smearing_mode) {
1552  logg[LExperiment].info() << "Smearing type: Covariant Gaussian";
1553  break;
1555  logg[LExperiment].info() << "Smearing type: Discrete with weight = "
1556  << parameters_.discrete_weight;
1557  break;
1559  logg[LExperiment].info() << "Smearing type: Triangular with range = "
1560  << parameters_.triangular_range;
1561  break;
1562  }
1563 
1564  // Create lattices
1565  if (config.has_value({"Lattice"})) {
1566  bool automatic = config.take({"Lattice", "Automatic"}, false);
1567  bool all_geometrical_properties_specified =
1568  config.has_value({"Lattice", "Cell_Number"}) &&
1569  config.has_value({"Lattice", "Origin"}) &&
1570  config.has_value({"Lattice", "Sizes"});
1571  if (!automatic && !all_geometrical_properties_specified) {
1572  throw std::invalid_argument(
1573  "The lattice was requested to be manually generated, but some\n"
1574  "lattice geometrical property was not specified. Be sure to provide\n"
1575  "both \"Cell_Number\" and \"Origin\" and \"Sizes\".");
1576  }
1577  if (automatic && all_geometrical_properties_specified) {
1578  throw std::invalid_argument(
1579  "The lattice was requested to be automatically generated, but all\n"
1580  "lattice geometrical properties were specified. In this case you\n"
1581  "need to set \"Automatic: False\".");
1582  }
1583  bool periodic = config.take({"Lattice", "Periodic"}, modus_.is_box());
1584  const auto [l, n, origin] = [&config, automatic, this]() {
1585  if (!automatic) {
1586  return std::make_tuple<std::array<double, 3>, std::array<int, 3>,
1587  std::array<double, 3>>(
1588  config.take({"Lattice", "Sizes"}),
1589  config.take({"Lattice", "Cell_Number"}),
1590  config.take({"Lattice", "Origin"}));
1591  } else {
1592  std::array<double, 3> l_default{20., 20., 20.};
1593  std::array<int, 3> n_default{10, 10, 10};
1594  std::array<double, 3> origin_default{-20., -20., -20.};
1595  if (modus_.is_collider() || (modus_.is_list() && !modus_.is_box())) {
1596  // Estimates on how far particles could get in x, y, z. The
1597  // default lattice is currently not contracted for afterburner runs
1598  const double gam = modus_.is_collider()
1599  ? modus_.sqrt_s_NN() / (2.0 * nucleon_mass)
1600  : 1.0;
1601  const double max_z = 5.0 / gam + end_time_;
1602  const double estimated_max_transverse_velocity = 0.7;
1603  const double max_xy =
1604  5.0 + estimated_max_transverse_velocity * end_time_;
1605  origin_default = {-max_xy, -max_xy, -max_z};
1606  l_default = {2 * max_xy, 2 * max_xy, 2 * max_z};
1607  // Go for approximately 0.8 fm cell size and contract
1608  // lattice in z by gamma factor
1609  const int n_xy = std::ceil(2 * max_xy / 0.8);
1610  int nz = std::ceil(2 * max_z / 0.8);
1611  // Contract lattice by gamma factor in case of smearing where
1612  // smearing length is bound to the lattice cell length
1613  if (parameters_.smearing_mode == SmearingMode::Discrete ||
1614  parameters_.smearing_mode == SmearingMode::Triangular) {
1615  nz = static_cast<int>(std::ceil(2 * max_z / 0.8 * gam));
1616  }
1617  n_default = {n_xy, n_xy, nz};
1618  } else if (modus_.is_box()) {
1619  origin_default = {0., 0., 0.};
1620  const double bl = modus_.length();
1621  l_default = {bl, bl, bl};
1622  const int n_xyz = std::ceil(bl / 0.5);
1623  n_default = {n_xyz, n_xyz, n_xyz};
1624  } else if (modus_.is_sphere()) {
1625  // Maximal distance from (0, 0, 0) at which a particle
1626  // may be found at the end of the simulation
1627  const double max_d = modus_.radius() + end_time_;
1628  origin_default = {-max_d, -max_d, -max_d};
1629  l_default = {2 * max_d, 2 * max_d, 2 * max_d};
1630  // Go for approximately 0.8 fm cell size
1631  const int n_xyz = std::ceil(2 * max_d / 0.8);
1632  n_default = {n_xyz, n_xyz, n_xyz};
1633  }
1634  // Take lattice properties from config to assign them to all lattices
1635  return std::make_tuple<std::array<double, 3>, std::array<int, 3>,
1636  std::array<double, 3>>(
1637  config.take({"Lattice", "Sizes"}, l_default),
1638  config.take({"Lattice", "Cell_Number"}, n_default),
1639  config.take({"Lattice", "Origin"}, origin_default));
1640  }
1641  }();
1642 
1643  logg[LExperiment].info()
1644  << "Lattice is ON. Origin = (" << origin[0] << "," << origin[1] << ","
1645  << origin[2] << "), sizes = (" << l[0] << "," << l[1] << "," << l[2]
1646  << "), number of cells = (" << n[0] << "," << n[1] << "," << n[2]
1647  << "), periodic = " << std::boolalpha << periodic;
1648 
1649  if (printout_lattice_td_ || printout_full_lattice_any_td_) {
1650  dens_type_lattice_printout_ = output_parameters.td_dens_type;
1651  printout_tmn_ = output_parameters.td_tmn;
1652  printout_tmn_landau_ = output_parameters.td_tmn_landau;
1653  printout_v_landau_ = output_parameters.td_v_landau;
1654  printout_j_QBS_ = output_parameters.td_jQBS;
1655  }
1656  if (printout_tmn_ || printout_tmn_landau_ || printout_v_landau_) {
1657  Tmn_ = std::make_unique<RectangularLattice<EnergyMomentumTensor>>(
1658  l, n, origin, periodic, LatticeUpdate::AtOutput);
1659  }
1660  if (printout_j_QBS_) {
1661  j_QBS_lat_ = std::make_unique<DensityLattice>(l, n, origin, periodic,
1662  LatticeUpdate::AtOutput);
1663  }
1664  /* Create baryon and isospin density lattices regardless of config
1665  if potentials are on. This is because they allow to compute
1666  potentials faster */
1667  if (potentials_) {
1668  // Create auxiliary lattices for baryon four-current calculation
1669  old_jmu_auxiliary_ = std::make_unique<RectangularLattice<FourVector>>(
1670  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1671  new_jmu_auxiliary_ = std::make_unique<RectangularLattice<FourVector>>(
1672  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1673  four_gradient_auxiliary_ =
1674  std::make_unique<RectangularLattice<std::array<FourVector, 4>>>(
1675  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1676 
1677  if (potentials_->use_skyrme()) {
1678  jmu_B_lat_ = std::make_unique<DensityLattice>(
1679  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1680  UB_lat_ = std::make_unique<RectangularLattice<FourVector>>(
1681  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1682  FB_lat_ = std::make_unique<
1683  RectangularLattice<std::pair<ThreeVector, ThreeVector>>>(
1684  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1685  }
1686  if (potentials_->use_symmetry()) {
1687  jmu_I3_lat_ = std::make_unique<DensityLattice>(
1688  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1689  UI3_lat_ = std::make_unique<RectangularLattice<FourVector>>(
1690  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1691  FI3_lat_ = std::make_unique<
1692  RectangularLattice<std::pair<ThreeVector, ThreeVector>>>(
1693  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1694  }
1695  if (potentials_->use_coulomb()) {
1696  jmu_el_lat_ = std::make_unique<DensityLattice>(
1697  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1698  EM_lat_ = std::make_unique<
1699  RectangularLattice<std::pair<ThreeVector, ThreeVector>>>(
1700  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1701  }
1702  if (potentials_->use_vdf()) {
1703  jmu_B_lat_ = std::make_unique<DensityLattice>(
1704  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1705  UB_lat_ = std::make_unique<RectangularLattice<FourVector>>(
1706  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1707  FB_lat_ = std::make_unique<
1708  RectangularLattice<std::pair<ThreeVector, ThreeVector>>>(
1709  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1710  }
1711  if (parameters_.field_derivatives_mode == FieldDerivativesMode::Direct) {
1712  // Create auxiliary lattices for field calculation
1713  old_fields_auxiliary_ =
1714  std::make_unique<RectangularLattice<FourVector>>(
1715  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1716  new_fields_auxiliary_ =
1717  std::make_unique<RectangularLattice<FourVector>>(
1718  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1719  fields_four_gradient_auxiliary_ =
1720  std::make_unique<RectangularLattice<std::array<FourVector, 4>>>(
1721  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1722 
1723  // Create the fields lattice
1724  fields_lat_ = std::make_unique<FieldsLattice>(
1725  l, n, origin, periodic, LatticeUpdate::EveryTimestep);
1726  }
1727  } else {
1728  if (dens_type_lattice_printout_ == DensityType::Baryon) {
1729  jmu_B_lat_ = std::make_unique<DensityLattice>(l, n, origin, periodic,
1730  LatticeUpdate::AtOutput);
1731  }
1732  if (dens_type_lattice_printout_ == DensityType::BaryonicIsospin) {
1733  jmu_I3_lat_ = std::make_unique<DensityLattice>(l, n, origin, periodic,
1734  LatticeUpdate::AtOutput);
1735  }
1736  }
1737  if (dens_type_lattice_printout_ != DensityType::None &&
1738  dens_type_lattice_printout_ != DensityType::BaryonicIsospin &&
1739  dens_type_lattice_printout_ != DensityType::Baryon) {
1740  jmu_custom_lat_ = std::make_unique<DensityLattice>(
1741  l, n, origin, periodic, LatticeUpdate::AtOutput);
1742  }
1743  } else if (printout_lattice_td_ || printout_full_lattice_any_td_) {
1744  logg[LExperiment].error(
1745  "If you want Therm. VTK or Lattice output, configure a lattice for "
1746  "it.");
1747  } else if (potentials_ && potentials_->use_coulomb()) {
1748  logg[LExperiment].error(
1749  "Coulomb potential requires a lattice. Please add one to the "
1750  "configuration");
1751  }
1752 
1753  // Warning for the mean field calculation if lattice is not on.
1754  if ((potentials_ != nullptr) && (jmu_B_lat_ == nullptr)) {
1755  logg[LExperiment].warn() << "Lattice is NOT used. Mean-field energy is "
1756  << "not going to be calculated.";
1757  }
1758 
1759  // Store pointers to potential and lattice accessible for Action
1760  if (parameters_.potential_affect_threshold) {
1761  UB_lat_pointer = UB_lat_.get();
1762  UI3_lat_pointer = UI3_lat_.get();
1763  pot_pointer = potentials_.get();
1764  }
1765 
1766  // Throw fatal if DerivativesMode == FiniteDifference and lattice is not on.
1767  if ((parameters_.derivatives_mode == DerivativesMode::FiniteDifference) &&
1768  (jmu_B_lat_ == nullptr)) {
1769  throw std::runtime_error(
1770  "Lattice is necessary to calculate finite difference gradients.");
1771  }
1772 
1773  // Create forced thermalizer
1774  if (config.has_value({"Forced_Thermalization"})) {
1775  Configuration th_conf =
1776  config.extract_sub_configuration({"Forced_Thermalization"});
1777  thermalizer_ = modus_.create_grandcan_thermalizer(th_conf);
1778  }
1779 
1780  /* Take the seed setting only after the configuration was stored to a file
1781  * in smash.cc */
1782  seed_ = config.take({"General", "Randomseed"});
1783 }
1784 
1786 const std::string hline(113, '-');
1787 
1812 std::string format_measurements(const std::vector<Particles> &ensembles,
1813  uint64_t scatterings_this_interval,
1814  const QuantumNumbers &conserved_initial,
1815  SystemTimePoint time_start, double time,
1816  double E_mean_field,
1817  double E_mean_field_initial);
1832  const Potentials &potentials,
1834  RectangularLattice<std::pair<ThreeVector, ThreeVector>> *em_lattice,
1835  const ExperimentParameters &parameters);
1836 
1853 EventInfo fill_event_info(const std::vector<Particles> &ensembles,
1854  double E_mean_field, double modus_impact_parameter,
1855  const ExperimentParameters &parameters,
1856  bool projectile_target_interact,
1857  bool kinematic_cut_for_SMASH_IC);
1858 
1859 template <typename Modus>
1861  random::set_seed(seed_);
1862  logg[LExperiment].info() << "random number seed: " << seed_;
1863  /* Set seed for the next event. It has to be positive, so it can be entered
1864  * in the config.
1865  *
1866  * We have to be careful about the minimal integer, whose absolute value
1867  * cannot be represented. */
1868  int64_t r = random::advance();
1869  while (r == INT64_MIN) {
1870  r = random::advance();
1871  }
1872  seed_ = std::abs(r);
1873  /* Set the random seed used in PYTHIA hadronization
1874  * to be same with the SMASH one.
1875  * In this way we ensure that the results are reproducible
1876  * for every event if one knows SMASH random seed. */
1877  if (process_string_ptr_ != NULL) {
1878  process_string_ptr_->init_pythia_hadron_rndm();
1879  }
1880 
1881  for (Particles &particles : ensembles_) {
1882  particles.reset();
1883  }
1884 
1885  // Sample particles according to the initial conditions
1886  double start_time = -1.0;
1887 
1888  // Sample impact parameter only once per all ensembles
1889  // It should be the same for all ensembles
1890  if (modus_.is_collider()) {
1891  modus_.sample_impact();
1892  logg[LExperiment].info("Impact parameter = ", modus_.impact_parameter(),
1893  " fm");
1894  }
1895  for (Particles &particles : ensembles_) {
1896  start_time = modus_.initial_conditions(&particles, parameters_);
1897  }
1898  /* For box modus make sure that particles are in the box. In principle, after
1899  * a correct initialization they should be, so this is just playing it safe.
1900  */
1901  for (Particles &particles : ensembles_) {
1902  modus_.impose_boundary_conditions(&particles, outputs_);
1903  }
1904  // Reset the simulation clock
1905  double timestep = delta_time_startup_;
1906 
1907  switch (time_step_mode_) {
1908  case TimeStepMode::Fixed:
1909  break;
1910  case TimeStepMode::None:
1911  timestep = end_time_ - start_time;
1912  // Take care of the box modus + timestepless propagation
1913  const double max_dt = modus_.max_timestep(max_transverse_distance_sqr_);
1914  if (max_dt > 0. && max_dt < timestep) {
1915  timestep = max_dt;
1916  }
1917  break;
1918  }
1919  std::unique_ptr<UniformClock> clock_for_this_event;
1920  if (modus_.is_list() && (timestep < 0.0)) {
1921  throw std::runtime_error(
1922  "Timestep for the given event is negative. \n"
1923  "This might happen if the formation times of the input particles are "
1924  "larger than the specified end time of the simulation.");
1925  }
1926  clock_for_this_event =
1927  std::make_unique<UniformClock>(start_time, timestep, end_time_);
1928  parameters_.labclock = std::move(clock_for_this_event);
1929 
1930  // Reset the output clock
1931  parameters_.outputclock->reset(start_time, true);
1932  // remove time before starting time in case of custom output times.
1933  parameters_.outputclock->remove_times_in_past(start_time);
1934 
1935  logg[LExperiment].debug(
1936  "Lab clock: t_start = ", parameters_.labclock->current_time(),
1937  ", dt = ", parameters_.labclock->timestep_duration());
1938 
1939  /* Save the initial conserved quantum numbers and total momentum in
1940  * the system for conservation checks */
1941  conserved_initial_ = QuantumNumbers(ensembles_);
1942  wall_actions_total_ = 0;
1943  previous_wall_actions_total_ = 0;
1944  interactions_total_ = 0;
1945  previous_interactions_total_ = 0;
1946  discarded_interactions_total_ = 0;
1947  total_pauli_blocked_ = 0;
1948  projectile_target_interact_.assign(parameters_.n_ensembles, false);
1949  total_hypersurface_crossing_actions_ = 0;
1950  total_energy_removed_ = 0.0;
1951  total_energy_violated_by_Pythia_ = 0.0;
1952  // Print output headers
1953  logg[LExperiment].info() << hline;
1954  logg[LExperiment].info() << "Time[fm] Ekin[GeV] E_MF[GeV] ETotal[GeV] "
1955  << "ETot/N[GeV] D(ETot/N)[GeV] Scatt&Decays "
1956  << "Particles Comp.Time";
1957  logg[LExperiment].info() << hline;
1958  double E_mean_field = 0.0;
1959  if (potentials_) {
1960  // update_potentials();
1961  // if (parameters.outputclock->current_time() == 0.0 )
1962  // using the lattice is necessary
1963  if ((jmu_B_lat_ != nullptr)) {
1964  update_lattice(jmu_B_lat_.get(), old_jmu_auxiliary_.get(),
1965  new_jmu_auxiliary_.get(), four_gradient_auxiliary_.get(),
1966  LatticeUpdate::EveryTimestep, DensityType::Baryon,
1967  density_param_, ensembles_,
1968  parameters_.labclock->timestep_duration(), true);
1969  // Because there was no lattice at t=-Delta_t, the time derivatives
1970  // drho_dt and dj^mu/dt at t=0 are huge, while they shouldn't be; we
1971  // overwrite the time derivative to zero by hand.
1972  for (auto &node : *jmu_B_lat_) {
1973  node.overwrite_drho_dt_to_zero();
1974  node.overwrite_djmu_dt_to_zero();
1975  }
1976  E_mean_field = calculate_mean_field_energy(*potentials_, *jmu_B_lat_,
1977  EM_lat_.get(), parameters_);
1978  }
1979  }
1980  initial_mean_field_energy_ = E_mean_field;
1982  ensembles_, 0u, conserved_initial_, time_start_,
1983  parameters_.labclock->current_time(), E_mean_field,
1984  initial_mean_field_energy_);
1985 
1986  // Output at event start
1987  for (const auto &output : outputs_) {
1988  for (int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
1989  auto event_info = fill_event_info(
1990  ensembles_, E_mean_field, modus_.impact_parameter(), parameters_,
1991  projectile_target_interact_[i_ens], kinematic_cuts_for_IC_output_);
1992  output->at_eventstart(ensembles_[i_ens],
1993  // Pretend each ensemble is an independent event
1994  event_ * parameters_.n_ensembles + i_ens,
1995  event_info);
1996  }
1997  // For thermodynamic output
1998  output->at_eventstart(ensembles_, event_);
1999  // For thermodynamic lattice output
2000  if (printout_full_lattice_any_td_) {
2001  switch (dens_type_lattice_printout_) {
2002  case DensityType::Baryon:
2003  output->at_eventstart(event_, ThermodynamicQuantity::EckartDensity,
2004  DensityType::Baryon, *jmu_B_lat_);
2005  break;
2006  case DensityType::BaryonicIsospin:
2007  output->at_eventstart(event_, ThermodynamicQuantity::EckartDensity,
2008  DensityType::BaryonicIsospin, *jmu_I3_lat_);
2009  break;
2010  case DensityType::None:
2011  break;
2012  default:
2013  output->at_eventstart(event_, ThermodynamicQuantity::EckartDensity,
2014  DensityType::BaryonicIsospin, *jmu_custom_lat_);
2015  }
2016  if (printout_tmn_) {
2017  output->at_eventstart(event_, ThermodynamicQuantity::Tmn,
2018  dens_type_lattice_printout_, *Tmn_);
2019  }
2020  if (printout_tmn_landau_) {
2021  output->at_eventstart(event_, ThermodynamicQuantity::TmnLandau,
2022  dens_type_lattice_printout_, *Tmn_);
2023  }
2024  if (printout_v_landau_) {
2025  output->at_eventstart(event_, ThermodynamicQuantity::LandauVelocity,
2026  dens_type_lattice_printout_, *Tmn_);
2027  }
2028  if (printout_j_QBS_) {
2029  output->at_eventstart(event_, ThermodynamicQuantity::j_QBS,
2030  dens_type_lattice_printout_, *j_QBS_lat_);
2031  }
2032  }
2033  }
2034 
2035  /* In the ColliderModus, if Fermi motion is frozen, assign the beam momenta
2036  * to the nucleons in both the projectile and the target. Every ensemble
2037  * gets the same beam momenta, so no need to create beam_momenta_ vector
2038  * for every ensemble.
2039  */
2040  if (modus_.is_collider() && modus_.fermi_motion() == FermiMotion::Frozen) {
2041  for (ParticleData &particle : ensembles_[0]) {
2042  const double m = particle.effective_mass();
2043  double v_beam = 0.0;
2044  if (particle.belongs_to() == BelongsTo::Projectile) {
2045  v_beam = modus_.velocity_projectile();
2046  } else if (particle.belongs_to() == BelongsTo::Target) {
2047  v_beam = modus_.velocity_target();
2048  }
2049  const double gamma = 1.0 / std::sqrt(1.0 - v_beam * v_beam);
2050  beam_momentum_.emplace_back(
2051  FourVector(gamma * m, 0.0, 0.0, gamma * v_beam * m));
2052  } // loop over particles
2053  }
2054 }
2055 
2056 template <typename Modus>
2057 bool Experiment<Modus>::perform_action(Action &action, int i_ensemble,
2058  bool include_pauli_blocking) {
2059  Particles &particles = ensembles_[i_ensemble];
2060  // Make sure to skip invalid and Pauli-blocked actions.
2061  if (!action.is_valid(particles)) {
2062  discarded_interactions_total_++;
2063  logg[LExperiment].debug(~einhard::DRed(), "✘ ", action,
2064  " (discarded: invalid)");
2065  return false;
2066  }
2067  action.generate_final_state();
2068  logg[LExperiment].debug("Process Type is: ", action.get_type());
2069  if (include_pauli_blocking && pauli_blocker_ &&
2070  action.is_pauli_blocked(ensembles_, *pauli_blocker_)) {
2071  total_pauli_blocked_++;
2072  return false;
2073  }
2074 
2075  // Prepare projectile_target_interact_, it's used for output
2076  // to signal that there was some interaction in this event
2077  if (modus_.is_collider()) {
2078  int count_target = 0, count_projectile = 0;
2079  for (const auto &incoming : action.incoming_particles()) {
2080  if (incoming.belongs_to() == BelongsTo::Projectile) {
2081  count_projectile++;
2082  } else if (incoming.belongs_to() == BelongsTo::Target) {
2083  count_target++;
2084  }
2085  }
2086  if (count_target > 0 && count_projectile > 0) {
2087  projectile_target_interact_[i_ensemble] = true;
2088  }
2089  }
2090 
2091  /* Make sure to pick a non-zero integer, because 0 is reserved for "no
2092  * interaction yet". */
2093  const auto id_process = static_cast<uint32_t>(interactions_total_ + 1);
2094  // we perform the action and collect possible energy violations by Pythia
2095  total_energy_violated_by_Pythia_ += action.perform(&particles, id_process);
2096 
2097  interactions_total_++;
2098  if (action.get_type() == ProcessType::Wall) {
2099  wall_actions_total_++;
2100  }
2101  if (action.get_type() == ProcessType::HyperSurfaceCrossing) {
2102  total_hypersurface_crossing_actions_++;
2103  total_energy_removed_ += action.incoming_particles()[0].momentum().x0();
2104  }
2105  // Calculate Eckart rest frame density at the interaction point
2106  double rho = 0.0;
2107  if (dens_type_ != DensityType::None) {
2108  const FourVector r_interaction = action.get_interaction_point();
2109  constexpr bool compute_grad = false;
2110  const bool smearing = true;
2111  // todo(oliiny): it's a rough density estimate from a single ensemble.
2112  // It might actually be appropriate for output. Discuss.
2113  rho = std::get<0>(current_eckart(r_interaction.threevec(), particles,
2114  density_param_, dens_type_, compute_grad,
2115  smearing));
2116  }
2132  for (const auto &output : outputs_) {
2133  if (!output->is_dilepton_output() && !output->is_photon_output()) {
2134  if (output->is_IC_output() &&
2135  action.get_type() == ProcessType::HyperSurfaceCrossing) {
2136  output->at_interaction(action, rho);
2137  } else if (!output->is_IC_output()) {
2138  output->at_interaction(action, rho);
2139  }
2140  }
2141  }
2142 
2143  // At every collision photons can be produced.
2144  // Note: We rely here on the lazy evaluation of the arguments to if.
2145  // It may happen that in a wall-crossing-action sqrt_s raises an exception.
2146  // Therefore we first have to check if the incoming particles can undergo
2147  // an em-interaction.
2148  if (photons_switch_ &&
2149  ScatterActionPhoton::is_photon_reaction(action.incoming_particles()) &&
2150  ScatterActionPhoton::is_kinematically_possible(
2151  action.sqrt_s(), action.incoming_particles())) {
2152  /* Time in the action constructor is relative to
2153  * current time of incoming */
2154  constexpr double action_time = 0.;
2155  ScatterActionPhoton photon_act(action.incoming_particles(), action_time,
2156  n_fractional_photons_,
2157  action.get_total_weight());
2158 
2169  photon_act.add_dummy_hadronic_process(action.get_total_weight());
2170 
2171  // Now add the actual photon reaction channel.
2172  photon_act.add_single_process();
2173 
2174  photon_act.perform_photons(outputs_);
2175  }
2176 
2177  if (bremsstrahlung_switch_ &&
2178  BremsstrahlungAction::is_bremsstrahlung_reaction(
2179  action.incoming_particles())) {
2180  /* Time in the action constructor is relative to
2181  * current time of incoming */
2182  constexpr double action_time = 0.;
2183 
2184  BremsstrahlungAction brems_act(action.incoming_particles(), action_time,
2185  n_fractional_photons_,
2186  action.get_total_weight());
2187 
2199  brems_act.add_dummy_hadronic_process(action.get_total_weight());
2200 
2201  // Now add the actual bremsstrahlung reaction channel.
2202  brems_act.add_single_process();
2203 
2204  brems_act.perform_bremsstrahlung(outputs_);
2205  }
2206 
2207  logg[LExperiment].debug(~einhard::Green(), "✔ ", action);
2208  return true;
2209 }
2210 
2221 void validate_and_adjust_particle_list(ParticleList &particle_list);
2222 
2223 template <typename Modus>
2225  ParticleList &&add_plist,
2226  ParticleList &&remove_plist) {
2227  if (!add_plist.empty() || !remove_plist.empty()) {
2228  if (ensembles_.size() > 1) {
2229  throw std::runtime_error(
2230  "Adding or removing particles from SMASH is only possible when one "
2231  "ensemble is used.");
2232  }
2233  const double action_time = parameters_.labclock->current_time();
2234  if (!add_plist.empty()) {
2236  }
2237  if (!add_plist.empty()) {
2238  // Create and perform action to add particle(s)
2239  auto action_add_particles = std::make_unique<FreeforallAction>(
2240  ParticleList{}, add_plist, action_time);
2241  perform_action(*action_add_particles, 0);
2242  }
2243  if (!remove_plist.empty()) {
2244  validate_and_adjust_particle_list(remove_plist);
2245  const auto number_of_particles_to_be_removed = remove_plist.size();
2246  remove_plist.erase(
2247  std::remove_if(
2248  remove_plist.begin(), remove_plist.end(),
2249  [this, &action_time](const ParticleData &particle_to_remove) {
2250  return std::find_if(
2251  ensembles_[0].begin(), ensembles_[0].end(),
2252  [&particle_to_remove,
2253  &action_time](const ParticleData &p) {
2254  return are_particles_identical_at_given_time(
2255  particle_to_remove, p, action_time);
2256  }) == ensembles_[0].end();
2257  }),
2258  remove_plist.end());
2259  if (auto delta = number_of_particles_to_be_removed - remove_plist.size();
2260  delta > 0) {
2261  logg[LExperiment].warn(
2262  "When trying to remove particle(s) at the beginning ",
2263  "of the system evolution,\n", delta,
2264  " particle(s) could not be found and will be ignored.");
2265  }
2266  }
2267  if (!remove_plist.empty()) {
2268  // Create and perform action to remove particles
2269  auto action_remove_particles = std::make_unique<FreeforallAction>(
2270  remove_plist, ParticleList{}, action_time);
2271  perform_action(*action_remove_particles, 0);
2272  }
2273  }
2274 
2275  while (parameters_.labclock->current_time() < t_end) {
2276  const double dt = parameters_.labclock->timestep_duration();
2277  logg[LExperiment].debug("Timestepless propagation for next ", dt, " fm.");
2278 
2279  // Perform forced thermalization if required
2280  if (thermalizer_ &&
2281  thermalizer_->is_time_to_thermalize(parameters_.labclock)) {
2282  const bool ignore_cells_under_treshold = true;
2283  // Thermodynamics in thermalizer is computed from all ensembles,
2284  // but thermalization actions act on each ensemble independently
2285  thermalizer_->update_thermalizer_lattice(ensembles_, density_param_,
2286  ignore_cells_under_treshold);
2287  const double current_t = parameters_.labclock->current_time();
2288  for (int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
2289  thermalizer_->thermalize(ensembles_[i_ens], current_t,
2290  parameters_.testparticles);
2291  ThermalizationAction th_act(*thermalizer_, current_t);
2292  if (th_act.any_particles_thermalized()) {
2293  perform_action(th_act, i_ens);
2294  }
2295  }
2296  }
2297 
2298  std::vector<Actions> actions(parameters_.n_ensembles);
2299  for (int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
2300  actions[i_ens].clear();
2301  if (ensembles_[i_ens].size() > 0 && action_finders_.size() > 0) {
2302  /* (1.a) Create grid. */
2303  const double min_cell_length = compute_min_cell_length(dt);
2304  logg[LExperiment].debug("Creating grid with minimal cell length ",
2305  min_cell_length);
2306  /* For the hyper-surface-crossing actions also unformed particles are
2307  * searched and therefore needed on the grid. */
2308  const bool include_unformed_particles = IC_output_switch_;
2309  const auto &grid =
2310  use_grid_ ? modus_.create_grid(ensembles_[i_ens], min_cell_length,
2311  dt, parameters_.coll_crit,
2312  include_unformed_particles)
2313  : modus_.create_grid(ensembles_[i_ens], min_cell_length,
2314  dt, parameters_.coll_crit,
2315  include_unformed_particles,
2316  CellSizeStrategy::Largest);
2317 
2318  const double gcell_vol = grid.cell_volume();
2319  /* (1.b) Iterate over cells and find actions. */
2320  grid.iterate_cells(
2321  [&](const ParticleList &search_list) {
2322  for (const auto &finder : action_finders_) {
2323  actions[i_ens].insert(finder->find_actions_in_cell(
2324  search_list, dt, gcell_vol, beam_momentum_));
2325  }
2326  },
2327  [&](const ParticleList &search_list,
2328  const ParticleList &neighbors_list) {
2329  for (const auto &finder : action_finders_) {
2330  actions[i_ens].insert(finder->find_actions_with_neighbors(
2331  search_list, neighbors_list, dt, beam_momentum_));
2332  }
2333  });
2334  }
2335  }
2336 
2337  /* \todo (optimizations) Adapt timestep size here */
2338 
2339  /* (2) Propagate from action to action until next output or timestep end */
2340  const double end_timestep_time = parameters_.labclock->next_time();
2341  while (next_output_time() < end_timestep_time) {
2342  for (int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
2343  run_time_evolution_timestepless(actions[i_ens], i_ens,
2344  next_output_time());
2345  }
2346  ++(*parameters_.outputclock);
2347 
2348  intermediate_output();
2349  }
2350  for (int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
2351  run_time_evolution_timestepless(actions[i_ens], i_ens, end_timestep_time);
2352  }
2353 
2354  /* (3) Update potentials (if computed on the lattice) and
2355  * compute new momenta according to equations of motion */
2356  if (potentials_) {
2357  update_potentials();
2358  update_momenta(ensembles_, parameters_.labclock->timestep_duration(),
2359  *potentials_, FB_lat_.get(), FI3_lat_.get(),
2360  EM_lat_.get());
2361  }
2362 
2363  /* (4) Expand universe if non-minkowskian metric; updates
2364  * positions and momenta according to the selected expansion */
2365  if (metric_.mode_ != ExpansionMode::NoExpansion) {
2366  for (Particles &particles : ensembles_) {
2367  expand_space_time(&particles, parameters_, metric_);
2368  }
2369  }
2370 
2371  ++(*parameters_.labclock);
2372 
2373  /* (5) Check conservation laws.
2374  *
2375  * Check conservation of conserved quantities if potentials and string
2376  * fragmentation are off. If potentials are on then momentum is conserved
2377  * only in average. If string fragmentation is on, then energy and
2378  * momentum are only very roughly conserved in high-energy collisions. */
2379  if (!potentials_ && !parameters_.strings_switch &&
2380  metric_.mode_ == ExpansionMode::NoExpansion && !IC_output_switch_) {
2381  std::string err_msg = conserved_initial_.report_deviations(ensembles_);
2382  if (!err_msg.empty()) {
2383  logg[LExperiment].error() << err_msg;
2384  throw std::runtime_error("Violation of conserved quantities!");
2385  }
2386  }
2387  }
2388 
2389  if (pauli_blocker_) {
2390  logg[LExperiment].info(
2391  "Interactions: Pauli-blocked/performed = ", total_pauli_blocked_, "/",
2392  interactions_total_ - wall_actions_total_);
2393  }
2394 }
2395 
2396 template <typename Modus>
2398  Particles &particles) {
2399  const double dt =
2400  propagate_straight_line(&particles, to_time, beam_momentum_);
2401  if (dilepton_finder_ != nullptr) {
2402  for (const auto &output : outputs_) {
2403  dilepton_finder_->shine(particles, output.get(), dt);
2404  }
2405  }
2406 }
2407 
2415 inline void check_interactions_total(uint64_t interactions_total) {
2416  constexpr uint64_t max_uint32 = std::numeric_limits<uint32_t>::max();
2417  if (interactions_total >= max_uint32) {
2418  throw std::runtime_error("Integer overflow in total interaction number!");
2419  }
2420 }
2421 
2422 template <typename Modus>
2424  Actions &actions, int i_ensemble, const double end_time_propagation) {
2425  Particles &particles = ensembles_[i_ensemble];
2426  logg[LExperiment].debug(
2427  "Timestepless propagation: ", "Actions size = ", actions.size(),
2428  ", end time = ", end_time_propagation);
2429 
2430  // iterate over all actions
2431  while (!actions.is_empty()) {
2432  if (actions.earliest_time() > end_time_propagation) {
2433  break;
2434  }
2435  // get next action
2436  ActionPtr act = actions.pop();
2437  if (!act->is_valid(particles)) {
2438  discarded_interactions_total_++;
2439  logg[LExperiment].debug(~einhard::DRed(), "✘ ", act,
2440  " (discarded: invalid)");
2441  continue;
2442  }
2443  logg[LExperiment].debug(~einhard::Green(), "✔ ", act,
2444  ", action time = ", act->time_of_execution());
2445 
2446  /* (1) Propagate to the next action. */
2447  propagate_and_shine(act->time_of_execution(), particles);
2448 
2449  /* (2) Perform action.
2450  *
2451  * Update the positions of the incoming particles, because the information
2452  * in the action object will be outdated as the particles have been
2453  * propagated since the construction of the action. */
2454  act->update_incoming(particles);
2455  const bool performed = perform_action(*act, i_ensemble);
2456 
2457  /* No need to update actions for outgoing particles
2458  * if the action is not performed. */
2459  if (!performed) {
2460  continue;
2461  }
2462 
2463  /* (3) Update actions for newly-produced particles. */
2464 
2465  const double end_time_timestep = parameters_.labclock->next_time();
2466  // New actions are always search until the end of the current timestep
2467  const double time_left = end_time_timestep - act->time_of_execution();
2468  const ParticleList &outgoing_particles = act->outgoing_particles();
2469  // Grid cell volume set to zero, since there is no grid
2470  const double gcell_vol = 0.0;
2471  for (const auto &finder : action_finders_) {
2472  // Outgoing particles can still decay, cross walls...
2473  actions.insert(finder->find_actions_in_cell(outgoing_particles, time_left,
2474  gcell_vol, beam_momentum_));
2475  // ... and collide with other particles.
2476  actions.insert(finder->find_actions_with_surrounding_particles(
2477  outgoing_particles, particles, time_left, beam_momentum_));
2478  }
2479 
2480  check_interactions_total(interactions_total_);
2481  }
2482 
2483  propagate_and_shine(end_time_propagation, particles);
2484 }
2485 
2486 template <typename Modus>
2488  const uint64_t wall_actions_this_interval =
2489  wall_actions_total_ - previous_wall_actions_total_;
2490  previous_wall_actions_total_ = wall_actions_total_;
2491  const uint64_t interactions_this_interval = interactions_total_ -
2492  previous_interactions_total_ -
2493  wall_actions_this_interval;
2494  previous_interactions_total_ = interactions_total_;
2495  double E_mean_field = 0.0;
2498  double computational_frame_time = 0.0;
2499  if (potentials_) {
2500  // using the lattice is necessary
2501  if ((jmu_B_lat_ != nullptr)) {
2502  E_mean_field = calculate_mean_field_energy(*potentials_, *jmu_B_lat_,
2503  EM_lat_.get(), parameters_);
2504  /*
2505  * Mean field calculated in a box should remain approximately constant if
2506  * the system is in equilibrium, and so deviations from its original value
2507  * may signal a phase transition or other dynamical process. This
2508  * comparison only makes sense in the Box Modus, hence the condition.
2509  */
2510  if (modus_.is_box()) {
2511  double tmp = (E_mean_field - initial_mean_field_energy_) /
2512  (E_mean_field + initial_mean_field_energy_);
2513  /*
2514  * This is displayed when the system evolves away from its initial
2515  * configuration (which is when the total mean field energy in the box
2516  * deviates from its initial value).
2517  */
2518  if (std::abs(tmp) > 0.01) {
2519  logg[LExperiment].info()
2520  << "\n\n\n\t The mean field at t = "
2521  << parameters_.outputclock->current_time()
2522  << " [fm] differs from the mean field at t = 0:"
2523  << "\n\t\t initial_mean_field_energy_ = "
2524  << initial_mean_field_energy_ << " [GeV]"
2525  << "\n\t\t abs[(E_MF - E_MF(t=0))/(E_MF + E_MF(t=0))] = "
2526  << std::abs(tmp)
2527  << "\n\t\t E_MF/E_MF(t=0) = "
2528  << E_mean_field / initial_mean_field_energy_ << "\n\n";
2529  }
2530  }
2531  }
2532  }
2533 
2535  ensembles_, interactions_this_interval, conserved_initial_, time_start_,
2536  parameters_.outputclock->current_time(), E_mean_field,
2537  initial_mean_field_energy_);
2538  const LatticeUpdate lat_upd = LatticeUpdate::AtOutput;
2539 
2540  // save evolution data
2541  if (!(modus_.is_box() && parameters_.outputclock->current_time() <
2542  modus_.equilibration_time())) {
2543  for (const auto &output : outputs_) {
2544  if (output->is_dilepton_output() || output->is_photon_output() ||
2545  output->is_IC_output()) {
2546  continue;
2547  }
2548  for (int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
2549  auto event_info = fill_event_info(
2550  ensembles_, E_mean_field, modus_.impact_parameter(), parameters_,
2551  projectile_target_interact_[i_ens], kinematic_cuts_for_IC_output_);
2552 
2553  output->at_intermediate_time(ensembles_[i_ens], parameters_.outputclock,
2554  density_param_, event_info);
2555  computational_frame_time = event_info.current_time;
2556  }
2557  // For thermodynamic output
2558  output->at_intermediate_time(ensembles_, parameters_.outputclock,
2559  density_param_);
2560 
2561  // Thermodynamic output on the lattice versus time
2562  switch (dens_type_lattice_printout_) {
2563  case DensityType::Baryon:
2564  update_lattice(jmu_B_lat_.get(), lat_upd, DensityType::Baryon,
2565  density_param_, ensembles_, false);
2566  output->thermodynamics_output(ThermodynamicQuantity::EckartDensity,
2567  DensityType::Baryon, *jmu_B_lat_);
2568  output->thermodynamics_lattice_output(*jmu_B_lat_,
2569  computational_frame_time);
2570  break;
2571  case DensityType::BaryonicIsospin:
2572  update_lattice(jmu_I3_lat_.get(), lat_upd,
2573  DensityType::BaryonicIsospin, density_param_,
2574  ensembles_, false);
2575  output->thermodynamics_output(ThermodynamicQuantity::EckartDensity,
2576  DensityType::BaryonicIsospin,
2577  *jmu_I3_lat_);
2578  output->thermodynamics_lattice_output(*jmu_I3_lat_,
2579  computational_frame_time);
2580  break;
2581  case DensityType::None:
2582  break;
2583  default:
2584  update_lattice(jmu_custom_lat_.get(), lat_upd,
2585  dens_type_lattice_printout_, density_param_,
2586  ensembles_, false);
2587  output->thermodynamics_output(ThermodynamicQuantity::EckartDensity,
2588  dens_type_lattice_printout_,
2589  *jmu_custom_lat_);
2590  output->thermodynamics_lattice_output(*jmu_custom_lat_,
2591  computational_frame_time);
2592  }
2593  if (printout_tmn_ || printout_tmn_landau_ || printout_v_landau_) {
2594  update_lattice(Tmn_.get(), lat_upd, dens_type_lattice_printout_,
2595  density_param_, ensembles_, false);
2596  if (printout_tmn_) {
2597  output->thermodynamics_output(ThermodynamicQuantity::Tmn,
2598  dens_type_lattice_printout_, *Tmn_);
2599  output->thermodynamics_lattice_output(
2600  ThermodynamicQuantity::Tmn, *Tmn_, computational_frame_time);
2601  }
2602  if (printout_tmn_landau_) {
2603  output->thermodynamics_output(ThermodynamicQuantity::TmnLandau,
2604  dens_type_lattice_printout_, *Tmn_);
2605  output->thermodynamics_lattice_output(
2607  computational_frame_time);
2608  }
2609  if (printout_v_landau_) {
2610  output->thermodynamics_output(ThermodynamicQuantity::LandauVelocity,
2611  dens_type_lattice_printout_, *Tmn_);
2612  output->thermodynamics_lattice_output(
2614  computational_frame_time);
2615  }
2616  }
2617  if (EM_lat_) {
2618  output->fields_output("Efield", "Bfield", *EM_lat_);
2619  }
2620  if (printout_j_QBS_) {
2621  output->thermodynamics_lattice_output(
2622  *j_QBS_lat_, computational_frame_time, ensembles_, density_param_);
2623  }
2624 
2625  if (thermalizer_) {
2626  output->thermodynamics_output(*thermalizer_);
2627  }
2628  }
2629  }
2630 }
2631 
2632 template <typename Modus>
2634  if (potentials_) {
2635  if (potentials_->use_symmetry() && jmu_I3_lat_ != nullptr) {
2636  update_lattice(jmu_I3_lat_.get(), old_jmu_auxiliary_.get(),
2637  new_jmu_auxiliary_.get(), four_gradient_auxiliary_.get(),
2638  LatticeUpdate::EveryTimestep, DensityType::BaryonicIsospin,
2639  density_param_, ensembles_,
2640  parameters_.labclock->timestep_duration(), true);
2641  }
2642  if ((potentials_->use_skyrme() || potentials_->use_symmetry()) &&
2643  jmu_B_lat_ != nullptr) {
2644  update_lattice(jmu_B_lat_.get(), old_jmu_auxiliary_.get(),
2645  new_jmu_auxiliary_.get(), four_gradient_auxiliary_.get(),
2646  LatticeUpdate::EveryTimestep, DensityType::Baryon,
2647  density_param_, ensembles_,
2648  parameters_.labclock->timestep_duration(), true);
2649  const size_t UBlattice_size = UB_lat_->size();
2650  for (size_t i = 0; i < UBlattice_size; i++) {
2651  auto jB = (*jmu_B_lat_)[i];
2652  const FourVector flow_four_velocity_B =
2653  std::abs(jB.rho()) > very_small_double ? jB.jmu_net() / jB.rho()
2654  : FourVector();
2655  double baryon_density = jB.rho();
2656  ThreeVector baryon_grad_j0 = jB.grad_j0();
2657  ThreeVector baryon_dvecj_dt = jB.dvecj_dt();
2658  ThreeVector baryon_curl_vecj = jB.curl_vecj();
2659  if (potentials_->use_skyrme()) {
2660  (*UB_lat_)[i] =
2661  flow_four_velocity_B * potentials_->skyrme_pot(baryon_density);
2662  (*FB_lat_)[i] =
2663  potentials_->skyrme_force(baryon_density, baryon_grad_j0,
2664  baryon_dvecj_dt, baryon_curl_vecj);
2665  }
2666  if (potentials_->use_symmetry() && jmu_I3_lat_ != nullptr) {
2667  auto jI3 = (*jmu_I3_lat_)[i];
2668  const FourVector flow_four_velocity_I3 =
2669  std::abs(jI3.rho()) > very_small_double
2670  ? jI3.jmu_net() / jI3.rho()
2671  : FourVector();
2672  (*UI3_lat_)[i] = flow_four_velocity_I3 *
2673  potentials_->symmetry_pot(jI3.rho(), baryon_density);
2674  (*FI3_lat_)[i] = potentials_->symmetry_force(
2675  jI3.rho(), jI3.grad_j0(), jI3.dvecj_dt(), jI3.curl_vecj(),
2676  baryon_density, baryon_grad_j0, baryon_dvecj_dt,
2677  baryon_curl_vecj);
2678  }
2679  }
2680  }
2681  if (potentials_->use_coulomb()) {
2682  update_lattice(jmu_el_lat_.get(), LatticeUpdate::EveryTimestep,
2683  DensityType::Charge, density_param_, ensembles_, true);
2684  for (size_t i = 0; i < EM_lat_->size(); i++) {
2685  ThreeVector electric_field = {0., 0., 0.};
2686  ThreeVector position = jmu_el_lat_->cell_center(i);
2687  jmu_el_lat_->integrate_volume(electric_field,
2688  Potentials::E_field_integrand,
2689  potentials_->coulomb_r_cut(), position);
2690  ThreeVector magnetic_field = {0., 0., 0.};
2691  jmu_el_lat_->integrate_volume(magnetic_field,
2692  Potentials::B_field_integrand,
2693  potentials_->coulomb_r_cut(), position);
2694  (*EM_lat_)[i] = std::make_pair(electric_field, magnetic_field);
2695  }
2696  } // if ((potentials_->use_skyrme() || ...
2697  if (potentials_->use_vdf() && jmu_B_lat_ != nullptr) {
2698  update_lattice(jmu_B_lat_.get(), old_jmu_auxiliary_.get(),
2699  new_jmu_auxiliary_.get(), four_gradient_auxiliary_.get(),
2700  LatticeUpdate::EveryTimestep, DensityType::Baryon,
2701  density_param_, ensembles_,
2702  parameters_.labclock->timestep_duration(), true);
2703  if (parameters_.field_derivatives_mode == FieldDerivativesMode::Direct) {
2705  fields_lat_.get(), old_fields_auxiliary_.get(),
2706  new_fields_auxiliary_.get(), fields_four_gradient_auxiliary_.get(),
2707  jmu_B_lat_.get(), LatticeUpdate::EveryTimestep, *potentials_,
2708  parameters_.labclock->timestep_duration());
2709  }
2710  const size_t UBlattice_size = UB_lat_->size();
2711  for (size_t i = 0; i < UBlattice_size; i++) {
2712  auto jB = (*jmu_B_lat_)[i];
2713  (*UB_lat_)[i] = potentials_->vdf_pot(jB.rho(), jB.jmu_net());
2714  switch (parameters_.field_derivatives_mode) {
2716  (*FB_lat_)[i] = potentials_->vdf_force(
2717  jB.rho(), jB.drho_dxnu().x0(), jB.drho_dxnu().threevec(),
2718  jB.grad_rho_cross_vecj(), jB.jmu_net().x0(), jB.grad_j0(),
2719  jB.jmu_net().threevec(), jB.dvecj_dt(), jB.curl_vecj());
2720  break;
2722  auto Amu = (*fields_lat_)[i];
2723  (*FB_lat_)[i] = potentials_->vdf_force(
2724  Amu.grad_A0(), Amu.dvecA_dt(), Amu.curl_vecA());
2725  break;
2726  }
2727  } // for (size_t i = 0; i < UBlattice_size; i++)
2728  } // if potentials_->use_vdf()
2729  }
2730 }
2731 
2732 template <typename Modus>
2734  /* At end of time evolution: Force all resonances to decay. In order to handle
2735  * decay chains, we need to loop until no further actions occur. */
2736  bool actions_performed, decays_found;
2737  uint64_t interactions_old;
2738  do {
2739  decays_found = false;
2740  interactions_old = interactions_total_;
2741  for (int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
2742  Actions actions;
2743 
2744  // Dileptons: shining of remaining resonances
2745  if (dilepton_finder_ != nullptr) {
2746  for (const auto &output : outputs_) {
2747  dilepton_finder_->shine_final(ensembles_[i_ens], output.get(), true);
2748  }
2749  }
2750  // Find actions.
2751  for (const auto &finder : action_finders_) {
2752  auto found_actions = finder->find_final_actions(ensembles_[i_ens]);
2753  if (!found_actions.empty()) {
2754  actions.insert(std::move(found_actions));
2755  decays_found = true;
2756  }
2757  }
2758  // Perform actions.
2759  while (!actions.is_empty()) {
2760  perform_action(*actions.pop(), i_ens, false);
2761  }
2762  }
2763  actions_performed = interactions_total_ > interactions_old;
2764  // Throw an error if actions were found but not performed
2765  if (decays_found && !actions_performed) {
2766  throw std::runtime_error("Final decays were found but not performed.");
2767  }
2768  // loop until no more decays occur
2769  } while (actions_performed);
2770 
2771  // Dileptons: shining of stable particles at the end
2772  if (dilepton_finder_ != nullptr) {
2773  for (const auto &output : outputs_) {
2774  for (Particles &particles : ensembles_) {
2775  dilepton_finder_->shine_final(particles, output.get(), false);
2776  }
2777  }
2778  }
2779 }
2780 
2781 template <typename Modus>
2783  /* make sure the experiment actually ran (note: we should compare this
2784  * to the start time, but we don't know that. Therefore, we check that
2785  * the time is positive, which should heuristically be the same). */
2786  double E_mean_field = 0.0;
2787  if (likely(parameters_.labclock > 0)) {
2788  const uint64_t wall_actions_this_interval =
2789  wall_actions_total_ - previous_wall_actions_total_;
2790  const uint64_t interactions_this_interval = interactions_total_ -
2791  previous_interactions_total_ -
2792  wall_actions_this_interval;
2793  if (potentials_) {
2794  // using the lattice is necessary
2795  if ((jmu_B_lat_ != nullptr)) {
2796  E_mean_field = calculate_mean_field_energy(*potentials_, *jmu_B_lat_,
2797  EM_lat_.get(), parameters_);
2798  }
2799  }
2800  if (std::abs(parameters_.labclock->current_time() - end_time_) >
2801  really_small) {
2802  logg[LExperiment].warn()
2803  << "SMASH not propagated until configured end time. Current time = "
2804  << parameters_.labclock->current_time()
2805  << "fm. End time = " << end_time_ << "fm.";
2806  } else {
2808  ensembles_, interactions_this_interval, conserved_initial_,
2809  time_start_, end_time_, E_mean_field, initial_mean_field_energy_);
2810  }
2811  int total_particles = 0;
2812  for (const Particles &particles : ensembles_) {
2813  total_particles += particles.size();
2814  }
2815  if (IC_output_switch_ && (total_particles == 0)) {
2816  const double initial_system_energy_plus_Pythia_violations =
2817  conserved_initial_.momentum().x0() + total_energy_violated_by_Pythia_;
2818  const double fraction_of_total_system_energy_removed =
2819  initial_system_energy_plus_Pythia_violations / total_energy_removed_;
2820  // Verify there is no more energy in the system if all particles were
2821  // removed when crossing the hypersurface
2822  if (std::fabs(fraction_of_total_system_energy_removed - 1.) >
2823  really_small) {
2824  throw std::runtime_error(
2825  "There is remaining energy in the system although all particles "
2826  "were removed.\n"
2827  "E_remain = " +
2828  std::to_string((initial_system_energy_plus_Pythia_violations -
2829  total_energy_removed_)) +
2830  " [GeV]");
2831  } else {
2832  logg[LExperiment].info() << hline;
2833  logg[LExperiment].info()
2834  << "Time real: " << SystemClock::now() - time_start_;
2835  logg[LExperiment].info()
2836  << "Interactions before reaching hypersurface: "
2837  << interactions_total_ - wall_actions_total_ -
2838  total_hypersurface_crossing_actions_;
2839  logg[LExperiment].info()
2840  << "Total number of particles removed on hypersurface: "
2841  << total_hypersurface_crossing_actions_;
2842  }
2843  } else {
2844  const double precent_discarded =
2845  interactions_total_ > 0
2846  ? static_cast<double>(discarded_interactions_total_) * 100.0 /
2847  interactions_total_
2848  : 0.0;
2849  std::stringstream msg_discarded;
2850  msg_discarded
2851  << "Discarded interaction number: " << discarded_interactions_total_
2852  << " (" << precent_discarded
2853  << "% of the total interaction number including wall crossings)";
2854 
2855  logg[LExperiment].info() << hline;
2856  logg[LExperiment].info()
2857  << "Time real: " << SystemClock::now() - time_start_;
2858  logg[LExperiment].debug() << msg_discarded.str();
2859 
2860  if (parameters_.coll_crit == CollisionCriterion::Stochastic &&
2861  precent_discarded > 1.0) {
2862  // The choosen threshold of 1% is a heuristical value
2863  logg[LExperiment].warn()
2864  << msg_discarded.str()
2865  << "\nThe number of discarded interactions is large, which means "
2866  "the assumption for the stochastic criterion of\n1 interaction "
2867  "per particle per timestep is probably violated. Consider "
2868  "reducing the timestep size.";
2869  }
2870 
2871  logg[LExperiment].info() << "Final interaction number: "
2872  << interactions_total_ - wall_actions_total_;
2873  }
2874 
2875  // Check if there are unformed particles
2876  int unformed_particles_count = 0;
2877  for (const Particles &particles : ensembles_) {
2878  for (const ParticleData &particle : particles) {
2879  if (particle.formation_time() > end_time_) {
2880  unformed_particles_count++;
2881  }
2882  }
2883  }
2884  if (unformed_particles_count > 0) {
2885  logg[LExperiment].warn(
2886  "End time might be too small. ", unformed_particles_count,
2887  " unformed particles were found at the end of the evolution.");
2888  }
2889  }
2890 
2891  // Keep track of how many ensembles had interactions
2892  count_nonempty_ensembles();
2893 
2894  for (const auto &output : outputs_) {
2895  for (int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
2896  auto event_info = fill_event_info(
2897  ensembles_, E_mean_field, modus_.impact_parameter(), parameters_,
2898  projectile_target_interact_[i_ens], kinematic_cuts_for_IC_output_);
2899  output->at_eventend(ensembles_[i_ens],
2900  // Pretend each ensemble is an independent event
2901  event_ * parameters_.n_ensembles + i_ens, event_info);
2902  }
2903  // For thermodynamic output
2904  output->at_eventend(ensembles_, event_);
2905 
2906  // For thermodynamic lattice output
2907  if (dens_type_lattice_printout_ != DensityType::None) {
2908  output->at_eventend(event_, ThermodynamicQuantity::EckartDensity,
2909  dens_type_lattice_printout_);
2910  output->at_eventend(ThermodynamicQuantity::EckartDensity);
2911  }
2912  if (printout_tmn_) {
2913  output->at_eventend(event_, ThermodynamicQuantity::Tmn,
2914  dens_type_lattice_printout_);
2915  output->at_eventend(ThermodynamicQuantity::Tmn);
2916  }
2917  if (printout_tmn_landau_) {
2918  output->at_eventend(event_, ThermodynamicQuantity::TmnLandau,
2919  dens_type_lattice_printout_);
2920  output->at_eventend(ThermodynamicQuantity::TmnLandau);
2921  }
2922  if (printout_v_landau_) {
2923  output->at_eventend(event_, ThermodynamicQuantity::LandauVelocity,
2924  dens_type_lattice_printout_);
2925  output->at_eventend(ThermodynamicQuantity::LandauVelocity);
2926  }
2927  if (printout_j_QBS_) {
2928  output->at_eventend(ThermodynamicQuantity::j_QBS);
2929  }
2930  }
2931 }
2932 
2933 template <typename Modus>
2935  for (bool has_interaction : projectile_target_interact_) {
2936  if (has_interaction) {
2937  nonempty_ensembles_++;
2938  }
2939  }
2940 }
2941 
2942 template <typename Modus>
2944  if (event_counting_ == EventCounting::FixedNumber) {
2945  return event_ >= nevents_;
2946  }
2947  if (event_counting_ == EventCounting::MinimumNonEmpty) {
2948  if (nonempty_ensembles_ >= minimum_nonempty_ensembles_) {
2949  return true;
2950  }
2951  if (event_ >= max_events_) {
2952  logg[LExperiment].warn()
2953  << "Maximum number of events (" << max_events_
2954  << ") exceeded. Stopping calculation. "
2955  << "The fraction of empty ensembles is "
2956  << (1.0 - static_cast<double>(nonempty_ensembles_) /
2957  (event_ * parameters_.n_ensembles))
2958  << ". If this fraction is expected, try increasing the "
2959  "Maximum_Ensembles_Run.";
2960  return true;
2961  }
2962  return false;
2963  }
2964  throw std::runtime_error("Event counting option is invalid");
2965  return false;
2966 }
2967 
2968 template <typename Modus>
2970  event_++;
2971 }
2972 
2973 template <typename Modus>
2975  const auto &mainlog = logg[LMain];
2976  for (event_ = 0; !is_finished(); event_++) {
2977  mainlog.info() << "Event " << event_;
2978 
2979  // Sample initial particles, start clock, some printout and book-keeping
2980  initialize_new_event();
2981 
2982  run_time_evolution(end_time_);
2983 
2984  if (force_decays_) {
2985  do_final_decays();
2986  }
2987 
2988  // Output at event end
2989  final_output();
2990  }
2991 }
2992 
2993 } // namespace smash
2994 
2995 #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:143
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.
bool has_value(std::initializer_list< const char * > keys) const
Return whether there is a non-empty value behind the requested keys.
void set_value(std::initializer_list< const char * > keys, T &&value)
Overwrite the value of the specified YAML node.
Configuration extract_sub_configuration(std::initializer_list< const char * > keys, Configuration::GetEmpty empty_if_not_existing=Configuration::GetEmpty::No)
Create a new configuration from a then-removed section of the present object.
Value take(std::initializer_list< const char * > keys)
The default interface for SMASH to read configuration values.
A class to pre-calculate and store parameters relevant for density calculation.
Definition: density.h:108
Non-template interface to Experiment<Modus>.
Definition: experiment.h:98
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:187
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:2397
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:425
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:641
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:724
std::vector< std::unique_ptr< ActionFinderInterface > > action_finders_
The Action finder objects.
Definition: experiment.h:419
bool printout_tmn_
Whether to print the energy-momentum tensor.
Definition: experiment.h:504
QuantumNumbers conserved_initial_
The conserved quantities of the system.
Definition: experiment.h:635
DensityParameters density_param_
Structure to precalculate and hold parameters for density computations.
Definition: experiment.h:370
double next_output_time() const
Shortcut for next output time.
Definition: experiment.h:345
bool printout_full_lattice_ascii_td_
Whether to print the thermodynamics quantities evaluated on the lattices, point by point,...
Definition: experiment.h:520
double total_energy_removed_
Total energy removed from the system in hypersurface crossing actions.
Definition: experiment.h:694
DensityType dens_type_lattice_printout_
Type of density for lattice printout.
Definition: experiment.h:455
double max_transverse_distance_sqr_
Maximal distance at which particles can interact in case of the geometric criterion,...
Definition: experiment.h:626
bool printout_j_QBS_
Whether to print the Q, B, S 4-currents.
Definition: experiment.h:513
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:2934
const TimeStepMode time_step_mode_
This indicates whether to use time steps.
Definition: experiment.h:620
std::unique_ptr< DensityLattice > jmu_custom_lat_
Custom density on the lattices.
Definition: experiment.h:452
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:2969
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:2423
Particles * first_ensemble()
Provides external access to SMASH particles.
Definition: experiment.h:256
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:2224
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:434
DensityType dens_type_
Type of density to be written to collision headers.
Definition: experiment.h:647
const int nevents_
Number of events.
Definition: experiment.h:553
void intermediate_output()
Intermediate output during an event.
Definition: experiment.h:2487
std::vector< FourVector > beam_momentum_
The initial nucleons in the ColliderModus propagate with beam_momentum_, if Fermi motion is frozen.
Definition: experiment.h:416
std::unique_ptr< RectangularLattice< std::pair< ThreeVector, ThreeVector > > > FI3_lat_
Lattices for the electric and magnetic component of the symmetry force.
Definition: experiment.h:478
std::unique_ptr< RectangularLattice< FourVector > > new_jmu_auxiliary_
Auxiliary lattice for values of jmu at a time step t0 + dt.
Definition: experiment.h:490
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:337
void initialize_new_event()
This is called in the beginning of each event.
Definition: experiment.h:1860
bool is_finished()
Checks wether the desired number events have been calculated.
Definition: experiment.h:2943
int n_fractional_photons_
Number of fractional photons produced per single reaction.
Definition: experiment.h:428
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:485
std::vector< Particles > ensembles_
Complete particle list, all ensembles in one vector.
Definition: experiment.h:379
std::unique_ptr< RectangularLattice< FourVector > > new_fields_auxiliary_
Auxiliary lattice for values of Amu at a time step t0 + dt.
Definition: experiment.h:498
SystemTimePoint time_start_
system starting time of the simulation
Definition: experiment.h:644
uint64_t previous_wall_actions_total_
Total number of wall-crossings for previous timestep.
Definition: experiment.h:671
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:474
OutputsList outputs_
A list of output formaters.
Definition: experiment.h:397
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:501
void do_final_decays()
Performs the final decays of an event.
Definition: experiment.h:2733
std::unique_ptr< PauliBlocker > pauli_blocker_
An instance of PauliBlocker class that stores parameters needed for Pauli blocking calculations and c...
Definition: experiment.h:391
bool printout_v_landau_
Whether to print the 4-velocity in Landau frame.
Definition: experiment.h:510
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:467
bool printout_lattice_td_
Whether to print the thermodynamics quantities evaluated on the lattices.
Definition: experiment.h:516
std::unique_ptr< RectangularLattice< std::pair< ThreeVector, ThreeVector > > > EM_lat_
Lattices for electric and magnetic field in fm^-2.
Definition: experiment.h:482
std::unique_ptr< FieldsLattice > fields_lat_
Mean-field A^mu on the lattice.
Definition: experiment.h:443
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:403
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:376
std::unique_ptr< RectangularLattice< FourVector > > old_jmu_auxiliary_
Auxiliary lattice for values of jmu at a time step t0.
Definition: experiment.h:488
std::unique_ptr< DecayActionsFinderDilepton > dilepton_finder_
The Dilepton Action Finder.
Definition: experiment.h:422
std::unique_ptr< DensityLattice > j_QBS_lat_
4-current for j_QBS lattice output
Definition: experiment.h:431
bool printout_tmn_landau_
Whether to print the energy-momentum tensor in Landau frame.
Definition: experiment.h:507
bool printout_full_lattice_binary_td_
Whether to print the thermodynamics quantities evaluated on the lattices, point by point,...
Definition: experiment.h:524
std::unique_ptr< RectangularLattice< std::array< FourVector, 4 > > > four_gradient_auxiliary_
Auxiliary lattice for calculating the four-gradient of jmu.
Definition: experiment.h:493
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:461
ExperimentParameters parameters_
Struct of several member variables.
Definition: experiment.h:367
std::unique_ptr< RectangularLattice< FourVector > > old_fields_auxiliary_
Auxiliary lattice for values of Amu at a time step t0.
Definition: experiment.h:496
std::unique_ptr< DensityLattice > jmu_el_lat_
Electric charge density on the lattice.
Definition: experiment.h:440
std::unique_ptr< Potentials > potentials_
An instance of potentials class, that stores parameters of potentials, calculates them and their grad...
Definition: experiment.h:385
uint64_t wall_actions_total_
Total number of wall-crossings for current timestep.
Definition: experiment.h:665
uint64_t total_hypersurface_crossing_actions_
Total number of particles removed from the evolution in hypersurface crossing actions.
Definition: experiment.h:683
uint64_t interactions_total_
Total number of interactions for current timestep.
Definition: experiment.h:653
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:437
uint64_t total_pauli_blocked_
Total number of Pauli-blockings for current timestep.
Definition: experiment.h:677
void final_output()
Output at the end of an event.
Definition: experiment.h:2782
std::vector< bool > projectile_target_interact_
Whether the projectile and the target collided.
Definition: experiment.h:409
Modus * modus()
Provides external access to SMASH calculation modus.
Definition: experiment.h:264
void update_potentials()
Recompute potentials on lattices if necessary.
Definition: experiment.h:2633
const bool IC_output_switch_
This indicates whether the IC output is enabled.
Definition: experiment.h:617
OutputPtr dilepton_output_
The Dilepton output.
Definition: experiment.h:400
int64_t seed_
random seed for the next event.
Definition: experiment.h:705
std::vector< Particles > * all_ensembles()
Getter for all ensembles.
Definition: experiment.h:258
uint64_t previous_interactions_total_
Total number of interactions for previous timestep.
Definition: experiment.h:659
uint64_t discarded_interactions_total_
Total number of discarded interactions, because they were invalidated before they could be performed.
Definition: experiment.h:689
void run() override
Runs the experiment.
Definition: experiment.h:2974
bool kinematic_cuts_for_IC_output_
This indicates whether kinematic cuts are enabled for the IC output.
Definition: experiment.h:702
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:699
The FourVector class holds relevant values in Minkowski spacetime with (+, −, −, −) metric signature.
Definition: fourvector.h:33
ParticleData contains the dynamic information of a certain particle.
Definition: particledata.h:58
static double formation_power_
Power with which the cross section scaling factor grows in time.
Definition: particledata.h:389
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:32
A container for storing conserved values.
A container class to hold all the arrays on the lattice and access them.
Definition: lattice.h:47
String excitation processes used in SMASH.
Definition: stringprocess.h:45
ThermalizationAction implements forced thermalization as an Action class.
The ThreeVector class represents a physical three-vector with the components .
Definition: threevector.h:31
@ Frozen
Use fermi motion without potentials.
TimeStepMode
The time step mode.
@ Fixed
Use fixed time step.
@ None
Don't use time steps; propagate from action to action.
@ Stochastic
Stochastic Criteiron.
EventCounting
Defines how the number of events is determined.
#define SMASH_SOURCE_LOCATION
Hackery that is required to output the location in the source code where the log statement occurs.
Definition: logging.h:150
std::ostream & operator<<(std::ostream &out, const ActionPtr &action)
Convenience: dereferences the ActionPtr to Action.
Definition: action.h:537
std::array< einhard::Logger<>, std::tuple_size< LogArea::AreaTuple >::value > logg
An array that stores all pre-configured Logger objects.
Definition: logging.cc:39
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:214
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:810
#define likely(x)
Tell the branch predictor that this expression is likely true.
Definition: macros.h:14
constexpr int n
Neutron.
Engine::result_type advance()
Advance the engine's state and return the generated value.
Definition: random.h:78
void set_seed(T &&seed)
Sets the seed of the random number engine.
Definition: random.h:71
Definition: action.h:24
void 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)
Updates the momenta of all particles at the current time step according to the equations of motion:
Definition: propagation.cc:111
static constexpr int LInitialConditions
Definition: experiment.h:89
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:604
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:336
void expand_space_time(Particles *particles, const ExperimentParameters &parameters, const ExpansionProperties &metric)
Modifies positions and momentum of all particles to account for space-time deformation.
Definition: propagation.cc:86
void update_lattice(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.
Definition: density.h:542
void check_interactions_total(uint64_t interactions_total)
Make sure interactions_total can be represented as a 32-bit integer.
Definition: experiment.h:2415
std::tuple< double, FourVector, ThreeVector, ThreeVector, FourVector, FourVector, FourVector, FourVector > current_eckart(const ThreeVector &r, const ParticleList &plist, const DensityParameters &par, DensityType dens_type, bool compute_gradient, bool smearing)
Calculates Eckart rest frame density and 4-current of a given density type and optionally the gradien...
Definition: density.cc:171
double propagate_straight_line(Particles *particles, double to_time, const std::vector< FourVector > &beam_momentum)
Propagates the positions of all particles on a straight line to a given moment.
Definition: propagation.cc:44
ExperimentParameters create_experiment_parameters(Configuration &config)
Gathers all general Experiment parameters.
Definition: experiment.cc:132
constexpr double very_small_double
A very small double, used to avoid division by zero.
Definition: constants.h:40
double calculate_mean_field_energy(const Potentials &potentials, RectangularLattice< smash::DensityOnLattice > &jmu_B_lat, RectangularLattice< std::pair< ThreeVector, ThreeVector >> *em_lattice, const ExperimentParameters &parameters)
Calculate the total mean field energy of the system; this will be printed to the screen when SMASH is...
Definition: experiment.cc:384
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
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:626
constexpr double nucleon_mass
Nucleon mass in GeV.
Definition: constants.h:58
LatticeUpdate
Enumerator option for lattice updates.
Definition: lattice.h:36
std::ostream & operator<<(std::ostream &out, const Experiment< Modus > &e)
Creates a verbose textual description of the setup of the Experiment.
Definition: experiment.h:717
Potentials * pot_pointer
Pointer to a Potential class.
static constexpr int LMain
Definition: experiment.h:88
constexpr double really_small
Numerical error tolerance.
Definition: constants.h:37
DensityType
Allows to choose which kind of density to calculate.
Definition: density.h:36
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:28
Structure to contain custom data for output.
Struct containing the type of the metric and the expansion parameter of the metric.
Definition: propagation.h:26
Exception class that is thrown if an invalid modus is requested from the Experiment factory.
Definition: experiment.h:143
Exception class that is thrown if the requested output path in the Experiment factory is not existing...
Definition: experiment.h:152
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.
Helper structure for Experiment to hold output options and parameters.
RivetOutputParameters rivet_parameters
Rivet specfic parameters.