CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2_patch1/src/GeneratorInterface/Pythia6Interface/plugins/Pythia6PartonPtGun.cc

Go to the documentation of this file.
00001 
00002 #include <iostream>
00003 
00004 #include "Pythia6PartonPtGun.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 Pythia6PartonPtGun::Pythia6PartonPtGun( const ParameterSet& pset ) :
00019    Pythia6PartonGun(pset)
00020 {
00021    
00022    // ParameterSet defpset ;
00023    ParameterSet pgun_params = 
00024       pset.getParameter<ParameterSet>("PGunParameters"); //, defpset ) ;
00025    fMinEta     = pgun_params.getParameter<double>("MinEta"); // ,-2.2);
00026    fMaxEta     = pgun_params.getParameter<double>("MaxEta"); // , 2.2);
00027    fMinPt       = pgun_params.getParameter<double>("MinPt"); // ,  20.);
00028    fMaxPt       = pgun_params.getParameter<double>("MaxPt"); // , 420.);
00029 
00030 }
00031 
00032 Pythia6PartonPtGun::~Pythia6PartonPtGun()
00033 {
00034 }
00035 
00036 void Pythia6PartonPtGun::generateEvent()
00037 {
00038    
00039    Pythia6Service::InstanceWrapper guard(fPy6Service);  // grab Py6 instance
00040 
00041    // now actualy, start cooking up the event gun 
00042    //
00043 
00044    // 1st, primary vertex
00045    //
00046    HepMC::GenVertex* Vtx = new HepMC::GenVertex( HepMC::FourVector(0.,0.,0.));
00047 
00048    // here re-create fEvt (memory)
00049    //
00050    fEvt = new HepMC::GenEvent() ;
00051      
00052    int ip=1;
00053    
00054    int py6PID = HepPID::translatePDTtoPythia( fPartonID );
00055    int dum = 0;
00056    double pt=0, mom=0, ee=0, the=0, eta=0;
00057    double mass = pymass_(py6PID);
00058          
00059    // fill p(ip,5) (in PYJETS) with mass value right now, 
00060    // because the (hardcoded) mstu(10)=1 will make py1ent
00061    // pick the mass from there
00062    pyjets.p[4][ip-1]=mass; 
00063                  
00064    double phi = (fMaxPhi-fMinPhi)*pyr_(&dum)+fMinPhi;
00065 
00066    eta  = (fMaxEta-fMinEta)*pyr_(&dum)+fMinEta;                                                      
00067 
00068    the  = 2.*atan(exp(-eta));                                                                          
00069          
00070    pt = (fMaxPt-fMinPt)*pyr_(&dum)+fMinPt;
00071          
00072    mom = pt/sin(the);
00073    ee = sqrt(mom*mom+mass*mass);
00074 
00075    py1ent_(ip, py6PID, ee, the, phi);
00076          
00077    double px     = pyjets.p[0][ip-1]; // pt*cos(phi) ;
00078    double py     = pyjets.p[1][ip-1]; // pt*sin(phi) ;
00079    double pz     = pyjets.p[2][ip-1]; // mom*cos(the) ;
00080          
00081    HepMC::FourVector p(px,py,pz,ee) ;
00082    HepMC::GenParticle* Part = new HepMC::GenParticle(p,fPartonID,1);
00083    Part->suggest_barcode( ip ) ;
00084    Vtx->add_particle_out(Part);
00085          
00086    // now add anti-quark
00087    ip = ip + 1;   
00088    HepMC::GenParticle* APart = addAntiParticle( ip, fPartonID, ee, eta, phi );   
00089    if ( APart ) 
00090    {
00091       Vtx->add_particle_out(APart) ;
00092    }
00093    else
00094    {
00095       // otherwise it should throw !
00096    }        
00097 
00098    // this should probably be configurable...
00099    //
00100    double qmax = 2.*ee;
00101    
00102    joinPartons( qmax );
00103          
00104    fEvt->add_vertex(Vtx);
00105      
00106    // run pythia
00107    pyexec_();
00108    
00109    return;
00110 }
00111 
00112 DEFINE_FWK_MODULE(Pythia6PartonPtGun);