Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:04:58

0001 #include <memory>
0002 
0003 #include "GeneratorInterface/Core/interface/GeneratorFilter.h"
0004 #include "GeneratorInterface/ExternalDecays/interface/ExternalDecayDriver.h"
0005 
0006 #include "GeneratorInterface/Pythia8Interface/interface/Py8GunBase.h"
0007 
0008 namespace gen {
0009 
0010   class Py8PtotGun : public Py8GunBase {
0011   public:
0012     Py8PtotGun(edm::ParameterSet const&);
0013     ~Py8PtotGun() override {}
0014 
0015     bool generatePartonsAndHadronize() override;
0016     const char* classname() const override;
0017 
0018   private:
0019     // Ptot Gun particle(s) characteristics
0020     double fMinEta;
0021     double fMaxEta;
0022     double fMinPtot;
0023     double fMaxPtot;
0024     bool fAddAntiParticle;
0025   };
0026 
0027   // implementation
0028   //
0029   Py8PtotGun::Py8PtotGun(edm::ParameterSet const& ps) : Py8GunBase(ps) {
0030     // ParameterSet defpset ;
0031     edm::ParameterSet pgun_params = ps.getParameter<edm::ParameterSet>("PGunParameters");  // , defpset ) ;
0032     fMinEta = pgun_params.getParameter<double>("MinEta");                                  // ,-2.2);
0033     fMaxEta = pgun_params.getParameter<double>("MaxEta");                                  // , 2.2);
0034     fMinPtot = pgun_params.getParameter<double>("MinPtot");                                // ,  0.);
0035     fMaxPtot = pgun_params.getParameter<double>("MaxPtot");                                // ,  0.);
0036     fAddAntiParticle = pgun_params.getParameter<bool>("AddAntiParticle");                  //, false) ;
0037   }
0038 
0039   bool Py8PtotGun::generatePartonsAndHadronize() {
0040     fMasterGen->event.reset();
0041 
0042     int NTotParticles = fPartIDs.size();
0043     if (fAddAntiParticle)
0044       NTotParticles *= 2;
0045 
0046     // energy below is dummy, it is not used
0047     (fMasterGen->event).append(990, -11, 0, 0, 2, 1 + NTotParticles, 0, 0, 0., 0., 0., 15000., 15000.);
0048 
0049     for (size_t i = 0; i < fPartIDs.size(); i++) {
0050       int particleID = fPartIDs[i];  // this is PDG
0051 
0052       double phi = (fMaxPhi - fMinPhi) * randomEngine().flat() + fMinPhi;
0053       double eta = (fMaxEta - fMinEta) * randomEngine().flat() + fMinEta;
0054       double the = 2. * atan(exp(-eta));
0055 
0056       double pp = (fMaxPtot - fMinPtot) * randomEngine().flat() + fMinPtot;
0057 
0058       double mass = (fMasterGen->particleData).m0(particleID);
0059 
0060       double pt = pp * sin(the);
0061       double ee = sqrt(pp * pp + mass * mass);
0062       double px = pt * cos(phi);
0063       double py = pt * sin(phi);
0064       double pz = pp * cos(the);
0065 
0066       if (!((fMasterGen->particleData).isParticle(particleID))) {
0067         particleID = std::abs(particleID);
0068       }
0069       if (1 <= std::abs(particleID) && std::abs(particleID) <= 6)  // quarks
0070         (fMasterGen->event).append(particleID, 23, 1, 0, 0, 0, 101, 0, px, py, pz, ee, mass);
0071       else if (std::abs(particleID) == 21)  // gluons
0072         (fMasterGen->event).append(21, 23, 1, 0, 0, 0, 101, 102, px, py, pz, ee, mass);
0073       // other
0074       else {
0075         (fMasterGen->event).append(particleID, 1, 1, 0, 0, 0, 0, 0, px, py, pz, ee, mass);
0076         int eventSize = (fMasterGen->event).size() - 1;
0077         // -log(flat) = exponential distribution
0078         double tauTmp = -(fMasterGen->event)[eventSize].tau0() * log(randomEngine().flat());
0079         (fMasterGen->event)[eventSize].tau(tauTmp);
0080       }
0081 
0082       // Here also need to add anti-particle (if any)
0083       // otherwise just add a 2nd particle of the same type
0084       // (for example, gamma)
0085       //
0086       if (fAddAntiParticle) {
0087         if (1 <= std::abs(particleID) && std::abs(particleID) <= 6) {  // quarks
0088           (fMasterGen->event).append(-particleID, 23, 1, 0, 0, 0, 0, 101, -px, -py, -pz, ee, mass);
0089         } else if (std::abs(particleID) == 21) {  // gluons
0090           (fMasterGen->event).append(21, 23, 1, 0, 0, 0, 102, 101, -px, -py, -pz, ee, mass);
0091         } else {
0092           if ((fMasterGen->particleData).isParticle(-particleID)) {
0093             (fMasterGen->event).append(-particleID, 1, 1, 0, 0, 0, 0, 0, -px, -py, -pz, ee, mass);
0094           } else {
0095             (fMasterGen->event).append(particleID, 1, 1, 0, 0, 0, 0, 0, -px, -py, -pz, ee, mass);
0096           }
0097           int eventSize = (fMasterGen->event).size() - 1;
0098           // -log(flat) = exponential distribution
0099           double tauTmp = -(fMasterGen->event)[eventSize].tau0() * log(randomEngine().flat());
0100           (fMasterGen->event)[eventSize].tau(tauTmp);
0101         }
0102       }
0103     }
0104 
0105     if (!fMasterGen->next())
0106       return false;
0107     evtGenDecay();
0108 
0109     event() = std::make_unique<HepMC::GenEvent>();
0110     return toHepMC.fill_next_event(fMasterGen->event, event().get());
0111   }
0112 
0113   const char* Py8PtotGun::classname() const { return "Py8PtotGun"; }
0114 
0115   typedef edm::GeneratorFilter<gen::Py8PtotGun, gen::ExternalDecayDriver> Pythia8PtotGun;
0116 
0117 }  // namespace gen
0118 
0119 using gen::Pythia8PtotGun;
0120 DEFINE_FWK_MODULE(Pythia8PtotGun);