CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Static Public Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes
lhef::LHEEvent Class Reference

#include <LHEEvent.h>

Public Types

typedef LHEEventProduct::PDF PDF
 
typedef LHEEventProduct::WGT WGT
 

Public Member Functions

void addComment (const std::string &line)
 
void addWeight (const WGT &wgt)
 
std::auto_ptr< HepMC::GenEventasHepMCEvent () const
 
void attempted ()
 
void count (LHERunInfo::CountMode count, double weight=1.0, double matchWeight=1.0)
 
void fillEventInfo (HepMC::GenEvent *hepmc) const
 
void fillPdfInfo (HepMC::PdfInfo *info) const
 
const std::vector< std::string > & getComments () const
 
const HEPEUPgetHEPEUP () const
 
const HEPRUPgetHEPRUP () const
 
const PDFgetPDF () const
 
const int getReadAttempts ()
 
const boost::shared_ptr< LHERunInfo > & getRunInfo () const
 
 LHEEvent (const boost::shared_ptr< LHERunInfo > &runInfo, std::istream &in)
 
 LHEEvent (const boost::shared_ptr< LHERunInfo > &runInfo, const HEPEUP &hepeup)
 
 LHEEvent (const boost::shared_ptr< LHERunInfo > &runInfo, const HEPEUP &hepeup, const LHEEventProduct::PDF *pdf, const std::vector< std::string > &comments)
 
 LHEEvent (const boost::shared_ptr< LHERunInfo > &runInfo, const LHEEventProduct &product)
 
int npLO () const
 
int npNLO () const
 
double originalXWGTUP () const
 
void removeResonances (const std::vector< int > &ids)
 
const std::vector< float > & scales () const
 
void setNpLO (int n)
 
void setNpNLO (int n)
 
void setPDF (std::auto_ptr< PDF > pdf)
 
void setScales (const std::vector< float > &scales)
 
const std::vector< WGT > & weights () const
 
 ~LHEEvent ()
 

Static Public Member Functions

static const HepMC::GenVertex * findSignalVertex (const HepMC::GenEvent *event, bool status3=true)
 
static void fixHepMCEventTimeOrdering (HepMC::GenEvent *event)
 
static void removeParticle (lhef::HEPEUP &hepeup, int index)
 

Private Member Functions

HepMC::GenParticle * makeHepMCParticle (unsigned int i) const
 

Static Private Member Functions

static bool checkHepMCTree (const HepMC::GenEvent *event)
 

Private Attributes

std::vector< std::string > comments
 
bool counted
 
HEPEUP hepeup
 
int npLO_
 
int npNLO_
 
double originalXWGTUP_
 
std::auto_ptr< PDFpdf
 
int readAttemptCounter
 
const boost::shared_ptr< LHERunInforunInfo
 
std::vector< float > scales_
 
std::vector< WGTweights_
 

Detailed Description

Definition at line 25 of file LHEEvent.h.

Member Typedef Documentation

Definition at line 39 of file LHEEvent.h.

Definition at line 40 of file LHEEvent.h.

Constructor & Destructor Documentation

lhef::LHEEvent::LHEEvent ( const boost::shared_ptr< LHERunInfo > &  runInfo,
std::istream &  in 
)

Definition at line 36 of file LHEEvent.cc.

References funct::abs(), lhef::HEPEUP::AQCDUP, lhef::HEPEUP::AQEDUP, comments, hepeup, mps_fire::i, lhef::HEPEUP::ICOLUP, lhef::HEPEUP::IDPRUP, lhef::HEPEUP::IDUP, lhef::HEPEUP::ISTUP, geometryCSVtoXML::line, lhef::HEPEUP::MOTHUP, lhef::HEPEUP::NUP, originalXWGTUP_, pdf, lhef::HEPEUP::PUP, lhef::HEPEUP::resize(), lhef::HEPEUP::SCALUP, skipWhitespace(), lhef::HEPEUP::SPINUP, AlCaHLTBitMon_QueryRunRegistry::string, GlobalPosition_Frontier_DevDB_cff::tag, lhef::HEPEUP::VTIMUP, lhef::HEPEUP::XPDWUP, and lhef::HEPEUP::XWGTUP.

37  :
38  runInfo(runInfo), weights_(0), counted(false),
39  readAttemptCounter(0), npLO_(-99), npNLO_(-99)
40 
41 {
42  hepeup.NUP = 0;
43  hepeup.XPDWUP.first = hepeup.XPDWUP.second = 0.0;
44 
47  if (!in.good())
48  throw cms::Exception("InvalidFormat")
49  << "Les Houches file contained invalid"
50  " event header." << std::endl;
51 
52  // store the original value of XWGTUP for the user
54 
55  int idwtup = runInfo->getHEPRUP()->IDWTUP;
56  if (idwtup >= 0 && hepeup.XWGTUP < 0) {
57  edm::LogWarning("Generator|LHEInterface")
58  << "Non-allowed negative event weight encountered."
59  << std::endl;
61  }
62 
63  if (std::abs(idwtup) == 3 && std::abs(hepeup.XWGTUP) != 1.) {
64  edm::LogInfo("Generator|LHEInterface")
65  << "Event weight not set to one for abs(IDWTUP) == 3"
66  << std::endl;
67  hepeup.XWGTUP = hepeup.XWGTUP > 0. ? 1.0 : -1.0;
68  }
69 
70  hepeup.resize();
71 
72  for(int i = 0; i < hepeup.NUP; i++) {
73  in >> hepeup.IDUP[i] >> hepeup.ISTUP[i]
74  >> hepeup.MOTHUP[i].first >> hepeup.MOTHUP[i].second
75  >> hepeup.ICOLUP[i].first >> hepeup.ICOLUP[i].second
76  >> hepeup.PUP[i][0] >> hepeup.PUP[i][1] >> hepeup.PUP[i][2]
77  >> hepeup.PUP[i][3] >> hepeup.PUP[i][4]
78  >> hepeup.VTIMUP[i] >> hepeup.SPINUP[i];
79  if (!in.good())
80  throw cms::Exception("InvalidFormat")
81  << "Les Houches file contained invalid event"
82  " in particle line " << (i + 1)
83  << "." << std::endl;
84  }
85 
86  while(skipWhitespace(in) == '#') {
88  std::getline(in, line);
89  std::istringstream ss(line);
91  ss >> tag;
92  if (tag == "#pdf") {
93  pdf.reset(new PDF);
94  ss >> pdf->id.first >> pdf->id.second
95  >> pdf->x.first >> pdf->x.second
96  >> pdf->scalePDF
97  >> pdf->xPDF.first >> pdf->xPDF.second;
98  if (ss.bad()) {
99  edm::LogWarning("Generator|LHEInterface")
100  << "Les Houches event contained"
101  " unparseable PDF information."
102  << std::endl;
103  pdf.reset();
104  } else
105  continue;
106  }
107  comments.push_back(line + "\n");
108  }
109 
110  if (!in.eof())
111  edm::LogWarning("Generator|LHEInterface")
112  << "Les Houches file contained spurious"
113  " content after event data." << std::endl;
114 }
HEPEUP hepeup
Definition: LHEEvent.h:90
std::vector< std::string > comments
Definition: LHEEvent.h:93
std::vector< double > VTIMUP
Definition: LesHouches.h:254
std::pair< double, double > XPDWUP
Definition: LesHouches.h:204
std::auto_ptr< PDF > pdf
Definition: LHEEvent.h:91
int readAttemptCounter
Definition: LHEEvent.h:95
std::vector< std::pair< int, int > > MOTHUP
Definition: LesHouches.h:236
void resize(int nup)
Definition: LesHouches.h:161
const boost::shared_ptr< LHERunInfo > runInfo
Definition: LHEEvent.h:88
std::vector< FiveVector > PUP
Definition: LesHouches.h:248
std::vector< double > SPINUP
Definition: LesHouches.h:261
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< int > ISTUP
Definition: LesHouches.h:230
LHEEventProduct::PDF PDF
Definition: LHEEvent.h:39
std::vector< int > IDUP
Definition: LesHouches.h:225
double AQCDUP
Definition: LesHouches.h:220
bool counted
Definition: LHEEvent.h:94
double AQEDUP
Definition: LesHouches.h:215
static int skipWhitespace(std::istream &in)
Definition: LHEEvent.cc:23
double XWGTUP
Definition: LesHouches.h:196
std::vector< std::pair< int, int > > ICOLUP
Definition: LesHouches.h:242
std::vector< WGT > weights_
Definition: LHEEvent.h:92
double originalXWGTUP_
Definition: LHEEvent.h:96
double SCALUP
Definition: LesHouches.h:210
lhef::LHEEvent::LHEEvent ( const boost::shared_ptr< LHERunInfo > &  runInfo,
const HEPEUP hepeup 
)

Definition at line 116 of file LHEEvent.cc.

117  :
119  npLO_(-99), npNLO_(-99)
120 {
121 }
HEPEUP hepeup
Definition: LHEEvent.h:90
int readAttemptCounter
Definition: LHEEvent.h:95
const boost::shared_ptr< LHERunInfo > runInfo
Definition: LHEEvent.h:88
bool counted
Definition: LHEEvent.h:94
lhef::LHEEvent::LHEEvent ( const boost::shared_ptr< LHERunInfo > &  runInfo,
const HEPEUP hepeup,
const LHEEventProduct::PDF pdf,
const std::vector< std::string > &  comments 
)

Definition at line 123 of file LHEEvent.cc.

126  :
127  runInfo(runInfo), hepeup(hepeup), pdf(pdf ? new PDF(*pdf) : nullptr),
129  npLO_(-99), npNLO_(-99)
130 {
131 }
HEPEUP hepeup
Definition: LHEEvent.h:90
std::vector< std::string > comments
Definition: LHEEvent.h:93
std::auto_ptr< PDF > pdf
Definition: LHEEvent.h:91
int readAttemptCounter
Definition: LHEEvent.h:95
const boost::shared_ptr< LHERunInfo > runInfo
Definition: LHEEvent.h:88
LHEEventProduct::PDF PDF
Definition: LHEEvent.h:39
bool counted
Definition: LHEEvent.h:94
lhef::LHEEvent::LHEEvent ( const boost::shared_ptr< LHERunInfo > &  runInfo,
const LHEEventProduct product 
)

Definition at line 133 of file LHEEvent.cc.

134  :
135  runInfo(runInfo), hepeup(product.hepeup()),
136  pdf(product.pdf() ? new PDF(*product.pdf()) : nullptr),
137  weights_(product.weights()),
138  comments(product.comments_begin(), product.comments_end()),
139  counted(false), readAttemptCounter(0),
140  originalXWGTUP_(product.originalXWGTUP()), scales_(product.scales()),
141  npLO_(product.npLO()), npNLO_(product.npNLO())
142 {
143 }
double originalXWGTUP() const
HEPEUP hepeup
Definition: LHEEvent.h:90
const lhef::HEPEUP & hepeup() const
const PDF * pdf() const
std::vector< std::string > comments
Definition: LHEEvent.h:93
std::auto_ptr< PDF > pdf
Definition: LHEEvent.h:91
const std::vector< WGT > & weights() const
int readAttemptCounter
Definition: LHEEvent.h:95
comments_const_iterator comments_begin() const
const boost::shared_ptr< LHERunInfo > runInfo
Definition: LHEEvent.h:88
LHEEventProduct::PDF PDF
Definition: LHEEvent.h:39
int npLO() const
const std::vector< float > & scales() const
comments_const_iterator comments_end() const
std::vector< float > scales_
Definition: LHEEvent.h:97
bool counted
Definition: LHEEvent.h:94
int npNLO() const
std::vector< WGT > weights_
Definition: LHEEvent.h:92
double originalXWGTUP_
Definition: LHEEvent.h:96
lhef::LHEEvent::~LHEEvent ( )

Definition at line 145 of file LHEEvent.cc.

146 {
147 }

Member Function Documentation

void lhef::LHEEvent::addComment ( const std::string &  line)
inline

Definition at line 64 of file LHEEvent.h.

References count(), photons_cff::ids, removeParticle(), and removeResonances().

64 { comments.push_back(line); }
std::vector< std::string > comments
Definition: LHEEvent.h:93
void lhef::LHEEvent::addWeight ( const WGT wgt)
inline

Definition at line 49 of file LHEEvent.h.

References weights_.

49 { weights_.push_back(wgt); }
std::vector< WGT > weights_
Definition: LHEEvent.h:92
std::auto_ptr< HepMC::GenEvent > lhef::LHEEvent::asHepMCEvent ( ) const

Definition at line 257 of file LHEEvent.cc.

References lhef::HEPEUP::AQCDUP, lhef::HEPEUP::AQEDUP, checkHepMCTree(), findSignalVertex(), GenParticle::GenParticle, GenHFHadronMatcher_cfi::genParticles, getHEPRUP(), hepeup, mps_fire::i, lhef::HEPEUP::IDPRUP, makeHepMCParticle(), lhef::HEPEUP::MOTHUP, lhef::HEPEUP::NUP, lhef::HEPEUP::SCALUP, and lhef::HEPEUP::VTIMUP.

Referenced by attempted().

258 {
259  std::auto_ptr<HepMC::GenEvent> hepmc(new HepMC::GenEvent);
260 
261  hepmc->set_signal_process_id(hepeup.IDPRUP);
262  hepmc->set_event_scale(hepeup.SCALUP);
263  hepmc->set_alphaQED(hepeup.AQEDUP);
264  hepmc->set_alphaQCD(hepeup.AQCDUP);
265 
266  unsigned int nup = hepeup.NUP; // particles in event
267 
268  // any particles in HEPEUP block?
269  if (!nup) {
270  edm::LogWarning("Generator|LHEInterface")
271  << "Les Houches Event does not contain any partons. "
272  << "Not much to convert." ;
273  return hepmc;
274  }
275 
276  // stores (pointers to) converted particles
277  std::vector<HepMC::GenParticle*> genParticles;
278  std::vector<HepMC::GenVertex*> genVertices;
279 
280  // I. convert particles
281  for(unsigned int i = 0; i < nup; i++)
282  genParticles.push_back(makeHepMCParticle(i));
283 
284  // II. loop again to build vertices
285  for(unsigned int i = 0; i < nup; i++) {
286  unsigned int mother1 = hepeup.MOTHUP.at(i).first;
287  unsigned int mother2 = hepeup.MOTHUP.at(i).second;
288  double cTau = hepeup.VTIMUP.at(i); // decay time
289 
290  // current particle has a mother? --- Sorry, parent! We're PC.
291  if (mother1) {
292  mother1--; // FORTRAN notation!
293  if (mother2)
294  mother2--;
295  else
296  mother2 = mother1;
297 
298  HepMC::GenParticle *in_par = genParticles.at(mother1);
299  HepMC::GenVertex *current_vtx = in_par->end_vertex(); // vertex of first mother
300 
301  if (!current_vtx) {
302  current_vtx = new HepMC::GenVertex(
303  HepMC::FourVector(0, 0, 0, cTau));
304 
305  // add vertex to event
306  genVertices.push_back(current_vtx);
307  }
308 
309  for(unsigned int j = mother1; j <= mother2; j++) // set mother-daughter relations
310  if (!genParticles.at(j)->end_vertex())
311  current_vtx->add_particle_in(genParticles.at(j));
312 
313  // connect THIS outgoing particle to current vertex
314  current_vtx->add_particle_out(genParticles.at(i));
315  }
316  }
317 
318  checkHepMCTree(hepmc.get());
319 
320  // III. restore color flow
321  // ok, nobody knows how to do it so far...
322 
323  // IV. fill run information
324  const HEPRUP *heprup = getHEPRUP();
325 
326  // set beam particles
328  HepMC::FourVector(0.0, 0.0, +heprup->EBMUP.first,
329  heprup->EBMUP.first),
330  heprup->IDBMUP.first);
332  HepMC::FourVector(0.0, 0.0, -heprup->EBMUP.second,
333  heprup->EBMUP.second),
334  heprup->IDBMUP.second);
335  b1->set_status(3);
336  b2->set_status(3);
337 
338  HepMC::GenVertex *v1 = new HepMC::GenVertex();
339  HepMC::GenVertex *v2 = new HepMC::GenVertex();
340  v1->add_particle_in(b1);
341  v2->add_particle_in(b2);
342 
343  hepmc->add_vertex(v1);
344  hepmc->add_vertex(v2);
345  hepmc->set_beam_particles(b1, b2);
346 
347  // first two particles have to be the hard partons going into the interaction
348  if (genParticles.size() >= 2) {
349  if (!genParticles.at(0)->production_vertex() &&
350  !genParticles.at(1)->production_vertex()) {
351  v1->add_particle_out(genParticles.at(0));
352  v2->add_particle_out(genParticles.at(1));
353  } else
354  edm::LogWarning("Generator|LHEInterface")
355  << "Initial partons do already have a"
356  " production vertex. " << std::endl
357  << "Beam particles not connected.";
358  } else
359  edm::LogWarning("Generator|LHEInterface")
360  << "Can't find any initial partons to be"
361  " connected to the beam particles.";
362 
363  for(std::vector<HepMC::GenVertex*>::const_iterator iter = genVertices.begin();
364  iter != genVertices.end(); ++iter)
365  hepmc->add_vertex(*iter);
366 
367  // do some more consistency checks
368  for(unsigned int i = 0; i < nup; i++) {
369  if (!genParticles.at(i)->parent_event()) {
370  edm::LogWarning("Generator|LHEInterface")
371  << "Not all LHE particles could be stored"
372  " stored in the HepMC event. "
373  << std::endl
374  << "Check the mother-daughter relations"
375  " in the given LHE input file.";
376  break;
377  }
378  }
379 
380  hepmc->set_signal_process_vertex(
381  const_cast<HepMC::GenVertex*>(
382  findSignalVertex(hepmc.get(), false)));
383 
384  return hepmc;
385 }
HEPEUP hepeup
Definition: LHEEvent.h:90
const HEPRUP * getHEPRUP() const
Definition: LHEEvent.h:44
std::vector< double > VTIMUP
Definition: LesHouches.h:254
std::vector< std::pair< int, int > > MOTHUP
Definition: LesHouches.h:236
HepMC::GenParticle * makeHepMCParticle(unsigned int i) const
Definition: LHEEvent.cc:387
double AQCDUP
Definition: LesHouches.h:220
static const HepMC::GenVertex * findSignalVertex(const HepMC::GenEvent *event, bool status3=true)
Definition: LHEEvent.cc:448
double AQEDUP
Definition: LesHouches.h:215
static bool checkHepMCTree(const HepMC::GenEvent *event)
Definition: LHEEvent.cc:404
double SCALUP
Definition: LesHouches.h:210
void lhef::LHEEvent::attempted ( )
inline
bool lhef::LHEEvent::checkHepMCTree ( const HepMC::GenEvent event)
staticprivate

Definition at line 404 of file LHEEvent.cc.

References mps_update::status.

Referenced by asHepMCEvent(), and attempted().

405 {
406  double px = 0, py = 0, pz = 0, E = 0;
407 
408  for(HepMC::GenEvent::particle_const_iterator iter = event->particles_begin();
409  iter != event->particles_end(); iter++) {
410  int status = (*iter)->status();
411  HepMC::FourVector fv = (*iter)->momentum();
412 
413  // incoming particles
414  if (status == 3 &&
415  *iter != event->beam_particles().first &&
416  *iter != event->beam_particles().second) {
417  px -= fv.px();
418  py -= fv.py();
419  pz -= fv.pz();
420  E -= fv.e();
421  }
422 
423  // outgoing particles
424  if (status == 1) {
425  px += fv.px();
426  py += fv.py();
427  pz += fv.pz();
428  E += fv.e();
429  }
430  }
431 
432  if (px*px + py*py + pz*pz + E*E > 0.1) {
433  edm::LogWarning("Generator|LHEInterface")
434  << "Energy-momentum badly conserved. "
435  << std::setprecision(3)
436  << "sum p_i = ["
437  << std::setw(7) << E << ", "
438  << std::setw(7) << px << ", "
439  << std::setw(7) << py << ", "
440  << std::setw(7) << pz << "]";
441 
442  return false;
443  }
444 
445  return true;
446 }
void lhef::LHEEvent::count ( LHERunInfo::CountMode  count,
double  weight = 1.0,
double  matchWeight = 1.0 
)

Definition at line 207 of file LHEEvent.cc.

References counted, hepeup, lhef::HEPEUP::IDPRUP, runInfo, and lhef::HEPEUP::XWGTUP.

Referenced by addComment(), gen::Pythia6Hadronizer::hadronize(), Pythia8Hadronizer::hadronize(), and Herwig6Hadronizer::hadronize().

209 {
210  if (counted)
211  edm::LogWarning("Generator|LHEInterface")
212  << "LHEEvent::count called twice on same event!"
213  << std::endl;
214 
216  weight, matchWeight);
217 
218  counted = true;
219 }
HEPEUP hepeup
Definition: LHEEvent.h:90
Definition: weight.py:1
const boost::shared_ptr< LHERunInfo > runInfo
Definition: LHEEvent.h:88
bool counted
Definition: LHEEvent.h:94
double XWGTUP
Definition: LesHouches.h:196
void lhef::LHEEvent::fillEventInfo ( HepMC::GenEvent hepmc) const

Definition at line 249 of file LHEEvent.cc.

References lhef::HEPEUP::AQCDUP, lhef::HEPEUP::AQEDUP, hepeup, lhef::HEPEUP::IDPRUP, and lhef::HEPEUP::SCALUP.

Referenced by attempted(), gen::Pythia6Hadronizer::finalizeEvent(), and Herwig6Hadronizer::finalizeEvent().

250 {
251  event->set_signal_process_id(hepeup.IDPRUP);
252  event->set_event_scale(hepeup.SCALUP);
253  event->set_alphaQED(hepeup.AQEDUP);
254  event->set_alphaQCD(hepeup.AQCDUP);
255 }
HEPEUP hepeup
Definition: LHEEvent.h:90
double AQCDUP
Definition: LesHouches.h:220
double AQEDUP
Definition: LesHouches.h:215
double SCALUP
Definition: LesHouches.h:210
void lhef::LHEEvent::fillPdfInfo ( HepMC::PdfInfo *  info) const

Definition at line 221 of file LHEEvent.cc.

References funct::abs(), lhef::HEPRUP::EBMUP, getHEPRUP(), hepeup, lhef::HEPEUP::IDUP, lhef::HEPEUP::NUP, pdf, lhef::HEPEUP::PUP, and lhef::HEPEUP::SCALUP.

Referenced by attempted(), gen::Pythia6Hadronizer::finalizeEvent(), and Herwig6Hadronizer::finalizeEvent().

222 {
223  if (pdf.get()) {
224  info->set_id1(pdf->id.first);
225  info->set_id2(pdf->id.second);
226  info->set_x1(pdf->x.first);
227  info->set_x2(pdf->x.second);
228  info->set_pdf1(pdf->xPDF.first);
229  info->set_pdf2(pdf->xPDF.second);
230  info->set_scalePDF(pdf->scalePDF);
231  } else if (hepeup.NUP >= 2) {
232  const HEPRUP *heprup = getHEPRUP();
233  info->set_id1(hepeup.IDUP[0] == 21 ? 0 : hepeup.IDUP[0]);
234  info->set_id2(hepeup.IDUP[1] == 21 ? 0 : hepeup.IDUP[1]);
235  info->set_x1(std::abs(hepeup.PUP[0][2] / heprup->EBMUP.first));
236  info->set_x2(std::abs(hepeup.PUP[1][2] / heprup->EBMUP.second));
237  info->set_pdf1(-1.0);
238  info->set_pdf2(-1.0);
239  info->set_scalePDF(hepeup.SCALUP);
240  } else {
241  info->set_x1(-1.0);
242  info->set_x2(-1.0);
243  info->set_pdf1(-1.0);
244  info->set_pdf2(-1.0);
245  info->set_scalePDF(hepeup.SCALUP);
246  }
247 }
HEPEUP hepeup
Definition: LHEEvent.h:90
static const TGPicture * info(bool iBackgroundIsBlack)
const HEPRUP * getHEPRUP() const
Definition: LHEEvent.h:44
std::auto_ptr< PDF > pdf
Definition: LHEEvent.h:91
std::vector< FiveVector > PUP
Definition: LesHouches.h:248
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< int > IDUP
Definition: LesHouches.h:225
double SCALUP
Definition: LesHouches.h:210
const HepMC::GenVertex * lhef::LHEEvent::findSignalVertex ( const HepMC::GenEvent event,
bool  status3 = true 
)
static

Definition at line 448 of file LHEEvent.cc.

Referenced by asHepMCEvent(), and attempted().

450 {
451  double largestMass2 = -9.0e-30;
452  const HepMC::GenVertex *vertex = nullptr;
453  for(HepMC::GenEvent::vertex_const_iterator iter = event->vertices_begin();
454  iter != event->vertices_end(); ++iter) {
455  if ((*iter)->particles_in_size() < 2)
456  continue;
457  if ((*iter)->particles_out_size() < 1 ||
458  ((*iter)->particles_out_size() == 1 &&
459  (!(*(*iter)->particles_out_const_begin())->end_vertex() ||
460  !(*(*iter)->particles_out_const_begin())
461  ->end_vertex()->particles_out_size())))
462  continue;
463 
464  double px = 0.0, py = 0.0, pz = 0.0, E = 0.0;
465  bool hadStatus3 = false;
466  for(HepMC::GenVertex::particles_out_const_iterator iter2 =
467  (*iter)->particles_out_const_begin();
468  iter2 != (*iter)->particles_out_const_end(); ++iter2) {
469  hadStatus3 = hadStatus3 || (*iter2)->status() == 3;
470  px += (*iter2)->momentum().px();
471  py += (*iter2)->momentum().py();
472  pz += (*iter2)->momentum().pz();
473  E += (*iter2)->momentum().e();
474  }
475  if (status3 && !hadStatus3)
476  continue;
477 
478  double mass2 = E * E - (px * px + py * py + pz * pz);
479  if (mass2 > largestMass2) {
480  vertex = *iter;
481  largestMass2 = mass2;
482  }
483  }
484 
485  return vertex;
486 }
void lhef::LHEEvent::fixHepMCEventTimeOrdering ( HepMC::GenEvent event)
static

Definition at line 513 of file LHEEvent.cc.

References lhef::fixSubTree(), and class-composition::visited.

Referenced by attempted(), gen::PomwigHadronizer::finalizeEvent(), and Herwig6Hadronizer::finalizeEvent().

514 {
515  std::set<const HepMC::GenVertex*> visited;
516  HepMC::FourVector zeroTime;
517  for(HepMC::GenEvent::vertex_iterator iter = event->vertices_begin();
518  iter != event->vertices_end(); ++iter)
519  fixSubTree(*iter, zeroTime, visited);
520 }
static void fixSubTree(HepMC::GenVertex *vertex, HepMC::FourVector &_time, std::set< const HepMC::GenVertex * > &visited)
Definition: LHEEvent.cc:488
VertexRefVector::iterator vertex_iterator
iterator over a vector of references to Vertex objects in the same collection
Definition: VertexFwd.h:19
const std::vector<std::string>& lhef::LHEEvent::getComments ( ) const
inline

Definition at line 46 of file LHEEvent.h.

References comments.

Referenced by Herwig6Hadronizer::upEvnt().

46 { return comments; }
std::vector< std::string > comments
Definition: LHEEvent.h:93
const HEPEUP* lhef::LHEEvent::getHEPEUP ( ) const
inline
const HEPRUP* lhef::LHEEvent::getHEPRUP ( ) const
inline

Definition at line 44 of file LHEEvent.h.

Referenced by asHepMCEvent(), and fillPdfInfo().

44 { return runInfo->getHEPRUP(); }
const boost::shared_ptr< LHERunInfo > runInfo
Definition: LHEEvent.h:88
const PDF* lhef::LHEEvent::getPDF ( ) const
inline

Definition at line 45 of file LHEEvent.h.

45 { return pdf.get(); }
std::auto_ptr< PDF > pdf
Definition: LHEEvent.h:91
const int lhef::LHEEvent::getReadAttempts ( )
inline

Definition at line 47 of file LHEEvent.h.

References readAttemptCounter.

47 { return readAttemptCounter; }
int readAttemptCounter
Definition: LHEEvent.h:95
const boost::shared_ptr<LHERunInfo>& lhef::LHEEvent::getRunInfo ( ) const
inline

Definition at line 42 of file LHEEvent.h.

References runInfo.

42 { return runInfo; }
const boost::shared_ptr< LHERunInfo > runInfo
Definition: LHEEvent.h:88
HepMC::GenParticle * lhef::LHEEvent::makeHepMCParticle ( unsigned int  i) const
private

Definition at line 387 of file LHEEvent.cc.

References GenParticle::GenParticle, hepeup, lhef::HEPEUP::IDUP, lhef::HEPEUP::ISTUP, lhef::HEPEUP::PUP, and mps_update::status.

Referenced by asHepMCEvent(), and attempted().

388 {
389  HepMC::GenParticle *particle = new HepMC::GenParticle(
390  HepMC::FourVector(hepeup.PUP.at(i)[0],
391  hepeup.PUP.at(i)[1],
392  hepeup.PUP.at(i)[2],
393  hepeup.PUP.at(i)[3]),
394  hepeup.IDUP.at(i));
395 
396  int status = hepeup.ISTUP.at(i);
397 
398  particle->set_generated_mass(hepeup.PUP.at(i)[4]);
399  particle->set_status(status > 0 ? (status == 2 ? 3 : status) : 3);
400 
401  return particle;
402 }
HEPEUP hepeup
Definition: LHEEvent.h:90
std::vector< FiveVector > PUP
Definition: LesHouches.h:248
std::vector< int > ISTUP
Definition: LesHouches.h:230
std::vector< int > IDUP
Definition: LesHouches.h:225
int lhef::LHEEvent::npLO ( ) const
inline

Definition at line 58 of file LHEEvent.h.

References npLO_.

58 { return npLO_; }
int lhef::LHEEvent::npNLO ( ) const
inline

Definition at line 59 of file LHEEvent.h.

References npNLO_.

59 { return npNLO_; }
double lhef::LHEEvent::originalXWGTUP ( ) const
inline

Definition at line 52 of file LHEEvent.h.

References originalXWGTUP_.

52 { return originalXWGTUP_; }
double originalXWGTUP_
Definition: LHEEvent.h:96
void lhef::LHEEvent::removeParticle ( lhef::HEPEUP hepeup,
int  index 
)
static

Definition at line 156 of file LHEEvent.cc.

References mps_fire::i, lhef::HEPEUP::ICOLUP, lhef::HEPEUP::IDUP, lhef::HEPEUP::ISTUP, lhef::HEPEUP::MOTHUP, lhef::HEPEUP::NUP, lhef::pop(), lhef::HEPEUP::PUP, lhef::HEPEUP::resize(), lhef::HEPEUP::SPINUP, and lhef::HEPEUP::VTIMUP.

Referenced by addComment(), and removeResonances().

157 {
158  index--;
159  if (index < 0 || index >= hepeup.NUP) {
160  edm::LogError("Generator|LHEInterface")
161  << "removeParticle: Index " << (index + 1)
162  << " out of bounds." << std::endl;
163  return;
164  }
165 
166  std::pair<int, int> mo = hepeup.MOTHUP[index];
167 
168  pop(hepeup.IDUP, index);
169  pop(hepeup.ISTUP, index);
170  pop(hepeup.MOTHUP, index);
171  pop(hepeup.ICOLUP, index);
172  pop(hepeup.PUP, index);
173  pop(hepeup.VTIMUP, index);
174  pop(hepeup.SPINUP, index);
175  hepeup.NUP--;
176  hepeup.resize();
177 
178  index++;
179  for(int i = 0; i < hepeup.NUP; i++) {
180  if (hepeup.MOTHUP[i].first == index) {
181  if (hepeup.MOTHUP[i].second > 0)
182  edm::LogError("Generator|LHEInterface")
183  << "removeParticle: Particle "
184  << (i + 2) << " has two mothers."
185  << std::endl;
186  hepeup.MOTHUP[i] = mo;
187  }
188 
189  if (hepeup.MOTHUP[i].first > index)
190  hepeup.MOTHUP[i].first--;
191  if (hepeup.MOTHUP[i].second > index)
192  hepeup.MOTHUP[i].second--;
193  }
194 }
std::vector< double > VTIMUP
Definition: LesHouches.h:254
static void pop(std::vector< T > &vec, unsigned int index)
Definition: LHEEvent.cc:150
std::vector< std::pair< int, int > > MOTHUP
Definition: LesHouches.h:236
void resize(int nup)
Definition: LesHouches.h:161
std::vector< FiveVector > PUP
Definition: LesHouches.h:248
std::vector< double > SPINUP
Definition: LesHouches.h:261
std::vector< int > ISTUP
Definition: LesHouches.h:230
std::vector< int > IDUP
Definition: LesHouches.h:225
std::vector< std::pair< int, int > > ICOLUP
Definition: LesHouches.h:242
void lhef::LHEEvent::removeResonances ( const std::vector< int > &  ids)

Definition at line 196 of file LHEEvent.cc.

References funct::abs(), spr::find(), hepeup, mps_fire::i, triggerObjects_cff::id, lhef::HEPEUP::IDUP, lhef::HEPEUP::MOTHUP, lhef::HEPEUP::NUP, and removeParticle().

Referenced by addComment().

197 {
198  for(int i = 1; i <= hepeup.NUP; i++) {
199  int id = std::abs(hepeup.IDUP[i - 1]);
200  if (hepeup.MOTHUP[i - 1].first > 0 &&
201  hepeup.MOTHUP[i - 1].second > 0 &&
202  std::find(ids.begin(), ids.end(), id) != ids.end())
203  removeParticle(hepeup, i--);
204  }
205 }
HEPEUP hepeup
Definition: LHEEvent.h:90
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
std::vector< std::pair< int, int > > MOTHUP
Definition: LesHouches.h:236
static void removeParticle(lhef::HEPEUP &hepeup, int index)
Definition: LHEEvent.cc:156
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< int > IDUP
Definition: LesHouches.h:225
const std::vector<float>& lhef::LHEEvent::scales ( ) const
inline

Definition at line 55 of file LHEEvent.h.

References scales_.

Referenced by setScales().

55 { return scales_; }
std::vector< float > scales_
Definition: LHEEvent.h:97
void lhef::LHEEvent::setNpLO ( int  n)
inline

Definition at line 61 of file LHEEvent.h.

References gen::n, and npLO_.

61 { npLO_ = n; }
void lhef::LHEEvent::setNpNLO ( int  n)
inline

Definition at line 62 of file LHEEvent.h.

References gen::n, and npNLO_.

62 { npNLO_ = n; }
void lhef::LHEEvent::setPDF ( std::auto_ptr< PDF pdf)
inline

Definition at line 50 of file LHEEvent.h.

References pdf.

50 { this->pdf = pdf; }
std::auto_ptr< PDF > pdf
Definition: LHEEvent.h:91
void lhef::LHEEvent::setScales ( const std::vector< float > &  scales)
inline

Definition at line 56 of file LHEEvent.h.

References scales(), and scales_.

56 { scales_ = scales; }
const std::vector< float > & scales() const
Definition: LHEEvent.h:55
std::vector< float > scales_
Definition: LHEEvent.h:97
const std::vector<WGT>& lhef::LHEEvent::weights ( ) const
inline

Definition at line 53 of file LHEEvent.h.

References weights_.

53 { return weights_; }
std::vector< WGT > weights_
Definition: LHEEvent.h:92

Member Data Documentation

std::vector<std::string> lhef::LHEEvent::comments
private

Definition at line 93 of file LHEEvent.h.

Referenced by getComments(), and LHEEvent().

bool lhef::LHEEvent::counted
private

Definition at line 94 of file LHEEvent.h.

Referenced by count().

HEPEUP lhef::LHEEvent::hepeup
private
int lhef::LHEEvent::npLO_
private

Definition at line 98 of file LHEEvent.h.

Referenced by npLO(), and setNpLO().

int lhef::LHEEvent::npNLO_
private

Definition at line 99 of file LHEEvent.h.

Referenced by npNLO(), and setNpNLO().

double lhef::LHEEvent::originalXWGTUP_
private

Definition at line 96 of file LHEEvent.h.

Referenced by LHEEvent(), and originalXWGTUP().

std::auto_ptr<PDF> lhef::LHEEvent::pdf
private

Definition at line 91 of file LHEEvent.h.

Referenced by fillPdfInfo(), LHEEvent(), and setPDF().

int lhef::LHEEvent::readAttemptCounter
private

Definition at line 95 of file LHEEvent.h.

Referenced by attempted(), and getReadAttempts().

const boost::shared_ptr<LHERunInfo> lhef::LHEEvent::runInfo
private

Definition at line 88 of file LHEEvent.h.

Referenced by count(), and getRunInfo().

std::vector<float> lhef::LHEEvent::scales_
private

Definition at line 97 of file LHEEvent.h.

Referenced by scales(), and setScales().

std::vector<WGT> lhef::LHEEvent::weights_
private

Definition at line 92 of file LHEEvent.h.

Referenced by addWeight(), and weights().