CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch1/src/IOMC/ParticleGuns/src/ExpoRandomPtGunSource.cc

Go to the documentation of this file.
00001 /*
00002  *  $Date: 2009/05/25 13:00:58 $
00003  *  $Revision: 1.4 $
00004  *  \author Jean-Roch Vlimant
00005  */
00006 
00007 #include <ostream>
00008 
00009 #include "IOMC/ParticleGuns/interface/ExpoRandomPtGunSource.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 // #include "FWCore/Utilities/interface/Exception.h"
00017 
00018 // #include "CLHEP/Random/RandExpo.h"
00019 
00020 using namespace edm;
00021 using namespace std;
00022 
00023 ExpoRandomPtGunSource::ExpoRandomPtGunSource(const ParameterSet& pset,
00024                                              const InputSourceDescription& desc) : 
00025    BaseFlatGunSource(pset, desc)
00026 {
00027 
00028 
00029    ParameterSet defpset ;
00030    ParameterSet pgun_params = 
00031       pset.getUntrackedParameter<ParameterSet>("PGunParameters",defpset) ;
00032   
00033    fMinPt = pgun_params.getUntrackedParameter<double>("MinPt",0.99);
00034    fMaxPt = pgun_params.getUntrackedParameter<double>("MaxPt",1.01);
00035 
00036    fMeanPt = pgun_params.getUntrackedParameter<double>("MeanPt",-1.);
00037 
00038    produces<HepMCProduct>();
00039 
00040    //the explonential generator
00041    fRandomExpoGenerator = new CLHEP::RandExponential(fRandomEngine,fMeanPt);
00042    
00043 }
00044 
00045 ExpoRandomPtGunSource::~ExpoRandomPtGunSource()
00046 {
00047    // no need to cleanup GenEvent memory - done in HepMCProduct
00048 }
00049 
00050 bool ExpoRandomPtGunSource::produce(Event &e) 
00051 {
00052 
00053    if ( fVerbosity > 0 )
00054    {
00055       cout << " ExpoRandomPtGunSource : Begin New Event Generation" << endl ; 
00056    }
00057    // event loop (well, another step in it...)
00058           
00059    // no need to clean up GenEvent memory - done in HepMCProduct
00060    // 
00061    
00062    // here re-create fEvt (memory)
00063    //
00064    fEvt = new HepMC::GenEvent() ;
00065    
00066    // now actualy, cook up the event from PDGTable and gun parameters
00067    //
00068    // 1st, primary vertex
00069    //
00070    //HepMC::GenVertex* Vtx = new HepMC::GenVertex(CLHEP::HepLorentzVector(0.,0.,0.));
00071    HepMC::GenVertex* Vtx = new HepMC::GenVertex(HepMC::FourVector(0.,0.,0.));
00072 
00073    // loop over particles
00074    //
00075    int barcode = 1 ;
00076    for (unsigned int ip=0; ip<fPartIDs.size(); ++ip)
00077    {
00078 
00079      //the max is to ensure you don't generate at 0
00080      //the 90% is to get rid of edge effect
00081      
00082      double pt     =  std::max(0.00001,0.90*fMinPt)+fRandomExpoGenerator->fire(fMeanPt);
00083      //shoot until in the designated range
00084      while (pt<fMinPt || pt>fMaxPt)
00085        {pt = std::max(0.00001,0.90*fMinPt) + fRandomExpoGenerator->fire(fMeanPt);}
00086      
00087      double eta    = fRandomGenerator->fire(fMinEta, fMaxEta) ;
00088        double phi    = fRandomGenerator->fire(fMinPhi, fMaxPhi) ;
00089        int PartID = fPartIDs[ip] ;
00090        const HepPDT::ParticleData* 
00091           PData = fPDGTable->particle(HepPDT::ParticleID(abs(PartID))) ;
00092        double mass   = PData->mass().value() ;
00093        double theta  = 2.*atan(exp(-eta)) ;
00094        double mom    = pt/sin(theta) ;
00095        double px     = pt*cos(phi) ;
00096        double py     = pt*sin(phi) ;
00097        double pz     = mom*cos(theta) ;
00098        double energy2= mom*mom + mass*mass ;
00099        double energy = sqrt(energy2) ; 
00100        //CLHEP::Hep3Vector p(px,py,pz) ;
00101        //HepMC::GenParticle* Part = 
00102        //    new HepMC::GenParticle(CLHEP::HepLorentzVector(p,energy),PartID,1);
00103        HepMC::FourVector p(px,py,pz,energy) ;
00104        HepMC::GenParticle* Part = 
00105            new HepMC::GenParticle(p,PartID,1);
00106        Part->suggest_barcode( barcode ) ;
00107        barcode++ ;
00108        Vtx->add_particle_out(Part);
00109 
00110        if ( fAddAntiParticle )
00111        {
00112           //CLHEP::Hep3Vector ap(-px,-py,-pz) ;
00113           HepMC::FourVector ap(-px,-py,-pz,energy) ;
00114           int APartID = -PartID ;
00115           if ( PartID == 22 || PartID == 23 )
00116           {
00117              APartID = PartID ;
00118           }       
00119           //HepMC::GenParticle* APart =
00120           //   new HepMC::GenParticle(CLHEP::HepLorentzVector(ap,energy),APartID,1);
00121           HepMC::GenParticle* APart =
00122              new HepMC::GenParticle(ap,APartID,1);
00123           APart->suggest_barcode( barcode ) ;
00124           barcode++ ;
00125           Vtx->add_particle_out(APart) ;
00126        }
00127 
00128    }
00129 
00130    fEvt->add_vertex(Vtx) ;
00131    fEvt->set_event_number(event()) ;
00132    fEvt->set_signal_process_id(20) ; 
00133         
00134    if ( fVerbosity > 0 )
00135    {
00136       fEvt->print() ;  
00137    }
00138 
00139    auto_ptr<HepMCProduct> BProduct(new HepMCProduct()) ;
00140    BProduct->addHepMCData( fEvt );
00141    e.put(BProduct);
00142     
00143    if ( fVerbosity > 0 )
00144    {
00145       // for testing purpose only
00146       // fEvt->print() ; // prints empty info after it's made into edm::Event
00147       cout << " FlatRandomPtGunSource : Event Generation Done " << endl;
00148    }
00149 
00150    return true;
00151 }
00152