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 Py8PtGun : public Py8GunBase {
0011   public:
0012     Py8PtGun(edm::ParameterSet const&);
0013     ~Py8PtGun() override {}
0014 
0015     bool generatePartonsAndHadronize() override;
0016     const char* classname() const override;
0017 
0018   private:
0019     // PtGun particle(s) characteristics
0020     double fMinEta;
0021     double fMaxEta;
0022     double fMinPt;
0023     double fMaxPt;
0024     bool fAddAntiParticle;
0025   };
0026 
0027   // implementation
0028   //
0029   Py8PtGun::Py8PtGun(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     fMinPt = pgun_params.getParameter<double>("MinPt");                                    // ,  0.);
0035     fMaxPt = pgun_params.getParameter<double>("MaxPt");                                    // ,  0.);
0036     fAddAntiParticle = pgun_params.getParameter<bool>("AddAntiParticle");                  //, false) ;
0037   }
0038 
0039   bool Py8PtGun::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 - need to convert to Py8 ???
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 pt = (fMaxPt - fMinPt) * randomEngine().flat() + fMinPt;
0057 
0058       double mass = (fMasterGen->particleData).m0(particleID);
0059 
0060       double pp = pt / sin(the);  // sqrt( ee*ee - mass*mass );
0061       double ee = sqrt(pp * pp + mass * mass);
0062 
0063       double px = pt * cos(phi);
0064       double py = pt * sin(phi);
0065       double pz = pp * cos(the);
0066 
0067       if (!((fMasterGen->particleData).isParticle(particleID))) {
0068         particleID = std::abs(particleID);
0069       }
0070       if (1 <= std::abs(particleID) && std::abs(particleID) <= 6)  // quarks
0071         (fMasterGen->event).append(particleID, 23, 1, 0, 0, 0, 101, 0, px, py, pz, ee, mass);
0072       else if (std::abs(particleID) == 21)  // gluons
0073         (fMasterGen->event).append(21, 23, 1, 0, 0, 0, 101, 102, px, py, pz, ee, mass);
0074       // other
0075       else {
0076         (fMasterGen->event).append(particleID, 1, 1, 0, 0, 0, 0, 0, px, py, pz, ee, mass);
0077         int eventSize = (fMasterGen->event).size() - 1;
0078         // -log(flat) = exponential distribution
0079         double tauTmp = -(fMasterGen->event)[eventSize].tau0() * log(randomEngine().flat());
0080         (fMasterGen->event)[eventSize].tau(tauTmp);
0081       }
0082 
0083       // Here also need to add anti-particle (if any)
0084       // otherwise just add a 2nd particle of the same type
0085       // (for example, gamma)
0086       //
0087       if (fAddAntiParticle) {
0088         if (1 <= std::abs(particleID) && std::abs(particleID) <= 6) {  // quarks
0089           (fMasterGen->event).append(-particleID, 23, 1, 0, 0, 0, 0, 101, -px, -py, -pz, ee, mass);
0090         } else if (std::abs(particleID) == 21) {  // gluons
0091           (fMasterGen->event).append(21, 23, 1, 0, 0, 0, 102, 101, -px, -py, -pz, ee, mass);
0092         } else {
0093           if ((fMasterGen->particleData).isParticle(-particleID)) {
0094             (fMasterGen->event).append(-particleID, 1, 1, 0, 0, 0, 0, 0, -px, -py, -pz, ee, mass);
0095           } else {
0096             (fMasterGen->event).append(particleID, 1, 1, 0, 0, 0, 0, 0, -px, -py, -pz, ee, mass);
0097           }
0098           int eventSize = (fMasterGen->event).size() - 1;
0099           // -log(flat) = exponential distribution
0100           double tauTmp = -(fMasterGen->event)[eventSize].tau0() * log(randomEngine().flat());
0101           (fMasterGen->event)[eventSize].tau(tauTmp);
0102         }
0103       }
0104     }
0105 
0106     if (!fMasterGen->next())
0107       return false;
0108     evtGenDecay();
0109 
0110     event() = std::make_unique<HepMC::GenEvent>();
0111     return toHepMC.fill_next_event(fMasterGen->event, event().get());
0112   }
0113 
0114   const char* Py8PtGun::classname() const { return "Py8PtGun"; }
0115 
0116   typedef edm::GeneratorFilter<gen::Py8PtGun, gen::ExternalDecayDriver> Pythia8PtGun;
0117 
0118 }  // namespace gen
0119 
0120 using gen::Pythia8PtGun;
0121 DEFINE_FWK_MODULE(Pythia8PtGun);