110ReadTree.py
Go to the documentation of this file.
1 #! /usr/bin/env python
2 #
3 # Read almost every field in the event tree.
4 #
5 
6 from __future__ import print_function
7 import glob, os, re, sys, getopt
8 from ROOT import *
9 gSystem.Load("libedepsim_io.so")
10 from ROOT import TG4Event
11 
12 # Print the fields in a TG4PrimaryParticle object
13 def printPrimaryParticle(depth, primaryParticle):
14  print(depth,"Class: ", primaryParticle.ClassName())
15  print(depth,"Track Id:", primaryParticle.GetTrackId())
16  print(depth,"Name:", primaryParticle.GetName())
17  print(depth,"PDG Code:",primaryParticle.GetPDGCode())
18  print(depth,"Momentum:",primaryParticle.GetMomentum().X(),
19  primaryParticle.GetMomentum().Y(),
20  primaryParticle.GetMomentum().Z(),
21  primaryParticle.GetMomentum().E(),
22  primaryParticle.GetMomentum().P(),
23  primaryParticle.GetMomentum().M())
24 
25 # Print the fields in an TG4PrimaryVertex object
26 def printPrimaryVertex(depth, primaryVertex):
27  print(depth,"Class: ", primaryVertex.ClassName())
28  print(depth,"Position:", primaryVertex.GetPosition().X(),
29  primaryVertex.GetPosition().Y(),
30  primaryVertex.GetPosition().Z(),
31  primaryVertex.GetPosition().T())
32  print(depth,"Generator:",primaryVertex.GetGeneratorName())
33  print(depth,"Reaction:",primaryVertex.GetReaction())
34  print(depth,"Filename:",primaryVertex.GetFilename())
35  print(depth,"InteractionNumber:",primaryVertex.GetInteractionNumber())
36  depth = depth + ".."
37  for infoVertex in primaryVertex.Informational:
38  printPrimaryVertex(depth,infoVertex)
39  for primaryParticle in primaryVertex.Particles:
40  printPrimaryParticle(depth,primaryParticle)
41 
42 # Print the fields in a TG4TrajectoryPoint object
43 def printTrajectoryPoint(depth, trajectoryPoint):
44  print(depth,"Class: ", trajectoryPoint.ClassName())
45  print(depth,"Position:", trajectoryPoint.GetPosition().X(),
46  trajectoryPoint.GetPosition().Y(),
47  trajectoryPoint.GetPosition().Z(),
48  trajectoryPoint.GetPosition().T())
49  print(depth,"Momentum:", trajectoryPoint.GetMomentum().X(),
50  trajectoryPoint.GetMomentum().Y(),
51  trajectoryPoint.GetMomentum().Z(),
52  trajectoryPoint.GetMomentum().Mag())
53  print(depth,"Process",trajectoryPoint.GetProcess())
54  print(depth,"Subprocess",trajectoryPoint.GetSubprocess())
55 
56 # Print the fields in a TG4Trajectory object
57 def printTrajectory(depth, trajectory):
58  print(depth,"Class: ", trajectory.ClassName())
59  depth = depth + ".."
60  print(depth,"Track Id/Parent Id:",
61  trajectory.GetTrackId(),
62  trajectory.GetParentId())
63  print(depth,"Name:",trajectory.GetName())
64  print(depth,"PDG Code",trajectory.GetPDGCode())
65  print(depth,"Initial Momentum:",trajectory.GetInitialMomentum().X(),
66  trajectory.GetInitialMomentum().Y(),
67  trajectory.GetInitialMomentum().Z(),
68  trajectory.GetInitialMomentum().E(),
69  trajectory.GetInitialMomentum().P(),
70  trajectory.GetInitialMomentum().M())
71  for trajectoryPoint in trajectory.Points:
72  printTrajectoryPoint(depth,trajectoryPoint)
73 
74 # Print the fields in a TG4HitSegment object
75 def printHitSegment(depth, hitSegment):
76  print(depth,"Class: ", hitSegment.ClassName())
77  print(depth,"Primary Id:", hitSegment.GetPrimaryId());
78  print(depth,"Energy Deposit:",hitSegment.GetEnergyDeposit())
79  print(depth,"Secondary Deposit:", hitSegment.GetSecondaryDeposit())
80  print(depth,"Track Length:",hitSegment.GetTrackLength())
81  print(depth,"Start:", hitSegment.GetStart().X(),
82  hitSegment.GetStart().Y(),
83  hitSegment.GetStart().Z(),
84  hitSegment.GetStart().T())
85  print(depth,"Stop:", hitSegment.GetStop().X(),
86  hitSegment.GetStop().Y(),
87  hitSegment.GetStop().Z(),
88  hitSegment.GetStop().T())
89  print(depth,"Contributor:", [contributor for contributor in hitSegment.Contrib])
90 
91 # Print the fields in a single element of the SegmentDetectors map.
92 # The container name is the key, and the hitSegments is the value (a
93 # vector of TG4HitSegment objects).
94 def printSegmentContainer(depth, containerName, hitSegments):
95  print(depth,"Detector: ", containerName, hitSegments.size())
96  depth = depth + ".."
97  for hitSegment in hitSegments: printHitSegment(depth, hitSegment)
98 
99 # Read a file and dump it.
100 def main(argv=None):
101  if argv is None:
102  argv = sys.argv
103 
104  # The input file is generated in a previous test (100TestTree.sh).
105  inputFile=TFile("100TestTree.root")
106 
107  # Get the input tree out of the file.
108  inputTree=inputFile.Get("EDepSimEvents")
109  print("Class:",inputTree.ClassName())
110 
111  # Attach a brach to the events.
112  event = TG4Event()
113  inputTree.SetBranchAddress("Event",event)
114 
115  # Read all of the events.
116  entries=inputTree.GetEntriesFast()
117  for jentry in xrange(entries):
118  nb = inputTree.GetEntry(jentry)
119  if nb<=0: continue
120  print("Class: ", event.ClassName())
121  print("Event number:", event.EventId)
122  # Dump the primary vertices
123  for primaryVertex in event.Primaries:
124  printPrimaryVertex("PP", primaryVertex)
125  # Dump the trajectories
126  for trajectory in event.Trajectories: printTrajectory("TT",trajectory)
127  # Dump the segment containers
128  print("Number of segment containers:", event.SegmentDetectors.size())
129  for containerName, hitSegments in event.SegmentDetectors:
130  printSegmentContainer("HH", containerName, hitSegments)
131 
132 if __name__ == "__main__":
133  sys.exit(main())
134 
def main(argv=None)
Definition: 110ReadTree.py:100
std::pair< float, std::string > P
def printTrajectoryPoint(depth, trajectoryPoint)
Definition: 110ReadTree.py:43
def printHitSegment(depth, hitSegment)
Definition: 110ReadTree.py:75
def printSegmentContainer(depth, containerName, hitSegments)
Definition: 110ReadTree.py:94
def printPrimaryParticle(depth, primaryParticle)
Definition: 110ReadTree.py:13
def printPrimaryVertex(depth, primaryVertex)
Definition: 110ReadTree.py:26
def printTrajectory(depth, trajectory)
Definition: 110ReadTree.py:57