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