CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/GeneratorInterface/Pythia6Interface/plugins/Pythia6JetGun.cc

Go to the documentation of this file.
00001 
00002 #include <iostream>
00003 
00004 #include "Pythia6JetGun.h"
00005 
00006 #include "FWCore/Utilities/interface/Exception.h"
00007 
00008 #include "FWCore/Framework/interface/EDProducer.h"
00009 #include "FWCore/Framework/interface/EventSetup.h"
00010 
00011 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00012 
00013 #include "FWCore/Framework/interface/MakerMacros.h"
00014 
00015 using namespace edm;
00016 using namespace gen;
00017 
00018 Pythia6JetGun::Pythia6JetGun( const ParameterSet& pset ) :
00019    Pythia6ParticleGun(pset),
00020    fMinEta(0.), fMaxEta(0.),
00021    fMinE(0.), fMaxE(0.),
00022    fMinP(0.), fMaxP(0.)
00023 {
00024    
00025    ParameterSet pgun_params = 
00026       pset.getParameter<ParameterSet>("PGunParameters"); 
00027    fMinEta     = pgun_params.getParameter<double>("MinEta"); 
00028    fMaxEta     = pgun_params.getParameter<double>("MaxEta");  
00029    fMinE = pgun_params.getParameter<double>("MinE");     
00030    fMaxE = pgun_params.getParameter<double>("MaxE");     
00031    fMinP = pgun_params.getParameter<double>("MinP"); 
00032    fMaxP = pgun_params.getParameter<double>("MaxP"); 
00033 
00034 }
00035 
00036 Pythia6JetGun::~Pythia6JetGun()
00037 {
00038 }
00039 
00040 void Pythia6JetGun::generateEvent()
00041 {
00042    Pythia6Service::InstanceWrapper guard(fPy6Service);  // grab Py6 instance
00043 
00044    // now actualy, start cooking up the event gun 
00045    //
00046 
00047    // 1st, primary vertex
00048    //
00049    HepMC::GenVertex* Vtx = new HepMC::GenVertex( HepMC::FourVector(0.,0.,0.));
00050 
00051    // here re-create fEvt (memory)
00052    //
00053    fEvt = new HepMC::GenEvent() ;
00054      
00055    int ip=1;
00056    double totPx = 0.;
00057    double totPy = 0.;
00058    double totPz = 0.;
00059    double totE  = 0.;
00060    double totM  = 0.;
00061    double phi, eta, the, ee, pp;
00062    int dum = 0;
00063    for ( size_t i=0; i<fPartIDs.size(); i++ )
00064    {
00065          int particleID = fPartIDs[i]; // this is PDG - need to convert to Py6 !!!
00066          int py6PID = HepPID::translatePDTtoPythia( particleID );
00067                  
00068          // internal numbers
00069          //
00070          phi = 2. * M_PI * pyr_(&dum);
00071          the = std::acos( -1. + 2.*pyr_(&dum) );
00072          
00073          // from input
00074          //
00075          ee   = (fMaxE-fMinE)*pyr_(&dum)+fMinE;
00076          
00077          // fill p(ip,5) (in PYJETS) with mass value right now,
00078          // because the (hardcoded) mstu(10)=1 will make py1ent
00079          // pick the mass from there
00080          double mass = pymass_(py6PID);
00081          pyjets.p[4][ip-1]=mass; 
00082 
00083          // add entry to py6
00084          //
00085          py1ent_(ip, py6PID, ee, the, phi);
00086 
00087          // values for computing total mass
00088          //
00089          totPx += pyjets.p[0][ip-1];
00090          totPy += pyjets.p[1][ip-1];
00091          totPz += pyjets.p[2][ip-1];
00092          totE  += pyjets.p[3][ip-1];
00093          
00094          ip++;
00095          
00096    } // end forming up py6 record of the jet
00097 
00098          
00099    // compute total mass
00100    //   
00101    totM = std::sqrt( totE*totE - (totPx*totPx+totPy*totPy+totPz*totPz) );
00102          
00103    //now the boost (from input params)
00104    //
00105    pp = (fMaxP-fMinP)*pyr_(&dum)+fMinP; 
00106    ee = std::sqrt( totM*totM + pp*pp );  
00107          
00108    //the boost direction (from input params)
00109    //
00110    phi = (fMaxPhi-fMinPhi)*pyr_(&dum)+fMinPhi;
00111    eta  = (fMaxEta-fMinEta)*pyr_(&dum)+fMinEta;                                                      
00112    the  = 2.*atan(exp(-eta));  
00113          
00114    double betaX = pp/ee * std::sin(the) * std::cos(phi);
00115    double betaY = pp/ee * std::sin(the) * std::sin(phi);
00116    double betaZ = pp/ee * std::cos(the);  
00117    
00118    // boost all particles
00119    // the first 2 params (-1) tell to boost all (fisrt-to-last),
00120    // and the next 2 params (0.) tell no rotation
00121    // 
00122    int first=-1, last=-1;
00123    double rothe=0, rophi=0.;
00124    
00125    pyrobo_( first, last, rothe, rophi, betaX, betaY, betaZ );                                                                      
00126                   
00127    // event should be formed from boosted record !!!
00128    // that's why additional loop
00129    //
00130    for ( int i=0; i<pyjets.n; i++ )
00131    {
00132       HepMC::FourVector p(pyjets.p[0][i],pyjets.p[1][i],pyjets.p[2][i],pyjets.p[3][i]) ;
00133       HepMC::GenParticle* Part = 
00134              new HepMC::GenParticle(p,
00135                                     HepPID::translatePythiatoPDT( pyjets.k[1][i] ),
00136                                     1); 
00137       Part->suggest_barcode( i+1 ) ;
00138       Vtx->add_particle_out(Part);
00139    }
00140    fEvt->add_vertex(Vtx);
00141      
00142    // run pythia
00143    pyexec_();
00144       
00145    return;
00146 }
00147 
00148 DEFINE_FWK_MODULE(Pythia6JetGun);