CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/IOMC/ParticleGuns/src/FlatRandomPtThetaGunSource.cc

Go to the documentation of this file.
00001 #include <ostream>
00002 
00003 #include "IOMC/ParticleGuns/interface/FlatRandomPtThetaGunSource.h"
00004 
00005 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00006 
00007 #include "FWCore/Framework/interface/Event.h"
00008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00010 
00011 using namespace edm;
00012 
00013 FlatRandomPtThetaGunSource::FlatRandomPtThetaGunSource(const edm::ParameterSet& pset,
00014                                                        const edm::InputSourceDescription& desc) : 
00015   FlatBaseThetaGunSource(pset, desc) {
00016 
00017   edm::ParameterSet defpset ;
00018   edm::ParameterSet pgun_params = 
00019     pset.getUntrackedParameter<edm::ParameterSet>("PGunParameters",defpset) ;
00020   
00021   fMinPt = pgun_params.getUntrackedParameter<double>("MinPt",0.99);
00022   fMaxPt = pgun_params.getUntrackedParameter<double>("MaxPt",1.01);
00023   
00024   produces<HepMCProduct>();
00025   edm::LogInfo("FlatThetaGun") << "Internal FlatRandomPtThetaGun is initialzed"
00026                                << "\nIt is going to generate " 
00027                                << remainingEvents() << " events";
00028 
00029 }
00030 
00031 FlatRandomPtThetaGunSource::~FlatRandomPtThetaGunSource() {}
00032 
00033 bool FlatRandomPtThetaGunSource::produce(edm::Event &e) {
00034 
00035   if ( fVerbosity > 0 ) {
00036     LogDebug("FlatThetaGun") << "FlatRandomPtThetaGunSource : Begin New Event Generation"; 
00037   }
00038 
00039   // event loop (well, another step in it...)
00040   
00041   // no need to clean up GenEvent memory - done in HepMCProduct
00042   // 
00043    
00044   // here re-create fEvt (memory)
00045   //
00046   fEvt = new HepMC::GenEvent() ;
00047    
00048   // now actualy, cook up the event from PDGTable and gun parameters
00049   //
00050   // 1st, primary vertex
00051   //
00052   HepMC::GenVertex* Vtx = new HepMC::GenVertex(HepMC::FourVector(0.,0.,0.));
00053 
00054   // loop over particles
00055   //
00056   int barcode = 1 ;
00057   for (unsigned int ip=0; ip<fPartIDs.size(); ++ip) {
00058 
00059     double pt     = fRandomGenerator->fire(fMinPt, fMaxPt) ;
00060     double theta  = fRandomGenerator->fire(fMinTheta, fMaxTheta) ;
00061     double phi    = fRandomGenerator->fire(fMinPhi, fMaxPhi) ;
00062     int PartID = fPartIDs[ip] ;
00063     const HepPDT::ParticleData* 
00064       PData = fPDGTable->particle(HepPDT::ParticleID(abs(PartID))) ;
00065     double mass   = PData->mass().value() ;
00066     double mom    = pt/sin(theta) ;
00067     double px     = pt*cos(phi) ;
00068     double py     = pt*sin(phi) ;
00069     double pz     = mom*cos(theta) ;
00070     double energy2= mom*mom + mass*mass ;
00071     double energy = sqrt(energy2) ; 
00072     HepMC::FourVector p(px,py,pz,energy) ;
00073     HepMC::GenParticle* Part = new HepMC::GenParticle(p,PartID,1);
00074     Part->suggest_barcode( barcode ) ;
00075     barcode++ ;
00076     Vtx->add_particle_out(Part);
00077 
00078     if ( fAddAntiParticle ) {
00079       HepMC::FourVector ap(-px,-py,-pz,energy) ;
00080       int APartID = -PartID ;
00081       if ( PartID == 22 || PartID == 23 ) {
00082         APartID = PartID ;
00083       }   
00084       HepMC::GenParticle* APart = new HepMC::GenParticle(ap,APartID,1);
00085       APart->suggest_barcode( barcode ) ;
00086       barcode++ ;
00087       Vtx->add_particle_out(APart) ;
00088     }
00089   }
00090 
00091   fEvt->add_vertex(Vtx) ;
00092   fEvt->set_event_number(event()) ;
00093   fEvt->set_signal_process_id(20) ; 
00094         
00095   if ( fVerbosity > 0 ) {
00096     fEvt->print() ;  
00097   }
00098 
00099   std::auto_ptr<HepMCProduct> BProduct(new HepMCProduct()) ;
00100   BProduct->addHepMCData( fEvt );
00101   e.put(BProduct);
00102     
00103   if ( fVerbosity > 0 ) {
00104     LogDebug("FlatThetaGun") << "FlatRandomPtThetaGunSource : Event Generation Done ";
00105   }
00106 
00107   return true;
00108 }
00109