CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_7/src/GeneratorInterface/Pythia8Interface/plugins/Py8JetGun.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 Py8JetGun : public Py8GunBase {
00011    
00012    public:
00013       
00014       Py8JetGun( edm::ParameterSet const& );
00015       ~Py8JetGun() {}
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  fMinP ;
00026       double  fMaxP ;
00027       double  fMinE ;
00028       double  fMaxE ;
00029 
00030 };
00031 
00032 // implementation 
00033 //
00034 Py8JetGun::Py8JetGun( edm::ParameterSet const& ps )
00035    : Py8GunBase(ps) 
00036 {
00037 
00038    // ParameterSet defpset ;
00039    edm::ParameterSet pgun_params = 
00040       ps.getParameter<edm::ParameterSet>("PGunParameters"); // , defpset ) ;
00041    fMinEta     = pgun_params.getParameter<double>("MinEta"); // ,-2.2);
00042    fMaxEta     = pgun_params.getParameter<double>("MaxEta"); // , 2.2);
00043    fMinP       = pgun_params.getParameter<double>("MinP"); // ,  0.);
00044    fMaxP       = pgun_params.getParameter<double>("MaxP"); // ,  0.);
00045    fMinE       = pgun_params.getParameter<double>("MinE"); // ,  0.);
00046    fMaxE       = pgun_params.getParameter<double>("MaxE"); // ,  0.);
00047 
00048 }
00049 
00050 bool Py8JetGun::generatePartonsAndHadronize()
00051 {
00052 
00053    fMasterGen->event.reset();
00054 
00055    double totPx = 0.;
00056    double totPy = 0.;
00057    double totPz = 0.;
00058    double totE  = 0.;
00059    double totM  = 0.;
00060    double phi, eta, the, ee, pp;
00061    
00062    for ( size_t i=0; i<fPartIDs.size(); i++ )
00063    {
00064 
00065       int particleID = fPartIDs[i]; // this is PDG - need to convert to Py8 ???
00066 
00067       // FIXME !!!
00068       // Ouch, it's using bare randomEngine pointer - that's NOT safe.
00069       // Need to hold a pointer somewhere properly !!!
00070       //
00071       phi = 2. * M_PI * randomEngine->flat() ;
00072       the = acos( -1. + 2.*randomEngine->flat() );
00073 
00074       // from input
00075       //
00076       ee   = (fMaxE-fMinE)*randomEngine->flat() + fMinE;
00077             
00078       double mass = (fMasterGen->particleData).mass( particleID );
00079 //      double mass = (pythia->particleData).m0( particleID );
00080 
00081       pp = sqrt( ee*ee - mass*mass );
00082       
00083       double px = pp * sin(the) * cos(phi);
00084       double py = pp * sin(the) * sin(phi);
00085       double pz = pp * cos(the);
00086 
00087       if ( !((fMasterGen->particleData).isParticle( particleID )) )
00088       {
00089          particleID = std::fabs(particleID) ;
00090       }
00091       (fMasterGen->event).append( particleID, 1, 0, 0, px, py, pz, ee, mass ); 
00092       
00093       // values for computing total mass
00094       //
00095       totPx += px;
00096       totPy += py;
00097       totPz += pz;
00098       totE  += ee;
00099 
00100    }
00101 
00102    totM = sqrt( totE*totE - (totPx*totPx+totPy*totPy+totPz*totPz) );
00103 
00104    //now the boost (from input params)
00105    //
00106    pp = (fMaxP-fMinP)*randomEngine->flat() + fMinP; 
00107    ee = sqrt( totM*totM + pp*pp );       
00108 
00109    //the boost direction (from input params)
00110    //
00111    phi = (fMaxPhi-fMinPhi)*randomEngine->flat() + fMinPhi;
00112    eta  = (fMaxEta-fMinEta)*randomEngine->flat() + fMinEta;                                                      
00113    the  = 2.*atan(exp(-eta));  
00114 
00115    double betaX = pp/ee * std::sin(the) * std::cos(phi);
00116    double betaY = pp/ee * std::sin(the) * std::sin(phi);
00117    double betaZ = pp/ee * std::cos(the);  
00118 
00119    // boost all particles
00120    //   
00121    (fMasterGen->event).bst( betaX, betaY, betaZ );
00122    
00123    if ( !fMasterGen->next() ) return false;
00124    
00125    event().reset(new HepMC::GenEvent);
00126    toHepMC.fill_next_event( fMasterGen->event, event().get() );
00127       
00128    return true;   
00129   
00130 }
00131 
00132 const char* Py8JetGun::classname() const
00133 {
00134    return "Py8JetGun"; 
00135 }
00136 
00137 typedef edm::GeneratorFilter<gen::Py8JetGun, gen::ExternalDecayDriver> Pythia8JetGun;
00138 
00139 } // end namespace
00140 
00141 using gen::Pythia8JetGun;
00142 DEFINE_FWK_MODULE(Pythia8JetGun);