CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/IOMC/ParticleGuns/src/FlatRandomOneOverPtGunProducer.cc

Go to the documentation of this file.
00001 #include <ostream>
00002 
00003 #include "IOMC/ParticleGuns/interface/FlatRandomOneOverPtGunProducer.h"
00004 
00005 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00006 #include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
00007 
00008 #include "FWCore/Framework/interface/Event.h"
00009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00011 
00012 #include "CLHEP/Units/GlobalSystemOfUnits.h"
00013 using namespace edm;
00014 
00015 FlatRandomOneOverPtGunProducer::FlatRandomOneOverPtGunProducer(const edm::ParameterSet& pset) : 
00016   BaseFlatGunProducer(pset) {
00017 
00018 
00019   edm::ParameterSet defpset ;
00020   edm::ParameterSet pgun_params = 
00021     pset.getParameter<ParameterSet>("PGunParameters") ;
00022   
00023   fMinOneOverPt = pgun_params.getParameter<double>("MinOneOverPt");
00024   fMaxOneOverPt = pgun_params.getParameter<double>("MaxOneOverPt");
00025   
00026   produces<HepMCProduct>();
00027   produces<GenEventInfoProduct>();
00028 
00029   edm::LogInfo("ParticleGun") << "FlatRandomOneOverPtGunProducer: initialized with minimum and maximum 1/pt " << fMinOneOverPt << ":" << fMaxOneOverPt;
00030 }
00031 
00032 FlatRandomOneOverPtGunProducer::~FlatRandomOneOverPtGunProducer() {
00033   // no need to cleanup GenEvent memory - done in HepMCProduct
00034 }
00035 
00036 void FlatRandomOneOverPtGunProducer::produce(Event &e, const EventSetup& es) {
00037 
00038   LogDebug("ParticleGun") << " FlatRandomOneOverPtGunProducer : Begin New Event Generation"; 
00039 
00040   // event loop (well, another step in it...)
00041           
00042   // no need to clean up GenEvent memory - done in HepMCProduct
00043   // 
00044    
00045   // here re-create fEvt (memory)
00046   //
00047   fEvt = new HepMC::GenEvent() ;
00048    
00049   // now actualy, cook up the event from PDGTable and gun parameters
00050   //
00051   // 1st, primary vertex
00052   //
00053   HepMC::GenVertex* Vtx = new HepMC::GenVertex(HepMC::FourVector(0.,0.,0.));
00054 
00055   // loop over particles
00056   //
00057   int barcode = 1 ;
00058   for (unsigned int ip=0; ip<fPartIDs.size(); ++ip) {
00059 
00060     double xx     = fRandomGenerator->fire(0.0,1.0);
00061     double pt     = std::exp((1.-xx)*std::log(fMinOneOverPt)+
00062                              xx*std::log(fMaxOneOverPt)) ;
00063     double eta    = fRandomGenerator->fire(fMinEta, fMaxEta) ;
00064     double phi    = fRandomGenerator->fire(fMinPhi, fMaxPhi) ;
00065     if (pt != 0) pt = 1./pt;
00066     int PartID = fPartIDs[ip] ;
00067     const HepPDT::ParticleData* 
00068       PData = fPDGTable->particle(HepPDT::ParticleID(abs(PartID))) ;
00069     double mass   = PData->mass().value() ;
00070     double theta  = 2.*atan(exp(-eta)) ;
00071     double mom    = pt/sin(theta) ;
00072     double px     = pt*cos(phi) ;
00073     double py     = pt*sin(phi) ;
00074     double pz     = mom*cos(theta) ;
00075     double energy2= mom*mom + mass*mass ;
00076     double energy = sqrt(energy2) ; 
00077     HepMC::FourVector p(px,py,pz,energy) ;
00078     HepMC::GenParticle* Part = new HepMC::GenParticle(p,PartID,1);
00079     Part->suggest_barcode( barcode ) ;
00080     barcode++ ;
00081     Vtx->add_particle_out(Part);
00082     LogDebug("ParticleGun") << "FlatRandomOneOverPtGunProducer: Event generated with pt:eta:phi " << pt << " " << eta << " " << phi << " (" << theta/CLHEP::deg << ":" << phi/CLHEP::deg << ")";
00083 
00084     if ( fAddAntiParticle ) {
00085       HepMC::FourVector ap(-px,-py,-pz,energy) ;
00086       int APartID = -PartID ;
00087       if ( PartID == 22 || PartID == 23 ) {
00088         APartID = PartID ;
00089       }   
00090       HepMC::GenParticle* APart = new HepMC::GenParticle(ap,APartID,1);
00091       APart->suggest_barcode( barcode ) ;
00092       barcode++ ;
00093       Vtx->add_particle_out(APart) ;
00094     }
00095 
00096   }
00097 
00098   fEvt->add_vertex(Vtx) ;
00099   fEvt->set_event_number(e.id().event()) ;
00100   fEvt->set_signal_process_id(20) ; 
00101         
00102   if ( fVerbosity > 0 ) {
00103     fEvt->print() ;  
00104   }
00105 
00106   std::auto_ptr<HepMCProduct> BProduct(new HepMCProduct()) ;
00107   BProduct->addHepMCData( fEvt );
00108   e.put(BProduct);
00109 
00110   std::auto_ptr<GenEventInfoProduct> genEventInfo(new GenEventInfoProduct(fEvt));
00111   e.put(genEventInfo);
00112     
00113   LogDebug("ParticleGun") << " FlatRandomOneOverPtGunProducer : Event Generation Done ";
00114 }