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 
8 #include <iostream>
9 #include <sstream>
10 #include <cstdlib>
11 #include <memory>
12 #include <cmath>
13 #include <stdlib.h>
14 
15 #include <boost/shared_ptr.hpp>
16 #include <boost/filesystem.hpp>
17 
18 #include <HepMC/GenEvent.h>
19 #include <HepMC/PdfInfo.h>
20 #include <HepMC/IO_GenEvent.h>
21 
22 #include <ThePEG/Utilities/DynamicLoader.h>
23 #include <ThePEG/Repository/Repository.h>
24 #include <ThePEG/Handlers/EventHandler.h>
25 #include <ThePEG/Handlers/XComb.h>
26 #include <ThePEG/EventRecord/Event.h>
27 #include <ThePEG/EventRecord/Particle.h>
28 #include <ThePEG/EventRecord/Collision.h>
29 #include <ThePEG/EventRecord/TmpTransform.h>
30 #include <ThePEG/Config/ThePEG.h>
31 #include <ThePEG/PDF/PartonExtractor.h>
32 #include <ThePEG/PDF/PDFBase.h>
33 #include <ThePEG/Utilities/UtilityBase.h>
34 
36 
38 
43 
44 using namespace std;
45 using namespace gen;
46 
48  randomEngineGlueProxy_(ThePEG::RandomEngineGlue::Proxy::create()),
49  dataLocation_(ParameterCollector::resolve(pset.getParameter<string>("dataLocation"))),
50  generator_(pset.getParameter<string>("generatorModule")),
51  run_(pset.getParameter<string>("run")),
52  dumpConfig_(pset.getUntrackedParameter<string>("dumpConfig", "")),
53  skipEvents_(pset.getUntrackedParameter<unsigned int>("skipEvents", 0))
54 {
55  // Write events in hepmc ascii format for debugging purposes
56  string dumpEvents = pset.getUntrackedParameter<string>("dumpEvents", "");
57  if (!dumpEvents.empty()) {
58  iobc_.reset(new HepMC::IO_GenEvent(dumpEvents.c_str(), ios::out));
59  edm::LogInfo("ThePEGSource") << "Event logging switched on (=> " << dumpEvents << ")";
60  }
61  // Clear dumpConfig target
62  if (!dumpConfig_.empty())
63  ofstream cfgDump(dumpConfig_.c_str(), ios_base::trunc);
64 }
65 
67 {
68  if (eg_)
69  eg_->finalize();
70  edm::LogInfo("ThePEGInterface") << "Event generator finalized";
71 }
72 
73 string ThePEGInterface::dataFile(const string &fileName) const
74 {
75  if (fileName.empty() || fileName[0] == '/')
76  return fileName;
77  else
78  return dataLocation_ + "/" + fileName;
79 }
80 
82  const string &paramName) const
83 {
84  const edm::Entry &entry = pset.retrieve(paramName);
85  if (entry.typeCode() == 'F')
86  return entry.getFileInPath().fullPath();
87  else
88  return dataFile(entry.getString());
89 }
90 
92 {
93  /* Initialize the repository from
94  * 1. the repository file (default values)
95  * 2. the ThePEG config files
96  * 3. the CMSSW config blocks
97  */
98  stringstream logstream;
99 
100  // Read the repository of serialized default values
101  string repository = dataFile(pset, "repository");
102  if (!repository.empty()) {
103  edm::LogInfo("ThePEGInterface") << "Loading repository (" << repository << ")";
104  ThePEG::Repository::load(repository);
105  }
106 
107  if (!getenv("ThePEG_INSTALL_PATH")) {
108  vector<string> libdirlist = ThePEG::DynamicLoader::allPaths();
109  for(vector<string>::const_iterator libdir = libdirlist.begin();
110  libdir < libdirlist.end(); ++libdir) {
111  if (libdir->empty() || (*libdir)[0] != '/')
112  continue;
113  if (boost::filesystem::exists(*libdir +
114  "/../../share/ThePEG/PDFsets.index")) {
115  setenv("ThePEG_INSTALL_PATH",
116  libdir->c_str(), 0);
117  break;
118  }
119  }
120  }
121 
122  // Read ThePEG config files to read
123  vector<string> configFiles = pset.getParameter<vector<string> >("configFiles");
124 
125  // Loop over the config files
126  for(vector<string>::const_iterator iter = configFiles.begin();
127  iter != configFiles.end(); ++iter) {
128  edm::LogInfo("ThePEGInterface") << "Reading config file (" << *iter << ")";
129  ThePEG::Repository::read(dataFile(*iter), logstream);
130  edm::LogInfo("ThePEGSource") << logstream.str();
131  }
132 
133  // Read CMSSW config file parameter sets starting from "parameterSets"
134 
135  ofstream cfgDump;
136  ParameterCollector collector(pset);
138  if (!dumpConfig_.empty()) {
139  cfgDump.open(dumpConfig_.c_str(), ios_base::app);
140  iter = collector.begin(cfgDump);
141  } else
142  iter = collector.begin();
143 
144  for(; iter != collector.end(); ++iter) {
145  string out = ThePEG::Repository::exec(*iter, logstream);
146  if (!out.empty()) {
147  edm::LogInfo("ThePEGInterface") << *iter << " => " << out;
148  cerr << "Error in ThePEG configuration!\n"
149  "\tLine: " << *iter << "\n" << out << endl;
150  }
151  }
152 
153  if (!dumpConfig_.empty()) {
154  cfgDump << "saverun " << run_ << " " << generator_ << endl;
155  cfgDump.close();
156  }
157 
158  // write the ProxyID for the RandomEngineGlue to fill its pointer in
159  ostringstream ss;
160  ss << randomEngineGlueProxy_->getID();
161  ThePEG::Repository::exec("set " + generator_ +
162  ":RandomNumberGenerator:ProxyID " + ss.str(),
163  logstream);
164 
165  // Print the directories where ThePEG looks for libs
166  vector<string> libdirlist = ThePEG::DynamicLoader::allPaths();
167  for(vector<string>::const_iterator libdir = libdirlist.begin();
168  libdir < libdirlist.end(); ++libdir)
169  edm::LogInfo("ThePEGInterface") << "DynamicLoader path = " << *libdir << endl;
170 
171  // Output status information about the repository
172  ThePEG::Repository::stats(logstream);
173  edm::LogInfo("ThePEGInterface") << logstream.str();
174 }
175 
177 {
178  // Get generator from the repository and initialize it
179  ThePEG::BaseRepository::CheckObjectDirectory(generator_);
180  ThePEG::EGPtr tmp = ThePEG::BaseRepository::GetObject<ThePEG::EGPtr>(generator_);
181  if (tmp) {
182  eg_ = ThePEG::Repository::makeRun(tmp, run_);
183  eg_->initialize();
184  edm::LogInfo("ThePEGInterface") << "EventGenerator initialized";
185  } else
186  throw cms::Exception("ThePEGInterface")
187  << "EventGenerator could not be initialized!" << endl;
188 
189  // Skip events
190  for (unsigned int i = 0; i < skipEvents_; i++) {
192  eg_->shoot();
193  edm::LogInfo("ThePEGInterface") << "Event discarded";
194  }
195 }
196 
198 {
199  ThePEG::RandomEngineGlue *rnd = randomEngineGlueProxy_->getInstance();
200 
201  if (!rnd)
202  edm::LogWarning("ProxyMissing")
203  << "ThePEG not initialised with RandomEngineGlue.";
204  else
205  rnd->flush();
206 }
207 
208 auto_ptr<HepMC::GenEvent> ThePEGInterface::convert(
209  const ThePEG::EventPtr &event)
210 {
211  return std::auto_ptr<HepMC::GenEvent>(
213 }
214 
215 void ThePEGInterface::clearAuxiliary(HepMC::GenEvent *hepmc,
216  HepMC::PdfInfo *pdf)
217 {
218  if (hepmc) {
219  hepmc->set_event_scale(-1.0);
220  hepmc->set_alphaQCD(-1.0);
221  hepmc->set_alphaQED(-1.0);
222  }
223  if (pdf) {
224  pdf->set_id1(-100);
225  pdf->set_id2(-100);
226  pdf->set_x1(-1.0);
227  pdf->set_x2(-1.0);
228  pdf->set_scalePDF(-1.0);
229  pdf->set_pdf1(-1.0);
230  pdf->set_pdf2(-1.0);
231  }
232 }
233 
234 void ThePEGInterface::fillAuxiliary(HepMC::GenEvent *hepmc,
235  HepMC::PdfInfo *pdf,
236  const ThePEG::EventPtr &event)
237 {
238  if (!event->primaryCollision())
239  return;
240 
241  ThePEG::tcEHPtr eh = ThePEG::dynamic_ptr_cast<ThePEG::tcEHPtr>(
242  event->primaryCollision()->handler());
243  double scale = eh->lastScale();
244 
245  if (hepmc) {
246  if (hepmc->event_scale() < 0 && scale > 0)
247  hepmc->set_event_scale(std::sqrt(scale) / ThePEG::GeV);
248 
249  if (hepmc->alphaQCD() < 0)
250  hepmc->set_alphaQCD(eh->lastAlphaS());
251  if (hepmc->alphaQED() < 0)
252  hepmc->set_alphaQED(eh->lastAlphaEM());
253  }
254 
255  if (pdf) {
256  const ThePEG::PPair &beams = eh->lastParticles();
257  const ThePEG::PPair &partons = eh->lastPartons();
258  ThePEG::tcPDFPtr pdf1 = eh->lastExtractor()->getPDF(
259  beams.first->dataPtr());
260  ThePEG::tcPDFPtr pdf2 = eh->lastExtractor()->getPDF(
261  beams.second->dataPtr());
262  double x1, x2;
263 
264  if (pdf->id1() == -100) {
265  int id = partons.first->id();
266  pdf->set_id1(id == 21 ? 0 : id);
267  }
268  if (pdf->id2() == -100) {
269  int id = partons.second->id();
270  pdf->set_id2(id == 21 ? 0 : id);
271  }
272 
273  if (pdf->scalePDF() < 0)
274  pdf->set_scalePDF(std::sqrt(scale) / ThePEG::GeV);
275  else
276  scale = ThePEG::sqr(pdf->scalePDF()) * ThePEG::GeV;
277 
278  if (pdf->x1() < 0) {
279  x1 = eh->lastX1();
280  pdf->set_x1(x1);
281  } else
282  x1 = pdf->x1();
283 
284  if (pdf->x2() < 0) {
285  x2 = eh->lastX2();
286  pdf->set_x2(x2);
287  } else
288  x2 = pdf->x2();
289 
290  if (pdf1 && pdf->pdf1() < 0) {
291  double v = pdf1->xfx(beams.first->dataPtr(),
292  partons.first->dataPtr(),
293  scale, x1);
294  if (v > 0 && x1 > 0)
295  pdf->set_pdf1(v / x1);
296  else
297  pdf->set_pdf2(-1.0);
298  }
299  if (pdf2 && pdf->pdf2() < 0) {
300  double v = pdf2->xfx(beams.first->dataPtr(),
301  partons.first->dataPtr(),
302  scale, x2);
303  if (v > 0 && x2 > 0)
304  pdf->set_pdf2(v / x2);
305  else
306  pdf->set_pdf2(-1.0);
307  }
308  }
309 
310 }
311 
313 {
314  using namespace ThePEG;
315 
316  if (!event->primaryCollision())
317  return -1.0;
318 
319  tSubProPtr sub = event->primaryCollision()->primarySubProcess();
320  TmpTransform<tSubProPtr> tmp(sub, Utilities::getBoostToCM(
321  sub->incoming()));
322 
323  double pthat = (*sub->outgoing().begin())->momentum().perp();
324  for(PVector::const_iterator it = sub->outgoing().begin();
325  it != sub->outgoing().end(); ++it)
326  pthat = std::min(pthat, (*it)->momentum().perp());
327 
328  return pthat / ThePEG::GeV;
329 }
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)
#define min(a, b)
Definition: mlp_lapack.h:161
boost::shared_ptr< boost::statechart::event_base > EventPtr
Definition: CommandQueue.h:23
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
std::pair< std::string, MonitorElement * > entry
Definition: ME_MAP.h:8
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)