HepMC3 event record library
full.md
1 # One page HepMC3 documentation
2 
3 ## 1 Build instructions
4 For the full list of the installation options including the description of the flags to build the
5 HepMC3 from the sources see at [HepMC3 page at CERN GitLab](https://gitlab.cern.ch/hepmc/HepMC3).
6 
7 A quick minimalist build that requires only C++11 compiler and a
8 recent version of CMake (3.15+) is described below.
9 To build HepMC3 is using CMake on the supported platforms, C++11 compiler
10 and a recent CMake is needed (3.15).
11 The commands executed in the unpacked source tarball of HepMC3
12 
13 ```
14  cmake -DCMAKE_INSTALL_PREFIX=desired_installation_path -DHEPMC3_ENABLE_ROOTIO=OFF -DHEPMC3_ENABLE_PYTHON=OFF CMakeLists.txt
15  cmake --build .
16  cmake --install .
17 ```
18 will configure the HepMC3 sources, compile them and install the library into "desired_installation_path".
19 
20 
21 ## 2 Differences between HepMC2 and HepMC3
22 The following is a list of main differences that should be taken into
23 account when transitioning from HepMC2 to HepMC3.
24 
25 ### 2.1 Structure change and header file organisation
26 Following changes in header files have been applied:
27 
28 ```
29  HepMC/HeavyIons.h -> HepMC3/GenHeavyIons.h (now a POD struct)
30  HepMC/PdfInfo.h -> HepMC3/GenPdfInfo.h (now a POD struct)
31  HepMC/SimpleVector.h -> HepMC3/FourVector.h (ThreeVector class removed)
32 ```
33 
34 The structure of GenCrossSection class has been changed to handle multiple
35 values of cross-sections. The cross-section values and errors (uncertainties)
36  can be accessed only by public functions, the corresponding data members are
37  private. By default, the number of cross-section values in every event is equal
38  to the number of event weights. Accordingly, each cross-section value can be
39  accessed using the corresponding event weight name (std::string) or event
40  weight index (int).
41 
42 Following header files are no longer available:
43 
44 ```
45  CompareGenEvent.h
46  Flow.h
47  GenRanges.h
48  HepMCDefs.h
49  HerwigWrapper.h
50  IteratorRange.h
51  Polarization.h
52  SearchVector.h
53  StreamHelpers.h
54  StreamInfo.h
55  TempParticleMap.h
56  WeightContainer.h
57  IO_GenEvent.h
58 
59  enable_if.h
60  is_arithmetic.h
61 
62  IO_AsciiParticles.h
63  IO_Exception.h
64  IO_HEPEVT.h
65  IO_HERWIG.h
66 
67  PythiaWrapper.h
68  PythiaWrapper6_4.h
69  PythiaWrapper6_4_WIN32.h
70 ```
71 
72 ### 2.2 Fortran generators
73 An example of interface to Pythia6 Fortran blocks is given in the examples.
74 Please note that the provided interface Pythia6ToHepMC3.cc and Pythia6ToHepMC3.inc
75 is an interface for HepMC3 from Pythia6 and not an interface to Pythia6 from HepMC,
76 as it was in the case of the HepMC2.
77 
78 ### 2.3 Changes to the I/O handling
79 Multiple file formats are supported. The implementation of reading and writing
80 is separated in HepMC3. All the reading operations are performed in
81 "reader" objects inherited from HepMC::Reader and the writing operations in the "writer"
82 objects inherited from HePMC::Writer. Therefore it is to use the desired headers explicitly, as needed.
83 
84 The IO_GenEvent.h header is not available anymore.
85 to write and/or read HepMC2 files the following includes
86 
87 ```
88  WriterAsciiHepMC2.h
89  ReaderAsciiHepMC2.h
90 ```
91 
92 should be used instead of
93 
94 ```
95  IO_GenEvent.h
96 ```
97 Please note that HepMC2 format is outdated and is not able to contain a lot
98  of information stored into event record by the modern Monte Carlo event generators.
99  It is recommended to use HepMC3 native event record in plain text or in ROOT TTree format.
100  The corresponding readers and writers are
101 
102 ```
103  WriterAscii.h
104  ReaderAscii.h
105 ```
106 and
107 
108 ```
109  WriterRootTree.h
110  ReaderRootTree.h
111 ```
112 
113 Implementation of custom Reader and Writer objects is possible as well.
114 
115 Please, note the difference in the behaviour of default Readers with respect to HepMC2.
116  when reading files with multiple headers. The ASCII files with multiple headers ( e.g. obtained with
117  cat 1.hepmc 2.hepmc > 12.hepmc) will be processed by the readers only till the
118  first occurrence of END_EVENT_LISTING.
119 
120  In addition to the standard readers, starting for the version 3.2.5 HepMC3 provides as set of
121  templates readers/writers to handle the zip-,lzma-,bz2-compressed files (ReaderGZ and WriterGZ) and to perform multithread
122  reading (ReaderMT).
123 
124 ### 2.4 Memory managed by shared pointers
125 Particles and vertices are managed using shared pointers, so they should
126 not be created through the call to 'new'.
127 
128 ```
129  GenEvent event;
130 
131  // Create particles
132  GenParticlePtr p1 = make_shared<GenParticle>();
133  GenParticlePtr p2 = make_shared<GenParticle>();
134  GenParticlePtr p3 = make_shared<GenParticle>();
135 
136  p1->set_pdg_id(23);
137  p2->set_pdg_id(15);
138  p3->set_pdg_id(-15);
139 
140  // Create vertex
141  GenVertexPtr v1 = make_shared<GenVertex>();
142  v1->add_particle_in(p1);
143  v1->add_particle_out(p2);
144  v1->add_particle_out(p3);
145 
146  event.add_vertex(v1);
147 ```
148 
149 ### 2.5 Iterators
150 The iterator-bases classes and access functions from HepMC2, e.g.
151 
152 ```
153  class particle_iterator;
154  class vertex_iterator;
155  class edge_iterator;
156  ...
157  inline int GenEvent::particles_size() const;
158  inline int GenEvent::vertices_size() const;
159  ...
160 
161 ```
162 were removed. The C++11 iterations should be used instead, e.g. instead of
163 
164 ```
165  for (GenEvent::particle_const_iterator p = evt->particles_begin(); p != evt->particles_end(); ++p ) {
166  ...
167  }
168 ```
169  one should use
170 ```
171  for (const auto& p: evt->particles()) {
172  ...
173  }
174 ```
175 or alternatively
176 
177 ```
178  for (size_t i=0;i<evt->particles().size();++i) {
179  ...
180  evt->particles().at(i)
181  ...
182  }
183 ```
184 
185 
186 
187 ### 2.6 Topological order
188 Particles and vertices in HepMC3 are stored in topological order. This means
189  that when creating vertices, incoming particles must have id lower than
190  any of the outgoing particles.
191 
192 This forces the tree structure to be constructed top-to-bottom
193  and disallows creating loops.
194 
195 ```
196  GenParticlePtr p1 = make_shared<GenParticle>();
197  GenParticlePtr p2 = make_shared<GenParticle>();
198  GenParticlePtr p3 = make_shared<GenParticle>();
199  GenParticlePtr p4 = make_shared<GenParticle>();
200 
201  GenVertexPtr v1 = make_shared<GenVertex>();
202  GenVertexPtr v2 = make_shared<GenVertex>();
203  GenVertexPtr v3 = make_shared<GenVertex>();
204 
205  event.add_particle(p1);
206  event.add_particle(p2);
207  event.add_particle(p3);
208  event.add_particle(p4);
209  event.add_vertex(v1);
210  event.add_vertex(v2);
211  event.add_vertex(v3);
212 
213  v1->add_particle_in (p2);
214  v1->add_particle_out(p3);
215  v1->add_particle_in (p4); // will cause error, because p3
216  // has higher index than p4
217 
218  v2->add_particle_in (p4);
219  v2->add_particle_out(p3); // will also cause error
220 
221  // Order of vertices does not matter. Index of end vertex
222  // can be lower than index of production vertex
223  v3->add_particle_in (p1);
224  v3->add_particle_out(p2);
225 ```
226 
227 
228 
229 ### 2.7 Deleting particles and vertices
230 Deleting a particle using GenEvent::remove_particle() will also remove
231  its end_vertex if this is the only particle that is on this vertex
232  particles_in() list.
233 
234 Deleting a vertex will delete all of its outgoing
235  particles. (and subsequently, all of their decays).
236 
237 
238 ### 2.8 Barcodes can no longer be se (Use constant ID instead)
239 The "barcode" integer in HepMC2 was an uncomfortable object, simultaneously
240  declared in the code documentation to be a meaningless unique identifier for
241  vertex and particle objects, and set to specific ranges by experiments'
242  production systems to encode information about a particle's origins. It proved
243  impossible to satisfactorily reconcile these twin uses, and experiments' demands
244  for particle provenance information have exceeded the capacity of an int (or
245  even a long int).
246 
247 Hence, barcodes are no longer available. Use attributes to provide additional
248  information that was previously encoded using barcodes
249  (see module @ref attributes).
250 
251 The unique identifier of particles and vertices is now called id() to
252  separate its role from barcodes. Id is set automatically and cannot
253  be changed. Id is not permanently attached to particle/vertex. When
254  a particle or vertex is removed from the event, id's of other particles
255  or vertices may change.
256 
257 ### 2.9 Flow is not a class on its own (it is an attribute).
258 The Flow class has been removed, since it was unused by any widespread event
259  generator, and to our knowledge the only active use-case is an abuse of it to
260  provide more ints in which to encode provenance information. As this is now done
261  via attributes, there is no case for Flow's continued existence. No backward
262  compatibility Flow class is provided since this usage is extremely localised in
263  one piece of user code and migration to the newer scheme should be simple.
264 
265 
266 ### 2.10 Units are no longer defined at compilation time
267 The default units are set to GEV and MM. They can be provided as
268  constructor parameters or changed later using HepMC::GenEvent::set_units
269 
270 ```
271  GenEvent event(Units::GEV,Units::CM);
272 
273  GenParticlePtr p = make_shared<GenParticle>();
274 
275  event.add_particle(p);
276  ...
277 
278  event.print(); // event printed in GEV/CM
279 
280  event.set_units(Units::MEV,Units::MM); // will trigger unit conversion for all particles and vertices
281 
282  event.print(); // event printed in MEV/MM
283 ```
284 
285 
286 ### 2.11 Deprecated code
287 A lot of HepMC2 functions has been declared obsolete and are marked as
288  deprecated. Warnings displayed at compilation time hint to what functions
289  or classes should be used instead.
290 
291 
292 ```
293  // HepMC2 code:
294  HepMC::FourVector position(pos[1],pos[2],pos[3],pos[0]);
295  vertex = new HepMC::GenVertex(position, id);
296 
297  // Replace with:
298  HepMC3::FourVector position(pos[1],pos[2],pos[3],pos[0]);
299  HepMC3::GenVertexPtr vertex = std::make_shared<HepMC3::GenVertex>(position);
300  vertex->set_status(1);
301 ```
302 
303 ```
304  // HepMC2 code:
305  std::vector<HepMC::GenParticle*> beamparticles
306  // ...
307  event.set_beam_particles(beamparticles[0],beamparticles[1]);
308 
309  // Replace with:
310  std::vector<HepMC3::GenParticlePtr> beamparticles;
311  // ...
312  event.add_particle(beamparticles[0]);
313  event.add_particle(beamparticles[1]);
314 ```
315 
316 ```
317  // HepMC2 code:
318  HepMC::GenVertex * vertex;
319  vertex->set_id(1);
320  vertex->id();
321 
322  // Replace with:
323  HepMC3::GenVertexPtr vertex = std::make_shared<HepMC3::GenVertex>();
324  vertex->set_status(1);
325  vertex->status();
326 ```
327 
328 ```
329  // HepMC2 code:
330  HepMC::GenVertex * vertex;
331  for (HepMC::GenVertex::particles_out_const_iterator pout
332  =v->particles_out_const_begin();
333  pout!=(*vit)->particles_out_const_end(); ++pout) { }
334 
335  // Replace with (similarly for particles_in):
336  HepMC3::GenVertexPtr vertex = std::make_shared<HepMC3::GenVertex>();
337  for (HepMC3::GenParticlePtr pout : vertex->particles_out() ) { }
338 ```
339 
340 ```
341  // HepMC2 code:
342  vertex->weights().push_back(1.);
343  // Replace with:
344  TODO
345 ```
346 
347 ```
348  GenEvent evt(Units::GEV,Units::MM);
349  // HepMC2 code:
350  evt.set_alphaQCD(m_alphas);
351  evt.set_alphaQED(m_alpha);
352  evt.set_event_scale(m_event_scale);
353  evt.set_mpi(m_mpi);
354  evt.set_signal_process_id(m_signal_process_id);
355  evt.set_random_states(m_random_states);
356  // Replace with:
357  evt.add_attribute("alphaQCD",
358  make_shared<DoubleAttribute>(m_alphas));
359  evt.add_attribute("alphaEM",
360  make_shared<DoubleAttribute>(m_alpha));
361  evt.add_attribute("event_scale",
362  make_shared<DoubleAttribute>(m_event_scale));
363  evt.add_attribute("mpi",
364  make_shared<IntAttribute>(m_mpi));
365  evt.add_attribute("signal_process_id",
366  make_shared<IntAttribute>(m_signal_process_id));
367  for ( size_t i=0;i<m_random_states.size();i++)
368  evt.add_attribute("random_states"+to_string(i),make_shared<IntAttribute>(m_random_states[i]));
369 ```
370 
371 ```
372  // HepMC2 code:
373  HepMC::GenVertex * vertex;
374  ...
375  vertex->weights().push_back(1.);
376  // Replace with:
377  vertex->add_attribute("weight1", make_shared<DoubleAttribute>(1.0));
378 
379 ```
380 
381 ```
382  // HepMC2 code:
383  HepMC::GenVertex * vertex;
384  vertex->check_momentum_conservation();
385  // Replace with:
386  TODO
387 ```
388 
389 ```
390  // HepMC2 code:
391  HepMC::GenParticle::set_flow(int code_index, int code = 0)
392  HepMC::GenParticle::set_polarization( Polarization(theta,phi))
393  // Replace with:
394  HepMC3::GenParticle::add_attribute("flow"+to_string(code_index),make_shared<IntAttribute>(code));
395  HepMC3::GenParticle::add_attribute("theta",make_shared<DoubleAttribute>(theta));
396  HepMC3::GenParticle::add_attribute("phi",make_shared<DoubleAttribute>(phi));
397 ```
398 
399 ```
400  // HepMC2 code:
401  GenCrossSection* XS;
402  double xs, xs_err;
403  XS->set_cross_section(xs,xs_err)
404  // Replace with:
405  GenCrossSectionPtr XS;
406  ....
407  std::vector<double> xs, xs_err;
408  ....
409  XS->set_cross_section(xs,xs_err)
410 ```
411 
412 
413 ## 3 Standard attributes
414 For the user convenience and backward compatibility the following standard attributes are
415  supported for the
416 
417 GenEvent
418 ```
419  double alphaQCD
420  double alphaEM
421  double event_scale
422  int mpi
423  int signal_process_id
424  int signal_vertex_id
425  int random_states1... random_statesN
426 ```
427 
428 GenVertex
429 ```
430  double weight1... weightN
431 ```
432 
433 GenParticle
434 ```
435  int flow
436  double theta
437  double phi
438 ```
439 
440 The presence of cycles in the event structure is indicated with an attribute
441 ```
442  int cycles
443 ```
444 
445  Note that attributes belong to the event, therefore these can be set only for particles and vertices
446  that belong to a GenEvent object.
447 
448 
449 ## 4 Interface to HEPEVT block
450 The most recent versions of HepMC3 has multiple implementations of the interfaces to HEPEVT Fortran common
451  block. These are
452 
453 - include/HepMC3/HEPEVT_Wrapper.h -- the default implementation. The size of common block is defined in
454  compile time via appropriate @c define directive. The block can hold float/double precision momenta.
455  This implementation is not compiled into any library. All functions and variables are static,
456  so only one instance of the interface can exists.
457 
458 - include/HepMC3/HEPEVT_Wrapper_Runtime.h -- The size of common block is defined in runtime.
459  The block can be held in the object. Multiple instances can exists.
460  The interface is compiled into the library. This interface is also available in Python.
461 
462 - include/HepMC3/HEPEVT_Wrapper_Runtime_Static.h -- The size of common block is defined in runtime.
463  All functions and variables are static,
464  so only one instance of the interface can exists. The interface is compiled into the library.
465 
466 - include/HepMC3/HEPEVT_Wrapper_Template.h -- The size of common block is defined in compile
467  time as a parameter of template. The block can hold float/double precision momenta.
468  The block can be held in the object. Multiple instances can exists.
469 
470 
471 ## 5 GenRunInfo class
472 A new class has been provided to store run-level information, such as
473  weight names, names and description of tools used to generate the event,
474  global attributes such as LHE run information or any other run information
475  provided by user. See HepMC::GenRunInfo class description for details.
476 
477 
478 
479 ## 6 Attributes
480 Attributes can be attached to GenEvent, GenParticle or GenVertex
481  and they can have any format defined by the user
482  (see @ref writing_attributes). An attribute is accessed through
483  a shared pointer and identified by its name.
484 
485 Example of reading an attribute from the event:
486 
487 ```
488  shared_ptr<GenPdfInfo> pdf_info = event.attribute<GenPdfInfo>("GenPdfInfo");
489 
490  if( pdf_info ) pdf_info->print();
491 ```
492 
493 Example of adding an attribute to the event:
494 
495 ```
496  shared_ptr<GenPdfInfo> pdf_info = make_shared<GenPdfInfo>();
497  evt.add_attribute("GenPdfInfo",pdf_info);
498 
499  // Setting values can be done before or after adding it to the event
500  pdf_info->set(1,2,3.4,5.6,7.8,9.0,1.2,3,4);
501 ```
502 
503 Adding and getting attributes of a vertex or particle uses the same
504  principles.
505 
506 Note: An event (or particle or vertex) can have more than one attribute
507  of the same type distinguished by different names. This might be
508  useful in some applications, however, we strongly encourage
509  to use just one instance named by its class name, as in these
510  examples.
511 
512 ### 6.1 Writing custom attributes
513 
514 Any class that derives from HepMC::Attribute class can be used as
515  an attribute that can be attached to the event, vertex or particle.
516 
517 User has to provide two abstract methods from HepMC::Attribute used
518  to parse the class content from/to string.
519 
520 Example:
521 
522 ```
523  #include "HepMC3/Attribute.h"
524 
525  struct MyAttribute : public HepMC::Attribute {
526 
527  double val1; /// First value
528  int val2; /// Second value
529 
530  public:
531  /// Implementation of Attribute::from_string
532  bool from_string(const string &att) {
533  val1 = stof( att );
534  val2 = stol( att.substr( att.find(' ')+1 ) );
535 
536  return true;
537  }
538 
539  /// Implementation of Attribute::to_string
540  bool to_string(string &att) const {
541  char buf[64];
542 
543  sprintf(buf,"%.8e %i",val1,val2);
544 
545  att = buf;
546 
547  return true;
548  }
549  };
550 ```
551 
552  For other examples see attributes provided in the HepMC3 package.
553 
554 
555 
556 ## 7 IO-related classes and interfaces
557 The main HepMC3 library contains the classes for the I/O of multiple event formats.
558 
559 Optionally the I/O capabilities can be implemented as plugin Reader/Writer classes compiled
560  separately into dynamically loadable libraries and used via RearedPlugin and WriterPlugin classes.
561  Please note that all required libraries/dlls should be loadable.
562  See examples for details.
563 
564 In some cases the fine tuning of the Reader/Writer classes behavior can be done using a
565  map of string "options"
566 
567 ```
568  void Reader::set_options(const std::map<std::string, std::string>& options)
569 
570  std::map<std::string, std::string> Reader::get_options() const
571 ```
572 
573 The options for ReaderAsciiHepMC2
574 
575  - "disable_pad_cross_sections", "pad_cross_section_value"/"pad_cross_section_error"
576  If "disable_pad_cross_sections" is present the reader will keep a single cross-section per event, just
577  in the HepMC2 style. This is pre-3.2.6 default behaviour.
578  Otherwise, the cross-section vector will be expanded to the size of event weights. This is 3.2.6+ default behaviour.
579  If present, "pad_cross_section_value"/"pad_cross_section_error" values will be inserted into the cross-section vector.
580  Otherwise, the cross-sections and errors will be filled with zeros.
581 
582 
583  - "particle_flows_are_separated" ,"event_random_states_are_separated"
584  "vertex_weights_are_separated"
585  "particle_flows_are_separated"
586  Regulate if the corresponding information from IO_GenEvent would be stored into multiple attributes as
587  individual numbers, i.e. "separated" or as a single std::vector. The former behavior is used if
588  the corresponding option name is present in the list of options, regardless of the option value.
589  The later behavior is the default one.
590 
591 The option for WriterAscii and WriterAsciiHepMC2
592 
593  - "float_printf_specifier"
594  specifies the float printf specifier used for the output format of the floats.
595  Two first characters from the option value are used.
596  The default behavior is equivalent to setting this option to "e" and results in the output formatted as
597  " %.*e". To save the disk space one can use the "g" option, e.g.
598 ```
599  WriterAscii outputA("someoutput.hepmc");
600  auto optionsA = outputA.get_options();
601  optionsA["float_printf_specifier"] = "g";
602  outputA.set_options(optionsA);
603 ```
604 This option will be the default on in the future.
605 
606 
607 
608 
609 ### 7.1 Links
610 The relations between vertices and particles in GenEventData are encoded via
611  members links1 and links2, wich are std::vector<int> containing object ids.
612  Direct manipulations with links1 and links2 can be useful. For instance,
613  when the events are saved in ROOT format, one can extract the information
614  from links1 and links2 without reading whole event.
615  In case links1[i] is particle, links2[i] is end vertex. In case links1[i] is
616  vertex, links2[i] is outgoing particle. An example of usage is given below.
617 ```
618  // Andrii Verbytskyi, 2017, MPI fuer Physik
619  // Here is a code to look for a scattered DIS electron in HepMC3 event record using links.
620  // The implementation is extended to provide example of links1, links2 usage.
621  // Dummy code.
622  GenEventData* A=...
623  ...
624  int i;
625  int j;
626  int current_l=0; // If the incoming electron is the first particle in the list
627  int vertex_l=0; // We set final vertex to some nonsense value.
628  bool found_next_vertex=true;
629  while(found_next_vertex) // Looking for the electron end vertex
630  {
631  found_next_vertex=false;
632  for (i=0; i<A->links1.size(); i++) // Points from to ...
633  if (A->links1[i]>0 && // The link1 should be particle, i.e. >0
634  A->links1[i]==current_l+1) // The link1 should be our electron
635  {
636  vertex_l=A->links2[i]; // The end vertex of this electron is found
637  found_next_vertex=true;
638  }
639  std::vector<int> out; // Here we will save the outgoing particles
640  if (found_next_vertex)
641  {
642  for (j=0; j<A->links1.size(); j++) // Points from to ...
643  if (A->links1[j]<0 && // The link1 should be a vertex, i.e. <0
644  A->links1[j]==vertex_l) // The link1 should be our end vertex
645  if (std::abs(A->particles_pid[A->links2[j]-1])==11) // If the outgoing particle is e+/e-.
646  out.push_back(A->links2[j]);
647  if (out.size()==0) {
648  printf("Warning: no electron in the new vertex.\n");
649  break;
650  }
651  else
652  {
653  if (out.size()>1) printf("Warning: more than one electron in the new vertex.\n");
654  current_l=out.at(0)-1; // Pick up the first electron in the list and use it as current electron.
655  }
656  }
657  if (A->particles_status[current_l]==1) break; // The end particle is stable. It is our scattered electron.
658  }
659  ...
660  // Here we have cuts on electron
661 ```
662 
663 ## 8 Search-related classes and interfaces
664 HepMC3 comes with an optional Search library for finding particles
665  related to other particles or vertices.
666  It provides a set of functions to perform simple search operations e.g.
667 
668 ```
669  std::vector<HepMC3::GenParticlePtr> children_particles(const HepMC3::GenVertexPtr& O); ///< Return children particles
670  std::vector<HepMC3::ConstGenVertexPtr> grandchildren_vertices(const HepMC3::ConstGenVertexPtr& O); ///< Return grandchildren vertices
671  std::vector<HepMC3::GenParticlePtr> parent_particles(const HepMC3::GenVertexPtr& O); ///< Return parent particles
672  std::vector<HepMC3::GenVertexPtr> ancestor_vertices(const HepMC3::GenVertexPtr& obj); ///< Return ancestor vertices
673 
674 ```
675 and interfaces for a more advanced usage. For the latter two main interfaces are defined:
676  Relatives, for finding a particular type of relative, and Feature, for
677  generating filters based on Features extracted from particles.
678  In addition, operator on Filters are also defined.
679 
680 
681 ### 8.1 Relatives Interface
682 The Relatives interface is defined within search/include/HepMC3/Relatives.h.
683  Classes that obey this interface must provide a set of operator functions
684  that take either a GenParticlePtr, ConstGenParticlePtr, GenVertexPtr or
685  ConstGenVertexPtr and return a vector of either GenParticlePtr or
686  ConstGenParticlePtr. Note that the versions of the operator functions that
687  receive a consted input parameter also return a vector<ConstGenParticlePtr>,
688  while the versions that take non-consted input return non-const output.
689  This ensures consistency with the rule that consted objects may only return
690  pointers to const objects.
691 
692 The Relatives base class is abstract, and has a concrete implementation in
693  the templated RelativesInterface class. The RelativesInterface uses type
694  erasure so that any class that obeys the defined Relatives interface
695  (i.e that has the necessary four operator functions) can be wrapped in the
696  RelativesInterface without needing to inherit from Relatives directly.
697 
698 For example, if class foo implements the four necessary functions then the
699  following will work
700 
701 ```
702  using FooRelatives = RelativesInterface<Foo>;
703  Relatives * relos = new FooRelatives();
704  GenParticlePtr someInput;
705  vector<GenParticlePtr> foos = (*relos)(someInput);
706 ```
707 
708 The purpose of Relatives is to be able to wrap any viable class in a common
709  interface for finding relatives from a particle or vertex. Examples are
710  provided in the form of the _parents and _children classes. These do not
711  inherit from Raltives, but they do implement the necessary functions.
712  The _parents and _children class are not intended to be used directly, but
713  they are aliased by wrapping in the RelativesInterface:
714 
715 ```
716  using Parents = RelativesInterface<_parents>;
717  using Children = RelativesInterface<_children>;
718 ```
719 
720 Note as well that the _parents and _children classes use some utility aliases
721  to help declare the appropriately consted return type. For example
722 
723 ```
724  template<typename GenObject_type>
725  GenParticles_type<GenObject_type> operator()(GenObject_type);
726 ```
727 
728 has a return type GenParticles_type that is a vector of GenParticlePtr that
729  is consted if GenObject_type is const, but is not consted if GenObject_type
730  is not consted. Note as well the use of enable_if so that a single
731  implementation can be used for both the const and non-const version of the
732  functions. For the simple case of _parents the four required funcs could
733  have been implemented directly without such templating, but for more
734  complicated relatives it avoids duplicated code.
735 
736 
737 ### 8.2 Recursive Relatives
738 In addition to the RelativesInterface wrapper, Relatives.h also contains a
739  Recursive class that can wrap the underlying relation in recursion. For
740  example, recursion applied to the parents relationship provides all of the
741  ancestors, i.e. parents repeatedly applied to the output of parents. The
742  only additional requirement to use the Recursive wrapper is that the
743  underlying class must implement a vertex(GenParticlePtr) method that returns
744  the appropriate vertex to follow from a given GenParticle. As long as a
745  class has such a method, it is possible to make a recursive version of it
746 
747 ```
748  using Ancestors = RelativesInterface<Recursive<_parents> >;
749 ```
750 
751 
752 ### 8.3 Existing Relatives
753 The Relatives class contains static implementations of the Parents,
754  Children, Ancestors and Descendants relatives, which can be accessed and
755  used as follows
756 
757 ```
758 
759  vector<const Relatives*> relos{&Relatives::PARENTS, &Relatives::ANCESTORS, &Relatives::CHILDREN, &Relatives::DESCENDANTS};
760  ConstGenVertexPtr startPosition;
761  // loop over different relationships.
762  for(const Relatives* r: relos){
763  for(auto p: r(startPosition)){
764  // Do something with search result p
765  }
766  }
767 
768 ```
769 
770 
771 ### 8.4 Filters
772 A Filter is any object that has an operator that takes as input a
773  ConstGenParticlePtr and returns a bool that reflects whether the input
774  particle passes the filter requirements or not. Filter is defined in
775  Filter.h as
776 
777 ```
778  using Filter = std::function<bool(ConstGenParticlePtr)>;
779 ```
780 
781  Filter.h also contains some logical operators that allow filters to be
782  combined to create new filters, for example
783 
784 ```
785  Filter filter1, filter2, filter3;
786  Filter filter4 = filter1 && filter2;
787  Filter filter5 = filter3 || filter4;
788  Filter filter6 = !filter1;
789 ```
790 
791 Filter.h additionally contains a dummy filter that always accepts every
792  possible particle. This may be needed in functions that require a default
793  filter. The dummy filter is accessed as
794 
795 ```
796  Filter dummy = ACCEPT_ALL;
797 ```
798 
799 It is possible to define a Filter by hand. However, there are some utility
800  classes to define Filters based on features that can be obtained from GenParticles
801 
802 
803 ### 8.5 Feature Interface
804 The Feature interface is defined in Feature.h. The interface is templated
805  on a Feature_type that is any type that can be extracted from a GenParticle.
806  This is very flexible, and the only criteria is that the Feature must have
807  the set of comparison operators. While the templated Feature is general
808  enough to be used with any type of Feature, there are specialisations
809  for both integral and floating point features. The specialisations will
810  cover the vast majority of Features that are likely to be useful, although
811  note that Attributes may be a source of more complicated Features.
812 
813 To create a Feature, one need only wrap a lambda expression in the Feature
814  interface. For example, to create a Feature based on particle status or pT:
815 
816 ```
817  Feature<int> status([](ConstGenParticlePtr p)->int{return p->status();});
818  Feature<double> pT([](ConstGenParticlePtr p)->double{return p->momentum().pt()});
819 ```
820 
821 The more general form for any type of Feature would be
822 
823 ```
824  Feature<type> foo([](ConstGenParticlePtr p)->type{return p->foo();});
825 ```
826 
827 Having created a Feature, it can be used to create Filters for particle
828  selection. Applying operators to Features creates the Filter, which is
829  a functor that evaluates on a particle. For example
830 
831 ```
832  ConstGenParticlePtr p;
833  Filter is_stable = (status == 1);
834  bool stable = is_stable(p);
835 
836  // this evaluates true if p has pT above 100.
837  bool passPTCut = (pT > 100.)(p);
838 
839  // The Features can be combined
840  bool combined = ((pT > 100.) && (status == 1))(p);
841 
842 ```
843 
844  It is also possible to make a new Feature from the absolute value of
845  a previous Feature, e.g.
846 
847 ```
848  Feature<double> rapidity([](ConstGenParticlePtr p)->double{return p->momentum().rapidity()});
849  bool passes_rapCut = (abs(rapidity) < 2.5)(p);
850 ```
851 
852 Some standard features are contained within the non-templated Selector class
853 
854 
855 ## 9 Selectors and SelectorWrapper
856 Selector is a simplified interface that contains some predefined Features
857  that can be used to search. Selector defines comparisons operators for
858  both integral and floating point types, as well as the following selection
859  Features:
860 
861 ```
862  Selector::STATUS
863  Selector::PDG_ID
864  Selector::PT
865  Selector::ENERGY
866  Selector::RAPIDITY
867  Selector::ETA
868  Selector::PHI
869  Selector::ET
870  Selector::MASS
871  Selector::ATTRIBUTE(const std::string)
872 ```
873 
874 So, for example, a filter can be defined as follows
875 
876 ```
877  Filter f = (Selector::STATUS == 1 && Selector::PT > 60.) || (Selector::MASS > 70. && Selector::MASS < 110.);
878  GenParticlePtr p;
879  bool passesCuts = f(p);
880 ```
881 
882 As with Feature, it is possible to take tbe absolute value of a Selector.
883  However, note that while Featue is templated, Selector is abstract and so
884  it is not possible for abs() to return a Selector object directly, only a
885  pointer
886 
887 ```
888  Filter f = *abs(Selector::RAPIDITY) < 2.5;
889  bool passRapidity = f(p);
890 ```
891 
892 Note that the ATTRIBUTE selection is different from the others and does not
893  have the full set of comparison operators. This is a current limitation of
894  the Attributes, which are not guaranteed to offer all comparisons.
895  ATTRIBUTE takes a string, which is the name of the attribute, and permits
896  the equality operator and the method exists, which checks if the attribute is
897  even present
898 
899 ```
900  string name = "MyAttribute";
901  Filter f = Selector::ATTRIBUTE(name).exists() && Selector::ATTRIBUTE(name) == "My Value";
902  bool passesAttribute = f(p);
903 ```
904 
905 
906 ### 9.1 Applying Filters
907 The function applyFilter is used to apply the Filter to a set of particles.
908  See for example examples/BasicExamples/basic_tree.cc
909 
910 ```
911  for(ConstGenParticlePtr p: applyFilter( *abs(Selector::PDG_ID) <= 6, someParticles)){
912  Print::line(p);
913  }
914 ```
915 
916 
917 ## 10 Python Bindings
918 HepMC3 includes Python bindings codes suitable for compilation of python
919  modules for Python3.
920 
921  The binding codes were generated automatically using the binder utility
922  version 1.4.0 created by Sergey Lyskov (Johns Hopkins University) et al..
923  See
924  - https://cppbinder.readthedocs.io/en/latest/ for details.
925 
926  The binding codes use the pybind11 library version 2.6.0 by Wenzel Jakob,
927  EPFL's School of Computer and Communication Sciences.
928  See
929  - https://pybind11.readthedocs.io/en/master/
930  - Wenzel Jakob and Jason Rhinelander and Dean Moldovan,
931  "pybind11 -- Seamless operability between C++11 and Python", 2017, https://github.com/pybind/pybind11
932 
933 ### 10.1 Installation from repositories
934 The Python bindings together with the HepMC3 itself can be installed from PyPy and multiple other repositories.
935  Please see [HepMC3 page](https://gitlab.cern.ch/hepmc/HepMC3) at CERN GitLab for details.
936 
937 ### 10.2 Installation from sources
938 To turn on the compilation of bindings use -DHEPMC3_ENABLE_PYTHON = ON.
939  By default the python modules will be generated for python3 if these are found in the system.
940  In case the test suite is enabled, tests of python bindings with all the enabled versions will run as well.
941 
942  Despite not recommended, it should be possible to compile the python bindings using the installed version of HepMC3.
943  To do this, copy the python directory outside of source tree, uncomment #project(pyHepMC3 CXX) in python/CMakeLists.txt and
944  run CMake inside python directory with -DUSE_INSTALLED_HEPMC3=ON option.
945 
946 ### 10.3 Selected aspects of Python bindings
947 In general, the syntax used in the python bindings is exactly the same as in the C++ code.
948  However, some C++ classes and routines are don't have their Python equivalents, namsly:
949 
950  - The constructors of Readers/Writers with ifstreams/ostreams were not binded.
951  - The multithread reader ReaderMT was not binded.
952  - The explicit readers and writers with compression support ReaderGZ/WriterGZ were not binded.
953  It is recommended to use build-in python compression modules in combination with the desired readers/writers.
954  - The `deduce_reader` was binded and uses the internal python compression libraries, i.e. it has no dependence on zlib, zstd etc.
955  The only requirement is that the corresponding module is available.
956  - A limited support for the ROOTTree format I/O with the uproot ( see [uproot](https://github.com/scikit-hep/uproot5)) is offered, namely, with uproot module installed
957  HepMC3 ROOTTree files can be read with ReaderuprootTree.
958 
959 
960 ## 11 Handling Les Houches Event Files
961 This module contains helper classes and Reader and Writer classes
962  for handling Les Houches event files - LHEF.
963 
964 ### 11.1 Introduction
965 The Les Houches accord on an event file format (LHEF) to be used
966  for passing events from a matrix element generator program (MEG)
967  to an event generator program (EG) implementing parton showers,
968  underlying event models, and hadronisation models etc., was not
969  originally included in the HepMC event record format. But as the
970  demand for more information to be included in HepMC, it was
971  decided to allow HepMC to include also the original information
972  from a MEG in the LHEF format (see the run attribute
973  HepMC3::HEPRUPAttribute and event attribute
974  HepMC3::HEPEUPAttribute). A separate /standard/ implementation in
975  C++ of the LHEF format had already been maintained by Leif
976  Lönnblad, and it was decided to include this (header only -
977  LHEF.h) as a part of HepMC3. This will both be used in above
978  mentioned HepMC3::Attribute classes and as a kind of definition of
979  the LHEF format, which so far has not been extremely well
980  documented. From now on these pages will serve as the defining
981  information about the LHEF format.
982 
983 ### 11.2 Background
984 The original Les Houches accord for communicating between MEGs and EGs
985  was agreed upon in
986  2001 [arXiv:hep-ph/0109068](http://archive.org/abs/hep-ph/0109068)
987  and consisted of two simple FORTRAN common blocks. In fact this
988  structure survived in the LHEF format, which was introduced in
989  2006
990  [arXiv:hep-ph/0609017](http://archive.org/abs/hep-ph/0609017), and is still there after the updated
991  versions 2 in 2009
992  [arXiv:1003.1643](http://archive.org/abs/1003.1643), and 3 in
993  2013
994  [arXiv:1405.1067](http://archive.org/abs/1405.1067), and in the current proposal developed at the
995  Les Houches workshop on TeV Colliders 2015.
996 
997 As the methods for combining MEGs and EGs has advanced since the first
998  accord, from the tree-level merging methods and NLO matching at the
999  turn of the millennium, to the multi-jet (N)NLO matching and merging
1000  methods being perfected to day, the LHEF format has developed and a lot
1001  of optional information can be passed beyond the original common block
1002  structures. In the following all features included will be described,
1003  also those that were added a bit prematurely and later became
1004  deprecated.
1005 
1006 ### 11.3 The basic structure
1007 The LHEF format is based on XML, but has some oddities that goes
1008  beyond pure XML. As the name indicates, XML is extensible, and anyone
1009  writing a LHEF file can add whatever information she or he wants,
1010  however the following basic structure must be observed.
1011 
1012 ```
1013  <LesHouchesEvents version="3.0">
1014  <!--
1015  # optional information in completely free format,
1016  # except for the reserved end tag (see next line)
1017  -->
1018  <header>
1019  <!-- individually designed XML tags, in fancy XML style -->
1020  </header>
1021  <init>
1022  compulsory initialization information
1023  # optional initialization information
1024  </init>
1025  <event>
1026  compulsory event information
1027  # optional event information
1028  </event>
1029  <event>
1030  compulsory event information
1031  <!-- more optional information -->
1032  </event>
1033  <!-- and as many events that you want, but ending with -->
1034  </LesHouchesEvents>
1035 ```
1036 
1037  This looks like fairly normal XML tags, and indeed they are. The only
1038  addition to the structure is that the <tt>init</tt> and
1039  <tt>event</tt> (and their respective end tags) are required to be
1040  alone on a line, and the content of these blocks are required to
1041  start with a number of lines on a specific format that follows
1042  exactly the structure of the fortran common block in original Les
1043  Houches Accord. This means that the first line in the
1044  <tt>init</tt> block must start with a line containing the numbers
1045 ```
1046  IDBMUP(1) IDBMUP(2) EBMUP(1) EBMUP(2) PDFGUP(1) PDFGUP(2) PDFSUP(1) PDFSUP(2) IDWTUP NPRUP
1047 ```
1048  and the following <tt>NPRUP</tt> lines should be numbers on the form
1049 ```
1050  XSECUP(IPR) XERRUP(IPR) XMAXUP(IPR) LPRUP(IPR)
1051 ```
1052  where the different variable names are defined as follows:
1053  <ul>
1054  <li>
1055  <tt>IDBMUP</tt>: the PDG numbers for the incoming beams;
1056  </li>
1057  <li>
1058  <tt>EBMUP</tt>: the energy (in GeV) of the incoming beams;
1059  </li>
1060  <li>
1061  <tt>PDFGUP</tt> and <tt>PDFSUP</tt>: the LHAPDF group and set
1062  numbers for the corresponding parton densities used;
1063  </li>
1064  <li>
1065  <tt>IDWTUP</tt>: the weight strategy used;
1066  </li>
1067  <li>
1068  <tt>NPRUP</tt> The number of different processes included in the file;
1069  </li>
1070  </ul>
1071  and for each process <tt>IPR</tt>:
1072  <ul>
1073  <li>
1074  <tt>XSECUP</tt>, <tt>XERRUP</tt>: the Monte Carlo integrated
1075  cross section and error estimate;
1076  </li>
1077  <li>
1078  <tt>XMAXUP</tt>: the overestimated cross section such that the
1079  sum over the individual weights in each <tt>&lt;event&gt;</tt> times
1080  this file times <tt>XMAXUP</tt> equals <tt>XSECUP</tt> times
1081  the number of events;
1082  </li>
1083  <li>
1084  <tt>LPRUP</tt>: a unique integer identifying the corresponding
1085  process.
1086  </li>
1087  </ul>
1088  In the LHEF::Reader and LHEF::Writer classes this information is
1089  available as the public <tt>heprup</tt> member of class
1090  LHEF::HEPRUP with public members mimicking the Fortran common
1091  block variables.
1092 
1093 
1094 
1095  Similarly, every <tt>event</tt> block must start with a line
1096  containing then numbers
1097 ```
1098  NUP IDPRUP XWGTUP SCALUP AQEDUP AQCDUP
1099 ```
1100  and the following <tt>NUP</tt> lines should be numbers describing
1101  each particle on the form
1102 ```
1103  IDUP(I) ISTUP(I) MOTHUP(1,I) MOTHUP(2,I) ICOLUP(1,I) ICOLUP(2,I) PUP(1,I) PUP(2,I) PUP(3,I) PUP(4,I) PUP(5,I) VTIMUP(I) SPINUP(I)
1104 ```
1105  where the different variable names are defined as follows:
1106  <ul>
1107  <li>
1108  <tt>NUP</tt>: the number of particles in the event;
1109  </li>
1110  <li>
1111  <tt>IDPRUP</tt>: the integer identifying the process;
1112  </li>
1113  <li>
1114  <tt>XWGTUP</tt>: the (default) weight of the event
1115  </li>
1116  <li>
1117  <tt>SCALUP</tt>: the scale of the hard process in GeV;
1118  </li>
1119  <li>
1120  <tt>AQEDUP</tt>: the value of &alpha;<sub>EM</sub> used;
1121  </li>
1122  <li>
1123  <tt>AQCDUP</tt>: the value of &alpha;<sub>S</sub> used;
1124  </li>
1125  </ul>
1126  and for each particle <tt>I</tt>:
1127  <ul>
1128  <li>
1129  <tt>IDUP</tt>: the PDG code for the particle;
1130  </li>
1131  <li>
1132  <tt>ISTUP</tt>: the status code for the particle;
1133  </li>
1134  <li>
1135  <tt>MOTHUP</tt>: the line numbers of two particles considered
1136  to be the mothers of this particle;
1137  </li>
1138  <li>
1139  <tt>ICOLUP</tt>: indices of the colour and anti-colour lines connected to this particle;
1140  </li>
1141  <li>
1142  <tt>PUP</tt>: the <i>x</i>, <i>y</i>, <i>z</i> and energy
1143  component of the 4-momentum and invariant masss (in units of
1144  GeV);
1145  </li>
1146  <li>
1147  <tt>VTIMUP</tt>; the proper lifetime of this particle;
1148  </li>
1149  <li>
1150  <tt>SPINUP</tt>: the spin.
1151  </li>
1152  </ul>
1153  In the LHEF::Reader and LHEF::Writer classes this information is
1154  available as the public <tt>hepeup</tt> member of class
1155  LHEF::HEPEUP with public members mimicking the Fortran common block
1156  variables.
1157 
1158 ### 11.4 LHEF-taglist Additional information
1159 
1160  Over the years several additional XML-tags has been formalised to
1161  specify information on top of what is given in the original Les
1162  Houches accord common block. These are listed below. In most cases
1163  the tag name corresponds to a class with a corresponding name
1164  available as suitably named public members in the LHEF::HEPRUP and
1165  LHEF::HEPEUP class respectively.
1166 
1167  Note that a tag may contain attributes in the following ways:
1168 ```
1169  <tag attribute1="value" attribute2="value">things marked by the tag</tag>
1170  <tag attribute1="value" attribute2="value" attribute3="value" />
1171 ```
1172  where the second case is a tag with only attributes an no contents.
1173  In the following an attribute may be described as required (R) or
1174  optional with a default value (D).
1175 
1176 ### 11.5 LHEF-init Standardised tags in the init block
1177 
1178  The <tt>&lt;init&gt;</tt> block contains information about the full run
1179  (similar to the information contained in HepMC3::GenRunInfo). The
1180  following tags are defined.
1181 
1182  <ul>
1183  <li><b><tt>&lt;generator&gt;</tt></b> (optional, multiple, see LHEF::HEPRUP::generators)<br> For easy
1184  access to the generator(s) used to generate this file. An
1185  optional attribute <tt>version</tt> can be given with a string
1186  containg a version string. The content of the tag can include
1187  any generator specific information.
1188  </li>
1189  <li><b><tt>&lt;xsecinfo&gt;</tt></b> (required, multiple, see LHEF::HEPRUP::xsecinfos)<br>
1190  The information in the HEPRUP common block is in principle
1191  sufficient to figure out the cross sections of the processes
1192  involved. However, the way things are specified is a bit
1193  confusing and complicated since it was assumed to be used to
1194  pass information between the MEG and PSG in both
1195  directions. For the event file, the communication is per
1196  definition one-way, and the information can be made more
1197  easily provided. The attributes for the xsecinfo tag are as
1198  follows.
1199  <ul>
1200  <li>
1201  <tt>weightname</tt>: in case of multiple weights for each
1202  event several <tt>xsecinfo</tt> tags can be given, in
1203  which case This should give the corresponding weight
1204  name. If missing This is assumed to describe the default
1205  weight.
1206  </li>
1207  <li>
1208  <tt>neve</tt> (R): the number of events in the file
1209  </li>
1210  <li>
1211  <tt>totxsec</tt> (R): the total cross section (in units of pb) of all processes in the file
1212  </li>
1213  <li>
1214  <tt>maxweight</tt> (D=1): the maximum weight of any event
1215  in the file (in an arbitrary unit)
1216  </li>
1217  <li>
1218  <tt>minweight</tt>: if the file contains negative weights,
1219  the minweight may contain the most negative of the
1220  negative weights in the file. If not given, minweight is
1221  assumed to be <tt>-maxweight</tt> .
1222  </li>
1223  <li>
1224  <tt>meanweight</tt> (D=1): The average weight of the
1225  events in the file (same unit as <tt>maxweight</tt> ).
1226  </li>
1227  <li>
1228  <tt>negweights</tt> (D=no): If yes, the file may contain
1229  events with negative weights.
1230  </li>
1231  <li>
1232  <tt>varweights</tt> (D=no): If yes, the file may contain
1233  varying weights (if no, all events are weighted with
1234  </li>
1235  <li>
1236  <tt>maxweight</tt> (or <tt>minweight</tt> )).
1237  </li>
1238  </ul>
1239  </li>
1240  <li><b><tt>&lt;cutsinfo&gt;</tt></b>(optional, multiple, see LHEF::HEPRUP::cuts)<br>
1241  Used by a generator to supply information about the cuts used
1242  in the generation. Inside this tag any number of <tt>cut</tt>
1243  and <tt>ptype</tt> tags can be supplied.
1244  </li>
1245  <li><b><tt>&lt;ptype&gt;</tt></b>(optional, multiple, see LHEF::HEPRUP::ptypes)<br>
1246  To be used inside the <tt>cutinfo</tt> tag to define a group
1247  of particles to which a cut can be applied. Its contents
1248  should be a white-space separated list of PDG particle
1249  codes. The only attribute is <tt>name</tt> (R): the string
1250  used to represent this <tt>ptype</tt> in a <tt>cut</tt> .
1251  </li>
1252  <li><b><tt>&lt;cut&gt;</tt></b>(optional, multiple, see LHEF::HEPRUP::cuts)<br>
1253 
1254  This tag has information of a particular cut used. The information
1255  needed is which particle(s) are affected, what variable is used
1256  and the maximum and/or minimum value of that parameter. The
1257  contents should be the minimum and maximum allowed values of this
1258  variable (in that order). If only one number is given, there is no
1259  maximum. If the minumum is larger or equal to the maximum, there
1260  is no minimum. The attributes are:
1261  <ul>
1262  <li>
1263  <tt>p1</tt> (D=0): The particles for which this cut
1264  applies. This can either be a number corresponding to a given
1265  particle PDG code or 0 for any particle or the name of a
1266  previously defined <tt>ptype</tt> tag.
1267  </li>
1268  <li>
1269  <tt>p2</tt> (D=0): If cut is between pairs of particles, this
1270  attribute should be non-zero. Allowed values are the same as
1271  for <tt>p1</tt>.
1272  </li>
1273  <li>
1274  <tt>type</tt> (R) This defines the variable which is cut. The
1275  following values are standardised (the lab frame is assumed in all
1276  cases):
1277  <ul>
1278  <li>"m" the invariant mass of two particles (or, if only one
1279  particle type is given, the invariant mass of that particle).
1280  </li>
1281  <li>"kt": the transverse momenta of a particle matching p1 (in GeV)
1282  </li>
1283  <li>"eta": the pseudorapidity of a particle matching p1
1284  </li>
1285  <li>"y": the true rapidity of a particle matching p1
1286  </li>
1287  <li>"deltaR": the pseudorapidity--azimuthal-angle difference between two particles matching p1 and p2 respectively.
1288  </li>
1289  <li>"E": the energy of a particle matching p1
1290  </li>
1291  <li>"ETmiss": the norm of the vectorial sum of the pt of particles matching p1 and not matching p2.
1292  </li>
1293  <li>"HT": the scalar sum of the transverse momentum of the particles matching p1 and not matching p2.
1294  </li>
1295  <li> other values are allowed but are not included in the standard.
1296  </li>
1297  </ul>
1298  </li>
1299  </ul>
1300  </li>
1301  <li><b><tt>&lt;procinfo&gt;</tt></b>(optional, multiple, see LHEF::HEPRUP::procinfo)<br>
1302  For each process number used in the LPRUP variable in the HEPEUP
1303  common block we can have additional information given in the
1304  following attributes:
1305  <ul>
1306  <li>
1307  <tt>iproc</tt> (D=0): The process number for which the information is given. 0 means all processes.
1308  </li>
1309  <li>
1310  <tt>loops</tt>: The number of loops used in calculating this process.
1311  </li>
1312  <li>
1313  <tt>qcdorder</tt>: The number of QCD vertices used in calculating this process.
1314  </li>
1315  <li>
1316  <tt>eworder</tt>: The number of electro-weak vertices used in calculating this process.
1317  </li>
1318  <li>
1319  <tt>rscheme</tt>: The renormalization scheme used, if applicable.
1320  </li>
1321  <li>
1322  <tt>fscheme</tt>: The factorization scheme used, if applicable.
1323  </li>
1324  <li>
1325  <tt>scheme</tt>: Information about the scheme used to
1326  calculate the matrix elements. If absent, a pure tree-level
1327  calculation is assumed. Other possible values could be
1328  CSdipole (NLO calculation with Catani-Seymour subtraction),
1329  MCatNLO, POWHEG and NLOexclusive (NLO calculation according to
1330  the exclusive cross section withing the given cuts).
1331  </li>
1332  </ul>
1333  The content of this tag is a string with arbitrary information about the process.
1334  </li>
1335  <li><b><tt>&lt;mergeinfo&gt;</tt></b>(DEPRECATED, multiple, see LHEF::HEPRUP::mergeinfo)<br>
1336  For some merging schemes (eg. for CKKW) it is possible to reweight the
1337  the events with Sudakov form factors already in the MEG. If this has
1338  been done the content of the mergetype tag for the corresponding
1339  process should give a name corresponding to the scheme used. The
1340  attributes are:
1341  <ul>
1342  <li>
1343  <tt>iproc</tt> (D=0): The process number for which the
1344  information is given. "0" means all processes.
1345  </li>
1346  <li>
1347  <tt>mergingscale</tt> (R): The value of the merging scale in GeV.
1348  </li>
1349  <li>
1350  <tt>maxmult</tt> (D=no): If yes the corresponding process is
1351  reweighted as if it is the maximum multiplicity process
1352  (ie. the Sudakov for the last step down to the merging scale
1353  is not included.
1354  </li>
1355  </ul>
1356  </li>
1357  <li><b><tt>&lt;eventfile&gt;</tt></b>(optional, multiple, see LHEF::HEPRUP::eventfiles)<br>
1358  Sometimes it is desirable to have the events in a separate file(s)
1359  from the <tt>init</tt> block, eg. when the event files become
1360  extremely long or when one wants to write the <tt>init</tt> block
1361  after having generated all events in order to collect the correct
1362  information from the run. The names of these files are then
1363  specified in <tt>eventfile</tt> tags with attributes:
1364  <ul>
1365  <li>
1366  <tt>name</tt> (R): the file name, interpreted as an absolute
1367  path if starting with "/", otherwise as a relative path to
1368  where the file with the <tt>init</tt> block is located;
1369  </li>
1370  <li>
1371  <tt>neve</tt>: the number of events in the file
1372  </li>
1373  <li>
1374  <tt>ntries</tt>: the number of attempts the generator needed
1375  to produce the <tt>neve</tt> events (for statistics purposes.
1376  </li>
1377  </ul>
1378  </li>
1379  <li><b><tt>&lt;weightinfo&gt;</tt></b>(optional, multiple, see LHEF::HEPRUP::weightinfo)<br>
1380  When using multiple weights for events, each weight is given an
1381  index and a name using the <tt>weightinfo</tt> tag. The default
1382  weight (as given by the LHEF:HEPEUP::XWGTUP variable) is always
1383  treated as index 0 and given the name "Default", while the
1384  additional weights are indexed by the order of
1385  the <tt>weightinfo</tt> tags. The attributes are:
1386  <ul>
1387  <li>
1388  <tt>name</tt> (R): the name of the weight (in the best of all
1389  worlds this would be standardised);
1390  </li>
1391  <li>
1392  <tt>muf</tt>: he factor multiplying the nominal factorisation
1393  scale of the event for the given weight;
1394  </li>
1395  <li>
1396  <tt>mur</tt>: he factor multiplying the nominal
1397  renormalisation scale of the event for the given weight;
1398  </li>
1399  <li>
1400  <tt>pdf</tt> (D=0): the LHAPDF code used for the given weight
1401  (where 0 corresponds to the default PDF of the run);
1402  </li>
1403  <li>
1404  <tt>pdf2</tt> (D=pdf): the LHAPDF code used for the second
1405  beam for the given weight;
1406  </li>
1407  </ul>
1408  </li>
1409  <li><b><tt>&lt;weightgroup&gt;</tt></b>(optional, multiple, see
1410  LHEF::HEPRUP::weightgroup)<br> This tag can be used to
1411  group <tt>weightinfo</tt> tags together. The only well defined
1412  attribute is
1413  <tt>name</tt> giving a string which will be combined with
1414  the <tt>name</tt> attribute of the included <tt>weightinfo</tt>
1415  tags to give the HepMC weight names. Also other attributes can be
1416  included to include information about the weight variation used
1417  and is available in LHEF::WeightGroup::attributes.
1418  </li>
1419  <li><b><tt>&lt;initrwgt&gt;</tt></b>(optional, multiple, see LHEF::HEPRUP::weightinfo)<br>
1420  Due to an initial confusion some generators uses this tag instead
1421  of the <tt>weightinfo</tt> tag. It accepts the same attributes
1422  as <tt>weightinfo</tt>, except that the name attribute is
1423  named <tt>id</tt>. Note that some generators puts these tags in
1424  the <tt>header</tt> block. The current implementation of
1425  LHEF::Reader cannot handle this. Note that this is handled by the
1426  same classes (LHEF::WeightInfo and LHEF::WeightGroup) in
1427  LHEF::Reader and LHEF::Writer.
1428  </li>
1429  </ul>
1430 
1431 ### 11.6 LHEF-events Standardised tags in the events block.
1432 
1433  After the <tt>init</tt> block any number of events can be given. In
1434  addition events can be given in separate files declared
1435  with <tt>eventfile</tt> tags in the <tt>init</tt> block.
1436 
1437  The main tag here is simply called <tt>event</tt> and can (for
1438  statistics purposes) have an attribute <tt>ntries</tt> specifying how
1439  many attempts the generator needed to produce the event. Also other
1440  attributes can be given (and will be stored in the
1441  LHEF::HEPEUP::attributes member variable).
1442 
1443  The <tt>event</tt> tags may be grouped together in
1444  a <tt>eventgroup</tt> tag. This is useful mainly for NLO generators
1445  that produces a (number) "real" event(s) accompanied with a number of
1446  "counter" events, where these events should be treated together for
1447  statistics purposes. For this reason the <tt>eventgroup</tt> tag can
1448  be provided with the optional tags <tt>nreal</tt>
1449  and <tt>ncounter</tt> to indicate the number of <tt>event</tt> tags
1450  included of each type.
1451 
1452  As indicated above the block must start with required information
1453  corresponding to the original Les Houches Accord Fortran common
1454  block. Here is a list of additional tags that may be provided.
1455 
1456  <ul>
1457  <li><b><tt>&lt;weights&gt;</tt></b> (optional, see LHEF::HEPEUP::weights)<br>
1458  This tag should contain a list of weights or the events given in
1459  the order they are defined by the <tt>weightinfo</tt> tags in the <tt>init</tt>
1460  block. If for some obscure reason a weight is absent, it should be
1461  anyway be included as a zero weight. Note that the weight for the
1462  nominal event is still given by <tt>XWGTUP</tt>. This tag has no attributes.
1463  </li>
1464  <li><b><tt>&lt;rwgt&gt;</tt></b> (optional, see LHEF::HEPEUP::weights)<br>
1465  For weights declared in <tt>initrwgt</tt> this should contain the
1466  weights for this event. The difference wrt. the <tt>weights</tt>
1467  tag is that each weight needs to be given in an <tt>wgt</tt> tag.
1468  </li>
1469  <li><b><tt>&lt;wgt&gt;</tt></b> (optional, see LHEF::HEPEUP::weights)<br>
1470  This tag contains a single weight correspoding to a weight
1471  declared in the <tt>initrwgt</tt> tag. The only attribute is
1472  the <tt>id</tt> of the weight defined in the
1473  corresponing <tt>weightinfo</tt> tag.
1474  </li>
1475  <li><b><tt>&lt;scales&gt;</tt></b> (optional, see LHEF::HEPEUP::scales)<br>
1476  Since it is often necessary to define more than one scale in an
1477  event, this tag is provided to specify these. The scales are given
1478  as the followin attributes:
1479  <ul>
1480  <li>
1481  <tt>muf</tt> (D=SCLAUP): the factorisation scale (in GeV);
1482  </li>
1483  <li>
1484  <tt>mur</tt> (D=SCLAUP): the renormalisation scale (in GeV);
1485  </li>
1486  <li>
1487  <tt>mups</tt> (D=SCLAUP): the suggested parton shower starting
1488  scale (in GeV).
1489  </li>
1490  </ul>
1491  Also other attributes are allowed and will be stored in
1492  LHEF::Scales::attributes.<br>
1493  </li>
1494  <li><b><tt>&lt;scale&gt;</tt></b> (optional inside a <tt>&lt;scales&gt;</tt>
1495  tag)<br> It is also possible to specify an individual scale for
1496  any particle in the event in a <tt>&lt;scale&gt;</tt> tag. This is
1497  typically used to define a starting scale for a parton shower in
1498  resonance decays. This tag can also be used to specify the scale
1499  of a set of particle types. The following attributes can be given:
1500  <ul>
1501  <li>
1502  <tt>pos</tt> (optional): the index of the particle for which this scale
1503  applies, optionally followed by a space-separated list of
1504  indices of recoilers involved when the particle was created.
1505  </li>
1506  <li>
1507  <tt>etype</tt> (optional): a space-separated list of PDG codes
1508  for which this scale applies. The short-hand notation "EW" can
1509  be used to specify all leptons and electro-weak bosons, and
1510  "QCD" can be used to specify the guon and all quarks.
1511  </li>
1512  </ul>
1513  The contents of the tag gives the scale in GeV.
1514  </li>
1515  </ul>
1516  <hr>
1517 
1518 ## 12 Handling of errors and warnings
1519 
1520 The errors and warnings since HepMC 3.3.0 have the following categories:
1521 
1522  - 100 Generic I/O problem, i.e. cannot open a file, cannot import a python module from C++, etc.
1523  - 200 Buffer overflows
1524  - 300 Problem serializing attributes, problems with the momentum/lengths units
1525  - 400 Math error
1526  - 500 Usupported expression in the input
1527  - 600 Warious event inconsistencies
1528  - 700 Various operational warnings, e.g. adding a null particle to vertex
1529  - 800 GenCrosssection warnings
1530  - 900 Outdate formats or objects
1531 
1532 The categories are numbered approximatelly according to their importance and if the
1533 current warning/error level is set below the warning level of the category, the warnigs/errors from the category
1534 will not be printed.
1535 
1536 The warning level can be obtained/set with ``int Setup::warnings_level()/void Setup::set_warnings_level(const int)`` functions.
1537 The functions ``bool Setup::print_warnings()/void Setup::set_print_errors(const bool flag)`` obtain/set the global flag for printing the warnings but don't change the warning levels.
1538 The default warning level is 750.
1539 
1540 
1541 
1542 The error level can be obtained/set with ``int Setup::errors_level()/void Setup::set_errors_level(const int)`` functions.
1543 The functions ``bool Setup::print_errors()/void Setup::set_print_errors(const bool flag)`` obtain/set the global flag for printing the errors but don't change the warning levels.
1544 The default error level is 1000.
1545 
1546 
1547 
1548 
1549 
1550