7 #ifndef SRC_INCLUDE_SMASH_EXPERIMENT_H_
8 #define SRC_INCLUDE_SMASH_EXPERIMENT_H_
41 #ifdef SMASH_USE_HEPMC
44 #ifdef SMASH_USE_RIVET
69 template <
typename T,
typename Ratio>
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]";
80 if (seconds < threshold_for_hours) {
81 return out << Minutes(seconds).count() <<
" [min]";
83 return out << Hours(seconds).count() <<
" [h]";
88 static constexpr
int LMain = LogArea::Main::id;
127 static std::unique_ptr<ExperimentBase>
create(
128 Configuration &config,
const std::filesystem::path &output_path);
144 using std::invalid_argument::invalid_argument;
153 using std::invalid_argument::invalid_argument;
157 template <
typename Modus>
159 template <
typename Modus>
186 template <
typename Modus>
213 const std::filesystem::path &output_path);
240 ParticleList &&remove_plist = {});
286 bool include_pauli_blocking =
true);
296 const std::filesystem::path &output_path,
321 const double end_time_propagation);
461 std::unique_ptr<RectangularLattice<FourVector>>
UB_lat_ =
nullptr;
467 std::unique_ptr<RectangularLattice<FourVector>>
UI3_lat_ =
nullptr;
473 std::unique_ptr<RectangularLattice<std::pair<ThreeVector, ThreeVector>>>
477 std::unique_ptr<RectangularLattice<std::pair<ThreeVector, ThreeVector>>>
481 std::unique_ptr<RectangularLattice<std::pair<ThreeVector, ThreeVector>>>
485 std::unique_ptr<RectangularLattice<EnergyMomentumTensor>>
Tmn_;
492 std::unique_ptr<RectangularLattice<std::array<FourVector, 4>>>
500 std::unique_ptr<RectangularLattice<std::array<FourVector, 4>>>
712 friend std::ostream &operator<<<>(std::ostream &out,
const Experiment &e);
716 template <
typename Modus>
718 out <<
"End time: " << e.
end_time_ <<
" fm\n";
723 template <
typename Modus>
725 const std::string &content,
726 const std::filesystem::path &output_path,
728 logg[
LExperiment].info() <<
"Adding output " << content <<
" of format "
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));
740 outputs_.emplace_back(
741 std::make_unique<RootOutput>(output_path, content, out_par));
745 "Root output requested, but Root support not compiled in");
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));
759 }
else if (
format ==
"Oscar1999" ||
format ==
"Oscar2013") {
760 outputs_.emplace_back(
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;
771 if (
format ==
"Lattice_Binary") {
772 printout_full_lattice_binary_td_ =
true;
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"));
797 "Requested HepMC_treeroot output not available, "
798 "ROOT or HepMC3 ROOTIO missing or not found by cmake.");
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"));
811 "Requested HepMC_treeroot output not available, "
812 "ROOT or HepMC3 ROOTIO missing or not found by cmake.");
817 "HepMC only available for Particles and "
818 "Collisions content. Requested for " +
823 "HepMC output requested, but HepMC support not compiled in");
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
831 static bool rivet_format_already_selected =
false;
833 if (rivet_format_already_selected) {
835 "Rivet output format can only be one, either YODA or YODA-full. "
836 "Only your first valid choice will be used.");
840 outputs_.emplace_back(std::make_unique<RivetOutput>(
842 rivet_format_already_selected =
true;
843 }
else if (
format ==
"YODA-full") {
844 outputs_.emplace_back(std::make_unique<RivetOutput>(
846 rivet_format_already_selected =
true;
849 "not one of YODA or YODA-full");
853 "Rivet output requested, but Rivet support not compiled in");
857 <<
"Unknown combination of format (" <<
format <<
") and content ("
858 << content <<
"). Fix the config.";
871 template <
typename Modus>
873 const std::filesystem::path &output_path)
876 modus_(std::invoke([&]() {
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);
891 config.
set_value(key_labels, temporary_taken_key);
893 return Modus{std::move(modus_config),
parameters_};
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()),
900 config.
take({
"Collision_Term",
"Force_Decays_At_End"},
true)),
901 use_grid_(config.
take({
"General",
"Use_Grid"},
true)),
904 config.
take({
"General",
"Expansion_Rate"}, 0.1)),
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"})),
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.");
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"});
927 std::ceil(
static_cast<double>(max_ensembles) / parameters_.n_ensembles);
935 throw std::invalid_argument(
936 "Covariant Gaussian derivatives only make sense for Covariant Gaussian "
944 parameters_.discrete_weight < (1. / 7.)) {
945 throw std::invalid_argument(
946 "The central weight for discrete smearing should be >= 1./7.");
951 throw std::invalid_argument(
952 "The stochastic criterion can only be employed for fixed time step "
953 "mode and with a grid!");
957 throw std::invalid_argument(
958 "The box modus can only be used with the fixed time step mode!");
962 " testparticles per particle.");
964 " parallel ensembles.");
967 if (dileptons_switch_) {
968 dilepton_finder_ = std::make_unique<DecayActionsFinderDilepton>();
970 if (photons_switch_ || bremsstrahlung_switch_) {
971 n_fractional_photons_ =
972 config.
take({
"Collision_Term",
"Photons",
"Fractional_Photons"}, 100);
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!");
981 "Resonance lifetime set to zero. Make sure resonances cannot "
983 "inelastically (e.g. resonance chains), else SMASH is known to "
986 action_finders_.emplace_back(std::make_unique<DecayActionsFinder>(
987 parameters_.res_lifetime_factor, parameters_.do_weak_decays));
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) &&
993 parameters_.use_monash_tune_default =
994 (modus_.is_collider() && modus_.sqrt_s_NN() >= 200.);
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));
1002 max_transverse_distance_sqr_ =
1003 parameters_.maximum_cross_section / M_PI *
fm2_mb;
1004 process_string_ptr_ = NULL;
1006 if (modus_.is_box()) {
1007 action_finders_.emplace_back(
1008 std::make_unique<WallCrossActionsFinder>(parameters_.box_length));
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.");
1016 if (config.
has_value({
"Output",
"Initial_Conditions",
"Proper_Time"})) {
1019 config.
take({
"Output",
"Initial_Conditions",
"Proper_Time"});
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;
1029 <<
"Nuclei passing time is too short, hypersurface proper time set "
1031 << lower_bound <<
" fm.";
1032 proper_time = lower_bound;
1036 double rapidity_cut = 0.0;
1037 if (config.
has_value({
"Output",
"Initial_Conditions",
"Rapidity_Cut"})) {
1039 config.
take({
"Output",
"Initial_Conditions",
"Rapidity_Cut"});
1040 if (rapidity_cut <= 0.0) {
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.");
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.");
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) {
1064 <<
"transverse momentum cut for initial conditions configured as "
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.");
1073 if (rapidity_cut > 0.0 || transverse_momentum_cut > 0.0) {
1074 kinematic_cuts_for_IC_output_ =
true;
1077 if (rapidity_cut > 0.0 && transverse_momentum_cut > 0.0) {
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) {
1084 <<
"Extracting initial conditions in kinematic range: "
1085 << -rapidity_cut <<
" <= y <= " << rapidity_cut <<
".";
1086 }
else if (transverse_momentum_cut > 0.0) {
1088 <<
"Extracting initial conditions in kinematic range: pT <= "
1089 << transverse_momentum_cut <<
" GeV.";
1092 <<
"Extracting initial conditions without kinematic cuts.";
1095 action_finders_.emplace_back(
1096 std::make_unique<HyperSurfaceCrossActionsFinder>(
1097 proper_time, rapidity_cut, transverse_momentum_cut));
1100 if (config.
has_value({
"Collision_Term",
"Pauli_Blocking"})) {
1102 pauli_blocker_ = std::make_unique<PauliBlocker>(
1108 {
"Collision_Term",
"String_Parameters",
"Power_Particle_Formation"},
1109 modus_.sqrt_s_NN() >= 200. ? -1. : 1.);
1352 " create OutputInterface objects");
1355 <<
"Density type printed to headers: " << dens_type_;
1361 if (output_path ==
"") {
1362 throw std::invalid_argument(
1363 "Invalid empty output path provided to Experiment constructor.");
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.";
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());
1373 output_contents.cbegin(), output_contents.cend(), list_of_formats.begin(),
1374 [&output_conf](std::string content) -> std::vector<std::string> {
1378 return output_conf.take({content.c_str(),
"Format"},
1379 std::vector<std::string>{});
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.");
1386 for (std::size_t i = 0; i < output_contents.size(); ++i) {
1387 if (list_of_formats[i].empty()) {
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) {
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();
1402 list_of_formats[i].clear();
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;
1410 if (outputs_.size() != total_number_of_requested_formats) {
1412 <<
"At least one invalid output format has been provided.";
1413 abort_because_of_invalid_input_file();
1424 throw std::invalid_argument(
"Can't use potentials without time steps!");
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!");
1435 << parameters_.labclock->timestep_duration();
1437 potentials_ = std::make_unique<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!");
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!");
1449 if (potentials_->use_skyrme()) {
1452 <<
"\t\tSkyrme_A [MeV] = " << potentials_->skyrme_a() <<
"\n";
1454 <<
"\t\tSkyrme_B [MeV] = " << potentials_->skyrme_b() <<
"\n";
1456 <<
"\t\t Skyrme_tau = " << potentials_->skyrme_tau() <<
"\n";
1458 if (potentials_->use_symmetry()) {
1460 <<
"Symmetry potential is:"
1461 <<
"\n S_pot [MeV] = " << potentials_->symmetry_S_pot() <<
"\n";
1463 if (potentials_->use_vdf()) {
1466 << potentials_->saturation_density() <<
"\n";
1467 for (
int i = 0; i < potentials_->number_of_terms(); i++) {
1469 <<
"\t\tCoefficient_" << i + 1 <<
" = "
1470 << 1000.0 * (potentials_->coeffs())[i] <<
" [MeV] \t Power_"
1471 << i + 1 <<
" = " << (potentials_->powers())[i] <<
"\n";
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\".");
1489 switch (parameters_.derivatives_mode) {
1500 switch (parameters_.rho_derivatives_mode) {
1509 if (potentials_->use_vdf()) {
1510 switch (parameters_.field_derivatives_mode) {
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!");
1534 throw std::runtime_error(
1535 "Can't use potentials without gradients of baryon current (Skyrme, "
1537 " or direct field derivatives (VDF)!");
1540 if (!(potentials_->use_vdf()) &&
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)");
1550 switch (parameters_.smearing_mode) {
1556 << parameters_.discrete_weight;
1560 << parameters_.triangular_range;
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"}) &&
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\".");
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\".");
1583 bool periodic = config.
take({
"Lattice",
"Periodic"}, modus_.is_box());
1584 const auto [l,
n, origin] = [&config, automatic,
this]() {
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"}));
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())) {
1598 const double gam = modus_.is_collider()
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};
1609 const int n_xy = std::ceil(2 * max_xy / 0.8);
1610 int nz = std::ceil(2 * max_z / 0.8);
1615 nz =
static_cast<int>(std::ceil(2 * max_z / 0.8 * gam));
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()) {
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};
1631 const int n_xyz = std::ceil(2 * max_d / 0.8);
1632 n_default = {n_xyz, n_xyz, n_xyz};
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));
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;
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;
1656 if (printout_tmn_ || printout_tmn_landau_ || printout_v_landau_) {
1657 Tmn_ = std::make_unique<RectangularLattice<EnergyMomentumTensor>>(
1658 l,
n, origin, periodic, LatticeUpdate::AtOutput);
1660 if (printout_j_QBS_) {
1661 j_QBS_lat_ = std::make_unique<DensityLattice>(l,
n, origin, periodic,
1662 LatticeUpdate::AtOutput);
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);
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);
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);
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);
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);
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);
1724 fields_lat_ = std::make_unique<FieldsLattice>(
1725 l,
n, origin, periodic, LatticeUpdate::EveryTimestep);
1728 if (dens_type_lattice_printout_ == DensityType::Baryon) {
1729 jmu_B_lat_ = std::make_unique<DensityLattice>(l,
n, origin, periodic,
1730 LatticeUpdate::AtOutput);
1732 if (dens_type_lattice_printout_ == DensityType::BaryonicIsospin) {
1733 jmu_I3_lat_ = std::make_unique<DensityLattice>(l,
n, origin, periodic,
1734 LatticeUpdate::AtOutput);
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);
1743 }
else if (printout_lattice_td_ || printout_full_lattice_any_td_) {
1745 "If you want Therm. VTK or Lattice output, configure a lattice for "
1747 }
else if (potentials_ && potentials_->use_coulomb()) {
1749 "Coulomb potential requires a lattice. Please add one to the "
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.";
1760 if (parameters_.potential_affect_threshold) {
1768 (jmu_B_lat_ ==
nullptr)) {
1769 throw std::runtime_error(
1770 "Lattice is necessary to calculate finite difference gradients.");
1774 if (config.
has_value({
"Forced_Thermalization"})) {
1775 Configuration th_conf =
1777 thermalizer_ = modus_.create_grandcan_thermalizer(th_conf);
1782 seed_ = config.
take({
"General",
"Randomseed"});
1813 uint64_t scatterings_this_interval,
1816 double E_mean_field,
1817 double E_mean_field_initial);
1854 double E_mean_field,
double modus_impact_parameter,
1856 bool projectile_target_interact,
1857 bool kinematic_cut_for_SMASH_IC);
1859 template <
typename Modus>
1869 while (r == INT64_MIN) {
1872 seed_ = std::abs(r);
1877 if (process_string_ptr_ != NULL) {
1878 process_string_ptr_->init_pythia_hadron_rndm();
1881 for (
Particles &particles : ensembles_) {
1886 double start_time = -1.0;
1890 if (modus_.is_collider()) {
1891 modus_.sample_impact();
1892 logg[
LExperiment].info(
"Impact parameter = ", modus_.impact_parameter(),
1895 for (
Particles &particles : ensembles_) {
1896 start_time = modus_.initial_conditions(&particles, parameters_);
1901 for (
Particles &particles : ensembles_) {
1902 modus_.impose_boundary_conditions(&particles, outputs_);
1905 double timestep = delta_time_startup_;
1907 switch (time_step_mode_) {
1911 timestep = end_time_ - start_time;
1913 const double max_dt = modus_.max_timestep(max_transverse_distance_sqr_);
1914 if (max_dt > 0. && max_dt < timestep) {
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.");
1926 clock_for_this_event =
1927 std::make_unique<UniformClock>(start_time, timestep, end_time_);
1928 parameters_.labclock = std::move(clock_for_this_event);
1931 parameters_.outputclock->reset(start_time,
true);
1933 parameters_.outputclock->remove_times_in_past(start_time);
1936 "Lab clock: t_start = ", parameters_.labclock->current_time(),
1937 ", dt = ", parameters_.labclock->timestep_duration());
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;
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";
1958 double E_mean_field = 0.0;
1963 if ((jmu_B_lat_ !=
nullptr)) {
1965 new_jmu_auxiliary_.get(), four_gradient_auxiliary_.get(),
1966 LatticeUpdate::EveryTimestep, DensityType::Baryon,
1967 density_param_, ensembles_,
1968 parameters_.labclock->timestep_duration(),
true);
1972 for (
auto &node : *jmu_B_lat_) {
1973 node.overwrite_drho_dt_to_zero();
1974 node.overwrite_djmu_dt_to_zero();
1977 EM_lat_.get(), parameters_);
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_);
1987 for (
const auto &output : outputs_) {
1988 for (
int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
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],
1994 event_ * parameters_.n_ensembles + i_ens,
1998 output->at_eventstart(ensembles_, event_);
2000 if (printout_full_lattice_any_td_) {
2001 switch (dens_type_lattice_printout_) {
2002 case DensityType::Baryon:
2004 DensityType::Baryon, *jmu_B_lat_);
2006 case DensityType::BaryonicIsospin:
2008 DensityType::BaryonicIsospin, *jmu_I3_lat_);
2010 case DensityType::None:
2014 DensityType::BaryonicIsospin, *jmu_custom_lat_);
2016 if (printout_tmn_) {
2018 dens_type_lattice_printout_, *Tmn_);
2020 if (printout_tmn_landau_) {
2022 dens_type_lattice_printout_, *Tmn_);
2024 if (printout_v_landau_) {
2026 dens_type_lattice_printout_, *Tmn_);
2028 if (printout_j_QBS_) {
2030 dens_type_lattice_printout_, *j_QBS_lat_);
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();
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));
2056 template <
typename Modus>
2058 bool include_pauli_blocking) {
2059 Particles &particles = ensembles_[i_ensemble];
2062 discarded_interactions_total_++;
2064 " (discarded: invalid)");
2069 if (include_pauli_blocking && pauli_blocker_ &&
2071 total_pauli_blocked_++;
2077 if (modus_.is_collider()) {
2078 int count_target = 0, count_projectile = 0;
2080 if (incoming.belongs_to() == BelongsTo::Projectile) {
2082 }
else if (incoming.belongs_to() == BelongsTo::Target) {
2086 if (count_target > 0 && count_projectile > 0) {
2087 projectile_target_interact_[i_ensemble] =
true;
2093 const auto id_process =
static_cast<uint32_t
>(interactions_total_ + 1);
2095 total_energy_violated_by_Pythia_ += action.
perform(&particles, id_process);
2097 interactions_total_++;
2098 if (action.
get_type() == ProcessType::Wall) {
2099 wall_actions_total_++;
2101 if (action.
get_type() == ProcessType::HyperSurfaceCrossing) {
2102 total_hypersurface_crossing_actions_++;
2107 if (dens_type_ != DensityType::None) {
2109 constexpr
bool compute_grad =
false;
2110 const bool smearing =
true;
2113 rho = std::get<0>(
current_eckart(r_interaction.threevec(), particles,
2114 density_param_, dens_type_, compute_grad,
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);
2148 if (photons_switch_ &&
2150 ScatterActionPhoton::is_kinematically_possible(
2154 constexpr
double action_time = 0.;
2156 n_fractional_photons_,
2172 photon_act.add_single_process();
2174 photon_act.perform_photons(outputs_);
2177 if (bremsstrahlung_switch_ &&
2178 BremsstrahlungAction::is_bremsstrahlung_reaction(
2182 constexpr
double action_time = 0.;
2185 n_fractional_photons_,
2202 brems_act.add_single_process();
2204 brems_act.perform_bremsstrahlung(outputs_);
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.");
2233 const double action_time = parameters_.labclock->current_time();
2234 if (!add_plist.empty()) {
2237 if (!add_plist.empty()) {
2239 auto action_add_particles = std::make_unique<FreeforallAction>(
2240 ParticleList{}, add_plist, action_time);
2241 perform_action(*action_add_particles, 0);
2243 if (!remove_plist.empty()) {
2245 const auto number_of_particles_to_be_removed = remove_plist.size();
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();
2258 remove_plist.end());
2259 if (
auto delta = number_of_particles_to_be_removed - remove_plist.size();
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.");
2267 if (!remove_plist.empty()) {
2269 auto action_remove_particles = std::make_unique<FreeforallAction>(
2270 remove_plist, ParticleList{}, action_time);
2271 perform_action(*action_remove_particles, 0);
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.");
2281 thermalizer_->is_time_to_thermalize(parameters_.labclock)) {
2282 const bool ignore_cells_under_treshold =
true;
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);
2292 if (th_act.any_particles_thermalized()) {
2293 perform_action(th_act, i_ens);
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) {
2303 const double min_cell_length = compute_min_cell_length(dt);
2308 const bool include_unformed_particles = IC_output_switch_;
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);
2318 const double gcell_vol = grid.cell_volume();
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_));
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_));
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());
2346 ++(*parameters_.outputclock);
2348 intermediate_output();
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);
2357 update_potentials();
2358 update_momenta(ensembles_, parameters_.labclock->timestep_duration(),
2359 *potentials_, FB_lat_.get(), FI3_lat_.get(),
2366 for (Particles &particles : ensembles_) {
2371 ++(*parameters_.labclock);
2379 if (!potentials_ && !parameters_.strings_switch &&
2381 std::string err_msg = conserved_initial_.report_deviations(ensembles_);
2382 if (!err_msg.empty()) {
2384 throw std::runtime_error(
"Violation of conserved quantities!");
2389 if (pauli_blocker_) {
2391 "Interactions: Pauli-blocked/performed = ", total_pauli_blocked_,
"/",
2392 interactions_total_ - wall_actions_total_);
2396 template <
typename Modus>
2401 if (dilepton_finder_ !=
nullptr) {
2402 for (
const auto &output : outputs_) {
2403 dilepton_finder_->shine(particles, output.get(), dt);
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!");
2422 template <
typename Modus>
2424 Actions &actions,
int i_ensemble,
const double end_time_propagation) {
2425 Particles &particles = ensembles_[i_ensemble];
2427 "Timestepless propagation: ",
"Actions size = ", actions.
size(),
2428 ", end time = ", end_time_propagation);
2436 ActionPtr act = actions.
pop();
2437 if (!act->is_valid(particles)) {
2438 discarded_interactions_total_++;
2440 " (discarded: invalid)");
2444 ", action time = ", act->time_of_execution());
2447 propagate_and_shine(act->time_of_execution(), particles);
2454 act->update_incoming(particles);
2455 const bool performed = perform_action(*act, i_ensemble);
2465 const double end_time_timestep = parameters_.labclock->next_time();
2467 const double time_left = end_time_timestep - act->time_of_execution();
2468 const ParticleList &outgoing_particles = act->outgoing_particles();
2470 const double gcell_vol = 0.0;
2471 for (
const auto &finder : action_finders_) {
2473 actions.
insert(finder->find_actions_in_cell(outgoing_particles, time_left,
2474 gcell_vol, beam_momentum_));
2476 actions.
insert(finder->find_actions_with_surrounding_particles(
2477 outgoing_particles, particles, time_left, beam_momentum_));
2483 propagate_and_shine(end_time_propagation, particles);
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;
2501 if ((jmu_B_lat_ !=
nullptr)) {
2503 EM_lat_.get(), parameters_);
2510 if (modus_.is_box()) {
2511 double tmp = (E_mean_field - initial_mean_field_energy_) /
2512 (E_mean_field + initial_mean_field_energy_);
2518 if (std::abs(tmp) > 0.01) {
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))] = "
2527 <<
"\n\t\t E_MF/E_MF(t=0) = "
2528 << E_mean_field / initial_mean_field_energy_ <<
"\n\n";
2535 ensembles_, interactions_this_interval, conserved_initial_, time_start_,
2536 parameters_.outputclock->current_time(), E_mean_field,
2537 initial_mean_field_energy_);
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()) {
2548 for (
int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
2550 ensembles_, E_mean_field, modus_.impact_parameter(), parameters_,
2551 projectile_target_interact_[i_ens], kinematic_cuts_for_IC_output_);
2553 output->at_intermediate_time(ensembles_[i_ens], parameters_.outputclock,
2554 density_param_, event_info);
2555 computational_frame_time = event_info.current_time;
2558 output->at_intermediate_time(ensembles_, parameters_.outputclock,
2562 switch (dens_type_lattice_printout_) {
2563 case DensityType::Baryon:
2565 density_param_, ensembles_,
false);
2567 DensityType::Baryon, *jmu_B_lat_);
2568 output->thermodynamics_lattice_output(*jmu_B_lat_,
2569 computational_frame_time);
2571 case DensityType::BaryonicIsospin:
2573 DensityType::BaryonicIsospin, density_param_,
2576 DensityType::BaryonicIsospin,
2578 output->thermodynamics_lattice_output(*jmu_I3_lat_,
2579 computational_frame_time);
2581 case DensityType::None:
2585 dens_type_lattice_printout_, density_param_,
2588 dens_type_lattice_printout_,
2590 output->thermodynamics_lattice_output(*jmu_custom_lat_,
2591 computational_frame_time);
2593 if (printout_tmn_ || printout_tmn_landau_ || printout_v_landau_) {
2595 density_param_, ensembles_,
false);
2596 if (printout_tmn_) {
2598 dens_type_lattice_printout_, *Tmn_);
2599 output->thermodynamics_lattice_output(
2602 if (printout_tmn_landau_) {
2604 dens_type_lattice_printout_, *Tmn_);
2605 output->thermodynamics_lattice_output(
2607 computational_frame_time);
2609 if (printout_v_landau_) {
2611 dens_type_lattice_printout_, *Tmn_);
2612 output->thermodynamics_lattice_output(
2614 computational_frame_time);
2618 output->fields_output(
"Efield",
"Bfield", *EM_lat_);
2620 if (printout_j_QBS_) {
2621 output->thermodynamics_lattice_output(
2622 *j_QBS_lat_, computational_frame_time, ensembles_, density_param_);
2626 output->thermodynamics_output(*thermalizer_);
2632 template <
typename Modus>
2635 if (potentials_->use_symmetry() && jmu_I3_lat_ !=
nullptr) {
2637 new_jmu_auxiliary_.get(), four_gradient_auxiliary_.get(),
2638 LatticeUpdate::EveryTimestep, DensityType::BaryonicIsospin,
2639 density_param_, ensembles_,
2640 parameters_.labclock->timestep_duration(),
true);
2642 if ((potentials_->use_skyrme() || potentials_->use_symmetry()) &&
2643 jmu_B_lat_ !=
nullptr) {
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];
2655 double baryon_density = jB.rho();
2659 if (potentials_->use_skyrme()) {
2661 flow_four_velocity_B * potentials_->skyrme_pot(baryon_density);
2663 potentials_->skyrme_force(baryon_density, baryon_grad_j0,
2664 baryon_dvecj_dt, baryon_curl_vecj);
2666 if (potentials_->use_symmetry() && jmu_I3_lat_ !=
nullptr) {
2667 auto jI3 = (*jmu_I3_lat_)[i];
2670 ? jI3.jmu_net() / jI3.rho()
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,
2681 if (potentials_->use_coulomb()) {
2683 DensityType::Charge, density_param_, ensembles_,
true);
2684 for (
size_t i = 0; i < EM_lat_->size(); i++) {
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);
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);
2697 if (potentials_->use_vdf() && jmu_B_lat_ !=
nullptr) {
2699 new_jmu_auxiliary_.get(), four_gradient_auxiliary_.get(),
2700 LatticeUpdate::EveryTimestep, DensityType::Baryon,
2701 density_param_, ensembles_,
2702 parameters_.labclock->timestep_duration(),
true);
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());
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());
2722 auto Amu = (*fields_lat_)[i];
2723 (*FB_lat_)[i] = potentials_->vdf_force(
2724 Amu.grad_A0(), Amu.dvecA_dt(), Amu.curl_vecA());
2732 template <
typename Modus>
2736 bool actions_performed, decays_found;
2737 uint64_t interactions_old;
2739 decays_found =
false;
2740 interactions_old = interactions_total_;
2741 for (
int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
2745 if (dilepton_finder_ !=
nullptr) {
2746 for (
const auto &output : outputs_) {
2747 dilepton_finder_->shine_final(ensembles_[i_ens], output.get(),
true);
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;
2760 perform_action(*actions.
pop(), i_ens,
false);
2763 actions_performed = interactions_total_ > interactions_old;
2765 if (decays_found && !actions_performed) {
2766 throw std::runtime_error(
"Final decays were found but not performed.");
2769 }
while (actions_performed);
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);
2781 template <
typename Modus>
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;
2795 if ((jmu_B_lat_ !=
nullptr)) {
2797 EM_lat_.get(), parameters_);
2800 if (std::abs(parameters_.labclock->current_time() - end_time_) >
2803 <<
"SMASH not propagated until configured end time. Current time = "
2804 << parameters_.labclock->current_time()
2805 <<
"fm. End time = " << end_time_ <<
"fm.";
2808 ensembles_, interactions_this_interval, conserved_initial_,
2809 time_start_, end_time_, E_mean_field, initial_mean_field_energy_);
2811 int total_particles = 0;
2812 for (
const Particles &particles : ensembles_) {
2813 total_particles += particles.
size();
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_;
2822 if (std::fabs(fraction_of_total_system_energy_removed - 1.) >
2824 throw std::runtime_error(
2825 "There is remaining energy in the system although all particles "
2828 std::to_string((initial_system_energy_plus_Pythia_violations -
2829 total_energy_removed_)) +
2834 <<
"Time real: " << SystemClock::now() - time_start_;
2836 <<
"Interactions before reaching hypersurface: "
2837 << interactions_total_ - wall_actions_total_ -
2838 total_hypersurface_crossing_actions_;
2840 <<
"Total number of particles removed on hypersurface: "
2841 << total_hypersurface_crossing_actions_;
2844 const double precent_discarded =
2845 interactions_total_ > 0
2846 ?
static_cast<double>(discarded_interactions_total_) * 100.0 /
2849 std::stringstream msg_discarded;
2851 <<
"Discarded interaction number: " << discarded_interactions_total_
2852 <<
" (" << precent_discarded
2853 <<
"% of the total interaction number including wall crossings)";
2857 <<
"Time real: " << SystemClock::now() - time_start_;
2861 precent_discarded > 1.0) {
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.";
2872 << interactions_total_ - wall_actions_total_;
2876 int unformed_particles_count = 0;
2877 for (
const Particles &particles : ensembles_) {
2879 if (particle.formation_time() > end_time_) {
2880 unformed_particles_count++;
2884 if (unformed_particles_count > 0) {
2886 "End time might be too small. ", unformed_particles_count,
2887 " unformed particles were found at the end of the evolution.");
2892 count_nonempty_ensembles();
2894 for (
const auto &output : outputs_) {
2895 for (
int i_ens = 0; i_ens < parameters_.n_ensembles; i_ens++) {
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],
2901 event_ * parameters_.n_ensembles + i_ens, event_info);
2904 output->at_eventend(ensembles_, event_);
2907 if (dens_type_lattice_printout_ != DensityType::None) {
2909 dens_type_lattice_printout_);
2912 if (printout_tmn_) {
2914 dens_type_lattice_printout_);
2917 if (printout_tmn_landau_) {
2919 dens_type_lattice_printout_);
2922 if (printout_v_landau_) {
2924 dens_type_lattice_printout_);
2927 if (printout_j_QBS_) {
2933 template <
typename Modus>
2935 for (
bool has_interaction : projectile_target_interact_) {
2936 if (has_interaction) {
2937 nonempty_ensembles_++;
2942 template <
typename Modus>
2945 return event_ >= nevents_;
2948 if (nonempty_ensembles_ >= minimum_nonempty_ensembles_) {
2951 if (event_ >= max_events_) {
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.";
2964 throw std::runtime_error(
"Event counting option is invalid");
2968 template <
typename Modus>
2973 template <
typename Modus>
2976 for (event_ = 0; !is_finished(); event_++) {
2977 mainlog.info() <<
"Event " << event_;
2980 initialize_new_event();
2982 run_time_evolution(end_time_);
2984 if (force_decays_) {
Collection of useful type aliases to measure and output the (real) runtime.
A stream modifier that allows to colorize the log output.
Action is the base class for a generic process that takes a number of incoming particles and transfor...
virtual ProcessType get_type() const
Get the process type.
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.
const ParticleList & incoming_particles() const
Get the list of particles that go into the action.
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].
FourVector get_interaction_point() const
Get the interaction point.
bool is_valid(const Particles &particles) const
Check whether the action still applies.
bool is_pauli_blocked(const std::vector< Particles > &ensembles, const PauliBlocker &p_bl) const
Check if the action is Pauli-blocked.
The Actions class abstracts the storage and manipulation of actions.
ActionPtr pop()
Return the first action in the list and removes it from the list.
double earliest_time() const
Return time of execution of earliest action.
ActionList::size_type size() const
void insert(ActionList &&new_acts)
Insert a list of actions into this object.
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.
Non-template interface to Experiment<Modus>.
static std::unique_ptr< ExperimentBase > create(Configuration &config, const std::filesystem::path &output_path)
Factory method that creates and initializes a new Experiment<Modus>.
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.
const bool bremsstrahlung_switch_
This indicates whether bremsstrahlung is switched on.
void propagate_and_shine(double to_time, Particles &particles)
Propagate all particles until time to_time without any interactions and shine dileptons.
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.
std::unique_ptr< ActionFinderInterface > photon_finder_
The (Scatter) Actions Finder for Direct Photons.
const bool force_decays_
This indicates whether we force all resonances to decay in the last timestep.
double initial_mean_field_energy_
The initial total mean field energy in the system.
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.
std::vector< std::unique_ptr< ActionFinderInterface > > action_finders_
The Action finder objects.
bool printout_tmn_
Whether to print the energy-momentum tensor.
QuantumNumbers conserved_initial_
The conserved quantities of the system.
DensityParameters density_param_
Structure to precalculate and hold parameters for density computations.
double next_output_time() const
Shortcut for next output time.
bool printout_full_lattice_ascii_td_
Whether to print the thermodynamics quantities evaluated on the lattices, point by point,...
double total_energy_removed_
Total energy removed from the system in hypersurface crossing actions.
DensityType dens_type_lattice_printout_
Type of density for lattice printout.
double max_transverse_distance_sqr_
Maximal distance at which particles can interact in case of the geometric criterion,...
bool printout_j_QBS_
Whether to print the Q, B, S 4-currents.
std::unique_ptr< GrandCanThermalizer > thermalizer_
Instance of class used for forced thermalization.
void count_nonempty_ensembles()
Counts the number of ensembles in wich interactions took place at the end of an event.
const TimeStepMode time_step_mode_
This indicates whether to use time steps.
std::unique_ptr< DensityLattice > jmu_custom_lat_
Custom density on the lattices.
bool printout_full_lattice_any_td_
Whether to print the thermodynamics quantities evaluated on the lattices, point by point,...
void increase_event_number()
Increases the event number by one.
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 ...
Particles * first_ensemble()
Provides external access to SMASH particles.
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,...
int minimum_nonempty_ensembles_
The number of ensembles, in which interactions take place, to be calculated.
std::unique_ptr< DensityLattice > jmu_B_lat_
Baryon density on the lattice.
DensityType dens_type_
Type of density to be written to collision headers.
const int nevents_
Number of events.
void intermediate_output()
Intermediate output during an event.
std::vector< FourVector > beam_momentum_
The initial nucleons in the ColliderModus propagate with beam_momentum_, if Fermi motion is frozen.
std::unique_ptr< RectangularLattice< std::pair< ThreeVector, ThreeVector > > > FI3_lat_
Lattices for the electric and magnetic component of the symmetry force.
std::unique_ptr< RectangularLattice< FourVector > > new_jmu_auxiliary_
Auxiliary lattice for values of jmu at a time step t0 + dt.
double compute_min_cell_length(double dt) const
Calculate the minimal size for the grid cells such that the ScatterActionsFinder will find all collis...
void initialize_new_event()
This is called in the beginning of each event.
bool is_finished()
Checks wether the desired number events have been calculated.
int n_fractional_photons_
Number of fractional photons produced per single reaction.
const double delta_time_startup_
The clock's timestep size at start up.
std::unique_ptr< RectangularLattice< EnergyMomentumTensor > > Tmn_
Lattices of energy-momentum tensors for printout.
std::vector< Particles > ensembles_
Complete particle list, all ensembles in one vector.
std::unique_ptr< RectangularLattice< FourVector > > new_fields_auxiliary_
Auxiliary lattice for values of Amu at a time step t0 + dt.
SystemTimePoint time_start_
system starting time of the simulation
uint64_t previous_wall_actions_total_
Total number of wall-crossings for previous timestep.
std::unique_ptr< RectangularLattice< std::pair< ThreeVector, ThreeVector > > > FB_lat_
Lattices for the electric and magnetic components of the Skyrme or VDF force.
OutputsList outputs_
A list of output formaters.
std::unique_ptr< RectangularLattice< std::array< FourVector, 4 > > > fields_four_gradient_auxiliary_
Auxiliary lattice for calculating the four-gradient of Amu.
void do_final_decays()
Performs the final decays of an event.
std::unique_ptr< PauliBlocker > pauli_blocker_
An instance of PauliBlocker class that stores parameters needed for Pauli blocking calculations and c...
bool printout_v_landau_
Whether to print the 4-velocity in Landau frame.
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...
bool printout_lattice_td_
Whether to print the thermodynamics quantities evaluated on the lattices.
std::unique_ptr< RectangularLattice< std::pair< ThreeVector, ThreeVector > > > EM_lat_
Lattices for electric and magnetic field in fm^-2.
std::unique_ptr< FieldsLattice > fields_lat_
Mean-field A^mu on the lattice.
int max_events_
Maximum number of events to be calculated in order obtain the desired number of non-empty events usin...
int nonempty_ensembles_
Number of ensembles containing an interaction.
OutputPtr photon_output_
The Photon output.
EventCounting event_counting_
The way in which the number of calculated events is specified.
const bool dileptons_switch_
This indicates whether dileptons are switched on.
Modus modus_
Instance of the Modus template parameter.
std::unique_ptr< RectangularLattice< FourVector > > old_jmu_auxiliary_
Auxiliary lattice for values of jmu at a time step t0.
std::unique_ptr< DecayActionsFinderDilepton > dilepton_finder_
The Dilepton Action Finder.
std::unique_ptr< DensityLattice > j_QBS_lat_
4-current for j_QBS lattice output
bool printout_tmn_landau_
Whether to print the energy-momentum tensor in Landau frame.
bool printout_full_lattice_binary_td_
Whether to print the thermodynamics quantities evaluated on the lattices, point by point,...
std::unique_ptr< RectangularLattice< std::array< FourVector, 4 > > > four_gradient_auxiliary_
Auxiliary lattice for calculating the four-gradient of jmu.
const bool photons_switch_
This indicates whether photons are switched on.
StringProcess * process_string_ptr_
Pointer to the string process class object, which is used to set the random seed for PYTHIA objects i...
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...
ExperimentParameters parameters_
Struct of several member variables.
std::unique_ptr< RectangularLattice< FourVector > > old_fields_auxiliary_
Auxiliary lattice for values of Amu at a time step t0.
std::unique_ptr< DensityLattice > jmu_el_lat_
Electric charge density on the lattice.
std::unique_ptr< Potentials > potentials_
An instance of potentials class, that stores parameters of potentials, calculates them and their grad...
uint64_t wall_actions_total_
Total number of wall-crossings for current timestep.
uint64_t total_hypersurface_crossing_actions_
Total number of particles removed from the evolution in hypersurface crossing actions.
uint64_t interactions_total_
Total number of interactions for current timestep.
const bool use_grid_
This indicates whether to use the grid.
std::unique_ptr< DensityLattice > jmu_I3_lat_
Isospin projection density on the lattice.
uint64_t total_pauli_blocked_
Total number of Pauli-blockings for current timestep.
void final_output()
Output at the end of an event.
std::vector< bool > projectile_target_interact_
Whether the projectile and the target collided.
Modus * modus()
Provides external access to SMASH calculation modus.
void update_potentials()
Recompute potentials on lattices if necessary.
const bool IC_output_switch_
This indicates whether the IC output is enabled.
OutputPtr dilepton_output_
The Dilepton output.
int64_t seed_
random seed for the next event.
std::vector< Particles > * all_ensembles()
Getter for all ensembles.
uint64_t previous_interactions_total_
Total number of interactions for previous timestep.
uint64_t discarded_interactions_total_
Total number of discarded interactions, because they were invalidated before they could be performed.
void run() override
Runs the experiment.
bool kinematic_cuts_for_IC_output_
This indicates whether kinematic cuts are enabled for the IC output.
const double end_time_
simulation time at which the evolution is stopped.
double total_energy_violated_by_Pythia_
Total energy violation introduced by Pythia.
The FourVector class holds relevant values in Minkowski spacetime with (+, −, −, −) metric signature.
ParticleData contains the dynamic information of a certain particle.
static double formation_power_
Power with which the cross section scaling factor grows in time.
The Particles class abstracts the storage and manipulation of particles.
A class that stores parameters of potentials, calculates potentials and their gradients.
A container for storing conserved values.
A container class to hold all the arrays on the lattice and access them.
String excitation processes used in SMASH.
ThermalizationAction implements forced thermalization as an Action class.
The ThreeVector class represents a physical three-vector with the components .
@ 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.
std::ostream & operator<<(std::ostream &out, const ActionPtr &action)
Convenience: dereferences the ActionPtr to Action.
std::array< einhard::Logger<>, std::tuple_size< LogArea::AreaTuple >::value > logg
An array that stores all pre-configured Logger objects.
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...
std::unique_ptr< OutputInterface > create_oscar_output(const std::string &format, const std::string &content, const std::filesystem::path &path, const OutputParameters &out_par)
#define likely(x)
Tell the branch predictor that this expression is likely true.
Engine::result_type advance()
Advance the engine's state and return the generated value.
void set_seed(T &&seed)
Sets the seed of the random number engine.
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:
static constexpr int LInitialConditions
EventInfo fill_event_info(const std::vector< Particles > &ensembles, double E_mean_field, double modus_impact_parameter, const ExperimentParameters ¶meters, bool projectile_target_interact, bool kinematic_cut_for_SMASH_IC)
Generate the EventInfo object which is passed to outputs_.
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.
void expand_space_time(Particles *particles, const ExperimentParameters ¶meters, const ExpansionProperties &metric)
Modifies positions and momentum of all particles to account for space-time deformation.
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.
void check_interactions_total(uint64_t interactions_total)
Make sure interactions_total can be represented as a 32-bit integer.
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...
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.
ExperimentParameters create_experiment_parameters(Configuration &config)
Gathers all general Experiment parameters.
constexpr double very_small_double
A very small double, used to avoid division by zero.
double calculate_mean_field_energy(const Potentials &potentials, RectangularLattice< smash::DensityOnLattice > &jmu_B_lat, RectangularLattice< std::pair< ThreeVector, ThreeVector >> *em_lattice, const ExperimentParameters ¶meters)
Calculate the total mean field energy of the system; this will be printed to the screen when SMASH is...
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.
void validate_and_adjust_particle_list(ParticleList &particle_list)
Validate a particle list adjusting each particle to be a valid SMASH particle.
constexpr double nucleon_mass
Nucleon mass in GeV.
LatticeUpdate
Enumerator option for lattice updates.
std::ostream & operator<<(std::ostream &out, const Experiment< Modus > &e)
Creates a verbose textual description of the setup of the Experiment.
Potentials * pot_pointer
Pointer to a Potential class.
static constexpr int LMain
constexpr double really_small
Numerical error tolerance.
DensityType
Allows to choose which kind of density to calculate.
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.
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.
Structure to contain custom data for output.
Struct containing the type of the metric and the expansion parameter of the metric.
Exception class that is thrown if an invalid modus is requested from the Experiment factory.
Exception class that is thrown if the requested output path in the Experiment factory is not existing...
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.