ATLAS Offline Software
Loading...
Searching...
No Matches
Pythia8B_i.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5// ======================================================================
6// Pythia8B_i
7// James.Catmore@cern.ch, Eva.Bouhova@cern.ch, Maria.Smizanska@cern.ch (July 2011)
8// Inherits from Pythia8_i (jmonk@hep.ucl.ac.uk)
9// Implementation of repeated hadronization of b-quark events in Pythia8
10// Based on pythia8150/examples/main15.cc in the Pythia8 distribution,
11// T. Sjostrand et al,
12// http://home.thep.lu.se/~torbjorn/pythiaaux/present.html
13// ======================================================================
17#include "AtlasHepMC/GenEvent.h"
18
19// calls to fortran routines
21
22// Pointer to random engine
23CLHEP::HepRandomEngine* Pythia8B_i::p_rndmEngine = nullptr;
24
26// User properties not defined by Pythia8_i
28(const std::string &name, ISvcLocator *pSvcLocator): Pythia8_i(name,pSvcLocator) {
29 declareProperty("NHadronizationLoops", m_had=1);
30 declareProperty("NDecayLoops", m_dec=1);
31 declareProperty("SelectBQuarks",m_selectBQuarks=true);
32 declareProperty("SelectCQuarks",m_selectCQuarks=false);
33 declareProperty("QuarkPtCut", m_qPtCut= 0.0);
34 declareProperty("AntiQuarkPtCut", m_aqPtCut= 0.0);
35 declareProperty("QuarkEtaCut", m_qEtaCut=102.5);
36 declareProperty("AntiQuarkEtaCut", m_aqEtaCut=102.5);
37 declareProperty("RequireBothQuarksPassCuts", m_and=false);
38 declareProperty("VetoDoubleBEvents", m_vetoDoubleB=true);
39 declareProperty("VetoDoubleCEvents", m_vetoDoubleC=true);
40 declareProperty("BPDGCodes", m_bcodes);
41 declareProperty("TriggerPDGCode", m_trigCode=0); // 0 = no selection applied
42 declareProperty("TriggerStatePtCut", m_trigPtCut);
43 declareProperty("TriggerStateEtaCut", m_trigEtaCut=2.5);
44 declareProperty("MinimumCountPerCut", m_cutCount);
45 declareProperty("TriggerInvariantMassCut",m_invMass=0.0);
46 declareProperty("TriggerOppositeCharges",m_oppCharges=false);
47 declareProperty("SignalPDGCodes", m_sigCodes);
48 declareProperty("SignalPtCuts", m_sigPtCuts);
49 declareProperty("SignalEtaCuts", m_sigEtaCuts);
50 declareProperty("NumberOfSignalsRequiredPerEvent", m_nSignalRequired=1);
51 declareProperty("UserSelection", m_userString="NONE");
52 declareProperty("UserSelectionVariables", m_userVar);
53 declareProperty("SuppressSmallPT", m_doSuppressSmallPT=false);
54
55 m_totalBQuark = 0;
57 m_totalCQuark = 0;
60 m_totalHard = 0;
61 m_totalClone = 0;
65 m_speciesCount.clear();
68
69 for (std::vector<int>::iterator iit=m_bcodes.begin(); iit!=m_bcodes.end(); ++iit) {
70 m_speciesCount[*iit] = 0;
71 }
72#ifdef HEPMC3
73 m_runinfo = std::make_shared<HepMC3::GenRunInfo>();
74 std::vector<std::string> names = {"Default"};
75 m_runinfo->set_weight_names(names);
76#endif
77}
78
82
83
85
86 // Logic checks
87 unsigned int trigPtCutsSize = m_trigPtCut.size();
88 if (trigPtCutsSize!=m_cutCount.size()) {
89 ATH_MSG_ERROR("You are requesting " << trigPtCutsSize << " trigger-like pt cuts but are providing required counts for " << m_cutCount.size() << " cuts. This doesn't make sense.");
90 return StatusCode::FAILURE;
91 }
92
93 // This over-rides the genInitialize in the base class Pythia8_i, but then calls it
94 // Sets the built-in UserHook called SuppressLowPT
95 // FOR ONIA USE ONLY
96 ATH_MSG_INFO("genInitialize() from Pythia8B_i");
97
98 bool canSetHook=true;
100 Pythia8_i::m_userHooks=std::vector<std::string>(1, "SuppressSmallPT");
101 }
102
103 if (!canSetHook) {
104 ATH_MSG_ERROR(" *** Unable to initialise PythiaB !! ***");
105 return StatusCode::FAILURE;
106 }
107
108 if(m_userString == "NONE") m_userString.clear();
109 // Call the base class genInitialize()
111
112 ATH_MSG_DEBUG("... seeding Athena random number generator");
113 if (m_useRndmGenSvc){
114 p_rndmEngine = m_atlasRndmEngine->getEngine(); // NOT THREAD-SAFE
115 if (!p_rndmEngine) {
116 ATH_MSG_FATAL("Unable to retrieve HepRandomEngine. Bailing out.");
117 return StatusCode::FAILURE;
118 }
119 }
120 return StatusCode::SUCCESS;
121}
122
123
126
127 // Just counter setting
128 ATH_MSG_INFO("Pythia8B_i from genuserInitialize()");
129
130 // Initialize global counters
131 m_totalBQuark = 0;
133 m_totalCQuark = 0;
136 m_totalHard = 0;
137 m_totalClone = 0;
138 m_atLeastOneAcc = 0;
141 m_speciesCount.clear();
142 for (std::vector<int>::iterator iit=m_bcodes.begin(); iit!=m_bcodes.end(); ++iit) {
143 m_speciesCount[*iit] = 0;
144 }
145
146 return StatusCode::SUCCESS;
147}
148
151
152 ATH_MSG_DEBUG(">>> Pythia8B_i from callGenerator");
153
154 // Check the buffers are empty - if not return so that the GenModule
155 // calls fillEvt again and takes another event from the buffer
156 if (m_BEventBuffer.size()!=0) {
157 ATH_MSG_DEBUG("BEventBuffer not empty; skipping callGenerator");
158 return StatusCode::SUCCESS;
159 }
160
162 //Re-seed the random number stream
163 long seeds[7];
164 const EventContext& ctx = Gaudi::Hive::currentContext();
165 ATHRNG::calculateSeedsMC21(seeds, "PYTHIA8B", ctx.eventID().event_number(), m_dsid, m_randomSeed);
166 m_atlasRndmEngine->getEngine()->setSeeds(seeds, 0); // NOT THREAD-SAFE
167 m_seeds.clear();
168 m_seeds.push_back(seeds[0]);
169 m_seeds.push_back(seeds[1]);
170 }
171
172 // Generator. Shorthand for event.
173 Pythia8::Event& event = Pythia8_i::m_pythia->event;
174
175 // Begin event loop.
176 int iEvent(0);
177 while (iEvent < 1) { // keep going until an event is accepted
178
179 ATH_MSG_DEBUG("Throwing the dice....");
180// if (!Pythia8_i::m_pythia->next()) continue;
181 if ( !Pythia8_i::m_pythia->next() ) {
182 // First check if it failed because it ran out of events.
183 if ( Pythia8_i::m_pythia->info.atEndOfFile() ) {
184 return StatusCode::FAILURE;
185 }
186
187 // Otherwise, just make sure that it's not failing too many times in a row.
189 if ( m_failureCount >= m_maxFailures ) {
190 ATH_MSG_ERROR("Exceeded the max number of consecutive event failures.");
191 return StatusCode::FAILURE;
192 } else {
193 ATH_MSG_INFO("Event generation failed - re-trying.");
194 continue;
195 }
196 }
197 // Reset failure counter and increment event counter after a successful event
198 m_failureCount = 0;
200
201 // Find b(c)/antib(c) quarks and enforce cuts as required
202 int nbBeforeSelection(0);
203 int nbbarBeforeSelection(0);
204 int ncBeforeSelection(0);
205 int ncbarBeforeSelection(0);
206 int nbQuark(0);
207 int nbbarQuark(0);
208 int ncQuark(0);
209 int ncbarQuark(0);
210 int stat;
211 for (int i = 0; i < event.size(); ++i) {
212 stat = event[i].statusAbs();
213 bool isBQuark(false);
214 bool isCQuark(false);
215 bool isAntiBQuark(false);
216 bool isAntiCQuark(false);
217 // b-quarks (=anti-B hadrons)
218 if (event[i].id() == 5 && (stat == 62 || stat == 63)) {isBQuark=true; ++nbBeforeSelection;}
219 if (event[i].id() == 4 && (stat == 62 || stat == 63)) {isCQuark=true; ++ncBeforeSelection;}
220 if ( (isBQuark && m_selectBQuarks) || ((isCQuark && m_selectCQuarks)) ) {
221 bool passesPtCut(false); bool passesEtaCut(false);
222 std::string accString = " : REJECTED";
223 double qpt = event[i].pT(); double qeta = std::abs(event[i].eta());
224 if (qpt>m_qPtCut) passesPtCut=true;
225 if (qeta<m_qEtaCut) passesEtaCut=true;
226 if (passesPtCut && passesEtaCut) {
227 if (isBQuark) ++nbQuark;
228 if (isCQuark) ++ncQuark;
229 accString = " : ACCEPTED";
230 }
231 if (isBQuark) ATH_MSG_DEBUG("bQuark pt/eta = " << qpt << "/" << qeta << accString);
232 if (isCQuark) ATH_MSG_DEBUG("cQuark pt/eta = " << qpt << "/" << qeta << accString);
233 }
234 // anti-b quarks (=B-hadrons)
235 if (event[i].id() == -5 && (stat == 62 || stat == 63)) {isAntiBQuark=true; ++nbbarBeforeSelection;}
236 if (event[i].id() == -4 && (stat == 62 || stat == 63)) {isAntiCQuark=true; ++ncbarBeforeSelection;}
237 if ( (isAntiBQuark && m_selectBQuarks) || ((isAntiCQuark && m_selectCQuarks)) ) {
238 bool passesPtCut(false); bool passesEtaCut(false);
239 std::string accString = " : REJECTED";
240 double aqpt = event[i].pT(); double aqeta = std::abs(event[i].eta());
241 if (aqpt>m_aqPtCut) passesPtCut=true;
242 if (aqeta<m_aqEtaCut) passesEtaCut=true;
243 if (passesPtCut && passesEtaCut) {
244 if (isAntiBQuark) ++nbbarQuark;
245 if (isAntiCQuark) ++ncbarQuark;
246 accString = " : ACCEPTED";
247 }
248 if (isAntiBQuark) ATH_MSG_DEBUG("bbarQuark pt/eta = " << aqpt << "/" << aqeta << accString);
249 if (isAntiCQuark) ATH_MSG_DEBUG("ccbarQuark pt/eta = " << aqpt << "/" << aqeta << accString);
250 }
251 }
252 if (nbBeforeSelection+nbbarBeforeSelection>=4 && m_vetoDoubleB) {
253 ATH_MSG_DEBUG("At user request, rejecting double b-bbar event; throwing dice again");
254 continue;
255 }
256 if (ncBeforeSelection+ncbarBeforeSelection>=4 && m_vetoDoubleC) {
257 ATH_MSG_DEBUG("At user request, rejecting double c-cbar event; throwing dice again");
258 continue;
259 }
260 bool rejectBBbar(false);
261 bool rejectCCbar(false);
262 // Enforce user selections
263 if (!m_and && (nbbarQuark==0 && nbQuark==0) && m_selectBQuarks) rejectBBbar=true;
264 if (m_and && (nbbarQuark==0 || nbQuark==0) && m_selectBQuarks) rejectBBbar=true;
265 if (!m_and && (ncbarQuark==0 && ncQuark==0) && m_selectCQuarks) rejectCCbar=true;
266 if (m_and && (ncbarQuark==0 || ncQuark==0) && m_selectCQuarks) rejectCCbar=true;
267 if (m_selectBQuarks && m_selectCQuarks && rejectBBbar && rejectCCbar) {
268 ATH_MSG_DEBUG("No b- or c- quarks accepted; throwing the dice again");
269 continue;
270 }
271 if (m_selectBQuarks && !m_selectCQuarks && rejectBBbar) {
272 ATH_MSG_DEBUG("No b-quarks accepted; throwing the dice again");
273 continue;
274 }
275 if (!m_selectBQuarks && m_selectCQuarks && rejectCCbar) {
276 ATH_MSG_DEBUG("No c-quarks accepted; throwing the dice again");
277 continue;
278 }
279 m_totalBQuark += nbQuark;
280 m_totalBBarQuark += nbbarQuark;
281 m_totalCQuark += ncQuark;
282 m_totalCBarQuark += ncbarQuark;
283 ++m_totalHard;
284
285 // Switch off decays of species of interest if requested (only for repeated decays)
286 bool doRepeatedDecays(false); if (m_dec>1) doRepeatedDecays=true;
287 if (doRepeatedDecays) {
288 for (unsigned int iC = 0; iC < m_bcodes.size(); ++iC) Pythia8_i::m_pythia->particleData.mayDecay( m_bcodes[iC], false);
289 }
290
291 // REPEATED HADRONIZATION LOOP
292 // Make a spare copy of event
293 Pythia8::Event eventCopy;
294 std::vector<Pythia8::Event> repeatHadronizedEvents;
295 std::vector<Pythia8::Event>::iterator eventItr;
296 eventCopy = event;
297 // Begin loop over rounds of hadronization for same event.
298 for (unsigned int iRepeat = 0; iRepeat < m_had; ++iRepeat) {
299 ++m_totalClone;
300 // Repeated hadronization: restore saved event record.
301 if (iRepeat > 0) event = eventCopy;
302 // Repeated hadronization: do HadronLevel (repeatedly).
303 if (!Pythia8_i::m_pythia->forceHadronLevel()) continue;
304 repeatHadronizedEvents.push_back(event); // save the events for selections
305 }
306
307 // REPEATED DECAY LOOP
308 std::vector<Pythia8::Event> savedEvents;
309 if (doRepeatedDecays) {
310 // switch decays back on
311 for (unsigned int iC = 0; iC < m_bcodes.size(); ++iC) Pythia8_i::m_pythia->particleData.mayDecay( m_bcodes[iC], true);
312 // Loop over repeatedly hadronized events
313 for (eventItr=repeatHadronizedEvents.begin(); eventItr!=repeatHadronizedEvents.end(); ++eventItr) {
314 for (unsigned int iRepeat = 0; iRepeat < m_dec; ++iRepeat) {
315 event = (*eventItr);
316 if (!Pythia8_i::m_pythia->moreDecays()) break;
317 savedEvents.push_back(event);
318 }
319 }
320 }
321 if (!doRepeatedDecays) savedEvents = std::move(repeatHadronizedEvents);
322
323 //event.list();
324
325 // Make final selections on savedEvents
326 std::vector<Pythia8::Event> finalEvents;
327 bool signalSelect(false);
328 if (m_sigCodes.size()>0) signalSelect=true;
329 for (eventItr=savedEvents.begin(); eventItr!=savedEvents.end(); ++eventItr) {
330 // Does event contain basic leptons needed to fire triggers?
333 // Does the events contain a signal decay as required by the user
334 if (signalSelect) {
335 if (!cleanUndecayed(*eventItr,m_bcodes)) continue;
336 if (!signalAccept(*eventItr,m_sigCodes,m_nSignalRequired)) continue;
337 }
338 // User-defined selections if any
339 if (!m_userString.empty()) {
340 if (!userSelection(*eventItr,m_userString,m_userVar)) continue;
341 }
342 finalEvents.push_back(*eventItr);
343 ++iEvent;
344 }
345 savedEvents.clear();
346 if (finalEvents.size()>0) ++m_atLeastOneAcc;
347
348 // Species counting and preparing vector for filling
349 for (eventItr=finalEvents.begin(); eventItr!=finalEvents.end(); ++eventItr) {
350 for (int i = 0; i < (*eventItr).size(); ++i) {
351 int id = (*eventItr)[i].id();
352 std::map<int,int>::iterator it;
353 it = m_speciesCount.find(id);
354 if (it != m_speciesCount.end()) {
355 it->second++;
356 }
357 }
358 m_BEventBuffer.push_back(*eventItr);
361 //(*eventItr).list();
362 }
363
364 // End of event loop.
365 }
366
367 StatusCode returnCode = StatusCode::SUCCESS;
368
369 return returnCode;
370}
371
373StatusCode Pythia8B_i::fillEvt(HepMC::GenEvent *evt){
374
375 ATH_MSG_DEBUG(">>> Pythia8B_i from fillEvt");
376 ATH_MSG_DEBUG("BEventBuffer contains " << m_BEventBuffer.size() << " events ");
377 if (m_BEventBuffer.size()==0 ) {
378 ATH_MSG_DEBUG("BEventBuffer now empty - going to next Pythia event");
379 return StatusCode::SUCCESS;
380 }
381
382 Pythia8::Event &pyev = *(m_BEventBuffer.begin());
383 evt->set_event_number(*(m_internalEventNumbers.begin()));
384 m_pythiaToHepMC.fill_next_event(pyev, evt, 1, &Pythia8_i::m_pythia->info);
385
386
387 // set the randomseeds
390
391 // set the event weight
392#ifdef HEPMC3
393 if (!evt->run_info()) evt->set_run_info(m_runinfo);
394 evt->set_units(HepMC3::Units::MEV, HepMC3::Units::MM);
395#else
396 evt->weights().push_back(m_pythia->info.weight());
397#endif
398
399//uncomment to list HepMC events
400//#ifdef HEPMC3
401// std::cout << " print::listing Pythia8B " << std::endl;
402// HepMC3::Print::listing(std::cout, *evt);
403//#else
404// std::cout << " print::printing Pythia8B " << std::endl;
405// evt->print();
406//#endif
407
408 // Remove event/number from buffer
409 m_BEventBuffer.erase(m_BEventBuffer.begin());
411
412
413 return StatusCode::SUCCESS;
414}
415
416
419
420 ATH_MSG_INFO(">>> Pythia8B_i from genFinalize");
421
422 Pythia8_i::m_pythia->stat();
423
424 Pythia8::Info info = Pythia8_i::m_pythia->info;
425 double xs = info.sigmaGen();// in mb
426 xs *= 1000. * 1000.;//convert to nb
427
428 std::cout << "MetaData: cross-section (nb)= " << xs << std::endl;
429 std::cout << "Total events from Pythia = " << m_totalPythiaCalls<< std::endl;
430 std::cout << "Number of accepted b quarks = " << m_totalBQuark<< std::endl;
431 std::cout << "Number of accepted bbar quarks = " << m_totalBBarQuark<< std::endl;
432 std::cout << "Number of accepted c quarks = " << m_totalCQuark<< std::endl;
433 std::cout << "Number of accepted cbar quarks = " << m_totalCBarQuark<< std::endl;
434 std::cout << "Number of accepted b/c events before hadronization = " << m_totalHard<< std::endl;
435 std::cout << "Number of hadronization loops per b/c-event = " << m_had<< std::endl;
436 std::cout << "Total number of hadronization loops = " << m_totalClone<< std::endl;
437 std::cout << "Number of accepted b/c events yielding at least one finally accepted event = " << m_atLeastOneAcc<< std::endl;
438 std::cout << "Average number of accepted events per hard event = " << static_cast<float>(m_internal_event_number)/static_cast<float>(m_atLeastOneAcc)<< std::endl;
439 std::cout << "Number of events passing trigger cuts = " << m_passingTriggerCuts<< std::endl;
440 std::cout << "Number of accepted events = " << m_internal_event_number<< std::endl;
441 std::cout << "Summary of cuts: " << std::endl;
443 if (m_selectBQuarks) std::cout << "Quarks cuts apply to b-quarks" << std::endl;
444 if (m_selectCQuarks) std::cout << "Quarks cuts apply to c-quarks" << std::endl;
445 std::cout << "Quark pt > " << m_qPtCut << std::endl;
446 std::cout << "Antiquark pt > " << m_aqPtCut << std::endl;
447 std::cout << "Quark eta < " << m_qEtaCut << std::endl;
448 std::cout << "Antiquark eta < " << m_aqEtaCut << std::endl;
449 if (m_and) std::cout << "Require quark and anti-quark pass cuts" << std::endl;
450 }
451 std::cout << "Trigger lepton type = " << m_trigCode << std::endl;
452 std::cout << "Trigger lepton pt cuts: ";
453 for (unsigned int prCntr=0; prCntr<m_trigPtCut.size(); ++prCntr) {std::cout << m_trigPtCut[prCntr] << " ";}
454 std::cout << std::endl;
455 std::cout << "Trigger lepton eta cut: " << m_trigEtaCut << std::endl;
456 std::cout << "Required number of leptons passing each trigger cut = ";
457 for (unsigned int prCntr=0; prCntr<m_cutCount.size(); ++prCntr) {std::cout << m_cutCount[prCntr] << " ";}
458 std::cout << std::endl;
459 std::cout << "Invariant mass of trigger leptons > " << m_invMass << std::endl;
460 if (m_oppCharges) std::cout << "Trigger leptons required to have opposite charges" << std::endl;
462 if (!m_userString.empty()) std::cout << "User selection: " << m_userString << std::endl;
463
464 if (m_speciesCount.size()>0) {
465 std::cout <<"\nSpecies\tCount\n"<< std::endl;
466 for ( std::map<int,int>::iterator iter = m_speciesCount.begin();
467 iter != m_speciesCount.end(); ++iter )
468 std::cout << iter->first << '\t' << iter->second << '\n'<< std::endl;
469 }
470 if (m_dec>1) std::cout << "Number of times each of these states were copied and decayed: " << m_dec << std::endl;
471
472 float cloningFactor = (float)(m_internal_event_number)/(float)(m_atLeastOneAcc);
473 float finalXS = ((float)xs*((float)m_internal_event_number/(float)m_totalPythiaCalls))/(float)m_had;
474
475 std::cout << " I===================================================================================== " << std::endl;
476 std::cout << " I CROSSSECTION OF YOUR B-CHANNEL IS I " << std::endl;
477 std::cout << " I BX= PX*NB/AC/MHAD= I " << finalXS << " nb" << std::endl;
478 std::cout << " I I " << std::endl;
479 std::cout << " I IN CASE YOU FORCED ANY DECAY YOU SHOULD I " << std::endl;
480 std::cout << " I CORRECT CROSS SECTION BX FURTHER, MULTIPLYING I " << std::endl;
481 std::cout << " I BX BY BRANCHING RATIO(S) OF YOUR FORCED I " << std::endl;
482 std::cout << " I DECAY(S) AND BY A FACTOR OF 2 FOR SYMMETRY I " << std::endl;
483 std::cout << " I I " << std::endl;
484 std::cout << " I MORE DETAILS ON CROSS SECTION I " << std::endl;
485 std::cout << " I PYTHIA CROSS SECTION IS PX= I " << xs <<" nb" << std::endl;
486 std::cout << " I NUMBER OF ACCEPTED HARD EVENTS AC= I " << m_totalPythiaCalls << std::endl;
487 std::cout << " I NUMBER OF ACCEPTED B-EVENTS IS NB= I " << m_internal_event_number<< std::endl;
488 std::cout << " I REPEATED HADRONIZATIONS IN EACH EVENT MHAD= I " << m_had << std::endl;
489 std::cout << " I AVERAGE NUM OF ACCEPTED EVTS IN HADRONIZATION LOOP I " << cloningFactor << std::endl;
490 std::cout << " I IN CASE YOU FORCED ANY DECAY YOU SHOULD I " << std::endl;
491 std::cout << " I CORRECT CROSS SECTION BX FURTHER, MULTIPLYING I " << std::endl;
492 std::cout << " I BX BY BRANCHING RATIO(S) OF YOUR FORCED I " << std::endl;
493 std::cout << " I DECAY(S) AND BY A FACTOR OF 2 FOR SYMMETRY I " << std::endl;
494 std::cout << " I I " << std::endl;
495 std::cout << " I===================================================================================== " << std::endl;
496 std::cout << "" << std::endl;
497 std::cout << "MetaData: cross-section (nb)= " << finalXS << std::endl;
498
499 return StatusCode::SUCCESS;
500}
501
503// Enforce trigger-like selections on final state leptons
504// User can require pt and eta cuts on a given number of leptons.
505// User can enforce at least one oppositely charged pair
506// User can enforce a minimum invariant mass on a pair
507
508bool Pythia8B_i::leptonSelect(Pythia8::Event &theEvent, const std::vector<double> &ptCut, double etaCut, const std::vector<int> &counts, int type_id, double massCut, bool opposite) {
509
510 if (type_id==0) return true;
511 bool passed(false);
512 std::string accString(" : REJECTED");
513
514 std::vector<int> leptonIDs;
515 int nCuts=ptCut.size();
516 std::vector<int> countGood(nCuts, 0);
517
518 for (int i = 0; i<theEvent.size(); ++i) {
519 const Pythia8::Particle &theParticle = theEvent[i];
520 int id = theParticle.idAbs();
521 if ( id == type_id ) { // find lepton flavour requested by user
522 double pt = theParticle.pT();
523 double eta = theParticle.eta();
524 ATH_MSG_DEBUG("Lepton of type " << id << " with pt/eta " << pt << "/" << eta);
525 for (int cutCntr=0; cutCntr<nCuts; ++cutCntr) {
526 if ( (pt>ptCut[cutCntr]) && (std::abs(eta)<etaCut) ) {
527 countGood[cutCntr] += 1;
528 leptonIDs.push_back(i); // remember leptons
529 }
530 }
531 }
532 // can the loop stop?
533 int nPassed(0);
534 for (int cutCntr=0; cutCntr<nCuts; ++cutCntr) {
535 if (countGood[cutCntr] >= counts[cutCntr]) ++nPassed;
536 }
537 if (nPassed==nCuts) break;
538 } // end of loop over particles
539
540 // Check the accumulated counts
541 int nPassed(0);
542 for (int cutCntr=0; cutCntr<nCuts; ++cutCntr) {
543 if (countGood[cutCntr] >= counts[cutCntr]) ++nPassed;
544 }
545 if (nPassed >= nCuts && nCuts==1) {accString=" : ACCEPTED"; passed = true;} // only 1 lepton required so no pair properties
546 if (nPassed >= nCuts && nCuts>1 && pairProperties(theEvent,leptonIDs,massCut,opposite)) {accString=" : ACCEPTED"; passed = true;}
547 ATH_MSG_DEBUG("Trigger-like selection" << accString);
548 return passed;
549}
550
552// Require that an event contains no undecayed particles of the types listed by
553// user
554bool Pythia8B_i::cleanUndecayed(Pythia8::Event &theEvent, const std::vector<int> &bCodes) {
555
556 bool cleanEvent(true);
557 std::string accString(" : ACCEPTED");
558 for (int i = 0; i<theEvent.size(); ++i) {
559 const Pythia8::Particle &theParticle = theEvent[i];
560 int id = theParticle.id();
561 int status = theParticle.status();
562 for (auto iItr = bCodes.begin(); iItr!=bCodes.end(); ++iItr) {
563 if ( (id == *iItr) && (status>0) ) {accString=" : REJECTED"; cleanEvent = false;}
564 }
565 }
566 ATH_MSG_DEBUG("Check event for undecayed signal particles" << accString);
567
568 return cleanEvent;
569
570}
571
573// Requires that a set of leptons represented by (index numbers in event record)
574// contains at least one pair with
575// - invariant mass above massCut
576// - opposite charges if user requires
577
578bool Pythia8B_i::pairProperties(Pythia8::Event &theEvent, const std::vector<int> &leptonIDs, double massCut, bool opposite) {
579
580 bool passesCuts(false);
581 std::string accString(" : REJECTED");
582 for (auto iit = leptonIDs.begin(); iit!=leptonIDs.end(); ++iit) {
583 for (auto iit2 = iit+1; iit2!=leptonIDs.end(); ++iit2) {
584 int q1=theEvent[*iit].charge();
585 int q2=theEvent[*iit2].charge();
586 if (opposite && (q1*q2>0)) continue;
587 double px1=theEvent[*iit].px();
588 double py1=theEvent[*iit].py();
589 double pz1=theEvent[*iit].pz();
590 double mass1=theEvent[*iit].mSel();
591 double e1=std::sqrt(px1*px1+py1*py1+pz1*pz1+mass1*mass1);
592 double px2=theEvent[*iit2].px();
593 double py2=theEvent[*iit2].py();
594 double pz2=theEvent[*iit2].pz();
595 double mass2=theEvent[*iit2].mSel();
596 double e2=std::sqrt(px2*px2+py2*py2+pz2*pz2+mass2*mass2);
597 double eSum=e1+e2;
598 double pxSum=px1+px2;
599 double pySum=py1+py2;
600 double pzSum=pz1+pz2;
601 double M=std::sqrt(eSum*eSum-pxSum*pxSum-pySum*pySum-pzSum*pzSum);
602 if (M>massCut) {
603 passesCuts=true;
604 ATH_MSG_DEBUG("Acceptable lepton pair with invariant mass : " << M);
605 break;
606 }
607 }
608 if (passesCuts) {accString=" : ACCEPTED"; break;}
609 }
610 ATH_MSG_DEBUG("Dilepton selection" << accString);
611 return passesCuts;
612}
613
615// Iteratively descends through a decay chain, recording all particles it comes
616// across
617void Pythia8B_i::descendThroughDecay(Pythia8::Event &theEvent, std::vector<Pythia8::Particle> &list, int i) const {
618
619 list.push_back(theEvent[i]);
620 std::vector<int> childrenIndices = theEvent.daughterList(i);
621 std::vector<int>::iterator iItr;
622 for (iItr = childrenIndices.begin(); iItr != childrenIndices.end(); ++iItr) {
623 descendThroughDecay(theEvent,list,*iItr);
624 }
625
626 return;
627
628}
629
631// Given a vector of Pythia8 particles, returns their PDG codes
632std::vector<int> Pythia8B_i::getCodes(const std::vector<Pythia8::Particle> &theParticles) const {
633
634 std::vector<int> codes;
635 codes.reserve(theParticles.size());
636 for (auto pItr = theParticles.begin(); pItr!=theParticles.end(); ++pItr ) {
637 codes.push_back( (*pItr).id() );
638 }
639 return codes;
640
641}
642
644// Compares two vectors of integers (any order) and if they are the same, returns
645// true
646bool Pythia8B_i::compare(std::vector<int> vect1, std::vector<int> vect2) const {
647
648 if (vect1.size()!=vect2.size()) return false;
649
650 bool isTheSame(true);
651 int size = vect1.size();
652 std::sort(vect1.rbegin(), vect1.rend());
653 std::sort(vect2.rbegin(), vect2.rend());
654 for (int i=0; i<size; ++i) {
655 if (vect1[i] != vect2[i]) {isTheSame = false; break;}
656 }
657 return isTheSame;
658
659}
660
662// Checks a vector of particles against a vector of pt and eta cuts and if all
663// particles do not pass the cuts, returns false
664// Last argument controls whether the cut is on pt or eta
665bool Pythia8B_i::passesPTCuts(const std::vector<Pythia8::Particle> &theParticles) const {
666
667 bool pass(true);
668 unsigned int i(0);
669 for (auto pItr=theParticles.cbegin(); pItr!=theParticles.cend(); ++pItr,++i) {
670 if ((*pItr).pT() < m_sigPtCuts[i]) pass = false;
671 if (!pass) break;
672 }
673 return pass;
674
675}
676
677bool Pythia8B_i::passesEtaCuts(const std::vector<Pythia8::Particle> &theParticles) const {
678
679 bool pass(true);
680 unsigned int i(0);
681 for (auto pItr=theParticles.cbegin(); pItr!=theParticles.cend(); ++pItr,++i) {
682 if (std::abs((*pItr).eta()) > m_sigEtaCuts[i]) pass = false;
683 if (!pass) break;
684 }
685
686 return pass;
687
688}
689
691// Given an event, checks whether it contains a decay process as defined by the
692// user and if so whether the final states pass user-defined cuts
693bool Pythia8B_i::signalAccept(Pythia8::Event &theEvent, const std::vector<int> &requiredDecay,
694 unsigned int nRequired) const {
695
696 bool acceptEvent(false);
697 std::vector<int> parentsIndices;
698 for (int i = 0; i<theEvent.size(); ++i) {
699 const Pythia8::Particle &theParticle = theEvent[i];
700 int id = theParticle.id();
701 if (id==requiredDecay[0]) parentsIndices.push_back(i);
702 }
703
704 // For the special case where the user merely requires a particular
705 // species to be in the event (mainly for EvtGen use)
706 if ( (requiredDecay.size()==1) && (parentsIndices.size()>0) ) {
707 acceptEvent = true;
708 return acceptEvent;
709 }
710
711 unsigned int goodDecayCounter(0);
712 for (std::vector<int>::iterator iItr = parentsIndices.begin(); iItr!=parentsIndices.end(); ++iItr) {
713 std::vector<Pythia8::Particle> decayMembers;
714 descendThroughDecay(theEvent,decayMembers,*iItr);
715 std::vector<int> pdgCodes = getCodes(decayMembers);
716 if (!compare(requiredDecay,std::move(pdgCodes))) {
717 ATH_MSG_DEBUG("Signal event REJECTED as does not contain required decay chain");
718 continue;
719 }
720 ATH_MSG_DEBUG("Event contains required signal decay chain");
721
722 if (decayMembers.size()==m_sigPtCuts.size()) {
723 if (!passesPTCuts(decayMembers)) {
724 ATH_MSG_DEBUG("Signal event REJECTED as signal chain does not pass pt cuts");
725 continue;
726 }
727 }
728
729 if (decayMembers.size()==m_sigEtaCuts.size()) {
730 if (!passesEtaCuts(decayMembers)) {
731 ATH_MSG_DEBUG("Signal event REJECTED as signal chain does not pass eta cuts");
732 continue;
733 }
734 }
735 if (nRequired==1) {
736 ATH_MSG_DEBUG("Signal decay good; event ACCEPTED");
737 acceptEvent = true;
738 break;
739 }
740 else {++goodDecayCounter;}
741 }
742 if (nRequired>1) {
743 if (goodDecayCounter>=nRequired) {
744 ATH_MSG_DEBUG(nRequired << "signal decays required; " << goodDecayCounter << " found; event ACCEPTED");
745 acceptEvent = true;
746 }
747 else if (goodDecayCounter<nRequired) {
748 ATH_MSG_DEBUG(nRequired << "signal decays required; " << goodDecayCounter << " found; event REJECTED");
749 acceptEvent = false;
750 }
751 }
752
753 return acceptEvent;
754
755}
756
759
760void Pythia8B_i::printSignalSelections(const std::vector<int> &signalProcess,const std::vector<double> &ptCuts, const std::vector<double> &etaCuts, unsigned int nRequired ) const {
761 std::cout << "Signal PDG codes required: ";
762 for (unsigned int k=0; k<m_sigCodes.size(); ++k) std::cout << signalProcess[k] << " ";
763 if (signalProcess.size()==ptCuts.size()) {
764 std::cout << "Cuts on the pt of the signal states: " << std::endl;
765 for (unsigned int l=0; l<ptCuts.size(); ++l) std::cout << ptCuts[l] << " ";
766 }
767 if (signalProcess.size()==etaCuts.size()) {
768 std::cout << "Cuts on the eta of the signal states: " << std::endl;
769 for (unsigned int l=0; l<etaCuts.size(); ++l) std::cout << etaCuts[l] << " ";
770 }
771 std::cout << "Number of such decays required per event: " << nRequired << std::endl;
772 std::cout << std::endl;
773}
774
775
Scalar eta() const
pseudorapidity method
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_DEBUG(x)
#define CHECK(...)
Evaluate an expression and check for errors.
bool passed(DecisionID id, const DecisionIDContainer &)
checks if required decision ID is in the set of IDs in the container
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
const HepPDT::ParticleData * particleData(int pid) const
Access an element in the particle data table.
Definition GenBase.h:126
IntegerProperty m_randomSeed
Seed for random number engine.
Definition GenModule.h:84
std::vector< long int > m_seeds
Definition Pythia8B_i.h:53
int m_totalClone
Definition Pythia8B_i.h:51
void printSignalSelections(const std::vector< int > &, const std::vector< double > &, const std::vector< double > &, unsigned int) const
std::vector< int > m_sigCodes
Definition Pythia8B_i.h:48
bool compare(std::vector< int >, std::vector< int >) const
static CLHEP::HepRandomEngine * p_rndmEngine
Definition Pythia8B_i.h:41
bool m_selectCQuarks
Definition Pythia8B_i.h:54
bool m_doSuppressSmallPT
Definition Pythia8B_i.h:60
int m_totalBQuark
Definition Pythia8B_i.h:51
int m_totalHard
Definition Pythia8B_i.h:51
std::vector< int > m_internalEventNumbers
Definition Pythia8B_i.h:59
double m_aqPtCut
Definition Pythia8B_i.h:55
std::vector< double > m_sigPtCuts
Definition Pythia8B_i.h:50
double m_aqEtaCut
Definition Pythia8B_i.h:55
unsigned int m_nSignalRequired
Definition Pythia8B_i.h:47
bool passesPTCuts(const std::vector< Pythia8::Particle > &) const
bool leptonSelect(Pythia8::Event &, const std::vector< double > &, double, const std::vector< int > &, int, double, bool)
int m_totalCQuark
Definition Pythia8B_i.h:51
Pythia8::SuppressSmallPT * m_SuppressSmallPT
Definition Pythia8B_i.h:61
unsigned int m_had
Definition Pythia8B_i.h:45
int m_atLeastOneAcc
Definition Pythia8B_i.h:51
int m_trigCode
Definition Pythia8B_i.h:46
virtual StatusCode callGenerator()
For calling the generator on each iteration of the event loop.
int m_totalBBarQuark
Definition Pythia8B_i.h:51
bool m_vetoDoubleC
Definition Pythia8B_i.h:54
void descendThroughDecay(Pythia8::Event &, std::vector< Pythia8::Particle > &, int) const
bool m_oppCharges
Definition Pythia8B_i.h:54
bool passesEtaCuts(const std::vector< Pythia8::Particle > &) const
unsigned int m_dec
Definition Pythia8B_i.h:45
std::vector< Pythia8::Event > m_BEventBuffer
Definition Pythia8B_i.h:58
std::vector< int > getCodes(const std::vector< Pythia8::Particle > &) const
std::vector< int > m_bcodes
Definition Pythia8B_i.h:48
double m_qEtaCut
Definition Pythia8B_i.h:55
std::map< int, int > m_speciesCount
Definition Pythia8B_i.h:52
double m_qPtCut
Definition Pythia8B_i.h:55
std::string m_userString
Definition Pythia8B_i.h:56
virtual StatusCode genFinalize()
For finalising the generator, if required.
bool m_vetoDoubleB
Definition Pythia8B_i.h:54
std::vector< int > m_cutCount
Definition Pythia8B_i.h:48
bool m_selectBQuarks
Definition Pythia8B_i.h:54
bool signalAccept(Pythia8::Event &, const std::vector< int > &, unsigned int) const
virtual StatusCode genuserInitialize()
For initialization of user code, if required. Called after genInitialize.
bool cleanUndecayed(Pythia8::Event &, const std::vector< int > &)
double m_invMass
Definition Pythia8B_i.h:55
unsigned int m_failureCount
Definition Pythia8B_i.h:62
std::vector< double > m_trigPtCut
Definition Pythia8B_i.h:49
virtual StatusCode fillEvt(HepMC::GenEvent *)
For filling the HepMC event object.
bool userSelection(Pythia8::Event &, std::string, std::vector< double >)
std::vector< double > m_userVar
Definition Pythia8B_i.h:57
bool pairProperties(Pythia8::Event &, const std::vector< int > &, double, bool)
virtual StatusCode genInitialize()
For initializing the generator, if required.
int m_totalPythiaCalls
Definition Pythia8B_i.h:51
std::vector< double > m_sigEtaCuts
Definition Pythia8B_i.h:50
int m_passingTriggerCuts
Definition Pythia8B_i.h:51
Pythia8B_i(const std::string &name, ISvcLocator *pSvcLocator)
int m_totalCBarQuark
Definition Pythia8B_i.h:51
int m_internal_event_number
Definition Pythia8B_i.h:51
double m_trigEtaCut
Definition Pythia8B_i.h:55
virtual StatusCode genInitialize()
For initializing the generator, if required.
Definition Pythia8_i.cxx:98
IntegerProperty m_dsid
Definition Pythia8_i.h:101
bool useRndmGenSvc() const
Definition Pythia8_i.h:89
HepMC::Pythia8ToHepMC m_pythiaToHepMC
Definition Pythia8_i.h:93
std::shared_ptr< customRndm > m_atlasRndmEngine
Definition Pythia8_i.h:97
StringArrayProperty m_userHooks
Definition Pythia8_i.h:102
std::unique_ptr< Pythia8::Pythia > m_pythia
Definition Pythia8_i.h:92
UnsignedIntegerProperty m_maxFailures
Definition Pythia8_i.h:94
BooleanProperty m_useRndmGenSvc
Definition Pythia8_i.h:96
Pythia8_i(const std::string &name, ISvcLocator *pSvcLocator)
Definition Pythia8_i.cxx:68
void calculateSeedsMC21(long *seeds, const std::string &algName, uint64_t ev, uint64_t run, uint64_t offset=0)
Set the random seed using a string (e.g.
void set_random_states(GenEvent *e, std::vector< T > a)
Definition GenEvent.h:649
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.