HackAnalysis  2
General_Pythia.h
1 #pragma once
2 
3 /*
4 
5 General routines for building your pythia event fill function
6 
7 
8 
9 */
10 
11 
12 HEP::P4 PVtoH(const Pythia8::Vec4 &pyv)
13 {
14  HEP::P4 newv(pyv.px(), pyv.py(), pyv.pz(), pyv.e());
15  return newv;
16  //return HEP::P4(double pyv.px(), double pyv.py(), double pyv.pz(), double pyv.E());
17 }
18 
19 
20 Pythia8::Vec4 getMissingMomentum(Pythia8::Event &event)
21 {
22  Pythia8::Vec4 missingETvec;
23  std::vector<int> metaStableMothers;
24  //First check for all final particles:
25  for (int i = 0; i < event.size(); ++i) {
26  Pythia8::Particle ptc = event[i];
27  if (!ptc.isFinal()) continue;
28  //If particle has been produced inside or after the calorimeter, skip it and include its parent
29  // These for ATLAS
30  //if (abs(ptc.vProd().pz()) > 6100. || ptc.vProd().pT() > 2300.){
31  // These for CMS
32  //if (abs(ptc.vProd().pz()) > 5500. || ptc.vProd().pT() > 2950.){
33  // CMS muon chamber. Assuming here something invisible produced by something charged
34  if (fabs(ptc.vProd().pz()) > 11000. || ptc.vProd().pT() > 7000.){
35  int mom = ptc.mother1();
36  //cout << "Adding mom " << mom << " from " << ptc.name() << " and index " << i << endl;
37  if(std::find(metaStableMothers.begin(), metaStableMothers.end(), mom) == metaStableMothers.end()) {
38  metaStableMothers.push_back(mom);
39  }
40  continue;
41  }
42 
43 
44  //Neutrinos are included in MET:
45  if (ptc.idAbs() == 12 || ptc.idAbs() == 14 || ptc.idAbs() == 16){
46  //cout << "Including: " << i << " " << ptc.name() << endl;
47  missingETvec += ptc.p();
48  continue;
49  }
50  if (ptc.isHadron()) continue; //Ignore usual hadrons
51  if (ptc.idAbs() == 22) continue; //Ignore photons
52  if (ptc.chargeType() !=0) continue; // Ignore anything charged!!
53  if (ptc.isVisible() && ptc.m() < 10.) continue; //Ignore light visible particles
54 
55  //Add to missing momentum
56  //cout << "Including: " << i << " " << ptc.name() << endl;
57  missingETvec += ptc.p();
58  }
59 
60  //Now check all meta-stable mothers:
61  for (int imom = 0; imom < metaStableMothers.size(); ++imom) {
62  Pythia8::Particle ptc = event[metaStableMothers[imom]];
63 
64 
65  //Neutrinos are included in MET:
66  if (ptc.idAbs() == 12 || ptc.idAbs() == 14 || ptc.idAbs() == 16){
67  //cout << "Including mom: " << metaStableMothers[imom] << " " << ptc.id() << endl;
68  missingETvec += ptc.p();
69  continue;
70  }
71  if (ptc.isHadron()) continue; //Ignore usual hadrons
72  if (ptc.idAbs() == 22) continue; //Ignore photons
73  if (ptc.chargeType() !=0) continue; // Ignore anything charged!!
74  if (ptc.isVisible() && ptc.m() < 10.) continue; //Ignore light visible particles
75 
76  //Add to missing momentum
77  //cout << "Including mom: " << metaStableMothers[imom] << " " << ptc.id() << endl;
78  missingETvec += ptc.p();
79  }
80 
81  return missingETvec;
82 
83 
84 }
85 
86 void simplecalcMissingMomentum(HEP::Event &OutEvt, std::mt19937 &engine)
87 {
88  // look at the invisibles. Bit stupid I know but whatever.
89  HEP::P4 pmiss;
90  pmiss.clear();
91  for(HEP::Particle* invisible : OutEvt.invisible_particles())
92  {
93  pmiss+=invisible->mom();
94  }
95 
96  OutEvt.set_missingmom(pmiss);
97 
98  // smear here ...
99 }
100 
101 
102 
103 bool isLLP(Pythia8::Particle *particle, const Pythia8::Event &event)
104 {
105  // if (particle->colType() != 0) return false;
106  if (particle->isHadron()) return false;
107 
108 
109  //Skip light particles
110  if (fabs(particle->m()) < 20.) return false;
111 
112  if (particle->isFinal()) return true;
113 
114  Pythia8::Vec4 decayVertex = particle->vDec();
115  if(decayVertex.pAbs() < 0.001) return false;
116  //Check if it is not an intermediate state (has itself as daughter)
117  std::vector<int> daughters = particle->daughterListRecursive();
118  for (int i = 0; i < daughters.size(); ++i){
119  int idaughter = daughters[i];
120  Pythia8::Particle daughterParticle = event[idaughter];
121  if (daughterParticle.id() == particle->id()) return false;
122  }
123 
124  return true;
125 
126 }
127 
128 bool isHSCP(Pythia8::Particle *particle, const Pythia8::Event &event)
129 //Return True (False) if the particle is (is not) a HSCP
130 {
131 
132  //Skip neutral particles
133  if (!particle->isCharged()) return false;
134  //Skip color charged particles or hadrons
135  if (particle->colType() != 0) return false;
136  if (particle->isHadron()) return false;
137  //Skip light particles
138  if (fabs(particle->m()) < 20.) return false;
139  //Check if it is not an intermediate state (has itself as daughter)
140  std::vector<int> daughters = particle->daughterListRecursive();
141  for (int i = 0; i < daughters.size(); ++i){
142  int idaughter = daughters[i];
143  Pythia8::Particle daughterParticle = event[idaughter];
144  if (daughterParticle.id() == particle->id()) return false;
145  }
146 
147  return true;
148 }
149 
150 
151 bool decayBeforeEndHcal(Pythia8::Particle &particle)
152 //Return True (False) if the particle decays before (after) the hadronic calorimeter
153 {
154  if (particle.isFinal()){return false;}
155  Pythia8::Vec4 decayVertex = particle.vDec();
156  if (decayVertex.pT() < 3.9e3 && abs(decayVertex.pz()) < 6.1e3){return true;}
157 
158  return false;
159 }
160 
161 bool decayInsideAtlas(Pythia8::Particle &particle)
162 //Return True (False) if the particle decays inside the detector
163 {
164  if (particle.isFinal()){return false;}
165  Pythia8::Vec4 decayVertex = particle.vDec();
166  if (decayVertex.pT() < 12e3 && abs(decayVertex.pz()) < 23e3){return true;}
167 
168  return false;
169 }
170 
171 bool ischargedTrackHadron(Pythia8::Particle* particle, const Pythia8::Event &event)
172 {
173  if(!particle->isHadron()) {return false;}
174  if(!particle->isCharged()) { return false;}
175 
176  Pythia8::Vec4 prodVertex = particle->vProd();
177  //if( (prodVertex.pT() > 1100.) || (abs(prodVertex.pz()) > 3100.)) {return false;}
178  if( (prodVertex.pT() > 200.) || (abs(prodVertex.pz()) > 270.)) {return false;}
179 
180  if(particle->isFinal())
181  {
182  return true;
183  }
184  else
185  {
186  std::vector<int> daughters = particle->daughterListRecursive();
187  for (int i = 0; i < daughters.size(); i++){
188  int idaughter = daughters[i];
189  const Pythia8::Particle* daughterParticle = &event[idaughter];
190  if (daughterParticle->id() == particle->id()) return false;
191  }
192 
193  if( (particle->vDec().pT() < 300.0) ){return false;}
194  }
195 
196  return true;
197 }
198 
199 
200 /*
201 
202 inline bool fromHadron(int n, const Pythia8::Event &evt) {
203  // Root particle is invalid
204  if (n == 0) return false;
205  const Pythia8::Particle p = evt[n];
206  if (p.isHadron()) return true;
207  if (p.isParton()) return false; // stop the walking at the end of the hadron level
208  for (int m : p.motherList()) {
209  if (fromHadron(m, evt)) return true;
210  }
211  return false;
212  }
213 */
214 
215 inline bool fromHadron(int n, const Pythia8::Event &evt) {
216  // Root particle is invalid
217  if (n == 0) return false;
218  const Pythia8::Particle* p = &evt[n];
219  if (p->isParton()) return false; // stop the walking at the end of the hadron level
220  if (p->isHadron()) return true;
221 
222 
223  for (int mm : p->motherList()) {
224  //cout << n << ": " << p->id() << " Mother!: " << mm << endl;
225  if (fromHadron(mm, evt)) return true;
226  }
227  return false;
228  }
229 
230 
231 /*
232 
233 bool registersinHCAL(Pythia8::Particle* prt)
234 {
235  Pythia8::Vec4 pos = prt->vProd();
236 
237 
238 }
239 
240 
241 */
242 
243 
244 // This function is not used.
245 HEP::PileupEvent* sortPileupEvent(Pythia8::Event &event)
246 {
247  HEP::PileupEvent* newPuE = new HEP::PileupEvent();
248 
249 
250  for(int npart = 0; npart< event.size(); npart++)
251  {
252  Pythia8::Particle* prt = &event[npart];
253  int pid=prt->id();
254  int apid=abs(pid);
255  // check if is final and produced inside tracker or if it decays outside enough (not yet)
256  //cout << pid << endl;
257  if(prt->isFinal())
258  {
259  HEP::Particle* newpart = new HEP::Particle(prt->px(),prt->py(),prt->pz(),prt->e(),pid);
260  // All particles other than invisibles and muons are jet constituents. Bollocks, really we only want photons
261  //if( pid == 12 || pid == 14 || pid == 16 || pid == 18 || pid == 13 || pid == 11)
262 
263 
264  if(prt->isHadron() && prt->isVisible())
265  {
266  if(prt->isCharged())
267  {
268  //OutEvt.add_charged_hadron(newpart->mom());
269  Pythia8::Vec4 pos = prt->vProd();
270  if (pos.pAbs() > 0.0001)
271  {
272  newpart->set_prod(HEP::P4(pos.px(),pos.py(),pos.pz(),pos.e()));
273  }
274  newpart->set_3Q(prt->chargeType());
275 
276  //newPuE->add_hadron_vertex(newpart);
277  if(newpart->pT() > 1.0)
278  {
279  newPuE->add_charged_hadron(newpart);
280  }
281  else
282  {
283  delete newpart;
284  }
285  }
286  else
287  {
288  newPuE->add_neutral_hadron(newpart->mom());
289  delete newpart;
290  }
291  }
292  else // not hadron
293  {
294 
295  // look at the vertices
296  //if ( (*p)->production_vertex() ) { // i.e. does it have one: we want the position
297  //HepMC::FourVector pos = (*p)->production_vertex()->position();
298  Pythia8::Vec4 pos = prt->vProd();
299 
300 
301  if (pos.pAbs() > 0.0001)
302  {
303  newpart->set_prod(HEP::P4(pos.px(),pos.py(),pos.pz(),pos.e()));
304  }
305  newpart->set_prompt();
306  newpart->set_3Q(prt->chargeType());
307  newPuE->add_particle(newpart);
308  // }
309 
310  }
311  // don't bother with jet constituents for this analysis
312  /*
313  if( pid == 12 || pid == 14 || pid == 16 || pid == 18 ) // include electrons and muons, not neutrinos
314  {
315  continue;
316  }
317  if((pid > 1000000) && (pid < 9000000))
318  {
319  continue;
320  }
321  //HepMC::FourVector v4((*p)->momentum());
322  Pythia8::Vec4 v4=prt.p();
323  fastjet::PseudoJet particleTemp(v4.px(),v4.py(),v4.pz(),v4.e());
324  fjInputs.push_back(particleTemp);
325  */
326  }
327  else
328  {
329  // cout << "is not final ";
330  if(ischargedTrackHadron(prt,event))
331  {
332 
333  Pythia8::Vec4 pos = prt->vProd();
334  Pythia8::Vec4 vdec = prt->vDec();
335  if(vdec.pT() < 33.0) continue; // check to make sure it doesn't decay too soon
336  HEP::Particle* newpart = new HEP::Particle(prt->px(),prt->py(),prt->pz(),prt->e(),pid);
337  if (pos.pAbs() > 0.0001)
338  {
339  newpart->set_prod(HEP::P4(pos.px(),pos.py(),pos.pz(),pos.e()));
340  }
341  newpart->set_decay(HEP::P4(vdec.px(),vdec.py(),vdec.pz(),vdec.e()));
342  newpart->set_3Q(prt->chargeType());
343  // cout << " adding charged hadron ";
344  //newPuE->add_hadron_vertex(newpart);
345  newPuE->add_charged_hadron(newpart);
346 
347 
348  }
349  }
350 
351 
352  }
353 
354  return newPuE;
355 
356 }
357 
358 /*
359 
360 NOONONONONO
361 
362 DEFINE HERE THE MAP FROM THE NAME OF FILL FUNCTIONS TO THE ACTUAL FUNCTION:
363 
364 std::map<std::string, std::function<int(int,int)>> funcMap =
365  {{ "add", add},
366  { "sub", sub}
367  };
368 
369 HEP::Event &OutEvt, Pythia8::Event &evt, double weight, int npileups, std::vector<HEP::PileupEvent*> &pileups, std::mt19937 &engine
370 
371 */
Simple event class, separating particles into classes.
Definition: heputil.h:227
void set_missingmom(const P4 &pmiss)
Set the missing momentum vector.
Definition: heputil.h:736
const std::vector< Particle * > & invisible_particles() const
Get invisible final state particles.
Definition: heputil.h:641
A 4-momentum class for vectors.
Definition: Vectors.h:45
void clear()
Set the components to zero.
Definition: Vectors.h:82
Definition: Particle.h:24
const P4 & mom() const
Get the 4 vector.
Definition: Particle.h:109
void set_prompt(bool isprompt=true)
Set promptness.
Definition: Particle.h:178
double pT() const
Get the squared transverse momentum.
Definition: Particle.h:165
Definition: heputil.h:98