CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
HepMC::Pythia8ToHepMCA Class Reference

#include <HepMCA2.h>

Inheritance diagram for HepMC::Pythia8ToHepMCA:

Public Member Functions

bool append_event (Pythia8::Event &pyev, GenEvent *evt, GenParticle *rootpart, int ibarcode=-1, Pythia8::Info *pyinfo=0, Pythia8::Settings *pyset=0)
 
bool convert_gluon_to_0 () const
 
bool crash_on_problem () const
 
bool free_parton_warnings () const
 
bool print_inconsistency () const
 
 Pythia8ToHepMCA ()
 
void set_convert_gluon_to_0 (bool b=false)
 
void set_crash_on_problem (bool b=false)
 
void set_free_parton_warnings (bool b=true)
 
void set_print_inconsistency (bool b=true)
 
virtual ~Pythia8ToHepMCA ()
 

Private Member Functions

virtual bool fill_next_event (GenEvent *)
 
 Pythia8ToHepMCA (const Pythia8ToHepMCA &)
 
virtual void write_event (const GenEvent *)
 

Private Attributes

bool m_convert_gluon_to_0
 
bool m_crash_on_problem
 
bool m_free_parton_warnings
 
int m_internal_event_number
 
bool m_print_inconsistency
 

Detailed Description

Definition at line 26 of file HepMCA2.h.

Constructor & Destructor Documentation

HepMC::Pythia8ToHepMCA::Pythia8ToHepMCA ( )
inline

Definition at line 31 of file HepMCA2.h.

virtual HepMC::Pythia8ToHepMCA::~Pythia8ToHepMCA ( )
inlinevirtual

Definition at line 34 of file HepMCA2.h.

34 {;}
HepMC::Pythia8ToHepMCA::Pythia8ToHepMCA ( const Pythia8ToHepMCA )
inlineprivate

Definition at line 59 of file HepMCA2.h.

59 : IO_BaseClass() {;}

Member Function Documentation

bool HepMC::Pythia8ToHepMCA::append_event ( Pythia8::Event &  pyev,
GenEvent *  evt,
GenParticle *  rootpart,
int  ibarcode = -1,
Pythia8::Info *  pyinfo = 0,
Pythia8::Settings *  pyset = 0 
)
inline

Definition at line 74 of file HepMCA2.h.

References funct::abs(), dtNoiseDBValidation_cfg::cerr, cuy::col, alignCSCRings::e, cmsRelvalreport::exit, GenParticle::GenParticle, i, contentValuesFiles::m, m_crash_on_problem, m_free_parton_warnings, m_print_inconsistency, HWWFunctions::MM, and AlCaHLTBitMon_ParallelJobs::p.

Referenced by gen::Py8GunBase::residualDecay(), and Pythia8Hadronizer::residualDecay().

75  {
76 
77  // 1. Error if no event passed.
78  if (!evt) {
79  std::cerr << "Pythia8ToHepMCA::fill_next_event error - passed null event."
80  << std::endl;
81  return 0;
82  }
83 
84  // Conversion factors from Pythia units GeV and mm to HepMC ones.
85  double momFac = HepMC::Units::conversion_factor(HepMC::Units::GEV,
86  evt->momentum_unit());
87  double lenFac = HepMC::Units::conversion_factor(HepMC::Units::MM,
88  evt->length_unit());
89 
90  int NewBarcode = evt->particles_size();
91  if (ibarcode > -1) NewBarcode = ibarcode;
92 
93  GenVertex* prod_vtx0 = new GenVertex();
94  prod_vtx0->add_particle_in( rootpart );
95  evt->add_vertex( prod_vtx0 );
96 
97  // 2. Create a particle instance for each entry and fill a map, and
98  // a vector which maps from the particle index to the GenParticle address.
99  std::vector<GenParticle*> hepevt_particles( pyev.size() );
100  for (int i = 2; i < pyev.size(); ++i) {
101 
102  // Fill the particle.
103 
104  hepevt_particles[i] = new GenParticle(
105  FourVector( momFac * pyev[i].px(), momFac * pyev[i].py(),
106  momFac * pyev[i].pz(), momFac * pyev[i].e() ),
107  pyev[i].id(), pyev[i].statusHepMC() );
108  if (ibarcode !=0) NewBarcode++;
109  hepevt_particles[i]->suggest_barcode(NewBarcode);
110  hepevt_particles[i]->set_generated_mass( momFac * pyev[i].m() );
111 
112  // Colour flow uses index 1 and 2.
113  int colType = pyev[i].colType();
114  if (colType == 1 || colType == 2)
115  hepevt_particles[i]->set_flow(1, pyev[i].col());
116  if (colType == -1 || colType == 2)
117  hepevt_particles[i]->set_flow(2, pyev[i].acol());
118  }
119 
120  // 3. Loop over particles AGAIN, this time creating vertices.
121  // We build the production vertex for each entry in hepevt.
122  // The HEPEVT pointers are bi-directional, so gives decay vertices as well.
123  for (int i = 2; i < pyev.size(); ++i) {
124  GenParticle* p = hepevt_particles[i];
125 
126  // 3a. Search to see if a production vertex already exists.
127  std::vector<int> mothers = pyev[i].motherList();
128 
129  unsigned int imother = 0;
130  int mother = -1; // note that in Pythia8 there is a particle number 0!
131  if ( !mothers.empty() ) mother = mothers[imother];
132  GenVertex* prod_vtx = p->production_vertex();
133  while ( !prod_vtx && mother > 0 ) {
134  if(mother == 1) {
135  prod_vtx = rootpart->end_vertex();
136  } else {
137  prod_vtx = hepevt_particles[mother]->end_vertex();
138  }
139  if ( prod_vtx ) prod_vtx->add_particle_out( p );
140  mother = ( ++imother < mothers.size() ) ? mothers[imother] : -1;
141  }
142 
143  // 3b. If no suitable production vertex exists - and the particle has
144  // at least one mother or position information to store - make one.
145  FourVector prod_pos( lenFac * pyev[i].xProd(), lenFac * pyev[i].yProd(),
146  lenFac * pyev[i].zProd(), lenFac * pyev[i].tProd() );
147  if ( !prod_vtx && ( mothers.size() > 0 || prod_pos != FourVector() ) ) {
148  prod_vtx = new GenVertex();
149  prod_vtx->add_particle_out( p );
150  evt->add_vertex( prod_vtx );
151  }
152 
153  // 3c. If prod_vtx doesn't already have position specified, fill it.
154  if ( prod_vtx && prod_vtx->position() == FourVector() )
155  prod_vtx->set_position( prod_pos );
156 
157  // 3d. loop over mothers to make sure their end_vertices are consistent.
158  imother = 0;
159  mother = -1;
160  if ( !mothers.empty() ) mother = mothers[imother];
161  while ( prod_vtx && mother > 0 ) {
162 
163  // If end vertex of the mother isn't specified, do it now.
164  GenParticle* ppp;
165  if(mother == 1) ppp = rootpart;
166  if(mother > 1) ppp = hepevt_particles[mother];
167 
168  if ( !ppp->end_vertex() ) {
169  prod_vtx->add_particle_in( ppp );
170 
171  // Problem scenario: the mother already has a decay vertex which
172  // differs from the daughter's production vertex. This means there is
173  // internal inconsistency in the HEPEVT event record. Print an error.
174  // Note: we could provide a fix by joining the two vertices with a
175  // dummy particle if the problem arises often.
176  } else if (ppp->end_vertex() != prod_vtx ) {
178  << "HepMC::Pythia8ToHepMC: inconsistent mother/daugher "
179  << "information in Pythia8 event " << std::endl
180  << "i = " << i << " mother = " << mother
181  << "\n This warning can be turned off with the "
182  << "Pythia8ToHepMC::print_inconsistency switch." << std::endl;
183  }
184 
185  // End of vertex-setting loops.
186  mother = ( ++imother < mothers.size() ) ? mothers[imother] : -1;
187  }
188  }
189 
190  // If hadronization switched on then no final coloured particles.
191  bool doHadr = (pyset == 0) ? m_free_parton_warnings
192  : pyset->flag("HadronLevel:all") && pyset->flag("HadronLevel:Hadronize");
193 
194  // 4. Check for particles which come from nowhere, i.e. are without
195  // mothers or daughters. These need to be attached to a vertex, or else
196  // they will never become part of the event.
197  for (int i = 2; i < pyev.size(); ++i) {
198  if ( !hepevt_particles[i]->end_vertex() &&
199  !hepevt_particles[i]->production_vertex() ) {
200  std::cerr << "hanging particle " << i << std::endl;
201  GenVertex* prod_vtx = new GenVertex();
202  prod_vtx->add_particle_out( hepevt_particles[i] );
203  evt->add_vertex( prod_vtx );
204  }
205 
206  // Also check for free partons (= gluons and quarks; not diquarks?).
207  if ( doHadr && m_free_parton_warnings ) {
208  if ( hepevt_particles[i]->pdg_id() == 21 &&
209  !hepevt_particles[i]->end_vertex() ) {
210  std::cerr << "gluon without end vertex " << i << std::endl;
211  if ( m_crash_on_problem ) exit(1);
212  }
213  if ( abs(hepevt_particles[i]->pdg_id()) <= 6 &&
214  !hepevt_particles[i]->end_vertex() ) {
215  std::cerr << "quark without end vertex " << i << std::endl;
216  if ( m_crash_on_problem ) exit(1);
217  }
218  }
219  }
220 
221  // Done.
222  return true;
223 
224 }
int i
Definition: DBlmapReader.cc:9
bool m_free_parton_warnings
Definition: HepMCA2.h:63
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool m_print_inconsistency
Definition: HepMCA2.h:63
int col
Definition: cuy.py:1008
bool HepMC::Pythia8ToHepMCA::convert_gluon_to_0 ( ) const
inline

Definition at line 44 of file HepMCA2.h.

References m_convert_gluon_to_0.

44 {return m_convert_gluon_to_0;}
bool HepMC::Pythia8ToHepMCA::crash_on_problem ( ) const
inline

Definition at line 43 of file HepMCA2.h.

References m_crash_on_problem.

43 {return m_crash_on_problem;}
virtual bool HepMC::Pythia8ToHepMCA::fill_next_event ( GenEvent *  )
inlineprivatevirtual

Definition at line 55 of file HepMCA2.h.

55 { return 0; }
bool HepMC::Pythia8ToHepMCA::free_parton_warnings ( ) const
inline

Definition at line 42 of file HepMCA2.h.

References m_free_parton_warnings.

42 {return m_free_parton_warnings;}
bool m_free_parton_warnings
Definition: HepMCA2.h:63
bool HepMC::Pythia8ToHepMCA::print_inconsistency ( ) const
inline

Definition at line 41 of file HepMCA2.h.

References m_print_inconsistency.

41 {return m_print_inconsistency;}
bool m_print_inconsistency
Definition: HepMCA2.h:63
void HepMC::Pythia8ToHepMCA::set_convert_gluon_to_0 ( bool  b = false)
inline

Definition at line 50 of file HepMCA2.h.

References b, and m_convert_gluon_to_0.

double b
Definition: hdecay.h:120
void HepMC::Pythia8ToHepMCA::set_crash_on_problem ( bool  b = false)
inline

Definition at line 49 of file HepMCA2.h.

References b, and m_crash_on_problem.

double b
Definition: hdecay.h:120
void HepMC::Pythia8ToHepMCA::set_free_parton_warnings ( bool  b = true)
inline

Definition at line 48 of file HepMCA2.h.

References b, and m_free_parton_warnings.

bool m_free_parton_warnings
Definition: HepMCA2.h:63
double b
Definition: hdecay.h:120
void HepMC::Pythia8ToHepMCA::set_print_inconsistency ( bool  b = true)
inline

Definition at line 47 of file HepMCA2.h.

References b, and m_print_inconsistency.

bool m_print_inconsistency
Definition: HepMCA2.h:63
double b
Definition: hdecay.h:120
virtual void HepMC::Pythia8ToHepMCA::write_event ( const GenEvent *  )
inlineprivatevirtual

Definition at line 56 of file HepMCA2.h.

56 {;}

Member Data Documentation

bool HepMC::Pythia8ToHepMCA::m_convert_gluon_to_0
private

Definition at line 63 of file HepMCA2.h.

Referenced by convert_gluon_to_0(), and set_convert_gluon_to_0().

bool HepMC::Pythia8ToHepMCA::m_crash_on_problem
private

Definition at line 63 of file HepMCA2.h.

Referenced by append_event(), crash_on_problem(), and set_crash_on_problem().

bool HepMC::Pythia8ToHepMCA::m_free_parton_warnings
private

Definition at line 63 of file HepMCA2.h.

Referenced by append_event(), free_parton_warnings(), and set_free_parton_warnings().

int HepMC::Pythia8ToHepMCA::m_internal_event_number
private

Definition at line 62 of file HepMCA2.h.

bool HepMC::Pythia8ToHepMCA::m_print_inconsistency
private

Definition at line 63 of file HepMCA2.h.

Referenced by append_event(), print_inconsistency(), and set_print_inconsistency().