CMS 3D CMS Logo

Public Member Functions | Private Attributes

edm::FlatRandomOneOverPtGunProducer Class Reference

#include <FlatRandomOneOverPtGunProducer.h>

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

List of all members.

Public Member Functions

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

Private Attributes

double fMaxOneOverPt
double fMinOneOverPt

Detailed Description

Definition at line 14 of file FlatRandomOneOverPtGunProducer.h.


Constructor & Destructor Documentation

FlatRandomOneOverPtGunProducer::FlatRandomOneOverPtGunProducer ( const ParameterSet pset)

Definition at line 15 of file FlatRandomOneOverPtGunProducer.cc.

References fMaxOneOverPt, fMinOneOverPt, and edm::ParameterSet::getParameter().

                                                                                          : 
  BaseFlatGunProducer(pset) {


  edm::ParameterSet defpset ;
  edm::ParameterSet pgun_params = 
    pset.getParameter<ParameterSet>("PGunParameters") ;
  
  fMinOneOverPt = pgun_params.getParameter<double>("MinOneOverPt");
  fMaxOneOverPt = pgun_params.getParameter<double>("MaxOneOverPt");
  
  produces<HepMCProduct>();
  produces<GenEventInfoProduct>();

  edm::LogInfo("ParticleGun") << "FlatRandomOneOverPtGunProducer: initialized with minimum and maximum 1/pt " << fMinOneOverPt << ":" << fMaxOneOverPt;
}
FlatRandomOneOverPtGunProducer::~FlatRandomOneOverPtGunProducer ( ) [virtual]

Definition at line 32 of file FlatRandomOneOverPtGunProducer.cc.

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

Member Function Documentation

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

Implements edm::EDProducer.

Definition at line 36 of file FlatRandomOneOverPtGunProducer.cc.

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

                                                                           {

  LogDebug("ParticleGun") << " FlatRandomOneOverPtGunProducer : Begin New Event Generation"; 

  // 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 xx     = fRandomGenerator->fire(0.0,1.0);
    double pt     = std::exp((1.-xx)*std::log(fMinOneOverPt)+
                             xx*std::log(fMaxOneOverPt)) ;
    double eta    = fRandomGenerator->fire(fMinEta, fMaxEta) ;
    double phi    = fRandomGenerator->fire(fMinPhi, fMaxPhi) ;
    if (pt != 0) pt = 1./pt;
    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);
    LogDebug("ParticleGun") << "FlatRandomOneOverPtGunProducer: Event generated with pt:eta:phi " << pt << " " << eta << " " << phi << " (" << theta/CLHEP::deg << ":" << phi/CLHEP::deg << ")";

    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() ;  
  }

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

  std::auto_ptr<GenEventInfoProduct> genEventInfo(new GenEventInfoProduct(fEvt));
  e.put(genEventInfo);
    
  LogDebug("ParticleGun") << " FlatRandomOneOverPtGunProducer : Event Generation Done ";
}

Member Data Documentation

Definition at line 28 of file FlatRandomOneOverPtGunProducer.h.

Referenced by FlatRandomOneOverPtGunProducer(), and produce().

Definition at line 27 of file FlatRandomOneOverPtGunProducer.h.

Referenced by FlatRandomOneOverPtGunProducer(), and produce().