Skip to content

Commit 552f72d

Browse files
authored
include tool in name, fix ratio for sibyll, attributes for toolset, particle optional (#29)
1 parent f6c77d4 commit 552f72d

File tree

5 files changed

+48
-27
lines changed

5 files changed

+48
-27
lines changed

src/core.cpp

Lines changed: 10 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -187,8 +187,11 @@ inline std::ostream& repr(std::ostream& os, const HepMC3::GenRunInfo::ToolInfo&
187187
}
188188

189189
inline std::ostream& repr(std::ostream& os, const HepMC3::FourVector& x) {
190+
const int saved = os.precision();
191+
os.precision(3);
190192
os << "FourVector(" << x.x() << ", " << x.y() << ", " << x.z() << ", " << x.t()
191193
<< ")";
194+
os.precision(saved);
192195
return os;
193196
}
194197

@@ -421,7 +424,13 @@ PYBIND11_MODULE(_core, m) {
421424
repr(os, self);
422425
return os.str();
423426
})
424-
.def(py::self == py::self);
427+
.def(py::self == py::self)
428+
// clang-format off
429+
ATTR(name, GenRunInfo::ToolInfo)
430+
ATTR(version, GenRunInfo::ToolInfo)
431+
ATTR(description, GenRunInfo::ToolInfo)
432+
// clang-format on
433+
;
425434

426435
py::implicitly_convertible<py::sequence, GenRunInfo::ToolInfo>();
427436

src/pybind.h

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -27,5 +27,6 @@ _T overload_cast(_T x) {
2727
.def_property(#name, overload_cast<rval, const cls>(&cls::name), &cls::set_##name)
2828
#define METH(name, cls) .def(#name, &cls::name)
2929
#define METH_OL(name, cls, rval, args) .def(#name, (rval(cls::*)(args)) & cls::name)
30+
#define ATTR(name, cls) .def_readwrite(#name, &cls::name)
3031

3132
#endif

src/pyhepmc/__init__.py

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -28,7 +28,7 @@
2828
WriterHEPEVT,
2929
pyhepmc_open as open,
3030
)
31-
from ._version import version as __version__ # noqa
31+
from ._version import __version__ # noqa
3232

3333
try:
3434
from .view import to_dot

src/pyhepmc/view.py

Lines changed: 35 additions & 24 deletions
Original file line numberDiff line numberDiff line change
@@ -1,26 +1,29 @@
11
from graphviz import Digraph
2-
from particle import Particle, ParticleNotFound, InvalidParticle
32
from pyhepmc._prettify import db as _prettify
43
from pyhepmc import Units
54
import numpy as np
65

76

87
def to_dot(evt, style=None):
9-
d = Digraph(name="event %i" % evt.event_number)
8+
if evt.run_info:
9+
name = "\n".join(f"{t.name}-{t.version}" for t in evt.run_info.tools)
10+
name += "\n"
11+
else:
12+
name = ""
13+
name += f"event number = {evt.event_number}"
14+
15+
d = Digraph(name=name)
1016
d.node_attr["shape"] = "point"
1117
d.graph_attr["rankdir"] = "LR"
1218
d.graph_attr["size"] = "8,100"
13-
d.graph_attr["ratio"] = "compress"
1419

1520
GeV = 1 if evt.momentum_unit == Units.GEV else 1e3
1621

1722
for v in evt.vertices:
18-
d.node(f"{v.id}", tooltip=f"{v.position} mm\nstatus={v.status}")
23+
d.node(f"{v.id}", tooltip=f"{v.position} mm\nstatus = {v.status}")
1924

2025
enmax = max(p.momentum.e / GeV for p in evt.particles)
2126
for p in evt.particles:
22-
style = "solid"
23-
2427
en = p.momentum.e / GeV
2528
penwidth = str(max(0.5, 4 * en / enmax))
2629

@@ -32,33 +35,41 @@ def to_dot(evt, style=None):
3235
en *= 1e3
3336
unit = "MeV"
3437

38+
color = "black"
39+
style = "solid"
40+
41+
pname = _prettify.get(p.pid, None)
42+
if pname is None:
43+
if p.pid == 0:
44+
pname = "Invalid(0)"
45+
color = "gainsboro"
46+
penwidth = "7"
47+
else:
48+
pname = f"Internal({p.pid})"
49+
color = "dodgerblue"
50+
penwidth = "7"
51+
52+
tooltip = f"{pname} [PDGID: {int(p.pid)}]"
53+
tooltip += f"\n{p.momentum} GeV"
54+
tooltip += f"\nm = {p.generated_mass:.4g} GeV"
55+
tooltip += f"\nstatus = {p.status}"
56+
3557
try:
58+
from particle import Particle, ParticleNotFound, InvalidParticle
59+
3660
pdb = Particle.from_pdgid(p.pid)
37-
pname = pdb.name
38-
pname = _prettify.get(pdb.pdgid, pname)
39-
tooltip = f"{pdb.name} [{int(pdb.pdgid)}]"
4061
quarks = pdb.quarks
4162
if quarks:
42-
tooltip += f"\n{quarks}"
63+
tooltip += f"\nquarks = {quarks}"
64+
else: # boson or lepton
65+
color = "goldenrod"
4366
tooltip += f"\nQ = {pdb.charge:.3g}"
44-
if pdb.mass:
45-
tooltip += f"\nm = {pdb.mass:.4g} MeV/c2"
4667
if pdb.ctau and np.isfinite(pdb.ctau):
4768
tooltip += f"\ncτ = {pdb.ctau:.3g} mm"
48-
color = "black" if quarks else "goldenrod"
49-
# if not quarks: # boson or lepton
5069
if pdb.charge == 0:
5170
style = "dashed"
52-
except ParticleNotFound:
53-
pname = f"Internal({p.pid})"
54-
tooltip = pname
55-
color = "dodgerblue"
56-
penwidth = "7"
57-
except InvalidParticle:
58-
pname = f"Invalid({p.pid})"
59-
tooltip = pname
60-
color = "gainsboro"
61-
penwidth = "7"
71+
except (ModuleNotFoundError, ParticleNotFound, InvalidParticle):
72+
pass
6273

6374
label = f"{pname} {en:.2g} {unit}"
6475

tests/test_view.py

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -9,7 +9,7 @@ def test_dot(evt): # noqa
99
d = view.to_dot(evt)
1010
s = str(d)
1111
print(s)
12-
assert s.startswith('digraph "event 1')
12+
assert s.startswith('digraph "dummy-1.0\nevent number = 1')
1313
assert "in_1" in s
1414
assert "in_2" in s
1515
assert "in_3" not in s

0 commit comments

Comments
 (0)