Belle II Software development
Ntuple Class Reference
Inheritance diagram for Ntuple:

Public Member Functions

def initialize (self)
 
def event (self)
 
def terminate (self)
 

Public Attributes

 file
 output file name
 
 tree
 output tree name
 
 data
 tree strruct
 

Static Public Attributes

dict t0const = {0: 0.0}
 t0 constant per channel per slot
 
ROOT f = ROOT.TFile.Open('t0const_slot' + str(args[1]) + '.root')
 input t0const root file
 
TH2F histTimeCh = TH2F('before T0Cal slot#' + str(args[1]), 'before T0Cal slot#' + str(args[1]), 512, 0, 511, 500, 50, 100)
 scatter plot before t0 calibartion
 
TH2F histCalTimeCh = TH2F('after T0Cal slot#' + str(args[1]), 'after T0Cal slot#' + str(args[1]), 512, 0, 511, 500, 50, 100)
 scatter plot after t0 calibartion
 

Detailed Description

 t0const ntpule infomation 

Definition at line 46 of file t0CalTime.py.

Member Function Documentation

◆ event()

def event (   self)
 Event processor: fill the tree and scatter plots 

Definition at line 85 of file t0CalTime.py.

85 def event(self):
86 ''' Event processor: fill the tree and scatter plots '''
87
88 digits = Belle2.PyStoreArray('TOPDigits')
89 for digit in digits:
90 if digit.getModuleID() == int(args[1]):
91 self.data.slot = digit.getModuleID()
92 self.data.channel = digit.getChannel()
93 self.data.pixel = digit.getPixelID()
94 self.data.time = digit.getTime()
95 self.data.caltime = digit.getTime() + t0const.get(digit.getChannel())
96 self.data.height = digit.getPulseHeight()
97 self.data.width = digit.getPulseWidth()
98 self.data.quality = digit.getHitQuality()
99 self.data.rawtime = digit.getRawTime()
100 self.data.firstwindow = digit.getFirstWindow()
101 self.data.pmtpixel = digit.getPMTPixel()
102 self.data.pmt = digit.getPMTNumber()
103
104 self.file.cd()
105 self.tree.Fill()
106
107 self.histTimeCh.Fill(digit.getChannel(), digit.getTime())
108 self.histCalTimeCh.Fill(digit.getChannel(), digit.getTime() - t0const.get(digit.getChannel()))
109
A (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:72

◆ initialize()

def initialize (   self)
 Initialize the Module: output root file 

Definition at line 68 of file t0CalTime.py.

68 def initialize(self):
69 ''' Initialize the Module: output root file '''
70
71
72 self.file = ROOT.TFile('t0CalTime_slot' + str(args[1]) + '.root', 'recreate')
73
74 self.tree = ROOT.TTree('laser', '')
75
76 self.data = TreeStruct()
77
78 for key in TreeStruct.__dict__.keys():
79 if '__' not in key:
80 formstring = '/F'
81 if isinstance(self.data.__getattribute__(key), int):
82 formstring = '/I'
83 self.tree.Branch(key, addressof(self.data, key), key + formstring)
84

◆ terminate()

def terminate (   self)
 Write the file 

Definition at line 110 of file t0CalTime.py.

110 def terminate(self):
111 ''' Write the file '''
112
113 self.file.cd()
114 self.file.Write()
115 self.histTimeCh.Write()
116 self.histCalTimeCh.Write()
117 self.file.Close()
118
119
120# Create path

Member Data Documentation

◆ data

data

tree strruct

Definition at line 76 of file t0CalTime.py.

◆ f

ROOT f = ROOT.TFile.Open('t0const_slot' + str(args[1]) + '.root')
static

input t0const root file

Definition at line 54 of file t0CalTime.py.

◆ file

file

output file name

Definition at line 72 of file t0CalTime.py.

◆ histCalTimeCh

TH2F histCalTimeCh = TH2F('after T0Cal slot#' + str(args[1]), 'after T0Cal slot#' + str(args[1]), 512, 0, 511, 500, 50, 100)
static

scatter plot after t0 calibartion

Definition at line 61 of file t0CalTime.py.

◆ histTimeCh

TH2F histTimeCh = TH2F('before T0Cal slot#' + str(args[1]), 'before T0Cal slot#' + str(args[1]), 512, 0, 511, 500, 50, 100)
static

scatter plot before t0 calibartion

Definition at line 59 of file t0CalTime.py.

◆ t0const

dict t0const = {0: 0.0}
static

t0 constant per channel per slot

Definition at line 51 of file t0CalTime.py.

◆ tree

tree

output tree name

Definition at line 74 of file t0CalTime.py.


The documentation for this class was generated from the following file: