CMS 3D CMS Logo

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

Go to the documentation of this file.
00001 #include <ostream>
00002 
00003 #include "IOMC/ParticleGuns/interface/FlatRandomEThetaGunSource.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 FlatRandomEThetaGunSource::FlatRandomEThetaGunSource(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   // doesn't seem necessary to check if pset is empty - if this
00022   // is the case, default values will be taken for params
00023   fMinE = pgun_params.getUntrackedParameter<double>("MinE",0.99);
00024   fMaxE = pgun_params.getUntrackedParameter<double>("MaxE",1.01); 
00025   
00026   produces<HepMCProduct>();
00027 
00028   edm::LogInfo("FlatThetaGun") << "Internal FlatRandomEThetaGun is initialzed"
00029                                << "\nIt is going to generate " 
00030                                << remainingEvents() << " events";
00031 }
00032 
00033 FlatRandomEThetaGunSource::~FlatRandomEThetaGunSource() {}
00034 
00035 bool FlatRandomEThetaGunSource::produce(edm::Event & e) {
00036 
00037   if ( fVerbosity > 0 ) {
00038     LogDebug("FlatThetaGun") << "FlatRandomEThetaGunSource : Begin New Event Generation"; 
00039   }
00040    
00041   // event loop (well, another step in it...)
00042           
00043   // no need to clean up GenEvent memory - done in HepMCProduct
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 
00052   // 1st, primary vertex
00053   //
00054   HepMC::GenVertex* Vtx = new HepMC::GenVertex( HepMC::FourVector(0.,0.,0.));
00055    
00056   // loop over particles
00057   //
00058   int barcode = 1;
00059   for (unsigned int ip=0; ip<fPartIDs.size(); ip++) {
00060     double energy = fRandomGenerator->fire(fMinE, fMaxE) ;
00061     double theta  = fRandomGenerator->fire(fMinTheta, fMaxTheta) ;
00062     double phi    = fRandomGenerator->fire(fMinPhi, fMaxPhi) ;
00063     int PartID = fPartIDs[ip] ;
00064     const HepPDT::ParticleData* 
00065       PData = fPDGTable->particle(HepPDT::ParticleID(abs(PartID))) ;
00066     double mass   = PData->mass().value() ;
00067     double mom2   = energy*energy - mass*mass ;
00068     double mom    = (mom2 > 0. ? std::sqrt(mom2) : 0.);
00069     double px     = mom*sin(theta)*cos(phi) ;
00070     double py     = mom*sin(theta)*sin(phi) ;
00071     double pz     = mom*cos(theta) ;
00072 
00073     HepMC::FourVector p(px,py,pz,energy) ;
00074     HepMC::GenParticle* Part =  new HepMC::GenParticle(p,PartID,1);
00075     Part->suggest_barcode( barcode ) ;
00076     barcode++ ;
00077     Vtx->add_particle_out(Part);
00078        
00079     if ( fAddAntiParticle ) {
00080       HepMC::FourVector ap(-px,-py,-pz,energy) ;
00081       int APartID = -PartID ;
00082       if ( PartID == 22 || PartID == 23 ) {
00083         APartID = PartID ;
00084       }
00085       HepMC::GenParticle* APart = new HepMC::GenParticle(ap,APartID,1);
00086       APart->suggest_barcode( barcode ) ;
00087       barcode++ ;
00088       Vtx->add_particle_out(APart) ;
00089     }
00090        
00091   }
00092   fEvt->add_vertex(Vtx) ;
00093   fEvt->set_event_number(event()) ;
00094   fEvt->set_signal_process_id(20) ;  
00095    
00096    
00097   if ( fVerbosity > 0 ) {
00098     fEvt->print() ;  
00099   }  
00100 
00101   std::auto_ptr<HepMCProduct> BProduct(new HepMCProduct()) ;
00102   BProduct->addHepMCData( fEvt );
00103   e.put(BProduct);
00104     
00105   if ( fVerbosity > 0 ) {
00106     LogDebug("FlatThetaGun") << "FlatRandomEThetaGunSource : Event Generation Done";
00107    }
00108   
00109   return true;
00110 }