CMS 3D CMS Logo

SuepShower.cc
Go to the documentation of this file.
1 // This is file contains the code to generate a dark sector shower in a strongly coupled, quasi-conformal hidden valley, often
2 // referred to as "soft unclustered energy patterns (SUEP)" or "softbomb" events. The shower is generated in its rest frame and
3 // for a realistic simulation this class needs to be interfaced with an event generator such as madgraph, pythia or herwig.
4 //
5 // The algorithm relies on arXiv:1305.5226. See arXiv:1612.00850 for a description of the model.
6 // Please cite both papers when using this code.
7 //
8 //
9 // Applicability:
10 // The ratio of the parameters m and T (m/T) should be an O(1) number. For m/T>>1 and m/T<<1 the theoretical description of the shower is not valid.
11 // The mass of the scalar which initiates the shower should be much larger than the mass of the mesons, in other words M>>m,T, by at least an order of magnitude.
12 //
13 // Written by Simon Knapen on 12/22/2019, following 1205.5226
14 // Adapted by Carlos Erice for CMSSW
16 
17 // constructor
18 SuepShower::SuepShower(double mass, double temperature, Pythia8::Rndm* rndmPtr) {
19  edm::LogInfo("Pythia8Interface") << "Creating a SuepShower module";
20  // Parameter setup
22  fRndmPtr_ = rndmPtr;
23  // No need to save temperature, as everything depends on it through m/T
25 
26  // 128 bit numerical precision, i.e. take machine precision
27  tolerance_ = boost::math::tools::eps_tolerance<double>(128);
28 
29  // Median momentum
31  double pmax = sqrt(2 + 2 * sqrt(1 + mass_over_T_ * mass_over_T_)) / mass_over_T_;
32  // Find the two cases where f(p)/f(pmax) = 1/e, given MB distribution, one appears at each side of pmax always
33  p_plus_ = (boost::math::tools::bisect(
34  boost::bind(&SuepShower::logTestFunction, this, boost::placeholders::_1), pmax, 50.0, tolerance_))
35  .first; // first root
36  p_minus_ = (boost::math::tools::bisect(
37  boost::bind(&SuepShower::logTestFunction, this, boost::placeholders::_1), 0.0, pmax, tolerance_))
38  .first; // second root
39  // Define the auxiliar quantities for the random generation of momenta
44  q_m_ = 1 - (q_plus_ + q_minus_);
45 }
46 
48 
49 // Generate a shower in the rest frame of the mediator
50 std::vector<Pythia8::Vec4> SuepShower::generateShower(double energy) {
52  std::vector<Pythia8::Vec4> shower;
53  double shower_energy = 0.0;
54 
55  // Fill up shower record
56  while (shower_energy < mediator_energy_) {
57  shower.push_back(generateFourVector());
58  shower_energy += (shower.back()).e();
59  }
60 
61  // Reballance momenta to ensure conservation
62  // Correction starts at 0
63  Pythia8::Vec4 correction = Pythia8::Vec4(0., 0., 0., 0.);
64  for (const auto& daughter : shower) {
65  correction = correction + daughter;
66  }
67  correction = correction / shower.size();
68  // We only want to correct momenta first
69  correction.e(0);
70  for (auto& daughter : shower) {
71  daughter = daughter - correction;
72  }
73 
74  // With momentum conserved, balance energy. scale is the multiplicative factor needed such that sum_daughters((scale*p)^2+m^2) = E_parent, i.e. energy is conserved
75  double scale;
76  scale =
77  (boost::math::tools::bisect(
78  boost::bind(&SuepShower::reballanceFunction, this, boost::placeholders::_1, shower), 0.0, 2.0, tolerance_))
79  .first;
80 
81  for (auto& daughter : shower) {
82  daughter.px(daughter.px() * scale);
83  daughter.py(daughter.py() * scale);
84  daughter.pz(daughter.pz() * scale);
85  // Force everything to be on-shell
86  daughter.e(sqrt(daughter.pAbs2() + darkmeson_mass_ * darkmeson_mass_));
87  }
88  return shower;
89 }
90 
92 // Maxwell-boltzman distribution, slightly massaged
93 const double SuepShower::fMaxwellBoltzmann(double p) {
94  return p * p * exp(-mass_over_T_ * p * p / (1 + sqrt(1 + p * p)));
95 }
96 
97 // Derivative of maxwell-boltzmann
98 const double SuepShower::fMaxwellBoltzmannPrime(double p) {
99  return exp(-mass_over_T_ * p * p / (1 + sqrt(1 + p * p))) * p * (2 - mass_over_T_ * p * p / sqrt(1 + p * p));
100 }
101 
102 // Test function to be solved for p_plus_ and p_minus_
103 const double SuepShower::logTestFunction(double p) { return log(fMaxwellBoltzmann(p) / fMaxwellBoltzmann(p_m_)) + 1.0; }
104 
105 // Generate one random 4 std::vector from the thermal distribution
107  double en, phi, theta, momentum; //kinematic variables of the 4 std::vector
108 
109  // First do momentum, following naming of arxiv:1305.5226
110  double U, V, X, Y, E;
111  int i = 0;
112  while (i < 100) {
113  // Very rarely (<1e-5 events) takes more than one loop to converge, set to 100 for safety
114  U = fRndmPtr_->flat();
115  V = fRndmPtr_->flat();
116 
117  if (U < q_m_) {
118  Y = U / q_m_;
119  X = (1 - Y) * (p_minus_ + lambda_minus_) + Y * (p_plus_ - lambda_plus_);
120  if (V < fMaxwellBoltzmann(X) / fMaxwellBoltzmann(p_m_) && X > 0) {
121  break;
122  }
123  } else {
124  if (U < q_m_ + q_plus_) {
125  E = -log((U - q_m_) / q_plus_);
126  X = p_plus_ - lambda_plus_ * (1 - E);
127  if (V < exp(E) * fMaxwellBoltzmann(X) / fMaxwellBoltzmann(p_m_) && X > 0) {
128  break;
129  }
130  } else {
131  E = -log((U - (q_m_ + q_plus_)) / q_minus_);
132  X = p_minus_ + lambda_minus_ * (1 - E);
133  if (V < exp(E) * fMaxwellBoltzmann(X) / fMaxwellBoltzmann(p_m_) && X > 0) {
134  break;
135  }
136  }
137  }
138  }
139  // X is the dimensionless momentum, p/m
140  momentum = X * darkmeson_mass_;
141 
142  // now do the angles, isotropic
143  phi = 2.0 * M_PI * (fRndmPtr_->flat());
144  theta = acos(2.0 * (fRndmPtr_->flat()) - 1.0);
145 
146  // compose the 4 std::vector
147  en = sqrt(momentum * momentum + darkmeson_mass_ * darkmeson_mass_);
148  Pythia8::Vec4 daughterFourMomentum =
149  Pythia8::Vec4(momentum * cos(phi) * sin(theta), momentum * sin(phi) * sin(theta), momentum * sin(theta), en);
150  return daughterFourMomentum;
151 }
152 
153 // Auxiliary function, to be solved in order to impose energy conservation
154 // To ballance energy, we solve for "scale" by demanding that this function vanishes, i.e. that shower energy and mediator energy are equal
155 // By rescaling the momentum rather than the energy, avoid having to do these annoying rotations from the previous version
156 const double SuepShower::reballanceFunction(double scale, const std::vector<Pythia8::Vec4>& shower) {
157  double showerEnergy = 0.0;
158  for (const auto& daughter : shower) {
159  showerEnergy += sqrt(scale * scale * daughter.pAbs2() + darkmeson_mass_ * darkmeson_mass_);
160  }
161  return showerEnergy - mediator_energy_;
162 }
double darkmeson_mass_
Definition: SuepShower.h:26
const double reballanceFunction(double scale, const std::vector< Pythia8::Vec4 > &shower)
Definition: SuepShower.cc:156
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
#define X(str)
Definition: MuonsGrabber.cc:38
Pythia8::Rndm * fRndmPtr_
Definition: SuepShower.h:46
double p_m_
Definition: SuepShower.h:37
double mass_over_T_
Definition: SuepShower.h:27
double lambda_plus_
Definition: SuepShower.h:41
boost::math::tools::eps_tolerance< double > tolerance_
Definition: SuepShower.h:33
uint32_t T const *__restrict__ uint32_t const *__restrict__ int32_t int Histo::index_type cudaStream_t V
T sqrt(T t)
Definition: SSEVec.h:19
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double q_plus_
Definition: SuepShower.h:43
double p_plus_
Definition: SuepShower.h:39
ExtVec< T, 4 > Vec4
Definition: ExtVec.h:64
const double fMaxwellBoltzmannPrime(double p)
Definition: SuepShower.cc:98
#define M_PI
Log< level::Info, false > LogInfo
const double logTestFunction(double p)
Definition: SuepShower.cc:103
const Pythia8::Vec4 generateFourVector()
Definition: SuepShower.cc:106
double mediator_energy_
Definition: SuepShower.h:30
double q_minus_
Definition: SuepShower.h:43
double p_minus_
Definition: SuepShower.h:39
double q_m_
Definition: SuepShower.h:43
double lambda_minus_
Definition: SuepShower.h:41
SuepShower(double mass, double temperature, Pythia8::Rndm *rndmPtr)
Definition: SuepShower.cc:18
Geom::Theta< T > theta() const
const double fMaxwellBoltzmann(double p)
Definition: SuepShower.cc:93
std::vector< Pythia8::Vec4 > generateShower(double energy)
Definition: SuepShower.cc:50