Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007 #include <ostream>
00008
00009 #include "IOMC/ParticleGuns/interface/FlatRandomPtGunSource.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
00017 using namespace edm;
00018 using namespace std;
00019
00020 FlatRandomPtGunSource::FlatRandomPtGunSource(const ParameterSet& pset,
00021 const InputSourceDescription& desc) :
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 }
00036
00037 FlatRandomPtGunSource::~FlatRandomPtGunSource()
00038 {
00039
00040 }
00041
00042 bool FlatRandomPtGunSource::produce(Event &e)
00043 {
00044
00045 if ( fVerbosity > 0 )
00046 {
00047 cout << " FlatRandomPtGunSource : Begin New Event Generation" << endl ;
00048 }
00049
00050
00051
00052
00053
00054
00055
00056 fEvt = new HepMC::GenEvent() ;
00057
00058
00059
00060
00061
00062 HepMC::GenVertex* Vtx = new HepMC::GenVertex(HepMC::FourVector(0.,0.,0.));
00063
00064
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
00124
00125 cout << " FlatRandomPtGunSource : Event Generation Done " << endl;
00126 }
00127
00128 return true;
00129 }
00130