CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/GeneratorInterface/Pythia8Interface/plugins/Py8PtGun.cc

Go to the documentation of this file.
00001 
00002 #include "GeneratorInterface/Core/interface/GeneratorFilter.h"
00003 #include "GeneratorInterface/ExternalDecays/interface/ExternalDecayDriver.h"
00004 
00005 #include "GeneratorInterface/Pythia8Interface/interface/Py8GunBase.h"
00006 #include "GeneratorInterface/Pythia8Interface/interface/RandomP8.h"
00007 
00008 namespace gen {
00009 
00010 class Py8PtGun : public Py8GunBase {
00011    
00012    public:
00013       
00014       Py8PtGun( edm::ParameterSet const& );
00015       ~Py8PtGun() {}
00016          
00017       bool generatePartonsAndHadronize();
00018       const char* classname() const;
00019          
00020    private:
00021       
00022       // PtGun particle(s) characteristics
00023       double  fMinEta;
00024       double  fMaxEta;
00025       double  fMinPt ;
00026       double  fMaxPt ;
00027       bool    fAddAntiParticle;
00028 
00029 };
00030 
00031 // implementation 
00032 //
00033 Py8PtGun::Py8PtGun( edm::ParameterSet const& ps )
00034    : Py8GunBase(ps) 
00035 {
00036 
00037    // ParameterSet defpset ;
00038    edm::ParameterSet pgun_params = 
00039       ps.getParameter<edm::ParameterSet>("PGunParameters"); // , defpset ) ;
00040    fMinEta     = pgun_params.getParameter<double>("MinEta"); // ,-2.2);
00041    fMaxEta     = pgun_params.getParameter<double>("MaxEta"); // , 2.2);
00042    fMinPt      = pgun_params.getParameter<double>("MinPt"); // ,  0.);
00043    fMaxPt      = pgun_params.getParameter<double>("MaxPt"); // ,  0.);
00044    fAddAntiParticle = pgun_params.getParameter<bool>("AddAntiParticle"); //, false) ;  
00045 
00046 }
00047 
00048 bool Py8PtGun::generatePartonsAndHadronize()
00049 {
00050 
00051    fMasterGen->event.reset();
00052    
00053    for ( size_t i=0; i<fPartIDs.size(); i++ )
00054    {
00055 
00056       int particleID = fPartIDs[i]; // this is PDG - need to convert to Py8 ???
00057 
00058       // FIXME !!!
00059       // Ouch, it's using bare randomEngine pointer - that's NOT safe.
00060       // Need to hold a pointer somewhere properly !!!
00061       //
00062       double phi = (fMaxPhi-fMinPhi) * randomEngine->flat() + fMinPhi;
00063       double eta  = (fMaxEta-fMinEta) * randomEngine->flat() + fMinEta;                                                      
00064       double the  = 2.*atan(exp(-eta));                                                                          
00065 
00066       double pt   = (fMaxPt-fMinPt) * randomEngine->flat() + fMinPt;
00067       
00068       double mass = (fMasterGen->particleData).mass( particleID );
00069 //      double mass = (pythia->particleData).m0( particleID );
00070 
00071       double pp = pt / sin(the); // sqrt( ee*ee - mass*mass );
00072       double ee = sqrt( pp*pp + mass*mass );
00073       
00074       double px = pt * cos(phi);
00075       double py = pt * sin(phi);
00076       double pz = pp * cos(the);
00077 
00078       if ( !((fMasterGen->particleData).isParticle( particleID )) )
00079       {
00080          particleID = std::fabs(particleID) ;
00081       }
00082       (fMasterGen->event).append( particleID, 1, 0, 0, px, py, pz, ee, mass ); 
00083 
00084 // Here also need to add anti-particle (if any)
00085 // otherwise just add a 2nd particle of the same type 
00086 // (for example, gamma)
00087 //
00088       if ( fAddAntiParticle )
00089       {
00090          if ( (fMasterGen->particleData).isParticle( -particleID ) )
00091          {
00092             (fMasterGen->event).append( -particleID, 1, 0, 0, px, py, pz, ee, mass );
00093          }
00094          else
00095          {
00096             (fMasterGen->event).append( particleID, 1, 0, 0, px, py, pz, ee, mass );
00097          }
00098       }
00099 
00100    }
00101    
00102    if ( !fMasterGen->next() ) return false;
00103    
00104    event().reset(new HepMC::GenEvent);
00105    toHepMC.fill_next_event( fMasterGen->event, event().get() );
00106       
00107    return true;   
00108   
00109 }
00110 
00111 const char* Py8PtGun::classname() const
00112 {
00113    return "Py8PtGun"; 
00114 }
00115 
00116 typedef edm::GeneratorFilter<gen::Py8PtGun, gen::ExternalDecayDriver> Pythia8PtGun;
00117 
00118 } // end namespace
00119 
00120 using gen::Pythia8PtGun;
00121 DEFINE_FWK_MODULE(Pythia8PtGun);