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::unique_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 std::shared_ptr< LHERunInfo > & getRunInfo () const
 
 LHEEvent (const std::shared_ptr< LHERunInfo > &runInfo, std::istream &in)
 
 LHEEvent (const std::shared_ptr< LHERunInfo > &runInfo, const HEPEUP &hepeup)
 
 LHEEvent (const std::shared_ptr< LHERunInfo > &runInfo, const HEPEUP &hepeup, const LHEEventProduct::PDF *pdf, const std::vector< std::string > &comments)
 
 LHEEvent (const std::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::unique_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::unique_ptr< PDFpdf
 
int readAttemptCounter
 
const std::shared_ptr< LHERunInforunInfo
 
std::vector< float > scales_
 
std::vector< WGTweights_
 

Detailed Description

Definition at line 23 of file LHEEvent.h.

Member Typedef Documentation

Definition at line 37 of file LHEEvent.h.

Definition at line 38 of file LHEEvent.h.

Constructor & Destructor Documentation

lhef::LHEEvent::LHEEvent ( const std::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, mps_splice::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:88
std::vector< std::string > comments
Definition: LHEEvent.h:91
const std::shared_ptr< LHERunInfo > runInfo
Definition: LHEEvent.h:86
std::vector< double > VTIMUP
Definition: LesHouches.h:267
std::pair< double, double > XPDWUP
Definition: LesHouches.h:217
int readAttemptCounter
Definition: LHEEvent.h:93
std::vector< std::pair< int, int > > MOTHUP
Definition: LesHouches.h:249
void resize(int nup)
Definition: LesHouches.h:174
std::vector< FiveVector > PUP
Definition: LesHouches.h:261
std::vector< double > SPINUP
Definition: LesHouches.h:274
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< int > ISTUP
Definition: LesHouches.h:243
LHEEventProduct::PDF PDF
Definition: LHEEvent.h:37
std::vector< int > IDUP
Definition: LesHouches.h:238
double AQCDUP
Definition: LesHouches.h:233
std::unique_ptr< PDF > pdf
Definition: LHEEvent.h:89
bool counted
Definition: LHEEvent.h:92
double AQEDUP
Definition: LesHouches.h:228
static int skipWhitespace(std::istream &in)
Definition: LHEEvent.cc:23
double XWGTUP
Definition: LesHouches.h:209
std::vector< std::pair< int, int > > ICOLUP
Definition: LesHouches.h:255
std::vector< WGT > weights_
Definition: LHEEvent.h:90
double originalXWGTUP_
Definition: LHEEvent.h:94
double SCALUP
Definition: LesHouches.h:223
lhef::LHEEvent::LHEEvent ( const std::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:88
const std::shared_ptr< LHERunInfo > runInfo
Definition: LHEEvent.h:86
int readAttemptCounter
Definition: LHEEvent.h:93
bool counted
Definition: LHEEvent.h:92
lhef::LHEEvent::LHEEvent ( const std::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:88
std::vector< std::string > comments
Definition: LHEEvent.h:91
const std::shared_ptr< LHERunInfo > runInfo
Definition: LHEEvent.h:86
int readAttemptCounter
Definition: LHEEvent.h:93
LHEEventProduct::PDF PDF
Definition: LHEEvent.h:37
std::unique_ptr< PDF > pdf
Definition: LHEEvent.h:89
bool counted
Definition: LHEEvent.h:92
lhef::LHEEvent::LHEEvent ( const std::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:88
const lhef::HEPEUP & hepeup() const
const PDF * pdf() const
std::vector< std::string > comments
Definition: LHEEvent.h:91
const std::shared_ptr< LHERunInfo > runInfo
Definition: LHEEvent.h:86
const std::vector< WGT > & weights() const
int readAttemptCounter
Definition: LHEEvent.h:93
comments_const_iterator comments_begin() const
LHEEventProduct::PDF PDF
Definition: LHEEvent.h:37
int npLO() const
const std::vector< float > & scales() const
comments_const_iterator comments_end() const
std::vector< float > scales_
Definition: LHEEvent.h:95
std::unique_ptr< PDF > pdf
Definition: LHEEvent.h:89
bool counted
Definition: LHEEvent.h:92
int npNLO() const
std::vector< WGT > weights_
Definition: LHEEvent.h:90
double originalXWGTUP_
Definition: LHEEvent.h:94
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 62 of file LHEEvent.h.

References count(), removeParticle(), and removeResonances().

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

Definition at line 47 of file LHEEvent.h.

References weights_.

47 { weights_.push_back(wgt); }
std::vector< WGT > weights_
Definition: LHEEvent.h:90
std::unique_ptr< HepMC::GenEvent > lhef::LHEEvent::asHepMCEvent ( ) const

Definition at line 250 of file LHEEvent.cc.

References lhef::HEPEUP::AQCDUP, lhef::HEPEUP::AQEDUP, checkHepMCTree(), pdg::cTau(), 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().

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

Definition at line 397 of file LHEEvent.cc.

References mps_update::status.

Referenced by asHepMCEvent(), and attempted().

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

Definition at line 200 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().

202 {
203  if (counted)
204  edm::LogWarning("Generator|LHEInterface")
205  << "LHEEvent::count called twice on same event!"
206  << std::endl;
207 
209  weight, matchWeight);
210 
211  counted = true;
212 }
HEPEUP hepeup
Definition: LHEEvent.h:88
const std::shared_ptr< LHERunInfo > runInfo
Definition: LHEEvent.h:86
Definition: weight.py:1
bool counted
Definition: LHEEvent.h:92
double XWGTUP
Definition: LesHouches.h:209
void lhef::LHEEvent::fillEventInfo ( HepMC::GenEvent hepmc) const

Definition at line 242 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().

243 {
244  event->set_signal_process_id(hepeup.IDPRUP);
245  event->set_event_scale(hepeup.SCALUP);
246  event->set_alphaQED(hepeup.AQEDUP);
247  event->set_alphaQCD(hepeup.AQCDUP);
248 }
HEPEUP hepeup
Definition: LHEEvent.h:88
double AQCDUP
Definition: LesHouches.h:233
double AQEDUP
Definition: LesHouches.h:228
double SCALUP
Definition: LesHouches.h:223
void lhef::LHEEvent::fillPdfInfo ( HepMC::PdfInfo *  info) const

Definition at line 214 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().

215 {
216  if (pdf.get()) {
217  info->set_id1(pdf->id.first);
218  info->set_id2(pdf->id.second);
219  info->set_x1(pdf->x.first);
220  info->set_x2(pdf->x.second);
221  info->set_pdf1(pdf->xPDF.first);
222  info->set_pdf2(pdf->xPDF.second);
223  info->set_scalePDF(pdf->scalePDF);
224  } else if (hepeup.NUP >= 2) {
225  const HEPRUP *heprup = getHEPRUP();
226  info->set_id1(hepeup.IDUP[0] == 21 ? 0 : hepeup.IDUP[0]);
227  info->set_id2(hepeup.IDUP[1] == 21 ? 0 : hepeup.IDUP[1]);
228  info->set_x1(std::abs(hepeup.PUP[0][2] / heprup->EBMUP.first));
229  info->set_x2(std::abs(hepeup.PUP[1][2] / heprup->EBMUP.second));
230  info->set_pdf1(-1.0);
231  info->set_pdf2(-1.0);
232  info->set_scalePDF(hepeup.SCALUP);
233  } else {
234  info->set_x1(-1.0);
235  info->set_x2(-1.0);
236  info->set_pdf1(-1.0);
237  info->set_pdf2(-1.0);
238  info->set_scalePDF(hepeup.SCALUP);
239  }
240 }
HEPEUP hepeup
Definition: LHEEvent.h:88
static const TGPicture * info(bool iBackgroundIsBlack)
const HEPRUP * getHEPRUP() const
Definition: LHEEvent.h:42
std::vector< FiveVector > PUP
Definition: LesHouches.h:261
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< int > IDUP
Definition: LesHouches.h:238
std::unique_ptr< PDF > pdf
Definition: LHEEvent.h:89
double SCALUP
Definition: LesHouches.h:223
const HepMC::GenVertex * lhef::LHEEvent::findSignalVertex ( const HepMC::GenEvent event,
bool  status3 = true 
)
static

Definition at line 441 of file LHEEvent.cc.

Referenced by asHepMCEvent(), and attempted().

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

Definition at line 506 of file LHEEvent.cc.

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

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

507 {
508  std::set<const HepMC::GenVertex*> visited;
509  HepMC::FourVector zeroTime;
510  for(HepMC::GenEvent::vertex_iterator iter = event->vertices_begin();
511  iter != event->vertices_end(); ++iter)
512  fixSubTree(*iter, zeroTime, visited);
513 }
static void fixSubTree(HepMC::GenVertex *vertex, HepMC::FourVector &_time, std::set< const HepMC::GenVertex * > &visited)
Definition: LHEEvent.cc:481
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 44 of file LHEEvent.h.

References comments.

Referenced by Herwig6Hadronizer::upEvnt().

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

Definition at line 42 of file LHEEvent.h.

Referenced by asHepMCEvent(), and fillPdfInfo().

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

Definition at line 43 of file LHEEvent.h.

43 { return pdf.get(); }
std::unique_ptr< PDF > pdf
Definition: LHEEvent.h:89
const int lhef::LHEEvent::getReadAttempts ( )
inline

Definition at line 45 of file LHEEvent.h.

References readAttemptCounter.

45 { return readAttemptCounter; }
int readAttemptCounter
Definition: LHEEvent.h:93
const std::shared_ptr<LHERunInfo>& lhef::LHEEvent::getRunInfo ( ) const
inline

Definition at line 40 of file LHEEvent.h.

References runInfo.

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

Definition at line 380 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().

381 {
382  HepMC::GenParticle *particle = new HepMC::GenParticle(
383  HepMC::FourVector(hepeup.PUP.at(i)[0],
384  hepeup.PUP.at(i)[1],
385  hepeup.PUP.at(i)[2],
386  hepeup.PUP.at(i)[3]),
387  hepeup.IDUP.at(i));
388 
389  int status = hepeup.ISTUP.at(i);
390 
391  particle->set_generated_mass(hepeup.PUP.at(i)[4]);
392  particle->set_status(status > 0 ? (status == 2 ? 3 : status) : 3);
393 
394  return particle;
395 }
HEPEUP hepeup
Definition: LHEEvent.h:88
std::vector< FiveVector > PUP
Definition: LesHouches.h:261
std::vector< int > ISTUP
Definition: LesHouches.h:243
std::vector< int > IDUP
Definition: LesHouches.h:238
int lhef::LHEEvent::npLO ( ) const
inline

Definition at line 56 of file LHEEvent.h.

References npLO_.

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

Definition at line 57 of file LHEEvent.h.

References npNLO_.

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

Definition at line 50 of file LHEEvent.h.

References originalXWGTUP_.

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

Definition at line 149 of file LHEEvent.cc.

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

Referenced by addComment(), and removeResonances().

150 {
151  index--;
152  if (index < 0 || index >= hepeup.NUP) {
153  edm::LogError("Generator|LHEInterface")
154  << "removeParticle: Index " << (index + 1)
155  << " out of bounds." << std::endl;
156  return;
157  }
158 
159  std::pair<int, int> mo = hepeup.MOTHUP[index];
160 
161  hepeup.IDUP.erase(hepeup.IDUP.begin() + index);
162  hepeup.ISTUP.erase(hepeup.ISTUP.begin() + index);
163  hepeup.MOTHUP.erase(hepeup.MOTHUP.begin() + index);
164  hepeup.ICOLUP.erase(hepeup.ICOLUP.begin() + index);
165  hepeup.PUP.erase(hepeup.PUP.begin() + index);
166  hepeup.VTIMUP.erase(hepeup.VTIMUP.begin() + index);
167  hepeup.SPINUP.erase(hepeup.SPINUP.begin() + index);
168 
169  hepeup.NUP--;
170 
171  index++;
172  for(int i = 0; i < hepeup.NUP; i++) {
173  if (hepeup.MOTHUP[i].first == index) {
174  if (hepeup.MOTHUP[i].second > 0)
175  edm::LogError("Generator|LHEInterface")
176  << "removeParticle: Particle "
177  << (i + 2) << " has two mothers."
178  << std::endl;
179  hepeup.MOTHUP[i] = mo;
180  }
181 
182  if (hepeup.MOTHUP[i].first > index)
183  hepeup.MOTHUP[i].first--;
184  if (hepeup.MOTHUP[i].second > index)
185  hepeup.MOTHUP[i].second--;
186  }
187 }
std::vector< double > VTIMUP
Definition: LesHouches.h:267
std::vector< std::pair< int, int > > MOTHUP
Definition: LesHouches.h:249
std::vector< FiveVector > PUP
Definition: LesHouches.h:261
std::vector< double > SPINUP
Definition: LesHouches.h:274
std::vector< int > ISTUP
Definition: LesHouches.h:243
std::vector< int > IDUP
Definition: LesHouches.h:238
std::vector< std::pair< int, int > > ICOLUP
Definition: LesHouches.h:255
void lhef::LHEEvent::removeResonances ( const std::vector< int > &  ids)

Definition at line 189 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().

190 {
191  for(int i = 1; i <= hepeup.NUP; i++) {
192  int id = std::abs(hepeup.IDUP[i - 1]);
193  if (hepeup.MOTHUP[i - 1].first > 0 &&
194  hepeup.MOTHUP[i - 1].second > 0 &&
195  std::find(ids.begin(), ids.end(), id) != ids.end())
196  removeParticle(hepeup, i--);
197  }
198 }
HEPEUP hepeup
Definition: LHEEvent.h:88
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:249
static void removeParticle(lhef::HEPEUP &hepeup, int index)
Definition: LHEEvent.cc:149
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< int > IDUP
Definition: LesHouches.h:238
const std::vector<float>& lhef::LHEEvent::scales ( ) const
inline

Definition at line 53 of file LHEEvent.h.

References scales_.

Referenced by setScales().

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

Definition at line 59 of file LHEEvent.h.

References gen::n, and npLO_.

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

Definition at line 60 of file LHEEvent.h.

References gen::n, and npNLO_.

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

Definition at line 48 of file LHEEvent.h.

References eostools::move().

48 { this->pdf = std::move(pdf); }
std::unique_ptr< PDF > pdf
Definition: LHEEvent.h:89
def move(src, dest)
Definition: eostools.py:511
void lhef::LHEEvent::setScales ( const std::vector< float > &  scales)
inline

Definition at line 54 of file LHEEvent.h.

References scales(), and scales_.

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

Definition at line 51 of file LHEEvent.h.

References weights_.

51 { return weights_; }
std::vector< WGT > weights_
Definition: LHEEvent.h:90

Member Data Documentation

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

Definition at line 91 of file LHEEvent.h.

Referenced by getComments(), and LHEEvent().

bool lhef::LHEEvent::counted
private

Definition at line 92 of file LHEEvent.h.

Referenced by count().

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

Definition at line 96 of file LHEEvent.h.

Referenced by npLO(), and setNpLO().

int lhef::LHEEvent::npNLO_
private

Definition at line 97 of file LHEEvent.h.

Referenced by npNLO(), and setNpNLO().

double lhef::LHEEvent::originalXWGTUP_
private

Definition at line 94 of file LHEEvent.h.

Referenced by LHEEvent(), and originalXWGTUP().

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

Definition at line 89 of file LHEEvent.h.

Referenced by fillPdfInfo(), and LHEEvent().

int lhef::LHEEvent::readAttemptCounter
private

Definition at line 93 of file LHEEvent.h.

Referenced by attempted(), and getReadAttempts().

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

Definition at line 86 of file LHEEvent.h.

Referenced by count(), and getRunInfo().

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

Definition at line 95 of file LHEEvent.h.

Referenced by scales(), and setScales().

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

Definition at line 90 of file LHEEvent.h.

Referenced by addWeight(), and weights().