CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_1/src/IOMC/ParticleGuns/src/ExpoRandomPtGunProducer.cc

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