CMS 3D CMS Logo

SuepDecay.cc
Go to the documentation of this file.
1 #include "SuepDecay.h"
2 
4  : idMediator_(iConfig.getParameter<int>("idMediator")),
5  idDark_(iConfig.getParameter<int>("idDark")),
6  temperature_(iConfig.getParameter<double>("temperature")) {}
7 
9  mDark_ = particleDataPtr->m0(idDark_);
10  bool medDecay = particleDataPtr->mayDecay(idMediator_);
11  if (!medDecay) {
12  loggerPtr->errorMsg("SuepDecay::initAfterBeams", "mediator decay should be enabled");
13  return false;
14  }
15 
16  //construct the shower helper
17  suep_shower_ = std::make_unique<SuepShower>(mDark_, temperature_, rndmPtr);
18 
19  return true;
20 }
21 
22 //based on https://gitlab.com/simonknapen/suep_generator/-/blob/master/suep_main.cc:AttachSuepShower
23 bool SuepDecay::doVetoProcessLevel(Pythia8::Event& event) {
24  Pythia8::Vec4 pMediator, pDark;
25 
26  // Find the mediator in the event
27  for (int i = 0; i < event.size(); ++i) {
28  //mediator w/ distinct daughters = last copy (decayed)
29  if (event[i].id() == idMediator_ && event[i].daughter1() != event[i].daughter2() && event[i].daughter1() > 0 &&
30  event[i].daughter2() > 0) {
31  pMediator = event[i].p();
32  // undo mediator decay
33  event[i].undoDecay();
34 
35  // Generate the shower, output are 4 vectors in the rest frame of the shower, adding energy here avoids issues if scalar is off-shell
36  std::vector<Pythia8::Vec4> suep_shower_fourmomenta = suep_shower_->generateShower(pMediator.mCalc());
37  // Loop over hidden sector mesons and append to the event
38  int firstDaughter = event.size();
39  for (auto& pDark : suep_shower_fourmomenta) {
40  // Boost to the lab frame, i.e. apply the mediator boost
41  pDark.bst(pMediator);
42  // Append particle to the event w/ hidden meson pdg code. Magic number 91 means it is produced as a normal decay product
43  event.append(idDark_, 91, i, 0, 0, 0, 0, 0, pDark.px(), pDark.py(), pDark.pz(), pDark.e(), mDark_);
44  }
45 
46  // Change the status code of the mediator to reflect that it has decayed.
47  event[i].statusNeg();
48 
49  //set daughters of the mediator: daughter1 < daughter2 > 0 -> the particle has a range of decay products from daughter1 to daughter2
50  event[i].daughters(firstDaughter, event.size() - 1);
51  break;
52  }
53  }
54 
55  //allow event to continue
56  return false;
57 }
bool initAfterBeams() override
Definition: SuepDecay.cc:8
SuepDecay(const edm::ParameterSet &iConfig)
Definition: SuepDecay.cc:3
int idDark_
Definition: SuepDecay.h:23
ExtVec< T, 4 > Vec4
Definition: ExtVec.h:60
float temperature_
Definition: SuepDecay.h:24
std::unique_ptr< SuepShower > suep_shower_
Definition: SuepDecay.h:25
float mDark_
Definition: SuepDecay.h:24
bool doVetoProcessLevel(Pythia8::Event &event) override
Definition: SuepDecay.cc:23
int idMediator_
Definition: SuepDecay.h:23
Definition: event.py:1