CMS 3D CMS Logo

FlatEGunASCIIWriter.cc
Go to the documentation of this file.
1 // ----------------------------------------------------------------------
2 // FlatEGunASCIIWriter.cc
3 // Author: Julia Yarba
4 //
5 // This code has been molded after examples in HepMC and HepPDT, and
6 // after single particle gun example (private contacts with Lynn Garren)
7 //
8 // Plus, it uses the ParameterSet funtionalities for "user interface"
9 //
10 // It creates a ParticleDataTable from PDG-2004 data file.
11 // It then creates one or several HepMC particle(s) and adds it(them)
12 // to the HepMC::GenEvent.
13 // For particle definition it uses flat random gun in a given eta, phi
14 // and energy ranges to calculate particle kinematics.
15 // Vertex smearing is currently off, can be easily added later.
16 // After all, it writes out the event into ASCII output file.
17 //
18 // ----------------------------------------------------------------------
19 
20 #include <iostream>
21 #include <fstream>
22 // header
24 
25 // essentials !!!
29 
31 
32 #include "CLHEP/Random/RandFlat.h"
33 
34 using namespace edm;
35 using namespace std;
36 
38  : fEvt(nullptr),
39  fOutFileName(pset.getUntrackedParameter<string>("OutFileName", "FlatEGunHepMC.dat")),
40  fCurrentEvent(0) {
41  ParameterSet pgun_params = pset.getParameter<ParameterSet>("PGunParameters");
42  fPartIDs = pgun_params.getParameter<vector<int> >("PartID");
43  fMinEta = pgun_params.getParameter<double>("MinEta");
44  fMaxEta = pgun_params.getParameter<double>("MaxEta");
45  fMinPhi = pgun_params.getParameter<double>("MinPhi");
46  fMaxPhi = pgun_params.getParameter<double>("MaxPhi");
47  fMinE = pgun_params.getParameter<double>("MinE");
48  fMaxE = pgun_params.getParameter<double>("MaxE");
49 }
50 
52  if (fEvt != nullptr)
53  delete fEvt;
54  if (fOutStream)
55  delete fOutStream;
56 }
57 
59 
61  fOutStream = new HepMC::IO_GenEvent(fOutFileName);
62  if (fOutStream->rdstate() == std::ios::failbit) {
63  throw cms::Exception("FileNotOpen", "FlatEGunASCIIWriter::beginJob()")
64  << "File " << fOutFileName << " was not open.\n";
65  }
66 
67  return;
68 }
69 
70 void FlatEGunASCIIWriter::analyze(const Event&, const EventSetup& /* not used so far, thus is "empty" */) {
71  // clean up GenEvent memory : also deletes all vtx/part in it
72  //
73  if (fEvt != nullptr)
74  delete fEvt;
75 
76  // here re-create fEvt (memory)
77  fEvt = new HepMC::GenEvent();
78 
79  // now actualy, cook up the event from PDGTable and gun parameters
80  //
81  //HepMC::GenVertex* Vtx = new HepMC::GenVertex( CLHEP::HepLorentzVector(0.,0.,0.) );
82  HepMC::GenVertex* Vtx = new HepMC::GenVertex(HepMC::FourVector(0., 0., 0.));
83 
84  // loop over particles
85  //
86  for (unsigned int ip = 0; ip < fPartIDs.size(); ip++) {
87  double energy = CLHEP::RandFlat::shoot(fMinE, fMaxE);
88  double eta = CLHEP::RandFlat::shoot(fMinEta, fMaxEta);
89  double phi = CLHEP::RandFlat::shoot(fMinPhi, fMaxPhi);
90  const HepPDT::ParticleData* PData = fPDGTable->particle(HepPDT::ParticleID(abs(fPartIDs[ip])));
91  double mass = PData->mass().value();
92  double mom2 = energy * energy - mass * mass;
93  double mom = 0.;
94  if (mom2 > 0.)
95  mom = sqrt(mom2);
96  double theta = 2. * atan(exp(-eta));
97  double px = mom * sin(theta) * cos(phi);
98  double py = mom * sin(theta) * sin(phi);
99  double pz = mom * cos(theta);
100  //CLHEP::Hep3Vector p(px,py,pz) ;
101  //HepMC::GenParticle* Part =
102  // new HepMC::GenParticle(CLHEP::HepLorentzVector(p,energy),fPartIDs[ip],1);
103  HepMC::FourVector p(px, py, pz, energy);
104  HepMC::GenParticle* Part = new HepMC::GenParticle(p, fPartIDs[ip], 1);
105  Vtx->add_particle_out(Part);
106  }
107  fEvt->add_vertex(Vtx);
108  fEvt->set_event_number(fCurrentEvent + 1);
109  fEvt->set_signal_process_id(20);
110 
111  // for testing purpose only
112  // fEvt->print() ;
113 
114  fOutStream->write_event(fEvt);
115 
116  fCurrentEvent++;
117 
118  return;
119 }
edm::FlatEGunASCIIWriter::~FlatEGunASCIIWriter
~FlatEGunASCIIWriter() override
Definition: FlatEGunASCIIWriter.cc:51
edm::FlatEGunASCIIWriter::fMinE
double fMinE
Definition: FlatEGunASCIIWriter.h:53
edm::Run
Definition: Run.h:45
multPhiCorr_741_25nsDY_cfi.py
py
Definition: multPhiCorr_741_25nsDY_cfi.py:12
edm
HLT enums.
Definition: AlignableModifier.h:19
edm::FlatEGunASCIIWriter::fOutStream
HepMC::IO_GenEvent * fOutStream
Definition: FlatEGunASCIIWriter.h:63
edm::FlatEGunASCIIWriter::fOutFileName
std::string fOutFileName
Definition: FlatEGunASCIIWriter.h:62
edm::FlatEGunASCIIWriter::analyze
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: FlatEGunASCIIWriter.cc:70
HepMC::GenEvent
Definition: hepmc_rootio.cc:9
ParticleData
HepPDT::ParticleData ParticleData
Definition: ParticleDataTable.h:9
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
MakerMacros.h
FlatEGunASCIIWriter.h
edm::FlatEGunASCIIWriter::FlatEGunASCIIWriter
FlatEGunASCIIWriter(const edm::ParameterSet &)
Definition: FlatEGunASCIIWriter.cc:37
funct::cos
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
PVValHelper::eta
Definition: PVValidationHelpers.h:70
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
theta
Geom::Theta< T > theta() const
Definition: Basic3DVectorLD.h:150
HCALHighEnergyHPDFilter_cfi.energy
energy
Definition: HCALHighEnergyHPDFilter_cfi.py:5
edm::FlatEGunASCIIWriter::fMinPhi
double fMinPhi
Definition: FlatEGunASCIIWriter.h:51
edm::FlatEGunASCIIWriter::fMaxPhi
double fMaxPhi
Definition: FlatEGunASCIIWriter.h:52
edm::ParameterSet
Definition: ParameterSet.h:47
AlCaHLTBitMon_ParallelJobs.p
def p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
Event.h
ParticleDataTable.h
edm::FlatEGunASCIIWriter::beginRun
void beginRun(const edm::Run &, const EventSetup &) override
Definition: FlatEGunASCIIWriter.cc:58
edm::FlatEGunASCIIWriter::beginJob
void beginJob() override
Definition: FlatEGunASCIIWriter.cc:60
edm::EventSetup
Definition: EventSetup.h:58
edm::FlatEGunASCIIWriter::fPartIDs
std::vector< int > fPartIDs
Definition: FlatEGunASCIIWriter.h:48
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
edm::FlatEGunASCIIWriter::fMinEta
double fMinEta
Definition: FlatEGunASCIIWriter.h:49
edm::EventSetup::getData
bool getData(T &iHolder) const
Definition: EventSetup.h:127
alignCSCRings.r
r
Definition: alignCSCRings.py:93
multPhiCorr_741_25nsDY_cfi.px
px
Definition: multPhiCorr_741_25nsDY_cfi.py:10
GenParticle.GenParticle
GenParticle
Definition: GenParticle.py:18
edm::FlatEGunASCIIWriter::fEvt
HepMC::GenEvent * fEvt
Definition: FlatEGunASCIIWriter.h:57
std
Definition: JetResolutionObject.h:76
edm::FlatEGunASCIIWriter::fPDGTable
ESHandle< HepPDT::ParticleDataTable > fPDGTable
Definition: FlatEGunASCIIWriter.h:60
Exception
Definition: hltDiff.cc:245
EgHLTOffHistBins_cfi.mass
mass
Definition: EgHLTOffHistBins_cfi.py:34
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
Exception.h
edm::FlatEGunASCIIWriter::fMaxE
double fMaxE
Definition: FlatEGunASCIIWriter.h:54
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
JetChargeProducer_cfi.exp
exp
Definition: JetChargeProducer_cfi.py:6
edm::Event
Definition: Event.h:73
LHEGenericFilter_cfi.ParticleID
ParticleID
Definition: LHEGenericFilter_cfi.py:6
edm::FlatEGunASCIIWriter::fMaxEta
double fMaxEta
Definition: FlatEGunASCIIWriter.h:50
muonDTDigis_cfi.pset
pset
Definition: muonDTDigis_cfi.py:27
edm::FlatEGunASCIIWriter::fCurrentEvent
int fCurrentEvent
Definition: FlatEGunASCIIWriter.h:65