#include <FlatEGunASCIIWriter.h>
Public Member Functions | |
virtual void | analyze (const edm::Event &, const edm::EventSetup &) override |
virtual void | beginJob () |
virtual void | beginRun (const edm::Run &, const EventSetup &) override |
FlatEGunASCIIWriter (const edm::ParameterSet &) | |
virtual | ~FlatEGunASCIIWriter () |
Private Attributes | |
int | fCurrentEvent |
HepMC::GenEvent * | fEvt |
double | fMaxE |
double | fMaxEta |
double | fMaxPhi |
double | fMinE |
double | fMinEta |
double | fMinPhi |
std::string | fOutFileName |
HepMC::IO_GenEvent * | fOutStream |
std::vector< int > | fPartIDs |
ESHandle < HepPDT::ParticleDataTable > | fPDGTable |
Definition at line 34 of file FlatEGunASCIIWriter.h.
FlatEGunASCIIWriter::FlatEGunASCIIWriter | ( | const edm::ParameterSet & | pset | ) | [explicit] |
Definition at line 37 of file FlatEGunASCIIWriter.cc.
References fMaxE, fMaxEta, fMaxPhi, fMinE, fMinEta, fMinPhi, fPartIDs, and edm::ParameterSet::getParameter().
: fEvt(0), fOutFileName( pset.getUntrackedParameter<string>("OutFileName","FlatEGunHepMC.dat") ), fCurrentEvent(0) { ParameterSet pgun_params = pset.getParameter<ParameterSet>("PGunParameters") ; fPartIDs = pgun_params.getParameter< vector<int> >("PartID"); fMinEta = pgun_params.getParameter<double>("MinEta"); fMaxEta = pgun_params.getParameter<double>("MaxEta"); fMinPhi = pgun_params.getParameter<double>("MinPhi"); fMaxPhi = pgun_params.getParameter<double>("MaxPhi"); fMinE = pgun_params.getParameter<double>("MinE"); fMaxE = pgun_params.getParameter<double>("MaxE"); }
FlatEGunASCIIWriter::~FlatEGunASCIIWriter | ( | ) | [virtual] |
Definition at line 54 of file FlatEGunASCIIWriter.cc.
References fEvt, fOutStream, and NULL.
{ if ( fEvt != NULL ) delete fEvt ; if( fOutStream) delete fOutStream; }
void FlatEGunASCIIWriter::analyze | ( | const edm::Event & | , |
const edm::EventSetup & | |||
) | [override, virtual] |
Implements edm::EDAnalyzer.
Definition at line 81 of file FlatEGunASCIIWriter.cc.
References abs, funct::cos(), relval_parameters_module::energy, eta(), create_public_lumi_plots::exp, fCurrentEvent, fEvt, fMaxE, fMaxEta, fMaxPhi, fMinE, fMinEta, fMinPhi, fOutStream, fPartIDs, fPDGTable, configurableAnalysis::GenParticle, NULL, AlCaHLTBitMon_ParallelJobs::p, RecoTau_DiTaus_pt_20-420_cfg::ParticleID, phi, funct::sin(), mathSSE::sqrt(), and theta().
{ // clean up GenEvent memory : also deletes all vtx/part in it // if ( fEvt != NULL ) delete fEvt ; // here re-create fEvt (memory) fEvt = new HepMC::GenEvent() ; // now actualy, cook up the event from PDGTable and gun parameters // //HepMC::GenVertex* Vtx = new HepMC::GenVertex( CLHEP::HepLorentzVector(0.,0.,0.) ); HepMC::GenVertex* Vtx = new HepMC::GenVertex( HepMC::FourVector(0.,0.,0.) ); // loop over particles // for ( unsigned int ip=0; ip<fPartIDs.size(); ip++ ) { double energy = CLHEP::RandFlat::shoot( fMinE, fMaxE ) ; double eta = CLHEP::RandFlat::shoot( fMinEta, fMaxEta ) ; double phi = CLHEP::RandFlat::shoot( fMinPhi, fMaxPhi ) ; const HepPDT::ParticleData* PData = fPDGTable->particle(HepPDT::ParticleID(abs(fPartIDs[ip]))) ; double mass = PData->mass().value() ; double mom2 = energy*energy - mass*mass ; double mom = 0. ; if ( mom2 > 0. ) mom = sqrt(mom2) ; double theta = 2.*atan(exp(-eta)) ; double px = mom*sin(theta)*cos(phi) ; double py = mom*sin(theta)*sin(phi) ; double pz = mom*cos(theta) ; //CLHEP::Hep3Vector p(px,py,pz) ; //HepMC::GenParticle* Part = // new HepMC::GenParticle(CLHEP::HepLorentzVector(p,energy),fPartIDs[ip],1); HepMC::FourVector p(px,py,pz,energy); HepMC::GenParticle* Part = new HepMC::GenParticle(p,fPartIDs[ip],1); Vtx->add_particle_out(Part); } fEvt->add_vertex( Vtx ) ; fEvt->set_event_number( fCurrentEvent+1 ) ; fEvt->set_signal_process_id(20) ; // for testing purpose only // fEvt->print() ; fOutStream->write_event( fEvt ); fCurrentEvent++ ; return ; }
void FlatEGunASCIIWriter::beginJob | ( | void | ) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 69 of file FlatEGunASCIIWriter.cc.
References Exception, fOutFileName, and fOutStream.
{ fOutStream = new HepMC::IO_GenEvent( fOutFileName.c_str() ); if ( fOutStream->rdstate() == std::ios::failbit ) { throw cms::Exception("FileNotOpen", "FlatEGunASCIIWriter::beginJob()") << "File " << fOutFileName << " was not open.\n"; } return ; }
void FlatEGunASCIIWriter::beginRun | ( | const edm::Run & | r, |
const EventSetup & | es | ||
) | [override, virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 62 of file FlatEGunASCIIWriter.cc.
References fPDGTable, and edm::EventSetup::getData().
int edm::FlatEGunASCIIWriter::fCurrentEvent [private] |
Definition at line 72 of file FlatEGunASCIIWriter.h.
Referenced by analyze().
HepMC::GenEvent* edm::FlatEGunASCIIWriter::fEvt [private] |
Definition at line 63 of file FlatEGunASCIIWriter.h.
Referenced by analyze(), and ~FlatEGunASCIIWriter().
double edm::FlatEGunASCIIWriter::fMaxE [private] |
Definition at line 60 of file FlatEGunASCIIWriter.h.
Referenced by analyze(), and FlatEGunASCIIWriter().
double edm::FlatEGunASCIIWriter::fMaxEta [private] |
Definition at line 56 of file FlatEGunASCIIWriter.h.
Referenced by analyze(), and FlatEGunASCIIWriter().
double edm::FlatEGunASCIIWriter::fMaxPhi [private] |
Definition at line 58 of file FlatEGunASCIIWriter.h.
Referenced by analyze(), and FlatEGunASCIIWriter().
double edm::FlatEGunASCIIWriter::fMinE [private] |
Definition at line 59 of file FlatEGunASCIIWriter.h.
Referenced by analyze(), and FlatEGunASCIIWriter().
double edm::FlatEGunASCIIWriter::fMinEta [private] |
Definition at line 55 of file FlatEGunASCIIWriter.h.
Referenced by analyze(), and FlatEGunASCIIWriter().
double edm::FlatEGunASCIIWriter::fMinPhi [private] |
Definition at line 57 of file FlatEGunASCIIWriter.h.
Referenced by analyze(), and FlatEGunASCIIWriter().
std::string edm::FlatEGunASCIIWriter::fOutFileName [private] |
Definition at line 69 of file FlatEGunASCIIWriter.h.
Referenced by beginJob().
HepMC::IO_GenEvent* edm::FlatEGunASCIIWriter::fOutStream [private] |
Definition at line 70 of file FlatEGunASCIIWriter.h.
Referenced by analyze(), beginJob(), and ~FlatEGunASCIIWriter().
std::vector<int> edm::FlatEGunASCIIWriter::fPartIDs [private] |
Definition at line 54 of file FlatEGunASCIIWriter.h.
Referenced by analyze(), and FlatEGunASCIIWriter().
Definition at line 67 of file FlatEGunASCIIWriter.h.
Referenced by analyze(), and beginRun().