Belle II Software development
EvtGenDecays Class Reference

A class to perform decays via the external EvtGen decay program, see http://evtgen.warwick.ac.uk/, the program manual provided with the EvtGen distribution, and D. More...

#include <EvtGenDecays.h>

Classes

struct  Signal
 Map of signal particle info. More...
 

Public Member Functions

 EvtGenDecays (Pythia8::Pythia *pythiaPtrIn, EvtGen *evtGen, bool limit=true)
 Expert constructor which takes an already initialized EvtGen instance.
 
 EvtGenDecays (Pythia8::Pythia *pythiaPtrIn, std::string decayFile, std::string particleDataFile, EvtExternalGenList *extPtrIn=0, EvtAbsRadCorr *fsrPtrIn=0, int mixing=1, bool xml=false, bool limit=true, bool extUse=true, bool fsrUse=true)
 Constructor.
 
 ~EvtGenDecays ()
 Destructor.
 
 EvtGenDecays (const EvtGenDecays &)=delete
 forbid copy
 
EvtGenDecaysoperator= (const EvtGenDecays &)=delete
 forbid assignment
 
void getDecayLimits (bool limit)
 Get the Decay limits from Pythia.
 
double decay ()
 Perform all decays and return the event weight.
 
void exclude (int id)
 Stop EvtGen decaying a particle.
 
void updatePythia ()
 Update the Pythia particle database from EvtGen.
 
void updateEvtGen ()
 Update the EvtGen particle database from Pythia.
 
void readDecayFile (std::string decayFile, bool xml=false)
 Read an EvtGen user decay file.
 

Public Attributes

bool extOwner
 bool has external model
 
bool fsrOwner
 bool has FSR model
 
EvtExternalGenList * extPtr
 External model pointer.
 
EvtAbsRadCorr * fsrPtr
 FSR model pointer.
 
std::list< EvtDecayBase * > models
 list of model pointers
 
std::map< int, Signalsignals
 signal particles and statuses
 
std::string signalSuffix
 The suffix indicating an EvtGen particle or alias is signal.
 

Protected Member Functions

void updateData (bool final=false)
 Update the particles to decay with EvtGen, and the selected signals.
 
void updateEvent (Pythia8::Particle *pyPro, EvtParticle *egPro, std::vector< int > *pySigs=0, std::vector< EvtParticle * > *egSigs=0, std::vector< double > *bfs=0, double *wgt=0)
 Update the Pythia event record with an EvtGen decay tree.
 
bool checkVertex (Pythia8::Particle *pyPro)
 Check if a particle can decay.
 
bool checkSignal (Pythia8::Particle *pyPro)
 Check if a particle is signal.
 
bool checkOsc (EvtParticle *egPro)
 Check if an EvtGen particle has oscillated.
 

Protected Attributes

Pythia8::Pythia * pythiaPtr
 The pointer to the associated Pythia object.
 
EvtGenRandom rndm
 Random number wrapper for EvtGen.
 
EvtGen * evtgen
 The EvtGen object.
 
std::set< int > incIds
 Set of particle IDs to include decays with EvtGen.
 
std::set< int > excIds
 Set of particle IDs to exclude decays with EvtGen.
 
bool updated
 Flag whether the final particle update has been performed.
 
std::map< int, Signal >::iterator signal
 The selected signal iterator.
 
double tau0Max
 Parameters used to check if a particle should decay (as set via Pythia).
 
double tauMax
 the pythia parameter ParticleDecays:tauMax
 
double rMax
 the pythia parameter ParticleDecays:rMax
 
double xyMax
 the pythia parameter ParticleDecays:xyMax
 
double zMax
 the pythia parameter ParticleDecays:zMax
 
bool limitTau0
 the pythia parameter ParticleDecays:limitTau0
 
bool limitTau
 the pythia parameter ParticleDecays:limitTau
 
bool limitRadius
 the pythia parameter ParticleDecays:limitRadius
 
bool limitCylinder
 the pythia parameter ParticleDecays:limitCylinder
 
bool limitDecay
 is the decay limited
 

Static Protected Attributes

static constexpr int NTRYDECAY = 1000
 Number of times to try a decay sampling (constant).
 

Detailed Description

A class to perform decays via the external EvtGen decay program, see http://evtgen.warwick.ac.uk/, the program manual provided with the EvtGen distribution, and D.

J. Lange, Nucl. Instrum. Meth. A462, 152 (2001) for details.

EvtGen performs a series of decays from some initial particle decay, rather than just a single decay, and so EvtGen cannot be interfaced through the standard external DecayHandler class without considerable complication. Consequently, EvtGen is called on the complete event record after all steps of Pythia are completed.

Oftentimes a specific "signal" decay is needed to occur once in an event, and all other decays performed normally. This is possible via reading in a user decay file (with readDecayFile) and creating aliased particles with names ending with signalSuffix. By default, this is "_SIGNAL". When decay() is called, all particles in the Pythia event record that are of the same types as the signal particles are collected. One is selected at random and decayed via the channel(s) defined for that aliased signal particle. All other particles are decayed normally. The weight for the event is calculated and returned.

It is also possible to specify a status needed to consider a particle as a signal candidate. This can be done by modifying the signals map, e.g. if the tau- is a signal candidate, then EvtGenDecays.signals[15].status = 201 will only only select as candidates any tau- with this status. This allows the event record to be changed before decays, so only certain particles are selected as possible signal candidates (e.g. passing kinematic requirements).

Please note that particles produced from a signal candidate decay are not searched for additional signal candidates. This means that if B0 and tau- have been designated as signal, then a tau- from a W- decay would be a signal candidate, while a tau- from a B0 decay would not. This restriction arises from the additional complexity of allowing recursive signal decays. The following statuses are used: 93 for particles decayed with EvtGen, 94 for particles oscillated with EvtGen, 95 for signal particles, and 96 for signal particles from an oscillation.

Definition at line 88 of file EvtGenDecays.h.

Constructor & Destructor Documentation

◆ EvtGenDecays() [1/2]

EvtGenDecays ( Pythia8::Pythia *  pythiaPtrIn,
EvtGen *  evtGen,
bool  limit = true 
)

Expert constructor which takes an already initialized EvtGen instance.

Definition at line 255 of file EvtGenDecays.h.

255 :
256 extOwner(false), fsrOwner(false), extPtr(nullptr), fsrPtr(nullptr),
257 signalSuffix("_SIGNAL"), pythiaPtr(pythiaPtrIn), rndm(nullptr),
258 evtgen(evtGen), updated(false)
259{
260 getDecayLimits(limit);
261}
EvtGen * evtgen
The EvtGen object.
Definition: EvtGenDecays.h:185
EvtAbsRadCorr * fsrPtr
FSR model pointer.
Definition: EvtGenDecays.h:139
bool fsrOwner
bool has FSR model
Definition: EvtGenDecays.h:137
bool updated
Flag whether the final particle update has been performed.
Definition: EvtGenDecays.h:191
EvtExternalGenList * extPtr
External model pointer.
Definition: EvtGenDecays.h:138
bool extOwner
bool has external model
Definition: EvtGenDecays.h:136
void getDecayLimits(bool limit)
Get the Decay limits from Pythia.
Definition: EvtGenDecays.h:286
std::string signalSuffix
The suffix indicating an EvtGen particle or alias is signal.
Definition: EvtGenDecays.h:154
EvtGenRandom rndm
Random number wrapper for EvtGen.
Definition: EvtGenDecays.h:182
Pythia8::Pythia * pythiaPtr
The pointer to the associated Pythia object.
Definition: EvtGenDecays.h:179

◆ EvtGenDecays() [2/2]

EvtGenDecays ( Pythia8::Pythia *  pythiaPtrIn,
std::string  decayFile,
std::string  particleDataFile,
EvtExternalGenList *  extPtrIn = 0,
EvtAbsRadCorr *  fsrPtrIn = 0,
int  mixing = 1,
bool  xml = false,
bool  limit = true,
bool  extUse = true,
bool  fsrUse = true 
)

Constructor.

Definition at line 263 of file EvtGenDecays.h.

266 : extPtr(extPtrIn), fsrPtr(fsrPtrIn),
267 signalSuffix("_SIGNAL"), pythiaPtr(pythiaPtrIn), rndm(&pythiaPtr->rndm),
268 updated(false)
269{
270
271 // Initialize EvtGen.
272 if (!extPtr && fsrPtr) {
273 B2ERROR("Error in EvtGenDecays::EvtGenDecays: extPtr is null but fsrPtr is provided");
274 return;
275 }
276 if (extPtr) extOwner = false;
277 else {extOwner = true; extPtr = new EvtExternalGenList();}
278 if (fsrPtr) fsrOwner = false;
279 else {fsrOwner = true; fsrPtr = extPtr->getPhotosModel();}
280 models = extPtr->getListOfModels();
281 evtgen = new EvtGen(decayFile.c_str(), particleDataFile.c_str(),
282 &rndm, fsrUse ? fsrPtr : 0, extUse ? &models : 0, mixing, xml);
283 getDecayLimits(limit);
284}
std::list< EvtDecayBase * > models
list of model pointers
Definition: EvtGenDecays.h:140

◆ ~EvtGenDecays()

~EvtGenDecays ( )
inline

Destructor.

Definition at line 102 of file EvtGenDecays.h.

103 {
104 if (evtgen) delete evtgen;
105 if (extOwner && extPtr) delete extPtr;
106 if (fsrOwner && fsrPtr) delete fsrPtr;
107 }

Member Function Documentation

◆ checkOsc()

bool checkOsc ( EvtParticle *  egPro)
protected

Check if an EvtGen particle has oscillated.

Definition at line 668 of file EvtGenDecays.h.

669{
670 if (egPro && egPro->getNDaug() == 1 &&
671 egPro->getPDGId() != egPro->getDaug(0)->getPDGId()) return true;
672 else return false;
673}

◆ checkSignal()

bool checkSignal ( Pythia8::Particle *  pyPro)
protected

Check if a particle is signal.

Definition at line 653 of file EvtGenDecays.h.

654{
655 signal = signals.find(pyPro->id());
656 if (signal != signals.end() && (signal->second.status < 0 ||
657 signal->second.status == pyPro->status())) return true;
658 else return false;
659}
std::map< int, Signal > signals
signal particles and statuses
Definition: EvtGenDecays.h:151
std::map< int, Signal >::iterator signal
The selected signal iterator.
Definition: EvtGenDecays.h:194

◆ checkVertex()

bool checkVertex ( Pythia8::Particle *  pyPro)
protected

Check if a particle can decay.

Definition at line 635 of file EvtGenDecays.h.

636{
637 using ::Pythia8::pow2;
638
639 if (!limitDecay) return true;
640 if (limitTau0 && pyPro->tau0() > tau0Max) return false;
641 if (limitTau && pyPro->tau() > tauMax) return false;
642 if (limitRadius && pow2(pyPro->xDec()) + pow2(pyPro->yDec())
643 + pow2(pyPro->zDec()) > pow2(rMax)) return false;
644 if (limitCylinder && (pow2(pyPro->xDec()) + pow2(pyPro->yDec())
645 > pow2(xyMax) || abs(pyPro->zDec()) > zMax)) return false;
646 return true;
647}
double tauMax
the pythia parameter ParticleDecays:tauMax
Definition: EvtGenDecays.h:201
bool limitCylinder
the pythia parameter ParticleDecays:limitCylinder
Definition: EvtGenDecays.h:208
double zMax
the pythia parameter ParticleDecays:zMax
Definition: EvtGenDecays.h:204
double xyMax
the pythia parameter ParticleDecays:xyMax
Definition: EvtGenDecays.h:203
bool limitTau
the pythia parameter ParticleDecays:limitTau
Definition: EvtGenDecays.h:206
double rMax
the pythia parameter ParticleDecays:rMax
Definition: EvtGenDecays.h:202
bool limitRadius
the pythia parameter ParticleDecays:limitRadius
Definition: EvtGenDecays.h:207
bool limitDecay
is the decay limited
Definition: EvtGenDecays.h:209
bool limitTau0
the pythia parameter ParticleDecays:limitTau0
Definition: EvtGenDecays.h:205
double tau0Max
Parameters used to check if a particle should decay (as set via Pythia).
Definition: EvtGenDecays.h:200

◆ decay()

double decay ( )

Perform all decays and return the event weight.

Definition at line 327 of file EvtGenDecays.h.

328{
329
330 // Reset the signal and signal counters.
331 if (!pythiaPtr || !evtgen) return -1;
332 if (!updated) updateData(true);
333
334 // Loop over all particles in the Pythia event.
335 Pythia8::Event& event = pythiaPtr->event;
336 std::vector<int> pySigs; std::vector<EvtParticle*> egSigs, egPrts;
337 std::vector<double> bfs; double wgt(1.);
338 for (int iPro = 0; iPro < event.size(); ++iPro) {
339
340 // Check particle is final and can be decayed by EvtGen.
341 Pythia8::Particle* pyPro = &event[iPro];
342 if (!pyPro->isFinal()) continue;
343 if (incIds.find(pyPro->id()) == incIds.end()) continue;
344 if (pyPro->status() == 93 || pyPro->status() == 94) continue;
345
346 // Decay the progenitor with EvtGen.
347 EvtParticle* egPro = EvtParticleFactory::particleFactory
348 (EvtPDL::evtIdFromStdHep(pyPro->id()),
349 EvtVector4R(pyPro->e(), pyPro->px(), pyPro->py(), pyPro->pz()));
350 egPrts.push_back(egPro);
351 egPro->setDiagonalSpinDensity();
352 evtgen->generateDecay(egPro);
353 pyPro->tau(egPro->getLifetime());
354 if (!checkVertex(pyPro)) continue;
355
356 // Add oscillations to event record.
357 if (checkOsc(egPro))
358 updateEvent(pyPro, egPro, &pySigs, &egSigs, &bfs, &wgt);
359
360 // Undo decay if signal (duplicate to stop oscillations).
361 else if (checkSignal(pyPro)) {
362 pySigs.push_back(pyPro->index());
363 egSigs.push_back(egPro);
364 bfs.push_back(signal->second.bfs[0]);
365 wgt *= 1 - bfs.back();
366 egPro->deleteDaughters();
367 // The following problematic part of the Pythia's code was removed, see
368 // https://gitlab.com/Pythia8/releases/
369 // -/blob/pythia8310/include/Pythia8Plugins/EvtGen.h?ref_type=tags#L321-325
370
371 // If not signal, add to event record.
372 } else updateEvent(pyPro, egPro, &pySigs, &egSigs, &bfs, &wgt);
373 }
374 if (pySigs.size() == 0) {
375 for (int iPrt = 0; iPrt < (int)egPrts.size(); ++iPrt)
376 egPrts[iPrt]->deleteTree();
377 return 0;
378 }
379
380 // Determine the decays of the signal particles (signal or background).
381 std::vector<int> modes; int force(-1);
382 for (int iTry = 1; iTry <= NTRYDECAY; ++iTry) {
383 int n(0);
384
385 modes.clear(); force = pythiaPtr->rndm.pick(bfs); n = 0;
386 for (int iSig = 0; iSig < (int)pySigs.size(); ++iSig) {
387 if (iSig == force) modes.push_back(0);
388 else modes.push_back(pythiaPtr->rndm.flat() > bfs[iSig]);
389 if (modes.back() == 0) ++n;
390 }
391 if (pythiaPtr->rndm.flat() <= 1.0 / n) break;
392 if (iTry == NTRYDECAY) {
393 wgt = 2.;
394 B2WARNING("Warning in EvtGenDecays::decay: maximum number of signal decay attempts exceeded"
395 << LogVar("maxTries", NTRYDECAY));
396 }
397 }
398
399 // Decay the signal particles and mark forced decay.
400 for (int iSig = 0; iSig < (int)pySigs.size(); ++iSig) {
401 EvtParticle* egSig = egSigs[iSig];
402 Pythia8::Particle* pySig = &event[pySigs[iSig]];
403 signal = signals.find(pySig->id());
404 if (egSig->getNDaug() > 0) egSig = egSig->getDaug(0);
405 if (modes[iSig] == 0) egSig->setId(signal->second.egId);
406 else {
407 signal->second.modes.getDecayModel(egSig);
408 egSig->setChannel(signal->second.map[egSig->getChannel()]);
409 }
410 if (iSig == force)
411 pySig->status(pySig->status() == 92 || pySig->status() == 94 ? 96 : 95);
412 evtgen->generateDecay(egSig);
413 updateEvent(pySig, egSigs[iSig]);
414 }
415
416 // Delete all EvtGen particles and return weight.
417 for (int iPrt = 0; iPrt < (int)egPrts.size(); ++iPrt)
418 egPrts[iPrt]->deleteTree();
419 return 1. - wgt;
420
421}
void updateData(bool final=false)
Update the particles to decay with EvtGen, and the selected signals.
Definition: EvtGenDecays.h:481
std::set< int > incIds
Set of particle IDs to include decays with EvtGen.
Definition: EvtGenDecays.h:187
bool checkSignal(Pythia8::Particle *pyPro)
Check if a particle is signal.
Definition: EvtGenDecays.h:653
bool checkVertex(Pythia8::Particle *pyPro)
Check if a particle can decay.
Definition: EvtGenDecays.h:635
bool checkOsc(EvtParticle *egPro)
Check if an EvtGen particle has oscillated.
Definition: EvtGenDecays.h:668
static constexpr int NTRYDECAY
Number of times to try a decay sampling (constant).
Definition: EvtGenDecays.h:176
void updateEvent(Pythia8::Particle *pyPro, EvtParticle *egPro, std::vector< int > *pySigs=0, std::vector< EvtParticle * > *egSigs=0, std::vector< double > *bfs=0, double *wgt=0)
Update the Pythia event record with an EvtGen decay tree.
Definition: EvtGenDecays.h:574
Class to store variables with their name which were sent to the logging service.

◆ exclude()

void exclude ( int  id)
inline

Stop EvtGen decaying a particle.

Definition at line 122 of file EvtGenDecays.h.

122{excIds.insert(id);}
std::set< int > excIds
Set of particle IDs to exclude decays with EvtGen.
Definition: EvtGenDecays.h:188

◆ getDecayLimits()

void getDecayLimits ( bool  limit)

Get the Decay limits from Pythia.

Definition at line 286 of file EvtGenDecays.h.

287{
288 // Get the Pythia decay limits.
289 if (!pythiaPtr) return;
290 limitTau0 = pythiaPtr->settings.flag("ParticleDecays:limitTau0");
291 tau0Max = pythiaPtr->settings.parm("ParticleDecays:tau0Max");
292 limitTau = pythiaPtr->settings.flag("ParticleDecays:limitTau");
293 tauMax = pythiaPtr->settings.parm("ParticleDecays:tauMax");
294 limitRadius = pythiaPtr->settings.flag("ParticleDecays:limitRadius");
295 rMax = pythiaPtr->settings.parm("ParticleDecays:rMax");
296 limitCylinder = pythiaPtr->settings.flag("ParticleDecays:limitCylinder");
297 xyMax = pythiaPtr->settings.parm("ParticleDecays:xyMax");
298 zMax = pythiaPtr->settings.parm("ParticleDecays:zMax");
299 limitDecay = limit && (limitTau0 || limitTau ||
301
302}

◆ readDecayFile()

void readDecayFile ( std::string  decayFile,
bool  xml = false 
)
inline

Read an EvtGen user decay file.

Definition at line 131 of file EvtGenDecays.h.

132 {
133 evtgen->readUDecay(decayFile.c_str(), xml); updateData();
134 }

◆ updateData()

void updateData ( bool  final = false)
protected

Update the particles to decay with EvtGen, and the selected signals.

Definition at line 481 of file EvtGenDecays.h.

482{
483
484 // Loop over the EvtGen entries.
485 if (!pythiaPtr) return;
486 EvtDecayTable* egTable = EvtDecayTable::getInstance();
487 if (!egTable) return;
488 for (int iEntry = 0; iEntry < (int)EvtPDL::entries(); ++iEntry) {
489 EvtId egId = EvtPDL::getEntry(iEntry);
490 int pyId = EvtPDL::getStdHep(egId);
491 if (egTable->getNModes(egId) == 0) continue;
492 if (excIds.find(pyId) != excIds.end()) continue;
493
494 // Stop Pythia from decaying the particle and include in decay set.
495 if (final) {
496 incIds.insert(pyId);
497 pythiaPtr->particleData.mayDecay(pyId, false);
498 }
499
500 // Check for signal.
501 std::string egName = EvtPDL::name(egId);
502 if (egName.size() <= signalSuffix.size() || egName.substr
503 (egName.size() - signalSuffix.size()) != signalSuffix) continue;
504 signal = signals.find(pyId);
505 if (signal == signals.end()) {
506 signals[pyId].status = -1;
507 signal = signals.find(pyId);
508 }
509 signal->second.egId = egId;
510 signal->second.bfs = std::vector<double>(2, 0);
511 if (!final) continue;
512
513 // Get the signal and background decay modes.
514 std::vector<EvtParticleDecayList> egList = egTable->getDecayTable();
515 int sigIdx = egId.getAlias();
516 int bkgIdx = EvtPDL::evtIdFromStdHep(pyId).getAlias();
517 if (sigIdx > (int)egList.size() || bkgIdx > (int)egList.size()) continue;
518 EvtParticleDecayList sigModes = egList[sigIdx];
519 EvtParticleDecayList bkgModes = egList[bkgIdx];
520 EvtParticleDecayList allModes = egList[bkgIdx];
521 double sum(0);
522
523 // Sum signal branching fractions.
524 for (int iMode = 0; iMode < sigModes.getNMode(); ++iMode) {
525 EvtDecayBase* mode = sigModes.getDecayModel(iMode);
526 if (!mode) continue;
527 signal->second.bfs[0] += mode->getBranchingFraction();
528 sum += mode->getBranchingFraction();
529 }
530
531 // Sum remaining background branching fractions.
532 for (int iMode = 0; iMode < allModes.getNMode(); ++iMode) {
533 EvtDecayBase* mode = allModes.getDecayModel(iMode);
534 if (!mode) continue;
535 bool match(false);
536 for (int jMode = 0; jMode < sigModes.getNMode(); ++jMode)
537 if (mode->matchingDecay(*(sigModes.getDecayModel(jMode)))) {
538 match = true; break;
539 }
540 if (match) bkgModes.removeMode(mode);
541 else {
542 signal->second.map.push_back(iMode);
543 signal->second.bfs[1] += mode->getBranchingFraction();
544 sum += mode->getBranchingFraction();
545 }
546 }
547 signal->second.modes = bkgModes;
548 for (int iBf = 0; iBf < (int)signal->second.bfs.size(); ++iBf)
549 signal->second.bfs[iBf] /= sum;
550 }
551 if (final) updated = true;
552
553}

◆ updateEvent()

void updateEvent ( Pythia8::Particle *  pyPro,
EvtParticle *  egPro,
std::vector< int > *  pySigs = 0,
std::vector< EvtParticle * > *  egSigs = 0,
std::vector< double > *  bfs = 0,
double *  wgt = 0 
)
protected

Update the Pythia event record with an EvtGen decay tree.

Definition at line 574 of file EvtGenDecays.h.

577{
578
579 // Set up the mother vector.
580 if (!pythiaPtr) return;
581 EvtParticle* egMom = egPro;
582 if (egPro->getNDaug() == 1 && egPro->getPDGId() ==
583 egPro->getDaug(0)->getPDGId()) egMom = egPro->getDaug(0);
584 Pythia8::Event& event = pythiaPtr->event;
585 std::vector< std::pair<EvtParticle*, int> >
586 moms(1, std::pair<EvtParticle*, int>(egMom, pyPro->index()));
587
588 // Loop over the mothers.
589 while (moms.size() != 0) {
590
591 // Check if particle can decay.
592 egMom = moms.back().first;
593 int iMom = moms.back().second;
594 Pythia8::Particle* pyMom = &event[iMom];
595 moms.pop_back();
596 if (!checkVertex(pyMom)) continue;
597 bool osc(checkOsc(egMom));
598
599 // Set the children of the mother.
600 pyMom->daughters(event.size(), event.size() + egMom->getNDaug() - 1);
601 pyMom->statusNeg();
602 Pythia8::Vec4 vProd = pyMom->vDec();
603 for (int iDtr = 0 ; iDtr < (int)egMom->getNDaug(); ++iDtr) {
604 EvtParticle* egDtr = egMom->getDaug(iDtr);
605 int id = egDtr->getPDGId();
606 EvtVector4R p = egDtr->getP4Lab();
607 int idx = event.append(id, 93, iMom, 0, 0, 0, 0, 0, p.get(1),
608 p.get(2), p.get(3), p.get(0), egDtr->mass());
609 Pythia8::Particle* pyDtr = &event.back();
610 pyDtr->vProd(vProd);
611 pyDtr->tau(egDtr->getLifetime());
612 if (pySigs && egSigs && bfs && wgt && checkSignal(pyDtr)) {
613 pySigs->push_back(pyDtr->index());
614 egSigs->push_back(egDtr);
615 bfs->push_back(signal->second.bfs[0]);
616 (*wgt) *= 1 - bfs->back();
617 egDtr->deleteDaughters();
618 }
619 if (osc) pyDtr->status(94);
620 if (egDtr->getAttribute("FSR"))
621 pyDtr->status(Belle2::GeneratorConst::FSR_STATUS_CODE);
622 if (egDtr->getNDaug() > 0)
623 moms.push_back(std::pair<EvtParticle*, int>(egDtr, idx));
624 }
625 }
626
627}

◆ updateEvtGen()

void updateEvtGen ( )

Update the EvtGen particle database from Pythia.

Definition at line 455 of file EvtGenDecays.h.

456{
457 if (!pythiaPtr || !evtgen) return;
458 int pyId = pythiaPtr->particleData.nextId(0);
459 while (pyId != 0) {
460 EvtId egId = EvtPDL::evtIdFromStdHep(pyId);
461 EvtPDL::reSetMass(egId, pythiaPtr->particleData.m0(pyId));
462 EvtPDL::reSetWidth(egId, pythiaPtr->particleData.mWidth(pyId));
463 EvtPDL::reSetMassMin(egId, pythiaPtr->particleData.mMin(pyId));
464 EvtPDL::reSetMassMax(egId, pythiaPtr->particleData.mMax(pyId));
465 pyId = pythiaPtr->particleData.nextId(pyId);
466 }
467}

◆ updatePythia()

void updatePythia ( )

Update the Pythia particle database from EvtGen.

Definition at line 431 of file EvtGenDecays.h.

432{
433 if (!pythiaPtr || !evtgen) return;
434 for (int entry = 0; entry < (int)EvtPDL::entries(); ++entry) {
435 EvtId egId = EvtPDL::getEntry(entry);
436 int pyId = EvtPDL::getStdHep(egId);
437 pythiaPtr->particleData.spinType(pyId, EvtPDL::getSpinType(egId));
438 pythiaPtr->particleData.chargeType(pyId, EvtPDL::chg3(egId));
439 pythiaPtr->particleData.m0(pyId, EvtPDL::getMass(egId));
440 pythiaPtr->particleData.mWidth(pyId, EvtPDL::getWidth(egId));
441 pythiaPtr->particleData.mMin(pyId, EvtPDL::getMinMass(egId));
442 pythiaPtr->particleData.mMax(pyId, EvtPDL::getMaxMass(egId));
443 pythiaPtr->particleData.tau0(pyId, EvtPDL::getctau(egId));
444 }
445}

Member Data Documentation

◆ evtgen

EvtGen* evtgen
protected

The EvtGen object.

Definition at line 185 of file EvtGenDecays.h.

◆ excIds

std::set<int> excIds
protected

Set of particle IDs to exclude decays with EvtGen.

Definition at line 188 of file EvtGenDecays.h.

◆ extOwner

bool extOwner

bool has external model

Definition at line 136 of file EvtGenDecays.h.

◆ extPtr

EvtExternalGenList* extPtr

External model pointer.

Definition at line 138 of file EvtGenDecays.h.

◆ fsrOwner

bool fsrOwner

bool has FSR model

Definition at line 137 of file EvtGenDecays.h.

◆ fsrPtr

EvtAbsRadCorr* fsrPtr

FSR model pointer.

Definition at line 139 of file EvtGenDecays.h.

◆ incIds

std::set<int> incIds
protected

Set of particle IDs to include decays with EvtGen.

Definition at line 187 of file EvtGenDecays.h.

◆ limitCylinder

bool limitCylinder
protected

the pythia parameter ParticleDecays:limitCylinder

Definition at line 208 of file EvtGenDecays.h.

◆ limitDecay

bool limitDecay
protected

is the decay limited

Definition at line 209 of file EvtGenDecays.h.

◆ limitRadius

bool limitRadius
protected

the pythia parameter ParticleDecays:limitRadius

Definition at line 207 of file EvtGenDecays.h.

◆ limitTau

bool limitTau
protected

the pythia parameter ParticleDecays:limitTau

Definition at line 206 of file EvtGenDecays.h.

◆ limitTau0

bool limitTau0
protected

the pythia parameter ParticleDecays:limitTau0

Definition at line 205 of file EvtGenDecays.h.

◆ models

std::list<EvtDecayBase*> models

list of model pointers

Definition at line 140 of file EvtGenDecays.h.

◆ NTRYDECAY

constexpr int NTRYDECAY = 1000
staticconstexprprotected

Number of times to try a decay sampling (constant).

Definition at line 176 of file EvtGenDecays.h.

◆ pythiaPtr

Pythia8::Pythia* pythiaPtr
protected

The pointer to the associated Pythia object.

Definition at line 179 of file EvtGenDecays.h.

◆ rMax

double rMax
protected

the pythia parameter ParticleDecays:rMax

Definition at line 202 of file EvtGenDecays.h.

◆ rndm

EvtGenRandom rndm
protected

Random number wrapper for EvtGen.

Definition at line 182 of file EvtGenDecays.h.

◆ signal

std::map<int,Signal>::iterator signal
protected

The selected signal iterator.

Definition at line 194 of file EvtGenDecays.h.

◆ signals

std::map<int, Signal> signals

signal particles and statuses

Definition at line 151 of file EvtGenDecays.h.

◆ signalSuffix

std::string signalSuffix

The suffix indicating an EvtGen particle or alias is signal.

Definition at line 154 of file EvtGenDecays.h.

◆ tau0Max

double tau0Max
protected

Parameters used to check if a particle should decay (as set via Pythia).

http://home.thep.lu.se/Pythia/pythia82html/ParticleDecays.html the pythia parameter ParticleDecays:tau0Max

Definition at line 200 of file EvtGenDecays.h.

◆ tauMax

double tauMax
protected

the pythia parameter ParticleDecays:tauMax

Definition at line 201 of file EvtGenDecays.h.

◆ updated

bool updated
protected

Flag whether the final particle update has been performed.

Definition at line 191 of file EvtGenDecays.h.

◆ xyMax

double xyMax
protected

the pythia parameter ParticleDecays:xyMax

Definition at line 203 of file EvtGenDecays.h.

◆ zMax

double zMax
protected

the pythia parameter ParticleDecays:zMax

Definition at line 204 of file EvtGenDecays.h.


The documentation for this class was generated from the following file: