Go to the documentation of this file.00001 #include <ostream>
00002
00003 #include "IOMC/ParticleGuns/interface/FlatRandomEThetaGunProducer.h"
00004
00005 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00006 #include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
00007
00008 #include "FWCore/Framework/interface/Event.h"
00009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00011
00012 using namespace edm;
00013
00014 FlatRandomEThetaGunProducer::FlatRandomEThetaGunProducer(const edm::ParameterSet& pset) :
00015 FlatBaseThetaGunProducer(pset) {
00016
00017 edm::ParameterSet defpset ;
00018 edm::ParameterSet pgun_params =
00019 pset.getParameter<edm::ParameterSet>("PGunParameters") ;
00020
00021
00022
00023 fMinE = pgun_params.getParameter<double>("MinE");
00024 fMaxE = pgun_params.getParameter<double>("MaxE");
00025
00026 produces<HepMCProduct>();
00027 produces<GenEventInfoProduct>();
00028
00029
00030
00031
00032 }
00033
00034 FlatRandomEThetaGunProducer::~FlatRandomEThetaGunProducer() {}
00035
00036 void FlatRandomEThetaGunProducer::produce(edm::Event & e, const edm::EventSetup& es) {
00037
00038 if ( fVerbosity > 0 ) {
00039 LogDebug("FlatThetaGun") << "FlatRandomEThetaGunProducer : Begin New Event Generation";
00040 }
00041
00042
00043
00044
00045
00046
00047
00048 fEvt = new HepMC::GenEvent() ;
00049
00050
00051
00052
00053
00054
00055 HepMC::GenVertex* Vtx = new HepMC::GenVertex( HepMC::FourVector(0.,0.,0.));
00056
00057
00058
00059 int barcode = 1;
00060 for (unsigned int ip=0; ip<fPartIDs.size(); ip++) {
00061 double energy = fRandomGenerator->fire(fMinE, fMaxE) ;
00062 double theta = fRandomGenerator->fire(fMinTheta, fMaxTheta) ;
00063 double phi = fRandomGenerator->fire(fMinPhi, fMaxPhi) ;
00064 int PartID = fPartIDs[ip] ;
00065 const HepPDT::ParticleData*
00066 PData = fPDGTable->particle(HepPDT::ParticleID(abs(PartID))) ;
00067 double mass = PData->mass().value() ;
00068 double mom2 = energy*energy - mass*mass ;
00069 double mom = (mom2 > 0. ? std::sqrt(mom2) : 0.);
00070 double px = mom*sin(theta)*cos(phi) ;
00071 double py = mom*sin(theta)*sin(phi) ;
00072 double pz = mom*cos(theta) ;
00073
00074 HepMC::FourVector p(px,py,pz,energy) ;
00075 HepMC::GenParticle* Part = new HepMC::GenParticle(p,PartID,1);
00076 Part->suggest_barcode( barcode ) ;
00077 barcode++ ;
00078 Vtx->add_particle_out(Part);
00079
00080 if ( fAddAntiParticle ) {
00081 HepMC::FourVector ap(-px,-py,-pz,energy) ;
00082 int APartID = -PartID ;
00083 if ( PartID == 22 || PartID == 23 ) {
00084 APartID = PartID ;
00085 }
00086 HepMC::GenParticle* APart = new HepMC::GenParticle(ap,APartID,1);
00087 APart->suggest_barcode( barcode ) ;
00088 barcode++ ;
00089 Vtx->add_particle_out(APart) ;
00090 }
00091
00092 }
00093 fEvt->add_vertex(Vtx) ;
00094 fEvt->set_event_number(e.id().event()) ;
00095 fEvt->set_signal_process_id(20) ;
00096
00097
00098 if ( fVerbosity > 0 ) {
00099 fEvt->print() ;
00100 }
00101
00102 std::auto_ptr<HepMCProduct> BProduct(new HepMCProduct()) ;
00103 BProduct->addHepMCData( fEvt );
00104 e.put(BProduct);
00105
00106 std::auto_ptr<GenEventInfoProduct> genEventInfo(new GenEventInfoProduct(fEvt));
00107 e.put(genEventInfo);
00108
00109 if ( fVerbosity > 0 ) {
00110 LogDebug("FlatThetaGun") << "FlatRandomEThetaGunProducer : Event Generation Done";
00111 }
00112 }