CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch2/src/IOMC/ParticleGuns/src/FlatRandomEGunSource.cc

Go to the documentation of this file.
00001 /*
00002  *  $Date: 2009/01/09 10:23:16 $
00003  *  $Revision: 1.17 $
00004  *  \author Julia Yarba
00005  */
00006 
00007 #include <ostream>
00008 
00009 #include "IOMC/ParticleGuns/interface/FlatRandomEGunSource.h"
00010 
00011 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00012 
00013 #include "FWCore/Framework/interface/Event.h"
00014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00015 
00016 using namespace edm;
00017 using namespace std;
00018 
00019 FlatRandomEGunSource::FlatRandomEGunSource(const ParameterSet& pset,
00020                                            const InputSourceDescription& desc) :
00021    BaseFlatGunSource(pset, desc)
00022 {
00023 
00024    ParameterSet defpset ;
00025    // ParameterSet pgun_params = pset.getParameter<ParameterSet>("PGunParameters") ;
00026    ParameterSet pgun_params = 
00027       pset.getUntrackedParameter<ParameterSet>("PGunParameters",defpset) ;
00028   
00029    // doesn't seem necessary to check if pset is empty - if this
00030    // is the case, default values will be taken for params
00031    fMinE = pgun_params.getUntrackedParameter<double>("MinE",0.99);
00032    fMaxE = pgun_params.getUntrackedParameter<double>("MaxE",1.01); 
00033   
00034    produces<HepMCProduct>();
00035 
00036    cout << "Internal FlatRandomEGun is initialzed" << endl ;
00037    cout << "It is going to generate " << remainingEvents() << "events" << endl ;
00038    
00039 }
00040 
00041 FlatRandomEGunSource::~FlatRandomEGunSource()
00042 {
00043    // no need to cleanup fEvt since it's done in HepMCProduct
00044 }
00045 
00046 bool FlatRandomEGunSource::produce(Event & e) 
00047 {
00048 
00049    if ( fVerbosity > 0 )
00050    {
00051       cout << " FlatRandomEGunSource : Begin New Event Generation" << endl ; 
00052    }
00053    
00054    // event loop (well, another step in it...)
00055           
00056    // no need to clean up GenEvent memory - done in HepMCProduct
00057 
00058    // here re-create fEvt (memory)
00059    //
00060    fEvt = new HepMC::GenEvent() ;
00061    
00062    // now actualy, cook up the event from PDGTable and gun parameters
00063    //
00064 
00065    // 1st, primary vertex
00066    //
00067    HepMC::GenVertex* Vtx = new HepMC::GenVertex( HepMC::FourVector(0.,0.,0.));
00068    
00069    // loop over particles
00070    //
00071    int barcode = 1;
00072    for (unsigned int ip=0; ip<fPartIDs.size(); ip++)
00073    {
00074        double energy = fRandomGenerator->fire(fMinE, fMaxE) ;
00075        double eta    = fRandomGenerator->fire(fMinEta, fMaxEta) ;
00076        double phi    = fRandomGenerator->fire(fMinPhi, fMaxPhi) ;
00077        int PartID = fPartIDs[ip] ;
00078        const HepPDT::ParticleData* 
00079           PData = fPDGTable->particle(HepPDT::ParticleID(abs(PartID))) ;
00080        double mass   = PData->mass().value() ;
00081        double mom2   = energy*energy - mass*mass ;
00082        double mom    = 0. ;
00083        if (mom2 > 0.) 
00084        {
00085           mom = sqrt(mom2) ;
00086        }
00087        else
00088        {
00089           mom = 0. ;
00090        }
00091        double theta  = 2.*atan(exp(-eta)) ;
00092        double px     = mom*sin(theta)*cos(phi) ;
00093        double py     = mom*sin(theta)*sin(phi) ;
00094        double pz     = mom*cos(theta) ;
00095 
00096        HepMC::FourVector p(px,py,pz,energy) ;
00097        HepMC::GenParticle* Part = 
00098            new HepMC::GenParticle(p,PartID,1);
00099        Part->suggest_barcode( barcode ) ;
00100        barcode++ ;
00101        Vtx->add_particle_out(Part);
00102        
00103        if ( fAddAntiParticle )
00104        {
00105           HepMC::FourVector ap(-px,-py,-pz,energy) ;
00106           int APartID = -PartID ;
00107           if ( PartID == 22 || PartID == 23 )
00108           {
00109              APartID = PartID ;
00110           }
00111           HepMC::GenParticle* APart =
00112              new HepMC::GenParticle(ap,APartID,1);
00113           APart->suggest_barcode( barcode ) ;
00114           barcode++ ;
00115           Vtx->add_particle_out(APart) ;
00116        }
00117        
00118    }
00119    fEvt->add_vertex(Vtx) ;
00120    fEvt->set_event_number(event()) ;
00121    fEvt->set_signal_process_id(20) ;  
00122    
00123    
00124    if ( fVerbosity > 0 )
00125    {
00126       fEvt->print() ;  
00127    }  
00128 
00129    auto_ptr<HepMCProduct> BProduct(new HepMCProduct()) ;
00130    BProduct->addHepMCData( fEvt );
00131    e.put(BProduct);
00132     
00133    if ( fVerbosity > 0 )
00134    {
00135       // for testing purpose only
00136       //fEvt->print() ;  // for some strange reason, it prints NO event info
00137       // after it's been put into edm::Event...
00138       cout << " FlatRandomEGunSource : Event Generation Done " << endl;
00139    }
00140       
00141    return true;
00142 }