Documentation ¶
Overview ¶
Package hepmc is a pure Go implementation of the C++ HepMC-2 library.
Index ¶
Examples ¶
Constants ¶
This section is empty.
Variables ¶
This section is empty.
Functions ¶
Types ¶
type CrossSection ¶
type CrossSection struct { Value float64 // value of the cross-section (in pb) Error float64 // error on the value of the cross-section (in pb) }
CrossSection is used to store the generated cross section. This type 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.
type Decoder ¶
type Decoder struct {
// contains filtered or unexported fields
}
Decoder decodes a hepmc Event from a stream
func NewDecoder ¶
NewDecoder returns a new hepmc Decoder that reads from the io.Reader.
type Encoder ¶
type Encoder struct {
// contains filtered or unexported fields
}
Encoder encodes a hepmc Event into a stream.
func NewEncoder ¶
NewEncoder returns a new hepmc Encoder that writes into the io.Writer.
type Event ¶
type Event struct { SignalProcessID int // id of the signal process EventNumber int // event number Mpi int // number of multi particle interactions Scale float64 // energy scale, AlphaQCD float64 // QCD coupling, see hep-ph/0109068 AlphaQED float64 // QED coupling, see hep-ph/0109068 SignalVertex *Vertex // signal vertex Beams [2]*Particle // incoming beams Weights Weights // weights for this event. first weight is used by default for hit and miss RandomStates []int64 // container of random number generator states Vertices map[int]*Vertex Particles map[int]*Particle CrossSection *CrossSection HeavyIon *HeavyIon PdfInfo *PdfInfo MomentumUnit MomentumUnit LengthUnit LengthUnit // contains filtered or unexported fields }
Event represents a record for MC generators (for use at any stage of generation)
This type is intended as both a "container class" ( to store a MC
event for interface between MC generators and detector simulation ) and also as a "work in progress class" ( that could be used inside a generator and modified as the event is built ).
Example (BuildFromScratch) ¶
This example will place the following event into HepMC "by hand"
name status pdg_id parent Px Py Pz Energy Mass 1 !p+! 3 2212 0,0 0.000 0.000 7000.000 7000.000 0.938 2 !p+! 3 2212 0,0 0.000 0.000-7000.000 7000.000 0.938
=========================================================================
3 !d! 3 1 1,1 0.750 -1.569 32.191 32.238 0.000 4 !u~! 3 -2 2,2 -3.047 -19.000 -54.629 57.920 0.000 5 !W-! 3 -24 1,2 1.517 -20.68 -20.605 85.925 80.799 6 !gamma! 1 22 1,2 -3.813 0.113 -1.833 4.233 0.000 7 !d! 1 1 5,5 -2.445 28.816 6.082 29.552 0.010 8 !u~! 1 -2 5,5 3.962 -49.498 -26.687 56.373 0.006
now we build the graph, which will look like
# p7 # # p1 / # # \v1__p3 p5---v4 # # \_v3_/ \ # # / \ p8 # # v2__p4 \ # # / p6 # # p2 # # #
package main import ( "bytes" "fmt" "log" "os" "go-hep.org/x/hep/fmom" "go-hep.org/x/hep/hepmc" ) func main() { var err error // first create the event container, with signal process 20, event number 1 evt := hepmc.Event{ SignalProcessID: 20, EventNumber: 1, Particles: make(map[int]*hepmc.Particle), Vertices: make(map[int]*hepmc.Vertex), } defer evt.Delete() // define the units evt.MomentumUnit = hepmc.GEV evt.LengthUnit = hepmc.MM // create vertex 1 and 2, together with their in-particles v1 := &hepmc.Vertex{} err = evt.AddVertex(v1) if err != nil { log.Fatal(err) } err = v1.AddParticleIn(&hepmc.Particle{ Momentum: fmom.PxPyPzE{0, 0, 7000, 7000}, PdgID: 2212, Status: 3, }) if err != nil { log.Fatal(err) } v2 := &hepmc.Vertex{} err = evt.AddVertex(v2) if err != nil { log.Fatal(err) } err = v2.AddParticleIn(&hepmc.Particle{ Momentum: fmom.PxPyPzE{0, 0, -7000, 7000}, PdgID: 2212, Status: 3, //Barcode: 2, }) if err != nil { log.Fatal(err) } // create the outgoing particles of v1 and v2 p3 := &hepmc.Particle{ Momentum: fmom.PxPyPzE{.750, -1.569, 32.191, 32.238}, PdgID: 1, Status: 3, // Barcode: 3, } err = v1.AddParticleOut(p3) if err != nil { log.Fatal(err) } p4 := &hepmc.Particle{ Momentum: fmom.PxPyPzE{-3.047, -19., -54.629, 57.920}, PdgID: -2, Status: 3, // Barcode: 4, } err = v2.AddParticleOut(p4) if err != nil { log.Fatal(err) } // create v3 v3 := &hepmc.Vertex{} err = evt.AddVertex(v3) if err != nil { log.Fatal(err) } err = v3.AddParticleIn(p3) if err != nil { log.Fatal(err) } err = v3.AddParticleIn(p4) if err != nil { log.Fatal(err) } err = v3.AddParticleOut(&hepmc.Particle{ Momentum: fmom.PxPyPzE{-3.813, 0.113, -1.833, 4.233}, PdgID: 22, Status: 1, }) if err != nil { log.Fatal(err) } p5 := &hepmc.Particle{ Momentum: fmom.PxPyPzE{1.517, -20.68, -20.605, 85.925}, PdgID: -24, Status: 3, } err = v3.AddParticleOut(p5) if err != nil { log.Fatal(err) } // create v4 v4 := &hepmc.Vertex{ Position: fmom.PxPyPzE{0.12, -0.3, 0.05, 0.004}, } err = evt.AddVertex(v4) if err != nil { log.Fatal(err) } err = v4.AddParticleIn(p5) if err != nil { log.Fatal(err) } err = v4.AddParticleOut(&hepmc.Particle{ Momentum: fmom.PxPyPzE{-2.445, 28.816, 6.082, 29.552}, PdgID: 1, Status: 1, }) if err != nil { log.Fatal(err) } err = v4.AddParticleOut(&hepmc.Particle{ Momentum: fmom.PxPyPzE{3.962, -49.498, -26.687, 56.373}, PdgID: -2, Status: 1, }) if err != nil { log.Fatal(err) } evt.SignalVertex = v3 err = evt.Print(os.Stdout) if err != nil { log.Fatal(err) } out := new(bytes.Buffer) enc := hepmc.NewEncoder(out) err = enc.Encode(&evt) if err != nil { log.Fatal(err) } fmt.Printf("%s\n", out.Bytes()) }
Output: ________________________________________________________________________________ GenEvent: #0001 ID= 20 SignalProcessGenVertex Barcode: -3 Momentum units: GEV Position units: MM Entries this event: 4 vertices, 8 particles. Beam Particles are not defined. RndmState(0)= Wgts(0)= EventScale 0.00000 [energy] alphaQCD=0.00000000 alphaQED=0.00000000 GenParticle Legend Barcode PDG ID ( Px, Py, Pz, E ) Stat DecayVtx ________________________________________________________________________________ GenVertex: -1 ID: 0 (X,cT):0 I: 1 1 2212 +0.00e+00,+0.00e+00,+7.00e+03,+7.00e+03 3 -1 O: 1 3 1 +7.50e-01,-1.57e+00,+3.22e+01,+3.22e+01 3 -3 GenVertex: -2 ID: 0 (X,cT):0 I: 1 2 2212 +0.00e+00,+0.00e+00,-7.00e+03,+7.00e+03 3 -2 O: 1 4 -2 -3.05e+00,-1.90e+01,-5.46e+01,+5.79e+01 3 -3 GenVertex: -3 ID: 0 (X,cT):0 I: 2 3 1 +7.50e-01,-1.57e+00,+3.22e+01,+3.22e+01 3 -3 4 -2 -3.05e+00,-1.90e+01,-5.46e+01,+5.79e+01 3 -3 O: 2 5 22 -3.81e+00,+1.13e-01,-1.83e+00,+4.23e+00 1 6 -24 +1.52e+00,-2.07e+01,-2.06e+01,+8.59e+01 3 -4 Vertex: -4 ID: 0 (X,cT)=+1.20e-01,-3.00e-01,+5.00e-02,+4.00e-03 I: 1 6 -24 +1.52e+00,-2.07e+01,-2.06e+01,+8.59e+01 3 -4 O: 2 7 1 -2.44e+00,+2.88e+01,+6.08e+00,+2.96e+01 1 8 -2 +3.96e+00,-4.95e+01,-2.67e+01,+5.64e+01 1 ________________________________________________________________________________ HepMC::Version 2.06.09 HepMC::IO_GenEvent-START_EVENT_LISTING E 1 0 0.0000000000000000e+00 0.0000000000000000e+00 0.0000000000000000e+00 20 -3 4 0 0 0 0 U GEV MM F 0 0 0.0000000000000000e+00 0.0000000000000000e+00 0.0000000000000000e+00 0.0000000000000000e+00 0.0000000000000000e+00 0 0 V -1 0 0.0000000000000000e+00 0.0000000000000000e+00 0.0000000000000000e+00 0.0000000000000000e+00 1 1 0 P 1 2212 0.0000000000000000e+00 0.0000000000000000e+00 7.0000000000000000e+03 7.0000000000000000e+03 0.0000000000000000e+00 3 0.0000000000000000e+00 0.0000000000000000e+00 -1 0 P 3 1 7.5000000000000000e-01 -1.5690000000000000e+00 3.2191000000000003e+01 3.2238000000000000e+01 0.0000000000000000e+00 3 0.0000000000000000e+00 0.0000000000000000e+00 -3 0 V -2 0 0.0000000000000000e+00 0.0000000000000000e+00 0.0000000000000000e+00 0.0000000000000000e+00 1 1 0 P 2 2212 0.0000000000000000e+00 0.0000000000000000e+00 -7.0000000000000000e+03 7.0000000000000000e+03 0.0000000000000000e+00 3 0.0000000000000000e+00 0.0000000000000000e+00 -2 0 P 4 -2 -3.0470000000000002e+00 -1.9000000000000000e+01 -5.4628999999999998e+01 5.7920000000000002e+01 0.0000000000000000e+00 3 0.0000000000000000e+00 0.0000000000000000e+00 -3 0 V -3 0 0.0000000000000000e+00 0.0000000000000000e+00 0.0000000000000000e+00 0.0000000000000000e+00 0 2 0 P 5 22 -3.8130000000000002e+00 1.1300000000000000e-01 -1.8330000000000000e+00 4.2329999999999997e+00 0.0000000000000000e+00 1 0.0000000000000000e+00 0.0000000000000000e+00 0 0 P 6 -24 1.5169999999999999e+00 -2.0680000000000000e+01 -2.0605000000000000e+01 8.5924999999999997e+01 0.0000000000000000e+00 3 0.0000000000000000e+00 0.0000000000000000e+00 -4 0 V -4 0 1.2000000000000000e-01 -2.9999999999999999e-01 5.0000000000000003e-02 4.0000000000000001e-03 0 2 0 P 7 1 -2.4449999999999998e+00 2.8815999999999999e+01 6.0819999999999999e+00 2.9552000000000000e+01 0.0000000000000000e+00 1 0.0000000000000000e+00 0.0000000000000000e+00 0 0 P 8 -2 3.9620000000000002e+00 -4.9497999999999998e+01 -2.6687000000000001e+01 5.6372999999999998e+01 0.0000000000000000e+00 1 0.0000000000000000e+00 0.0000000000000000e+00 0 0
type Flow ¶
type Flow struct { Particle *Particle // the particle this flow describes Icode map[int]int // flow patterns as (code_index, icode) }
Flow represents a particle's flow and keeps track of an arbitrary number of flow patterns within a graph (i.e. color flow, charge flow, lepton number flow,...)
Flow patterns are coded with an integer, in the same manner as in Herwig. Note: 0 is NOT allowed as code index nor as flow code since it
is used to indicate null.
This class can be used to keep track of flow patterns within
a graph. An example is color flow. If we have two quarks going through an s-channel gluon to form two more quarks: \q1 /q3 then we can keep track of the color flow with the \_______/ HepMC::Flow class as follows: / g \. /q2 \q4 lets say the color flows from q2-->g-->q3 and q1-->g-->q4 the individual colors are unimportant, but the flow pattern is. We can capture this flow by assigning the first pattern (q2-->g-->q3) a unique (arbitrary) flow code 678 and the second pattern (q1-->g-->q4) flow code 269 ( you can ask HepMC::Flow to choose a unique code for you using Flow::set_unique_icode() ). these codes with the particles as follows: q2->flow().set_icode(1,678); g->flow().set_icode(1,678); q3->flow().set_icode(1,678); q1->flow().set_icode(1,269); g->flow().set_icode(2,269); q4->flow().set_icode(1,269); later on if we wish to know the color partner of q1 we can ask for a list of all particles connected via this code to q1 which do have less than 2 color partners using: vector<GenParticle*> result=q1->dangling_connected_partners(q1->icode(1),1,2); this will return a list containing q1 and q4. vector<GenParticle*> result=q1->connected_partners(q1->icode(1),1,2); would return a list containing q1, g, and q4.
type HeavyIon ¶
type HeavyIon struct { NCollHard int // number of hard scatterings NPartProj int // number of projectile participants NPartTarg int // number of target participants NColl int // number of NN (nucleon-nucleon) collisions NNwColl int // Number of N-Nwounded collisions NwNColl int // Number of Nwounded-N collisons NwNwColl int // Number of Nwounded-Nwounded collisions SpectatorNeutrons int // Number of spectators neutrons SpectatorProtons int // Number of spectators protons ImpactParameter float32 // Impact Parameter(fm) of collision EventPlaneAngle float32 // Azimuthal angle of event plane Eccentricity float32 // eccentricity of participating nucleons in the transverse plane (as in phobos nucl-ex/0510031) SigmaInelNN float32 // nucleon-nucleon inelastic (including diffractive) cross-section }
HeavyIon holds additional information for heavy-ion collisions
type LengthUnit ¶
type LengthUnit int
LengthUnit describes the units of length quantities (mm or cm)
const ( // MM is a Length in mm (default) MM LengthUnit = iota // CM is a Length in cm CM )
func LengthUnitFromString ¶
func LengthUnitFromString(s string) (LengthUnit, error)
LengthUnitFromString creates a LengthUnit value from its string representation
func (LengthUnit) String ¶
func (lu LengthUnit) String() string
type MomentumUnit ¶
type MomentumUnit int
MomentumUnit describes the units of momentum quantities (MeV or GeV)
const ( // MEV is a Momentum in MeV (default) MEV MomentumUnit = iota // GEV is a Momentum in GeV GEV )
func MomentumUnitFromString ¶
func MomentumUnitFromString(s string) (MomentumUnit, error)
MomentumUnitFromString creates a MomentumUnit value from its string representation
func (MomentumUnit) String ¶
func (mu MomentumUnit) String() string
type Particle ¶
type Particle struct { Momentum fmom.PxPyPzE // momentum vector PdgID int64 // id according to PDG convention Status int // status code as defined for HEPEVT Flow Flow // flow of this particle Polarization Polarization // polarization of this particle ProdVertex *Vertex // pointer to production vertex (nil if vacuum or beam) EndVertex *Vertex // pointer to decay vertex (nil if not-decayed) Barcode int // unique identifier in the event GeneratedMass float64 // mass of this particle when it was generated }
Particle represents a generator particle within an event coming in/out of a vertex
Particle is the basic building block of the event record
type Particles ¶
type Particles []*Particle
Particles is a []*Particle sorted by increasing-barcodes
type PdfInfo ¶
type PdfInfo struct { ID1 int // flavour code of first parton ID2 int // flavour code of second parton LHAPdf1 int // LHA PDF id of first parton LHAPdf2 int // LHA PDF id of second parton X1 float64 // fraction of beam momentum carried by first parton ("beam side") X2 float64 // fraction of beam momentum carried by second parton ("target side") ScalePDF float64 // Q-scale used in evaluation of PDF's (in GeV) Pdf1 float64 // PDF (id1, x1, Q) Pdf2 float64 // PDF (id2, x2, Q) }
PdfInfo holds informations about the partons distribution functions
type Polarization ¶
type Polarization struct { Theta float64 // polar angle of polarization in radians [0, math.Pi) Phi float64 // azimuthal angle of polarization in radians [0, 2*math.Pi) }
Polarization holds informations about a particle's polarization
type Vertex ¶
type Vertex struct { Position fmom.PxPyPzE // 4-vector of vertex [mm] ParticlesIn []*Particle // all incoming particles ParticlesOut []*Particle // all outgoing particles ID int // vertex id Weights Weights // weights for this vertex Event *Event // pointer to event owning this vertex Barcode int // unique identifier in the event }
Vertex represents a generator vertex within an event A vertex is indirectly (via particle "edges") linked to other
vertices ("nodes") to form a composite "graph"
func (*Vertex) AddParticleIn ¶
AddParticleIn adds a particle to the list of in-coming particles to this vertex
func (*Vertex) AddParticleOut ¶
AddParticleOut adds a particle to the list of out-going particles to this vertex
Source Files ¶
Directories ¶
Path | Synopsis |
---|---|
go-hepmc-dump is a simple command to dump in an almost human-friendly format the content of a hepmc file.
|
go-hepmc-dump is a simple command to dump in an almost human-friendly format the content of a hepmc file. |