CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/GeneratorInterface/Pythia6Interface/plugins/Pythia6PartonEGun.cc

Go to the documentation of this file.
00001 
00002 #include <iostream>
00003 
00004 #include "Pythia6PartonEGun.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 Pythia6PartonEGun::Pythia6PartonEGun( 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    fMinE       = pgun_params.getParameter<double>("MinE"); // ,  20.);
00028    fMaxE       = pgun_params.getParameter<double>("MaxE"); // , 420.);
00029 
00030 }
00031 
00032 Pythia6PartonEGun::~Pythia6PartonEGun()
00033 {
00034 }
00035 
00036 void Pythia6PartonEGun::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 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    ee   = (fMaxE-fMinE)*pyr_(&dum)+fMinE;
00066    eta  = (fMaxEta-fMinEta)*pyr_(&dum)+fMinEta;                                                      
00067    the  = 2.*atan(exp(-eta));                                                                          
00068          
00069    py1ent_(ip, py6PID, ee, the, phi);
00070          
00071    double px     = pyjets.p[0][ip-1]; // pt*cos(phi) ;
00072    double py     = pyjets.p[1][ip-1]; // pt*sin(phi) ;
00073    double pz     = pyjets.p[2][ip-1]; // mom*cos(the) ;
00074          
00075    HepMC::FourVector p(px,py,pz,ee) ;
00076    HepMC::GenParticle* Part = new HepMC::GenParticle(p,fPartonID,1);
00077    Part->suggest_barcode( ip ) ;
00078    Vtx->add_particle_out(Part);
00079          
00080    // now add anti-quark
00081    ip = ip + 1;
00082    HepMC::GenParticle* APart = addAntiParticle( ip, fPartonID, ee, eta, phi );   
00083    if ( APart ) 
00084    {
00085       Vtx->add_particle_out(APart) ;
00086    }
00087    else
00088    {
00089       // otherwise it should throw !
00090    }        
00091 
00092    // this should probably be configurable...
00093    //
00094    double qmax = 2.*ee;
00095    
00096    joinPartons( qmax );
00097          
00098    fEvt->add_vertex(Vtx);
00099      
00100    // run pythia
00101    pyexec_();
00102    
00103    return;
00104 }
00105 
00106 DEFINE_FWK_MODULE(Pythia6PartonEGun);