Reference

pyhepmc

pyhepmc is a pythonic and Jupyter-friendly Python API for HepMC3.

Differences between HepMC3 C++ and pyhepmc

  • The pyhepmc API uses properties where the C++ API uses setters/getters (where possible).

  • Sequences with matching types and lengths are implicitly convertible to FourVector und ToolInfo.

  • In addition to the C++ Reader/Writer classes, we offer an easy to use open(). It can read and write any supported HepMC3 format, including compressed files (gzip, bzip2, lzma are supported).

  • Attributes for GenRunInfo, GenEvent, GenParticle, GenVertex can be accessed via a dict-like view returned by the attributes property. The view automatically converts between native C++ types to native Python types.

  • The Print class is missing, but listing() and content() are present as free functions.

  • The member functions delta_X of FourVector are free functions with two arguments.

  • HEPEVT_Wrapper and friends are missing, use GenEvent.from_hepevt() instead.

  • ReaderGZ and WriterGZ are missing, since open() offers this functionality.

  • API marked as deprecated in HepMC3 is not available in Python.

  • pyhepmc offers event visualization and renders in Jupyter notebooks if all required extra packages are installed, see pyhepmc.view.to_dot().

Missing functionality

  • ReaderMT will not be implemented. If you want to use multi-threaded IO, it is better to use the high-level threading API of Python to achieve this. You can read the next GenEvent in a child thread from a file while you are processing the current event in the main thread.

class pyhepmc.FourVector

Generic 4-vector.

Interpretation of its content depends on accessors used: it’s much simpler to do this than to distinguish between space and momentum vectors via the type system (especially given the need for backward compatibility with HepMC2). Be sensible and don’t call energy functions on spatial vectors! To avoid duplication, most definitions are only implemented on the spatial function names, with the energy-momentum functions as aliases.

This is not intended to be a fully featured 4-vector, but does contain the majority of common non-boosting functionality, as well as a few support operations on 4-vectors.

The implementations in this class are fully inlined.

abs_eta(self: pyhepmc._core.FourVector) float

Absolute pseudorapidity.

abs_rap(self: pyhepmc._core.FourVector) float

Absolute rapidity.

eta(self: pyhepmc._core.FourVector) float

Pseudorapidity.

interval(self: pyhepmc._core.FourVector) float

Spacetime invariant interval s^2 = t^2 - x^2 - y^2 - z^2.

is_zero(self: pyhepmc._core.FourVector) bool

Check if the length of this vertex is zero.

length(self: pyhepmc._core.FourVector) float

Magnitude of spatial (x, y, z) 3-vector.

length2(self: pyhepmc._core.FourVector) float

Squared magnitude of (x, y, z) 3-vector.

m(self: pyhepmc._core.FourVector) float

Invariant mass. Returns -sqrt(-m) if e^2 - P^2 is negative.

m2(self: pyhepmc._core.FourVector) float

Squared invariant mass m^2 = E^2 - px^2 - py^2 - pz^2.

p3mod(self: pyhepmc._core.FourVector) float

Magnitude of p3 = (px, py, pz) vector.

p3mod2(self: pyhepmc._core.FourVector) float

Squared magnitude of p3 = (px, py, pz) vector.

perp(self: pyhepmc._core.FourVector) float

Magnitude of (x, y) vector.

perp2(self: pyhepmc._core.FourVector) float

Squared magnitude of (x, y) vector.

phi(self: pyhepmc._core.FourVector) float

Azimuthal angle.

pt(self: pyhepmc._core.FourVector) float

Transverse momentum.

pt2(self: pyhepmc._core.FourVector) float

Squared transverse momentum px^2 + py^2.

rap(self: pyhepmc._core.FourVector) float

Rapidity.

rho(self: pyhepmc._core.FourVector) float

Magnitude of spatial (x, y, z) 3-vector, for HepMC2 compatibility.

theta(self: pyhepmc._core.FourVector) float

Polar angle w.r.t. z direction.

property e

Energy component of momentum.

property px

x-component of momentum.

property py

y-component of momentum.

property pz

z-component of momentum.

property t

Time component of position/displacement.

property x

x-component of position/displacement.

property y

y-component of position/displacement.

property z

z-component of position/displacement.

class pyhepmc.GenCrossSection

Stores additional information about cross-section.

This is an example of event attribute used to store cross-section information.

This class is meant to be used to pass, on an event by event basis, the current best guess of the total cross section. it is expected that the final cross section will be stored elsewhere.

Several cross sections and related info can be included in case of runs with multiple weights.

The units of cross-sections are expected to be pb.

set_cross_section(self: pyhepmc._core.GenCrossSection, cross_section: float, cross_section_error: float, accepted_events: int = -1, attempted_events: int = -1) None

Set all fields.

set_xsec(self: pyhepmc._core.GenCrossSection, index_or_name: object, value: float) None
set_xsec_err(self: pyhepmc._core.GenCrossSection, index_or_name: object, value: float) None
xsec(self: pyhepmc._core.GenCrossSection, index_or_name: object = 0) float

Access the cross section corresponding to the weight named wName.

Access the cross section corresponding to the weight with index indx.

Access the cross section corresponding to the weight named wName.

Access the cross section corresponding to the weight with index indx.

xsec_err(self: pyhepmc._core.GenCrossSection, index_or_name: object = 0) float

Access the cross section error corresponding to the weight named wName.

Access the cross section error corresponding to the weight with index indx.

Access the cross section error corresponding to the weight named wName.

Access the cross section error corresponding to the weight with index indx.

property accepted_events

Access the number of accepted events.

property attempted_events

Access the number of attempted events.

property is_valid

Verify that the instance contains non-zero information.

class pyhepmc.GenEvent

Stores event-related information.

Manages event-related information. Contains lists of GenParticle and GenVertex objects.

add_particle(self: pyhepmc._core.GenEvent, arg0: HepMC3::GenParticle) None

Add particle.

add_vertex(self: pyhepmc._core.GenEvent, arg0: HepMC3::GenVertex) None

Add vertex.

clear(self: pyhepmc._core.GenEvent) None

Remove contents of this event.

event_pos(self: pyhepmc._core.GenEvent) pyhepmc._core.FourVector

Vertex representing the overall event position.

from_hepevt(self: pyhepmc._core.GenEvent, event_number: int, px: numpy.ndarray[numpy.float64], py: numpy.ndarray[numpy.float64], pz: numpy.ndarray[numpy.float64], en: numpy.ndarray[numpy.float64], m: numpy.ndarray[numpy.float64], pid: numpy.ndarray[numpy.int32], status: numpy.ndarray[numpy.int32], parents: object = None, children: object = None, vx: object = None, vy: object = None, vz: object = None, vt: object = None, fortran: bool = True) None

Convert HEPEVT record to GenEvent.

GenEvent is cleared and then filled with the information provided in the argument arrays. Recovering the particle history only requires either parents or children. If both are not None, parents are used.

All particle arrays must have equal length.

Parameters:
  • event_number (int) – Event number, starting with zero.

  • px (array-like) – X-component of momentum of particles in GeV.

  • py (array-like) – Y-component of momentum of particles in GeV.

  • pz (array-like) – Z-component of momentum of particles in GeV.

  • en (array-like) – Energy of particles in GeV.

  • m (array-like) – Generated mass of particles in GeV.

  • pid (array-like) – PDG IDs of particles.

  • status (array-like) – Status codes of particles.

  • parents (array-like or None, optional) – Array of int with shape (N, 2). Start and stop index (inclusive) in the record for the parents of each particle. Indices must be zero-based (C) or one-based (Fortran). In the latter case, keyword subtract_one must be set to true. No parents for a particle are indicated by the pair (-1, -1). Single parents are indicated either by (N, -1) or (N, N) with N > 0.

  • children (array-like or None, optional) – Like parents but for the children of this particle.

  • vx (array-like or None, optional) – X-component of location of production vertex of particles in mm.

  • vy (array-like or None, optional) – Y-component of location of production vertex of particles in mm.

  • vz (array-like or None, optional) – Z-component of location of production vertex of particles in mm.

  • vt (array-like or None, optional) – Time (ct) of production vertex of particles in mm.

  • fortran (bool, optional) – If True (default), the source indices are 1-based (Fortran, Pythia8). Set this to False, if the indices are 0-based (C-style).

read_data(self: pyhepmc._core.GenEvent, data: pyhepmc._core.GenEventData) None

Fill GenEvent based on GenEventData.

remove_particle(self: pyhepmc._core.GenEvent, arg0: HepMC3::GenParticle) None

Remove particle from the event.

This function will remove whole sub-tree starting from this particle if it is the only incoming particle of this vertex. It will also production vertex of this particle if this vertex has no more outgoing particles.

remove_vertex(self: pyhepmc._core.GenEvent, arg0: HepMC3::GenVertex) None

Remove vertex from the event.

This will remove all sub-trees of all outgoing particles of this vertex.

reserve(self: pyhepmc._core.GenEvent, particles: int, vertices: int = 0) None

Reserve memory for particles and vertices.

Helps optimize event creation when size of the event is known beforehand.

set_beam_particles(self: pyhepmc._core.GenEvent, arg0: HepMC3::GenParticle, arg1: HepMC3::GenParticle) None
set_units(self: pyhepmc._core.GenEvent, arg0: pyhepmc._core.Units.MomentumUnit, arg1: pyhepmc._core.Units.LengthUnit) None

Change event units Converts event from current units to new ones.

set_weight(self: pyhepmc._core.GenEvent, name: str, value: float) None

Set event weight accessed by index (or the canonical/first one if there is no argument) or name.

Access by weight name requires a GenRunInfo attached to the event, otherwise this will throw an exception.

It’s the user’s responsibility to ensure that the index or name exists!

weight(self: pyhepmc._core.GenEvent, index_or_name: object = 0) float

Get event weight accessed by index (or the canonical/first one if there is no argument) or name.

Access by weight name requires a GenRunInfo attached to the event, otherwise this will throw an exception.

It’s the user’s responsibility to ensure that the index or name exists!

write_data(self: pyhepmc._core.GenEvent, data: pyhepmc._core.GenEventData) None

Fill GenEventData object.

property attributes

Access attributes with a dict-like view.

It is possible to read and write attributes. Primitive C++ types (and vectors therefore) are converted from/to native Python types.

property beams

Access beam particles as a sequence.

HepMC3 considers every particle as a beam particles, which does not have another ancestor particle. Therefore, this is a read-only property.

property cross_section

Access to the GenCrossSection.

If this attribute is not set, returns None.

property event_number

Access event number.

property heavy_ion

Access to the GenHeavyIon.

If this attribute is not set, returns None.

property length_unit

Get length unit.

property momentum_unit

Get momentum unit.

property numpy
property particles

Access list of particles.

property pdf_info

Access to the GenPdfInfo.

If this attribute is not set, returns None.

property run_info

Access to the GenRunInfo.

If this attribute is not set, returns None.

property vertices

Access list of vertices.

property weight_names

Get event weight names.

Access requires a GenRunInfo attached to the event, otherwise this will throw an exception.

property weights

Access event weight values as a sequence.

Can be used to set multiple weights by assigning a sequence.

class pyhepmc.GenEventData

Stores serializable event information.

property attribute_id

Attribute owner id.

property attribute_name

Attribute name.

property attribute_string

Attribute serialized as string.

property event_number

Event number.

property event_pos

Event position.

property length_unit

Length unit.

property links1

First id of the vertex links.

If this id is positive - it is the incoming particle id of a vertex which id is written in GenEventData::links2.

If this id is negative - it’s the id of a vertex which outgoing particle id is written in GenEventData::links2.

The links1[i] points to links2[i]. In case links1[i] is particle, links2[i] is end vertex. In case links2[i] is vertex, links2[i] is outgoing particle. An example of usage is given in documentation.

property links2

Second id of the vertex links.

property momentum_unit

Momentum unit.

property particles

Particles.

property vertices

Vertices.

property weights

Weights.

class pyhepmc.GenHeavyIon

Stores additional information about Heavy Ion generator.

This is an example of event attribute used to store Heavy Ion information.

property N_Nwounded_collisions

Collisions with a diffractively excited target nucleon.

The number of single diffractive nucleon-nucleon collisions where the target nucleon is excited. A negative value means that the information is not available.

property Ncoll

the number of inelastic nucleon-nucleon collisions.

Note that a one participating nucleon can be involved in many inelastic collisions, and that inelastic also includes diffractive excitation. A negative value means that the information is not available.

property Ncoll_hard

the number of hard nucleon-nucleon collisions.

Model-dependent. Usually the number of nucleon-nucleon collisions containing a special signal process. A negative value means that the information is not available.

property Npart_proj

the number of participating nucleons in the projectile.

The number of nucleons in the projectile participating in an inelastic collision (see Ncoll). A negative value means that the information is not available.

property Npart_targ

the number of participating nucleons in the target.

The number of nucleons in the target participating in an inelastic collision (see Ncoll). A negative value means that the information is not available.

property Nspec_proj_n

The number of spectator neutrons in the projectile.

ie. those that thave not participated in any inelastic nucleon-nucleon collision. A negative value indicatess that the information is not available.

property Nspec_proj_p

The number of spectator protons in the projectile.

ie. those that thave not participated in any inelastic nucleon-nucleon collision. A negative value indicatess that the information is not available.

property Nspec_targ_n

The number of spectator neutrons in the target.

ie. those that thave not participated in any inelastic nucleon-nucleon collision. A negative value indicatess that the information is not available.

property Nspec_targ_p

The number of spectator protons in the target.

ie. those that thave not participated in any inelastic nucleon-nucleon collision. A negative value indicatess that the information is not available.

property Nwounded_N_collisions

Collisions with a diffractively excited projectile nucleon.

The number of single diffractive nucleon-nucleon collisions where the projectile nucleon is excited. A negative value means that the information is not available.

property Nwounded_Nwounded_collisions

Non-diffractive or doubly diffractive collisions.

The number of nucleon-nucleon collisions where both projectile and target nucleons are wounded. A negative value means that the information is not available.

property centrality

The centrality.

The generated centrality in percentiles, where 0 is the maximally central and 100 is the minimally central. A negative value means that the information is not available.

property eccentricities

Eccentricities.

Calculated to different orders. The key of the map specifies the order, and the value gives the corresponding eccentricity.

property event_plane_angle

The event plane angle.

The angle wrt. the x-axix of the impact parameter vector (pointing frm the target to the projectile). A positive number between 0 and two pi. A negative value means that the information is not available.

property impact_parameter

The impact parameter.

The impact parameter given in units of femtometer. A negative value means that the information is not available.

property participant_plane_angles

Participant plane angles.

calculated to different orders. The key of the map specifies the order, and the value gives to the angle wrt. the event plane.

property sigma_inel_NN

The assumed inelastic nucleon-nucleon cross section.

in units of millibarn. As used in a Glauber calculation to simulate the distribution in Ncoll. A negative value means that the information is not available.

property user_cent_estimate

A user defined centrality estimator.

This variable may contain anything a generator feels is reasonable for estimating centrality. The value should be non-negative, and a low value corresponds to a low centrality. A negative value indicatess that the information is not available.

class pyhepmc.GenParticle

Stores particle-related information.

is_generated_mass_set(self: pyhepmc._core.GenParticle) bool

See generated_mass.

unset_generated_mass(self: pyhepmc._core.GenParticle) None

Declare that generated mass is not set. See generated_mass.

property abs_pid

Equivalent to abs(genparticle.pid)

property attributes

Access attributes with a dict-like view.

It is possible to read and write attributes. Primitive C++ types (and vectors therefore) are converted from/to native Python types.

property children

Convenience access to immediate outgoing particles via end vertex.

Less efficient than via the vertex since return must be by value (in case there is no vertex).

property end_vertex

Access end vertex.

property generated_mass

Get or set generated mass.

This function will return mass as set by a generator/tool. If not set, it will return momentum().m().

property id

Number which uniquely identifies this particle in the current event. This is not a particle ID, see pid for the PDG ID. The id is a persistent pointer to this particle.

property in_event

Check if this particle belongs to an event.

property momentum

Access momentum.

property parent_event

Access parent event.

property parents

Convenience access to immediate incoming particles via production vertex.

Less efficient than via the vertex since return must be by value (in case there is no vertex).

property pid

PDG ID of particle

property production_vertex

Access production vertex.

property status

Access status code.

class pyhepmc.GenPdfInfo

Stores additional information about PDFs.

This is an example of event attribute used to store PDF-related information.

Input parton flavour codes id1 & id2 are expected to obey the PDG code conventions, especially g = 21.

The contents of pdf1 and pdf2 are expected to be x*f(x). The LHAPDF set ids are the entries in the first column of http:///projects.hepforge.org/lhapdf/PDFsets.index.

property parton_id1

Parton PDG ID.

property parton_id2

Parton PDG ID.

property pdf_id1

LHAPDF ID code.

property pdf_id2

LHAPDF ID code.

property scale

Factorisation scale (in GEV).

property x1

Parton momentum fraction.

property x2

Parton momentum fraction.

property xf1

PDF value.

property xf2

PDF value.

class pyhepmc.GenRunInfo

Stores run-related information.

Manages run-related information. Contains run-wide attributes.

class ToolInfo

Interrnal struct for keeping track of tools.

property description
property name
property version
property attributes

Access attributes with a dict-like view.

It is possible to read and write attributes. Primitive C++ types (and vectors therefore) are converted from/to native Python types.

property tools
property weight_names

Access the vector of weight names.

Access the names of the weights in this run.

For consistency, the length of the vector should be the same as the number of weights in the events in the run.

class pyhepmc.GenVertex

Stores vertex-related information.

add_particle_in(self: pyhepmc._core.GenVertex, arg0: pyhepmc._core.GenParticle) None

Add incoming particle.

add_particle_out(self: pyhepmc._core.GenVertex, arg0: pyhepmc._core.GenParticle) None

Add outgoing particle.

has_set_position(self: pyhepmc._core.GenVertex) bool

Check if position of this vertex is set.

remove_particle_in(self: pyhepmc._core.GenVertex, arg0: pyhepmc._core.GenParticle) None

Remove incoming particle.

remove_particle_out(self: pyhepmc._core.GenVertex, arg0: pyhepmc._core.GenParticle) None

Remove outgoing particle.

property attributes
property id

Number which uniquely identifies this vertex in the current event. The id is a persistent pointer to this vertex.

property in_event

Check if this vertex belongs to an event.

property parent_event

Access parent event.

property particles_in

Access list of incoming particles.

property particles_out

Access list of outgoing particles.

property position

Access vertex position.

Returns the position of this vertex. If a position is not set on this vertex, the production vertices of ancestors are searched to find the inherited position. FourVector(0,0,0,0) is returned if no position information is found.

Access vertex position.

property status

Access vertex status code.

class pyhepmc.HEPEUPAttribute
momentum(self: pyhepmc._core.HEPEUPAttribute, arg0: int) pyhepmc._core.FourVector
property hepeup
class pyhepmc.HEPRUPAttribute
property heprup
class pyhepmc.Setup

Imitates the Setup namespace.

You can directly read and write to the attributes of this class without creating an instance. They manipulate the corresponding global values in the HepMC3 C++ library.

debug_level
print_errors
print_warnings
class pyhepmc.Units

Units in which HepMC3 stores momentum and position.

class LengthUnit

Members:

CM

MM

CM = <LengthUnit.CM: 1>
MM = <LengthUnit.MM: 0>
property name
property value
class MomentumUnit

Members:

MEV

GEV

GEV = <MomentumUnit.GEV: 1>
MEV = <MomentumUnit.MEV: 0>
property name
property value
CM = <LengthUnit.CM: 1>
GEV = <MomentumUnit.GEV: 1>
MEV = <MomentumUnit.MEV: 0>
MM = <LengthUnit.MM: 0>
pyhepmc.content(arg0: pyhepmc._core.GenEvent) str

Print content of all GenEvent containers.

pyhepmc.delta_eta(arg0: pyhepmc._core.FourVector, arg1: pyhepmc._core.FourVector) float

Pseudorapidity separation.

pyhepmc.delta_phi(arg0: pyhepmc._core.FourVector, arg1: pyhepmc._core.FourVector) float

Signed azimuthal angle separation in [-pi, pi].

pyhepmc.delta_r2_eta(arg0: pyhepmc._core.FourVector, arg1: pyhepmc._core.FourVector) float

R_eta^2-distance separation dR^2 = dphi^2 + deta^2.

pyhepmc.delta_r2_rap(arg0: pyhepmc._core.FourVector, arg1: pyhepmc._core.FourVector) float

R_rap^2-distance separation dR^2 = dphi^2 + drap^2.

pyhepmc.delta_r_eta(arg0: pyhepmc._core.FourVector, arg1: pyhepmc._core.FourVector) float

R_eta-distance separation dR = sqrt(dphi^2 + deta^2).

pyhepmc.delta_r_rap(arg0: pyhepmc._core.FourVector, arg1: pyhepmc._core.FourVector) float

R-rap-distance separation dR = sqrt(dphi^2 + drap^2).

pyhepmc.delta_rap(arg0: pyhepmc._core.FourVector, arg1: pyhepmc._core.FourVector) float

Rapidity separation.

pyhepmc.equal_particle_sets(arg0: List[pyhepmc._core.GenParticle], arg1: List[pyhepmc._core.GenParticle]) bool

Whether two sets of particles are equal. Algorithm is analog to equal_vertex_sets().

pyhepmc.equal_vertex_sets(arg0: List[pyhepmc._core.GenVertex], arg1: List[pyhepmc._core.GenVertex]) bool

Whether two sets of vertices are equal. The vertices can have arbitrary order and their id() is ignored in the comparison. The sets are considered equal, if each vertex in one set compares equal to another vertex in the other set.

pyhepmc.listing(arg0: pyhepmc._core.GenEvent, arg1: int) str

Print event in listing (HepMC2) format.

pyhepmc.open(fileobj: str | PurePath, mode: str = 'r', precision: int | None = None, format: str | None = None) Any

Open HepMC files for reading or writing.

See HepMCFile.

pyhepmc.io

HepMC3 IO classes.

This module contains various Reader and Writer classes, which have pythonic interfaces. They act as context managers and close the file automatically. The Readers can be iterated over and yield events.

The open() function is even easier to use. It can read any supported file and will auto-detect the format. It can be used for reading and writing.

class pyhepmc.io.ReaderAscii

Reader for HepMC3 ASCII files.

class pyhepmc.io.ReaderAsciiHepMC2

Reader for HepMC2 ASCII files.

class pyhepmc.io.ReaderHEPEVT

Reader for HEPEVT files.

class pyhepmc.io.ReaderLHEF

Reader for LHEF files.

read_event(self: pyhepmc._core.Reader, event: pyhepmc._core.GenEvent) bool
class pyhepmc.io.UnparsedAttribute

Unparsed attribute after deserialization.

HepMC3 does not serialize the type of attributes, therefore the correct type cannot be restored upon deserialization (this is a limition of the HepMC3 C++ library and its serialization format). Use the astype() method to parse the attribute into a concrete type; this has important side-effects, see method description.

astype(self: pyhepmc._core.UnparsedAttribute, pytype: object) object

Convert unparsed attribute to concrete type.

If the conversion is successful, the unparsed attribute is replaced with the parsed attribute, so this method has to be called only once. If the conversion fails, a TypeError is raised.

Parameters:

pytype (type) – Type of the attribute. Allowed values: bool, int, float, str, GenParticle, GenPdfInfo, GenHeavyIon, GenCrossSection, HEPRUPAttribute, HEPEUPAttribute, typing.List[int], typing.List[float], typing.List[str]. In Python-3.9+, typing.List can be replaced by list.

class pyhepmc.io.WriterAscii
write()

write_event(self: pyhepmc._core.Writer, event: pyhepmc._core.GenEvent) -> None

property precision

Return output precision.

Access output precision.

So far available range is [2,24]. Default is 16.

class pyhepmc.io.WriterAsciiHepMC2
write()

write_event(self: pyhepmc._core.Writer, event: pyhepmc._core.GenEvent) -> None

property precision

Return output precision.

Access output precision.

Available range is [2,24]. Default is 16.

class pyhepmc.io.WriterHEPEVT
write()

write_event(self: pyhepmc._core.Writer, event: pyhepmc._core.GenEvent) -> None

pyhepmc.io.open(fileobj: str | PurePath, mode: str = 'r', precision: int | None = None, format: str | None = None) Any

Open HepMC files for reading or writing.

See HepMCFile.

pyhepmc.view

Visualization for GenEvent.

pyhepmc.view.savefig(event: GenEvent | Digraph, fname: str | BinaryIO, *, format: str | None = None, **kwargs: Any) None

Save event as an image.

The SVG format is recommended, because it contains tooltips with extra information.

Supported formats: .

Parameters:
  • event (GenEvent or Digraph) – The event to be saved.

  • fname (path-like or file-like) – Path or file-like handle to which the image is written.

  • format (str, optional) – Output format.

  • **kwargs – Other arguments are forwarded to pyhepmc.view.to_dot.

pyhepmc.view.to_dot(evt: GenEvent, *, size: Tuple[int, int] | None = None, color_hadron: str = 'black', color_lepton_or_boson: str = 'goldenrod', color_quark_or_gluon: str = 'darkred', color_internal: str = 'dodgerblue', color_invalid: str = 'gainsboro', color_special: str = 'green') Digraph

Convert GenEvent to Digraph.

This converts the GenEvent into a graphviz Digraph in the DOT language, which in turn can be rendered in SVG format and displayed in Jupyter notebooks.

Parameters:
  • evt (GenEvent) – The event to convert.

  • size ((int, int) or None, optional) – Maximum size of the graph in inches.

  • color_hadron (str, optional) – Color (HTML color specification) for hadrons.

  • color_lepton_or_boson (str, optional) – Color (HTML color specification) for leptons or bosons (gamma, W+/-, Z, H).

  • color_quark_or_gluon (str, optional) – Color (HTML color specification) for quarks, diquarks, and gluons.

  • color_internal (str, optional) – Color (HTML color specification) for generator-internal particles (e.g. strings, clusters, …).

  • color_special (str, optional) – Color (HTML color specification) for any valid particle which does not fit into the other categories (e.g. BSM particles).

  • color_invalid (str, optional) – Color (HTML color specification) for invalid particles (e.g. PID==0).