CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
BeamMomentumGunProducer.cc
Go to the documentation of this file.
1 #include <ostream>
3 
6 
12 
13 #include "CLHEP/Random/RandFlat.h"
14 
15 #include "TFile.h"
16 
17 namespace CLHEP {
18  class HepRandomEngine;
19 }
20 
21 namespace edm {
22  BeamMomentumGunProducer::BeamMomentumGunProducer(const edm::ParameterSet& pset)
24  parPDGId_(nullptr),
25  parX_(nullptr),
26  parY_(nullptr),
27  parZ_(nullptr),
28  parPx_(nullptr),
29  parPy_(nullptr),
30  parPz_(nullptr),
31  b_npar_(nullptr),
32  b_eventId_(nullptr),
33  b_parPDGId_(nullptr),
34  b_parX_(nullptr),
35  b_parY_(nullptr),
36  b_parZ_(nullptr),
37  b_parPx_(nullptr),
38  b_parPy_(nullptr),
39  b_parPz_(nullptr) {
40  edm::ParameterSet pgun_params = pset.getParameter<edm::ParameterSet>("PGunParameters");
41 
42  // doesn't seem necessary to check if pset is empty
43  xoff_ = pgun_params.getParameter<double>("XOffset");
44  yoff_ = pgun_params.getParameter<double>("YOffset");
45  zpos_ = pgun_params.getParameter<double>("ZPosition");
46  if (fVerbosity > 0)
47  edm::LogVerbatim("BeamMomentumGun")
48  << "Beam vertex offset (cm) " << xoff_ << ":" << yoff_ << " and z position " << zpos_;
49 
50  edm::FileInPath fp = pgun_params.getParameter<edm::FileInPath>("FileName");
52 
53  fFile_ = new TFile(infileName.c_str());
54  fFile_->GetObject("EventTree", fTree_);
55  nentries_ = fTree_->GetEntriesFast();
56  if (fVerbosity > 0)
57  edm::LogVerbatim("BeamMomentumGun") << "Total Events: " << nentries_ << " in " << infileName;
58 
59  // Set branch addresses and branch pointers
60  int npart = fTree_->SetBranchAddress("npar", &npar_, &b_npar_);
61  int event = fTree_->SetBranchAddress("eventId", &eventId_, &b_eventId_);
62  int pdgid = fTree_->SetBranchAddress("parPDGId", &parPDGId_, &b_parPDGId_);
63  int parxx = fTree_->SetBranchAddress("parX", &parX_, &b_parX_);
64  int paryy = fTree_->SetBranchAddress("parY", &parY_, &b_parY_);
65  int parzz = fTree_->SetBranchAddress("parZ", &parZ_, &b_parZ_);
66  int parpx = fTree_->SetBranchAddress("parPx", &parPx_, &b_parPx_);
67  int parpy = fTree_->SetBranchAddress("parPy", &parPy_, &b_parPy_);
68  int parpz = fTree_->SetBranchAddress("parPz", &parPz_, &b_parPz_);
69  if ((npart != 0) || (event != 0) || (pdgid != 0) || (parxx != 0) || (paryy != 0) || (parzz != 0) || (parpx != 0) ||
70  (parpy != 0) || (parpz != 0))
71  throw cms::Exception("GenException") << "Branch address wrong in i/p file\n";
72 
73  produces<HepMCProduct>("unsmeared");
74  produces<GenEventInfoProduct>();
75 
76  if (fVerbosity > 0)
77  edm::LogVerbatim("BeamMomentumGun") << "BeamMonetumGun is initialzed";
78  }
79 
81  if (fVerbosity > 0)
82  edm::LogVerbatim("BeamMomentumGun") << "BeamMomentumGunProducer : Begin New Event Generation";
83 
85  CLHEP::HepRandomEngine* engine = &rng->getEngine(e.streamID());
86 
87  // event loop (well, another step in it...)
88  // no need to clean up GenEvent memory - done in HepMCProduct
89  // here re-create fEvt (memory)
90  //
91  fEvt = new HepMC::GenEvent();
92 
93  // random entry generation for peaking event randomly from tree
94  long int rjentry = static_cast<long int>(CLHEP::RandFlat::shoot(engine, 0, nentries_ - 1));
95  fTree_->GetEntry(rjentry);
96  if (fVerbosity > 0)
97  edm::LogVerbatim("BeamMomentumGun") << "Entry " << rjentry << " : " << eventId_ << " : " << npar_;
98 
99  // loop over particles
100  int barcode = 1;
101  for (unsigned int ip = 0; ip < parPDGId_->size(); ip++) {
102  int partID = parPDGId_->at(ip);
103  const HepPDT::ParticleData* pData = fPDGTable->particle(HepPDT::ParticleID(std::abs(partID)));
104  double mass = pData->mass().value();
105  if (fVerbosity > 0)
106  edm::LogVerbatim("BeamMomentumGun") << "PDGId: " << partID << " mass: " << mass;
107  double xp = (xoff_ + mm2cm_ * parX_->at(ip));
108  double yp = (yoff_ + mm2cm_ * parY_->at(ip));
109  HepMC::GenVertex* Vtx = new HepMC::GenVertex(HepMC::FourVector(xp, yp, zpos_));
110  double pxGeV = MeV2GeV_ * parPx_->at(ip);
111  double pyGeV = MeV2GeV_ * parPy_->at(ip);
112  double pzGeV = MeV2GeV_ * parPz_->at(ip);
113  double momRand2 = pxGeV * pxGeV + pyGeV * pyGeV + pzGeV * pzGeV;
114  double energy = std::sqrt(momRand2 + mass * mass);
115  double mom = std::sqrt(momRand2);
116  double theta = CLHEP::RandFlat::shoot(engine, fMinTheta, fMaxTheta);
117  double phi = CLHEP::RandFlat::shoot(engine, fMinPhi, fMaxPhi);
118  double px = mom * sin(theta) * cos(phi);
119  double py = mom * sin(theta) * sin(phi);
120  double pz = mom * cos(theta);
121 
122  if (fVerbosity > 0)
123  edm::LogVerbatim("BeamMomentumGun") << "px:py:pz " << px << ":" << py << ":" << pz;
124 
125  HepMC::FourVector p(px, py, pz, energy);
126  HepMC::GenParticle* part = new HepMC::GenParticle(p, partID, 1);
127  part->suggest_barcode(barcode);
128  barcode++;
129  Vtx->add_particle_out(part);
130 
131  if (fAddAntiParticle) {
132  HepMC::FourVector ap(-px, -py, -pz, energy);
133  int apartID = (partID == 22 || partID == 23) ? partID : -partID;
134  HepMC::GenParticle* apart = new HepMC::GenParticle(ap, apartID, 1);
135  apart->suggest_barcode(barcode);
136  if (fVerbosity > 0)
137  edm::LogVerbatim("BeamMomentumGun")
138  << "Add anti-particle " << apartID << ":" << -px << ":" << -py << ":" << -pz;
139  barcode++;
140  Vtx->add_particle_out(apart);
141  }
142 
143  fEvt->add_vertex(Vtx);
144  }
145 
146  fEvt->set_event_number(e.id().event());
147  fEvt->set_signal_process_id(20);
148 
149  if (fVerbosity > 0)
150  fEvt->print();
151 
152  std::unique_ptr<HepMCProduct> BProduct(new HepMCProduct());
153  BProduct->addHepMCData(fEvt);
154  e.put(std::move(BProduct), "unsmeared");
155 
156  std::unique_ptr<GenEventInfoProduct> genEventInfo(new GenEventInfoProduct(fEvt));
157  e.put(std::move(genEventInfo));
158 
159  if (fVerbosity > 0)
160  edm::LogVerbatim("BeamMomentumGun") << "BeamMomentumGunProducer : Event Generation Done";
161  }
162 } // namespace edm
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:40
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
#define nullptr
ESHandle< HepPDT::ParticleDataTable > fPDGTable
Geom::Theta< T > theta() const
double npart
Definition: HydjetWrapper.h:46
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
T sqrt(T t)
Definition: SSEVec.h:19
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void produce(Event &e, const EventSetup &es) override
HepPDT::ParticleData ParticleData
part
Definition: HCALResponse.h:20
edm::EventID id() const
Definition: EventBase.h:59
HLT enums.
StreamID streamID() const
Definition: Event.h:96
std::string fullPath() const
Definition: FileInPath.cc:163
def move(src, dest)
Definition: eostools.py:511
Definition: event.py:1