CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ThePEGInterface.cc
Go to the documentation of this file.
1 
7 #include <iostream>
8 #include <sstream>
9 #include <cstdlib>
10 #include <memory>
11 #include <cmath>
12 #include <stdlib.h>
13 
14 #include <boost/shared_ptr.hpp>
15 #include <boost/filesystem.hpp>
16 
17 #include <HepMC/GenEvent.h>
18 #include <HepMC/PdfInfo.h>
19 #include <HepMC/IO_GenEvent.h>
20 
21 #include <ThePEG/Utilities/DynamicLoader.h>
22 #include <ThePEG/Repository/Repository.h>
23 #include <ThePEG/Handlers/EventHandler.h>
24 #include <ThePEG/Handlers/XComb.h>
25 #include <ThePEG/EventRecord/Event.h>
26 #include <ThePEG/EventRecord/Particle.h>
27 #include <ThePEG/EventRecord/Collision.h>
28 #include <ThePEG/EventRecord/TmpTransform.h>
29 #include <ThePEG/Config/ThePEG.h>
30 #include <ThePEG/PDF/PartonExtractor.h>
31 #include <ThePEG/PDF/PDFBase.h>
32 #include <ThePEG/Utilities/UtilityBase.h>
33 
35 
37 
42 
43 using namespace std;
44 using namespace gen;
45 
47  randomEngineGlueProxy_(ThePEG::RandomEngineGlue::Proxy::create()),
48  dataLocation_(ParameterCollector::resolve(pset.getParameter<string>("dataLocation"))),
49  generator_(pset.getParameter<string>("generatorModule")),
50  run_(pset.getParameter<string>("run")),
51  dumpConfig_(pset.getUntrackedParameter<string>("dumpConfig", "")),
52  skipEvents_(pset.getUntrackedParameter<unsigned int>("skipEvents", 0))
53 {
54  // Write events in hepmc ascii format for debugging purposes
55  string dumpEvents = pset.getUntrackedParameter<string>("dumpEvents", "");
56  if (!dumpEvents.empty()) {
57  iobc_.reset(new HepMC::IO_GenEvent(dumpEvents.c_str(), ios::out));
58  edm::LogInfo("ThePEGSource") << "Event logging switched on (=> " << dumpEvents << ")";
59  }
60  // Clear dumpConfig target
61  if (!dumpConfig_.empty())
62  ofstream cfgDump(dumpConfig_.c_str(), ios_base::trunc);
63 }
64 
66 {
67  if (eg_)
68  eg_->finalize();
69  edm::LogInfo("ThePEGInterface") << "Event generator finalized";
70 }
71 
72 string ThePEGInterface::dataFile(const string &fileName) const
73 {
74  if (fileName.empty() || fileName[0] == '/')
75  return fileName;
76  else
77  return dataLocation_ + "/" + fileName;
78 }
79 
81  const string &paramName) const
82 {
83  const edm::Entry &entry = pset.retrieve(paramName);
84  if (entry.typeCode() == 'F')
85  return entry.getFileInPath().fullPath();
86  else
87  return dataFile(entry.getString());
88 }
89 
91 {
92  /* Initialize the repository from
93  * 1. the repository file (default values)
94  * 2. the ThePEG config files
95  * 3. the CMSSW config blocks
96  */
97  stringstream logstream;
98 
99  // Read the repository of serialized default values
100  string repository = dataFile(pset, "repository");
101  if (!repository.empty()) {
102  edm::LogInfo("ThePEGInterface") << "Loading repository (" << repository << ")";
103  ThePEG::Repository::load(repository);
104  }
105 
106  if (!getenv("ThePEG_INSTALL_PATH")) {
107  vector<string> libdirlist = ThePEG::DynamicLoader::allPaths();
108  for(vector<string>::const_iterator libdir = libdirlist.begin();
109  libdir < libdirlist.end(); ++libdir) {
110  if (libdir->empty() || (*libdir)[0] != '/')
111  continue;
112  if (boost::filesystem::exists(*libdir +
113  "/../../share/ThePEG/PDFsets.index")) {
114  setenv("ThePEG_INSTALL_PATH",
115  libdir->c_str(), 0);
116  break;
117  }
118  }
119  }
120 
121  // Read ThePEG config files to read
122  vector<string> configFiles = pset.getParameter<vector<string> >("configFiles");
123 
124  // Loop over the config files
125  for(vector<string>::const_iterator iter = configFiles.begin();
126  iter != configFiles.end(); ++iter) {
127  edm::LogInfo("ThePEGInterface") << "Reading config file (" << *iter << ")";
128  ThePEG::Repository::read(dataFile(*iter), logstream);
129  edm::LogInfo("ThePEGSource") << logstream.str();
130  }
131 
132  // Read CMSSW config file parameter sets starting from "parameterSets"
133 
134  ofstream cfgDump;
135  ParameterCollector collector(pset);
137  if (!dumpConfig_.empty()) {
138  cfgDump.open(dumpConfig_.c_str(), ios_base::app);
139  iter = collector.begin(cfgDump);
140  } else
141  iter = collector.begin();
142 
143  for(; iter != collector.end(); ++iter) {
144  string out = ThePEG::Repository::exec(*iter, logstream);
145  if (!out.empty()) {
146  edm::LogInfo("ThePEGInterface") << *iter << " => " << out;
147  cerr << "Error in ThePEG configuration!\n"
148  "\tLine: " << *iter << "\n" << out << endl;
149  }
150  }
151 
152  if (!dumpConfig_.empty()) {
153  cfgDump << "saverun " << run_ << " " << generator_ << endl;
154  cfgDump.close();
155  }
156 
157  // write the ProxyID for the RandomEngineGlue to fill its pointer in
158  ostringstream ss;
159  ss << randomEngineGlueProxy_->getID();
160  ThePEG::Repository::exec("set " + generator_ +
161  ":RandomNumberGenerator:ProxyID " + ss.str(),
162  logstream);
163 
164  // Print the directories where ThePEG looks for libs
165  vector<string> libdirlist = ThePEG::DynamicLoader::allPaths();
166  for(vector<string>::const_iterator libdir = libdirlist.begin();
167  libdir < libdirlist.end(); ++libdir)
168  edm::LogInfo("ThePEGInterface") << "DynamicLoader path = " << *libdir << endl;
169 
170  // Output status information about the repository
171  ThePEG::Repository::stats(logstream);
172  edm::LogInfo("ThePEGInterface") << logstream.str();
173 }
174 
176 {
177  // Get generator from the repository and initialize it
178  ThePEG::BaseRepository::CheckObjectDirectory(generator_);
179  ThePEG::EGPtr tmp = ThePEG::BaseRepository::GetObject<ThePEG::EGPtr>(generator_);
180  if (tmp) {
181  eg_ = ThePEG::Repository::makeRun(tmp, run_);
182  eg_->initialize();
183  edm::LogInfo("ThePEGInterface") << "EventGenerator initialized";
184  } else
185  throw cms::Exception("ThePEGInterface")
186  << "EventGenerator could not be initialized!" << endl;
187 
188  // Skip events
189  for (unsigned int i = 0; i < skipEvents_; i++) {
191  eg_->shoot();
192  edm::LogInfo("ThePEGInterface") << "Event discarded";
193  }
194 }
195 
197 {
198  ThePEG::RandomEngineGlue *rnd = randomEngineGlueProxy_->getInstance();
199 
200  if (!rnd)
201  edm::LogWarning("ProxyMissing")
202  << "ThePEG not initialised with RandomEngineGlue.";
203  else
204  rnd->flush();
205 }
206 
207 auto_ptr<HepMC::GenEvent> ThePEGInterface::convert(
208  const ThePEG::EventPtr &event)
209 {
210  return std::auto_ptr<HepMC::GenEvent>(
212 }
213 
214 void ThePEGInterface::clearAuxiliary(HepMC::GenEvent *hepmc,
215  HepMC::PdfInfo *pdf)
216 {
217  if (hepmc) {
218  hepmc->set_event_scale(-1.0);
219  hepmc->set_alphaQCD(-1.0);
220  hepmc->set_alphaQED(-1.0);
221  }
222  if (pdf) {
223  pdf->set_id1(-100);
224  pdf->set_id2(-100);
225  pdf->set_x1(-1.0);
226  pdf->set_x2(-1.0);
227  pdf->set_scalePDF(-1.0);
228  pdf->set_pdf1(-1.0);
229  pdf->set_pdf2(-1.0);
230  }
231 }
232 
233 void ThePEGInterface::fillAuxiliary(HepMC::GenEvent *hepmc,
234  HepMC::PdfInfo *pdf,
235  const ThePEG::EventPtr &event)
236 {
237  if (!event->primaryCollision())
238  return;
239 
240  ThePEG::tcEHPtr eh = ThePEG::dynamic_ptr_cast<ThePEG::tcEHPtr>(
241  event->primaryCollision()->handler());
242  double scale = eh->lastScale();
243 
244  if (hepmc) {
245  if (hepmc->event_scale() < 0 && scale > 0)
246  hepmc->set_event_scale(std::sqrt(scale) / ThePEG::GeV);
247 
248  if (hepmc->alphaQCD() < 0)
249  hepmc->set_alphaQCD(eh->lastAlphaS());
250  if (hepmc->alphaQED() < 0)
251  hepmc->set_alphaQED(eh->lastAlphaEM());
252  }
253 
254  if (pdf) {
255  const ThePEG::PPair &beams = eh->lastParticles();
256  const ThePEG::PPair &partons = eh->lastPartons();
257  ThePEG::tcPDFPtr pdf1 = eh->lastExtractor()->getPDF(
258  beams.first->dataPtr());
259  ThePEG::tcPDFPtr pdf2 = eh->lastExtractor()->getPDF(
260  beams.second->dataPtr());
261  double x1, x2;
262 
263  if (pdf->id1() == -100) {
264  int id = partons.first->id();
265  pdf->set_id1(id == 21 ? 0 : id);
266  }
267  if (pdf->id2() == -100) {
268  int id = partons.second->id();
269  pdf->set_id2(id == 21 ? 0 : id);
270  }
271 
272  if (pdf->scalePDF() < 0)
273  pdf->set_scalePDF(std::sqrt(scale) / ThePEG::GeV);
274  else
275  scale = ThePEG::sqr(pdf->scalePDF()) * ThePEG::GeV;
276 
277  if (pdf->x1() < 0) {
278  x1 = eh->lastX1();
279  pdf->set_x1(x1);
280  } else
281  x1 = pdf->x1();
282 
283  if (pdf->x2() < 0) {
284  x2 = eh->lastX2();
285  pdf->set_x2(x2);
286  } else
287  x2 = pdf->x2();
288 
289  if (pdf1 && pdf->pdf1() < 0) {
290  double v = pdf1->xfx(beams.first->dataPtr(),
291  partons.first->dataPtr(),
292  scale, x1);
293  if (v > 0 && x1 > 0)
294  pdf->set_pdf1(v / x1);
295  else
296  pdf->set_pdf2(-1.0);
297  }
298  if (pdf2 && pdf->pdf2() < 0) {
299  double v = pdf2->xfx(beams.first->dataPtr(),
300  partons.first->dataPtr(),
301  scale, x2);
302  if (v > 0 && x2 > 0)
303  pdf->set_pdf2(v / x2);
304  else
305  pdf->set_pdf2(-1.0);
306  }
307  }
308 
309 }
310 
311 double ThePEGInterface::pthat(const ThePEG::EventPtr &event)
312 {
313  using namespace ThePEG;
314 
315  if (!event->primaryCollision())
316  return -1.0;
317 
318  tSubProPtr sub = event->primaryCollision()->primarySubProcess();
319  TmpTransform<tSubProPtr> tmp(sub, Utilities::getBoostToCM(
320  sub->incoming()));
321 
322  double pthat = (*sub->outgoing().begin())->momentum().perp();
323  for(PVector::const_iterator it = sub->outgoing().begin();
324  it != sub->outgoing().end(); ++it)
325  pthat = std::min(pthat, (*it)->momentum().perp());
326 
327  return pthat / ThePEG::GeV;
328 }
T getParameter(std::string const &) const
Entry const & retrieve(char const *) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
void initRepository(const edm::ParameterSet &params) const
ThePEGInterface(const edm::ParameterSet &params)
static double pthat(const ThePEG::EventPtr &event)
std::string dataFile(const std::string &fileName) const
const std::string run_
T sqrt(T t)
Definition: SSEVec.h:48
def load
Definition: svgfig.py:546
void flushRandomNumberGenerator()
std::auto_ptr< HepMC::IO_BaseClass > iobc_
ThePEG::EGPtr eg_
const std::string dumpConfig_
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
static std::auto_ptr< HepMC::GenEvent > convert(const ThePEG::EventPtr &event)
static GenEvent * convert(const Event &ev, bool nocopies=false, Energy eunit=Traits::defaultEnergyUnit(), Length lunit=Traits::defaultLengthUnit())
tuple out
Definition: dbtoconf.py:99
const std::string generator_
FileInPath getFileInPath() const
Definition: Entry.cc:787
virtual ~ThePEGInterface()
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
const unsigned int skipEvents_
const std::string dataLocation_
Square< F >::type sqr(const F &f)
Definition: Square.h:13
const_iterator end() const
const_iterator begin() const
char typeCode() const
Definition: Entry.h:175
std::string fullPath() const
Definition: FileInPath.cc:171
boost::shared_ptr< ThePEG::RandomEngineGlue::Proxy > randomEngineGlueProxy_
SurfaceDeformation * create(int type, const std::vector< double > &params)
static void clearAuxiliary(HepMC::GenEvent *hepmc, HepMC::PdfInfo *pdf)
std::string getString() const
Definition: Entry.cc:764
static void fillAuxiliary(HepMC::GenEvent *hepmc, HepMC::PdfInfo *pdf, const ThePEG::EventPtr &event)