CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_7/src/IOMC/ParticleGuns/src/FlatRandomEGunProducer.cc

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