Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-05-02 05:09:38

0001 #include "FastSimulation/SimplifiedGeometryPropagator/interface/Decayer.h"
0002 #include "FastSimulation/SimplifiedGeometryPropagator/interface/Particle.h"
0003 #include "FastSimulation/SimplifiedGeometryPropagator/interface/ParticleManager.h"
0004 #include "FWCore/ServiceRegistry/interface/RandomEngineSentry.h"
0005 #include "GeneratorInterface/Pythia8Interface/interface/P8RndmEngine.h"
0006 
0007 #include <Pythia8/Pythia.h>
0008 #include "Pythia8Plugins/HepMC2.h"
0009 
0010 fastsim::Decayer::~Decayer() { ; }
0011 
0012 fastsim::Decayer::Decayer() : pythia_(new Pythia8::Pythia()), pythiaRandomEngine_(new gen::P8RndmEngine()) {
0013   pythia_->setRndmEnginePtr(pythiaRandomEngine_);
0014   pythia_->settings.flag("ProcessLevel:all", false);
0015   pythia_->settings.flag("PartonLevel:FSRinResonances", false);
0016   pythia_->settings.flag("ProcessLevel:resonanceDecays", false);
0017   pythia_->init();
0018 
0019   // forbid all decays
0020   // (decays are allowed selectively in the decay function)
0021   Pythia8::ParticleData& pdt = pythia_->particleData;
0022   int pid = 0;
0023   while (pdt.nextId(pid) > pid) {
0024     pid = pdt.nextId(pid);
0025     pdt.mayDecay(pid, false);
0026   }
0027 }
0028 
0029 void fastsim::Decayer::decay(const Particle& particle,
0030                              std::vector<std::unique_ptr<fastsim::Particle> >& secondaries,
0031                              CLHEP::HepRandomEngine& engine) const {
0032   // make sure pythia takes random numbers from the engine past through via the function arguments
0033   edm::RandomEngineSentry<gen::P8RndmEngine> sentry(pythiaRandomEngine_.get(), &engine);
0034 
0035   // inspired by method Pythia8Hadronizer::residualDecay() in GeneratorInterface/Pythia8Interface/src/Py8GunBase.cc
0036   int pid = particle.pdgId();
0037   // snip decay products of exotic particles or their children. These decay products are preserved from the event record.
0038   // limitation: if exotic incurs heavy energy loss during propagation, the saved decay products could be too hard.
0039 
0040   if (isExotic(pid) || isExotic(particle.getMotherPdgId())) {
0041     return;
0042   }
0043 
0044   pythia_->event.reset();
0045 
0046   // create a pythia particle which has the same properties as the FastSim particle
0047   Pythia8::Particle pythiaParticle(pid,
0048                                    93,
0049                                    0,
0050                                    0,
0051                                    0,
0052                                    0,
0053                                    0,
0054                                    0,
0055                                    particle.momentum().X(),
0056                                    particle.momentum().Y(),
0057                                    particle.momentum().Z(),
0058                                    particle.momentum().E(),
0059                                    particle.momentum().M());
0060   pythiaParticle.vProd(
0061       particle.position().X(), particle.position().Y(), particle.position().Z(), particle.position().T());
0062   pythia_->event.append(pythiaParticle);
0063 
0064   int nentries_before = pythia_->event.size();
0065   // switch on the decay of this and only this particle (avoid double decays)
0066   pythia_->particleData.mayDecay(pid, true);
0067   // do the decay
0068   pythia_->next();
0069   // switch it off again
0070   pythia_->particleData.mayDecay(pid, false);
0071   int nentries_after = pythia_->event.size();
0072 
0073   if (nentries_after <= nentries_before)
0074     return;
0075 
0076   // add decay products back to the event
0077   for (int ipart = nentries_before; ipart < nentries_after; ipart++) {
0078     Pythia8::Particle& daughter = pythia_->event[ipart];
0079 
0080     secondaries.emplace_back(new fastsim::Particle(
0081         daughter.id(),
0082         math::XYZTLorentzVector(daughter.xProd(), daughter.yProd(), daughter.zProd(), daughter.tProd()),
0083         math::XYZTLorentzVector(daughter.px(), daughter.py(), daughter.pz(), daughter.e())));
0084 
0085     // daughter can inherit the SimTrackIndex of mother (if both charged): necessary for FastSim (cheat) tracking
0086     if (particle.charge() != 0 && std::abs(particle.charge() - daughter.charge()) < 1E-10) {
0087       secondaries.back()->setMotherDeltaR(particle.momentum());
0088       secondaries.back()->setMotherPdgId(particle.getMotherDeltaR() == -1 ? particle.pdgId()
0089                                                                           : particle.getMotherPdgId());
0090       secondaries.back()->setMotherSimTrackIndex(particle.simTrackIndex());
0091     }
0092   }
0093 
0094   return;
0095 }