CMS 3D CMS Logo

edm::FlatRandomPtGunSource Class Reference

#include <IOMC/ParticleGuns/interface/FlatRandomPtGunSource.h>

Inheritance diagram for edm::FlatRandomPtGunSource:

edm::BaseFlatGunSource edm::GeneratedInputSource edm::ConfigurableInputSource edm::InputSource edm::ProductRegistryHelper

List of all members.

Public Member Functions

 FlatRandomPtGunSource (const ParameterSet &, const InputSourceDescription &)
virtual ~FlatRandomPtGunSource ()

Protected Attributes

double fMaxPt
double fMinPt

Private Member Functions

virtual bool produce (Event &e)


Detailed Description

Definition at line 15 of file FlatRandomPtGunSource.h.


Constructor & Destructor Documentation

FlatRandomPtGunSource::FlatRandomPtGunSource ( const ParameterSet pset,
const InputSourceDescription desc 
)

Definition at line 20 of file FlatRandomPtGunSource.cc.

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

00021                                                                                  : 
00022    BaseFlatGunSource(pset, desc)
00023 {
00024 
00025 
00026    ParameterSet defpset ;
00027    ParameterSet pgun_params = 
00028       pset.getUntrackedParameter<ParameterSet>("PGunParameters",defpset) ;
00029   
00030    fMinPt = pgun_params.getUntrackedParameter<double>("MinPt",0.99);
00031    fMaxPt = pgun_params.getUntrackedParameter<double>("MaxPt",1.01);
00032   
00033   produces<HepMCProduct>();
00034    
00035 }

FlatRandomPtGunSource::~FlatRandomPtGunSource (  )  [virtual]

Definition at line 37 of file FlatRandomPtGunSource.cc.

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


Member Function Documentation

bool FlatRandomPtGunSource::produce ( Event e  )  [private, virtual]

Implements edm::ConfigurableInputSource.

Definition at line 42 of file FlatRandomPtGunSource.cc.

References funct::abs(), funct::cos(), GenMuonPlsPt100GeV_cfg::cout, lat::endl(), relval_parameters_module::energy, eta, edm::ConfigurableInputSource::event(), funct::exp(), edm::BaseFlatGunSource::fAddAntiParticle, edm::BaseFlatGunSource::fEvt, edm::BaseFlatGunSource::fMaxEta, edm::BaseFlatGunSource::fMaxPhi, fMaxPt, edm::BaseFlatGunSource::fMinEta, edm::BaseFlatGunSource::fMinPhi, fMinPt, edm::BaseFlatGunSource::fPartIDs, edm::BaseFlatGunSource::fPDGTable, edm::BaseFlatGunSource::fRandomGenerator, edm::BaseFlatGunSource::fVerbosity, p, gen_jpsi2muons_cfg::ParticleID, phi, edm::Event::put(), funct::sin(), funct::sqrt(), and theta.

00043 {
00044 
00045    if ( fVerbosity > 0 )
00046    {
00047       cout << " FlatRandomPtGunSource : Begin New Event Generation" << endl ; 
00048    }
00049    // event loop (well, another step in it...)
00050           
00051    // no need to clean up GenEvent memory - done in HepMCProduct
00052    // 
00053    
00054    // here re-create fEvt (memory)
00055    //
00056    fEvt = new HepMC::GenEvent() ;
00057    
00058    // now actualy, cook up the event from PDGTable and gun parameters
00059    //
00060    // 1st, primary vertex
00061    //
00062    HepMC::GenVertex* Vtx = new HepMC::GenVertex(HepMC::FourVector(0.,0.,0.));
00063 
00064    // loop over particles
00065    //
00066    int barcode = 1 ;
00067    for (unsigned int ip=0; ip<fPartIDs.size(); ++ip)
00068    {
00069 
00070        double pt     = fRandomGenerator->fire(fMinPt, fMaxPt) ;
00071        double eta    = fRandomGenerator->fire(fMinEta, fMaxEta) ;
00072        double phi    = fRandomGenerator->fire(fMinPhi, fMaxPhi) ;
00073        int PartID = fPartIDs[ip] ;
00074        const HepPDT::ParticleData* 
00075           PData = fPDGTable->particle(HepPDT::ParticleID(abs(PartID))) ;
00076        double mass   = PData->mass().value() ;
00077        double theta  = 2.*atan(exp(-eta)) ;
00078        double mom    = pt/sin(theta) ;
00079        double px     = pt*cos(phi) ;
00080        double py     = pt*sin(phi) ;
00081        double pz     = mom*cos(theta) ;
00082        double energy2= mom*mom + mass*mass ;
00083        double energy = sqrt(energy2) ; 
00084        HepMC::FourVector p(px,py,pz,energy) ;
00085        HepMC::GenParticle* Part = 
00086            new HepMC::GenParticle(p,PartID,1);
00087        Part->suggest_barcode( barcode ) ;
00088        barcode++ ;
00089        Vtx->add_particle_out(Part);
00090 
00091        if ( fAddAntiParticle )
00092        {
00093           HepMC::FourVector ap(-px,-py,-pz,energy) ;
00094           int APartID = -PartID ;
00095           if ( PartID == 22 || PartID == 23 )
00096           {
00097              APartID = PartID ;
00098           }       
00099           HepMC::GenParticle* APart =
00100              new HepMC::GenParticle(ap,APartID,1);
00101           APart->suggest_barcode( barcode ) ;
00102           barcode++ ;
00103           Vtx->add_particle_out(APart) ;
00104        }
00105 
00106    }
00107 
00108    fEvt->add_vertex(Vtx) ;
00109    fEvt->set_event_number(event()) ;
00110    fEvt->set_signal_process_id(20) ; 
00111         
00112    if ( fVerbosity > 0 )
00113    {
00114       fEvt->print() ;  
00115    }
00116 
00117    auto_ptr<HepMCProduct> BProduct(new HepMCProduct()) ;
00118    BProduct->addHepMCData( fEvt );
00119    e.put(BProduct);
00120     
00121    if ( fVerbosity > 0 )
00122    {
00123       // for testing purpose only
00124       // fEvt->print() ; // prints empty info after it's made into edm::Event
00125       cout << " FlatRandomPtGunSource : Event Generation Done " << endl;
00126    }
00127 
00128    return true;
00129 }


Member Data Documentation

double edm::FlatRandomPtGunSource::fMaxPt [protected]

Definition at line 31 of file FlatRandomPtGunSource.h.

Referenced by FlatRandomPtGunSource(), and produce().

double edm::FlatRandomPtGunSource::fMinPt [protected]

Definition at line 30 of file FlatRandomPtGunSource.h.

Referenced by FlatRandomPtGunSource(), and produce().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:41:02 2009 for CMSSW by  doxygen 1.5.4