CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2/src/GeneratorInterface/Pythia6Interface/plugins/Pythia6PtYDistGun.cc

Go to the documentation of this file.
00001 
00002 #include <iostream>
00003 
00004 #include "Pythia6PtYDistGun.h"
00005 #include "GeneratorInterface/Pythia6Interface/interface/PtYDistributor.h"
00006 
00007 #include "FWCore/Utilities/interface/Exception.h"
00008 
00009 #include "FWCore/Framework/interface/EDProducer.h"
00010 #include "FWCore/Framework/interface/EventSetup.h"
00011 
00012 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00013 
00014 #include "GeneratorInterface/Core/interface/RNDMEngineAccess.h"
00015 
00016 #include "FWCore/Framework/interface/MakerMacros.h"
00017 
00018 using namespace edm;
00019 using namespace gen;
00020 
00021 Pythia6PtYDistGun::Pythia6PtYDistGun( const ParameterSet& pset ) :
00022    Pythia6ParticleGun(pset), fPtYGenerator(0)
00023 {
00024    
00025    ParameterSet pgun_params = 
00026       pset.getParameter<ParameterSet>("PGunParameters"); 
00027    
00028    fPtYGenerator = new PtYDistributor( pgun_params.getParameter<FileInPath>("kinematicsFile"), 
00029                                        getEngineReference(), 
00030                                        pgun_params.getParameter<double>("MaxPt"),
00031                                        pgun_params.getParameter<double>("MinPt"),
00032                                        pgun_params.getParameter<double>("MaxY"), 
00033                                        pgun_params.getParameter<double>("MinY"), 
00034                                        pgun_params.getParameter<int>("PtBinning"), 
00035                                        pgun_params.getParameter<int>("YBinning") );
00036 }
00037 
00038 Pythia6PtYDistGun::~Pythia6PtYDistGun()
00039 {
00040    if ( fPtYGenerator ) delete fPtYGenerator;
00041 }
00042 
00043 void Pythia6PtYDistGun::generateEvent()
00044 {
00045    Pythia6Service::InstanceWrapper guard(fPy6Service);  // grab Py6 instance
00046 
00047    // now actualy, start cooking up the event gun 
00048    //
00049 
00050    // 1st, primary vertex
00051    //
00052    HepMC::GenVertex* Vtx = new HepMC::GenVertex( HepMC::FourVector(0.,0.,0.) );
00053 
00054    // here re-create fEvt (memory)
00055    //
00056    fEvt = new HepMC::GenEvent() ;
00057      
00058    int ip=1;
00059 
00060    for ( size_t i=0; i<fPartIDs.size(); i++ )
00061    {
00062 
00063          int particleID = fPartIDs[i]; // this is PDG - need to convert to Py6 !!!
00064          int py6PID = HepPID::translatePDTtoPythia( particleID );
00065          
00066          int dum = 0;
00067          double pt=0, y=0, u=0, ee=0, the=0;
00068          
00069          pt = fPtYGenerator->firePt();
00070          y  = fPtYGenerator->fireY();
00071          u  = exp(y); 
00072          
00073          double mass = pymass_(py6PID);
00074          
00075          // fill p(ip,5) (in PYJETS) with mass value right now,
00076          // because the (hardcoded) mstu(10)=1 will make py1ent
00077          // pick the mass from there
00078          pyjets.p[4][ip-1]=mass; 
00079 
00080          ee = 0.5*std::sqrt(mass*mass+pt*pt)*(u*u+1)/u;
00081          
00082          double pz = std::sqrt(ee*ee-pt*pt-mass*mass);
00083          if ( y < 0. ) pz *= -1;
00084          
00085          double phi = (fMaxPhi-fMinPhi)*pyr_(&dum)+fMinPhi;
00086 
00087          the  = atan(pt/pz);
00088          if ( pz < 0. ) the += M_PI ;                                                                      
00089          
00090          py1ent_(ip, py6PID, ee, the, phi);
00091          
00092 /*
00093          double px     = pt*cos(phi) ;
00094          double py     = pt*sin(phi) ;
00095 */         
00096          double px     = pyjets.p[0][ip-1]; // pt*cos(phi) ;
00097          double py     = pyjets.p[1][ip-1]; // pt*sin(phi) ;
00098          pz     = pyjets.p[2][ip-1]; // mom*cos(the) ;
00099 
00100          HepMC::FourVector p(px,py,pz,ee) ;
00101          HepMC::GenParticle* Part = 
00102              new HepMC::GenParticle(p,particleID,1);
00103          Part->suggest_barcode( ip ) ;
00104          Vtx->add_particle_out(Part);
00105          
00106          ip++;
00107    
00108    }
00109   
00110    fEvt->add_vertex(Vtx);
00111      
00112    // run pythia
00113    pyexec_();
00114    
00115    return;
00116 }
00117 
00118 DEFINE_FWK_MODULE(Pythia6PtYDistGun);