CMS 3D CMS Logo

Public Member Functions | Private Attributes

edm::FlatRandomPtGunProducer Class Reference

#include <FlatRandomPtGunProducer.h>

Inheritance diagram for edm::FlatRandomPtGunProducer:
edm::BaseFlatGunProducer edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

 FlatRandomPtGunProducer (const ParameterSet &pset)
virtual void produce (Event &e, const EventSetup &es)
virtual ~FlatRandomPtGunProducer ()

Private Attributes

double fMaxPt
double fMinPt

Detailed Description

Definition at line 15 of file FlatRandomPtGunProducer.h.


Constructor & Destructor Documentation

FlatRandomPtGunProducer::FlatRandomPtGunProducer ( const ParameterSet pset)

Definition at line 21 of file FlatRandomPtGunProducer.cc.

References fMaxPt, fMinPt, and edm::ParameterSet::getParameter().

                                                                         : 
   BaseFlatGunProducer(pset)
{


   ParameterSet defpset ;
   ParameterSet pgun_params = 
      pset.getParameter<ParameterSet>("PGunParameters") ;
  
   fMinPt = pgun_params.getParameter<double>("MinPt");
   fMaxPt = pgun_params.getParameter<double>("MaxPt");
  
  produces<HepMCProduct>();
  produces<GenEventInfoProduct>();
}
FlatRandomPtGunProducer::~FlatRandomPtGunProducer ( ) [virtual]

Definition at line 37 of file FlatRandomPtGunProducer.cc.

{
   // no need to cleanup GenEvent memory - done in HepMCProduct
}

Member Function Documentation

void FlatRandomPtGunProducer::produce ( Event e,
const EventSetup es 
) [virtual]

Implements edm::EDProducer.

Definition at line 42 of file FlatRandomPtGunProducer.cc.

References abs, funct::cos(), gather_cfg::cout, relval_parameters_module::energy, eta(), edm::EventID::event(), funct::exp(), edm::BaseFlatGunProducer::fAddAntiParticle, edm::BaseFlatGunProducer::fEvt, edm::BaseFlatGunProducer::fMaxEta, edm::BaseFlatGunProducer::fMaxPhi, fMaxPt, edm::BaseFlatGunProducer::fMinEta, edm::BaseFlatGunProducer::fMinPhi, fMinPt, edm::BaseFlatGunProducer::fPartIDs, edm::BaseFlatGunProducer::fPDGTable, edm::BaseFlatGunProducer::fRandomGenerator, edm::BaseFlatGunProducer::fVerbosity, configurableAnalysis::GenParticle, edm::EventBase::id(), L1TEmulatorMonitor_cff::p, RecoTau_DiTaus_pt_20-420_cfg::ParticleID, phi, ExpressReco_HICollisions_FallBack::pt, edm::Event::put(), funct::sin(), mathSSE::sqrt(), and theta().

{

   if ( fVerbosity > 0 )
   {
      cout << " FlatRandomPtGunProducer : Begin New Event Generation" << endl ; 
   }
   // event loop (well, another step in it...)
          
   // no need to clean up GenEvent memory - done in HepMCProduct
   // 
   
   // here re-create fEvt (memory)
   //
   fEvt = new HepMC::GenEvent() ;
   
   // now actualy, cook up the event from PDGTable and gun parameters
   //
   // 1st, primary vertex
   //
   HepMC::GenVertex* Vtx = new HepMC::GenVertex(HepMC::FourVector(0.,0.,0.));

   // loop over particles
   //
   int barcode = 1 ;
   for (unsigned int ip=0; ip<fPartIDs.size(); ++ip)
   {

       double pt     = fRandomGenerator->fire(fMinPt, fMaxPt) ;
       double eta    = fRandomGenerator->fire(fMinEta, fMaxEta) ;
       double phi    = fRandomGenerator->fire(fMinPhi, fMaxPhi) ;
       int PartID = fPartIDs[ip] ;
       const HepPDT::ParticleData* 
          PData = fPDGTable->particle(HepPDT::ParticleID(abs(PartID))) ;
       double mass   = PData->mass().value() ;
       double theta  = 2.*atan(exp(-eta)) ;
       double mom    = pt/sin(theta) ;
       double px     = pt*cos(phi) ;
       double py     = pt*sin(phi) ;
       double pz     = mom*cos(theta) ;
       double energy2= mom*mom + mass*mass ;
       double energy = sqrt(energy2) ; 
       HepMC::FourVector p(px,py,pz,energy) ;
       HepMC::GenParticle* Part = 
           new HepMC::GenParticle(p,PartID,1);
       Part->suggest_barcode( barcode ) ;
       barcode++ ;
       Vtx->add_particle_out(Part);

       if ( fAddAntiParticle )
       {
          HepMC::FourVector ap(-px,-py,-pz,energy) ;
          int APartID = -PartID ;
          if ( PartID == 22 || PartID == 23 )
          {
             APartID = PartID ;
          }       
          HepMC::GenParticle* APart =
             new HepMC::GenParticle(ap,APartID,1);
          APart->suggest_barcode( barcode ) ;
          barcode++ ;
          Vtx->add_particle_out(APart) ;
       }

   }

   fEvt->add_vertex(Vtx) ;
   fEvt->set_event_number(e.id().event()) ;
   fEvt->set_signal_process_id(20) ; 
        
   if ( fVerbosity > 0 )
   {
      fEvt->print() ;  
   }

   auto_ptr<HepMCProduct> BProduct(new HepMCProduct()) ;
   BProduct->addHepMCData( fEvt );
   e.put(BProduct);

   auto_ptr<GenEventInfoProduct> genEventInfo(new GenEventInfoProduct(fEvt));
   e.put(genEventInfo);
    
   if ( fVerbosity > 0 )
   {
      // for testing purpose only
      // fEvt->print() ; // prints empty info after it's made into edm::Event
      cout << " FlatRandomPtGunProducer : Event Generation Done " << endl;
   }
}

Member Data Documentation

Definition at line 29 of file FlatRandomPtGunProducer.h.

Referenced by FlatRandomPtGunProducer(), and produce().

Definition at line 28 of file FlatRandomPtGunProducer.h.

Referenced by FlatRandomPtGunProducer(), and produce().