Belle II Software development
SVDValidationTTreeSimhit.py
1#!/usr/bin/env python3
2
3
10
11"""
12<header>
13 <description>
14 This module is used for the SVD validation.
15 It gets information about truehits and clusters, saving
16 in a ttree in a ROOT file.
17 </description>
18 <noexecute>SVD validation helper class</noexecute>
19</header>
20"""
21
22import basf2 as b2
23
24# Some ROOT tools
25import ROOT
26from ROOT import Belle2 # noqa: make Belle2 namespace available
27from ROOT import gROOT, addressof
28
29# Define a ROOT struct to hold output data in the TTree
30gROOT.ProcessLine('struct EventDataSimhit {\
31 int sensor_id;\
32 int layer;\
33 int ladder;\
34 int sensor;\
35 int sensor_type;\
36 float simhit_length;\
37 float simhit_energy;\
38 float simhit_dEdx;\
39 };')
40
41from ROOT import EventDataSimhit # noqa
42
43
44class SVDValidationTTreeSimhit(b2.Module):
45 '''class to create sim hit ttree'''
46
47 def __init__(self):
48 """Initialize the module"""
49
50 super().__init__()
51
52 self.file = ROOT.TFile('../SVDValidationTTreeSimhit.root', 'recreate')
53
54 self.tree = ROOT.TTree('tree', 'Event data of SVD validation events')
55
56 self.data = EventDataSimhit()
57
58 # Declare tree branches
59 for key in EventDataSimhit.__dict__:
60 if '__' not in key:
61 formstring = '/F'
62 if isinstance(self.data.__getattribute__(key), int):
63 formstring = '/I'
64 self.tree.Branch(key, addressof(self.data, key), key + formstring)
65
66 def event(self):
67 """Find simhits with a truehit and save needed information"""
68
69 # Start with truehits and use the relation to get the corresponding
70 # simhits
71 svd_truehits = Belle2.PyStoreArray('SVDTrueHits')
72 for truehit in svd_truehits:
73 simhits = truehit.getRelationsTo('SVDSimHits')
74 for simhit in simhits:
75 # Sensor identification
76 sensorID = simhit.getSensorID()
77 self.data.sensor_id = int(sensorID)
78 sensorNum = sensorID.getSensorNumber()
79 self.data.sensor = sensorNum
80 layerNum = sensorID.getLayerNumber()
81 self.data.layer = layerNum
82 ladderNum = sensorID.getLadderNumber()
83 self.data.ladder = ladderNum
84 if (layerNum == 3):
85 sensorType = 1
86 else:
87 if (sensorNum == 1):
88 sensorType = 0
89 else:
90 sensorType = 1
91 self.data.sensor_type = sensorType
92 # We are interested only in the primary particles, not
93 # delta electrons, etc.
94 particle = simhit.getRelatedFrom('MCParticles')
95 if not particle.hasStatus(Belle2.MCParticle.c_PrimaryParticle):
96 continue
97
98 length = (simhit.getPosOut() - simhit.getPosIn()).R()
99 # The deposited energy is the number of electrons multiplied
100 # by the energy required to create an electron-hole pair
101 energy = simhit.getElectrons() * Belle2.Const.ehEnergy
102 self.data.simhit_length = length
103 self.data.simhit_energy = energy
104 if (length == 0):
105 continue
106 self.data.simhit_dEdx = energy / length
107 # A reasonable cut to see a nice Landau distribution
108 # remove this cut to be sensitive to delta rays
109 # if self.data.simhit_dEdx > 0.015:
110 # continue
111 # Fill tree
112 self.file.cd()
113 self.tree.Fill()
114
115 def terminate(self):
116 """Close the output file. """
117 self.file.cd()
118 self.file.Write()
119 self.file.Close()
A (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:72