Belle II Software development
channelMasking_hltdqm.py
1#!/usr/bin/env python
2
3
10
11import sys
12import os
13import glob
14import logging
15from ROOT import TFile, TH1F
16
17
18# checks if the run can be analyzed based on the HLT histograms, returns -1 if not and the number of events otherwise
19def checkIfRunUsable(file):
20 # consider only runs identified as physics runs (possible since experiment 8)
21 h = file.FindObject('DQMInfo/rtype')
22 if not h:
23 logging.warning(fileName + ' ... histogram runtype not found')
24 return -1
25 if not h.GetTitle() == "physics":
26 logging.info(fileName + ' ... not a physics run, skipping')
27 return -1
28
29 # check if we have enough events for meaningful masking
30 h = file.FindObject('TOP/good_hits_per_event1')
31 if not h:
32 logging.warning(fileName + ' ... histogram good_hits_per_event1 not found')
33 return -1
34 nev = int(h.GetEntries())
35 if nev < 10000:
36 logging.warning(fileName + 'run =' + str(run) + 'events =' + str(nev) + ' ... skipped, not enough events')
37 return -1
38
39 # check if we have the necessary histogram to determine masking
40 h = file.FindObject('TOP/good_hits')
41 if not h:
42 logging.warning(fileName + ' ... histogram good_hits not found')
43 return -1
44 if h.GetEntries() == 0:
45 logging.warning(fileName + ' ... histogram good_hits has no entries ... skipped')
46 return -1
47 return nev
48
49
50# creates a file with histograms for masked channels and returns the total number of masked channels in a run
51def makeChannelMasks(file, outFileName):
52
53 masks = [TH1F('slot_' + str(slot), 'Channel mask for slot ' + str(slot),
54 512, 0.0, 512.0) for slot in range(1, 17)]
55 masked = 0
56
57 # dead/hot channels: <10% or >1000% of the average in this slot
58 for slot in range(1, 17):
59 h = file.FindObject('TOP/good_channel_hits_' + str(slot))
60 if not h:
61 logging.error('no good_channel_hits found for slot'+str(slot))
62 continue
63 # calculate average number of hits per channel
64 mean = 0
65 n = 0
66 for chan in range(h.GetNbinsX()):
67 y = h.GetBinContent(chan+1)
68 if y > 0:
69 mean += y
70 n += 1
71 if n > 0:
72 mean /= n
73 # define threshold for dead and hot
74 deadCut = mean / 10
75 hotCut = mean * 10
76 for chan in range(h.GetNbinsX()):
77 y = h.GetBinContent(chan+1)
78 if y <= deadCut:
79 masks[slot-1].SetBinContent(chan+1, 1) # 1 = dead
80 masked += 1
81 elif y > hotCut:
82 masks[slot-1].SetBinContent(chan+1, 2) # 2 = noisy
83 masked += 1
84
85 # asics with window corruption
86 for slot in range(1, 17):
87 h = file.FindObject('TOP/window_vs_asic_' + str(slot))
88 if not h:
89 logging.error('Error: no window_vs_asic found for slot' + str(slot))
90 continue
91 h0 = h.ProjectionX()
92 h1 = h.ProjectionX('_tmp', 222, 245)
93 for asic in range(h.GetNbinsX()):
94 if h0.GetBinContent(asic+1) > 0:
95 r = 1 - h1.GetBinContent(asic+1) / h0.GetBinContent(asic+1)
96 if r > 0.20:
97 for chan in range(8):
98 masks[slot-1].SetBinContent(asic*8+chan+1, 2) # mark as noisy
99 masked += 1
100
101 # save to file
102 if masked < 3000: # don't save for nonsense
103 outfile = TFile(outFileName, 'recreate')
104 for h in masks:
105 h.Write()
106 outfile.Close()
107 return masked
108
109
110experiment = 10
111experimentstring = f"{experiment:04d}"
112outdir = 'masks'
113if not os.path.exists(outdir):
114 os.makedirs(outdir)
115fileNames = sorted(glob.glob('/group/belle2/phase3/dqm/dqmsrv1/e'+experimentstring+'/dqmhisto/hltdqm*.root'))
116
117logging.basicConfig(level=logging.INFO, filename="channelmasking.log")
118logging.info("Starting channelmasking from HLT histograms")
119logging.info("Experiment: "+str(experiment))
120
121numFiles = len(fileNames)
122if numFiles == 0:
123 logging.error('No files found, exiting')
124 sys.exit()
125
126for fileName in fileNames:
127 run = ((fileName.split('/')[-1]).split('r')[1]).split('.')[0] # exp7-8
128 outFileName = outdir + '/channelMask_e' + experimentstring + '_r' + run + '.root'
129
130 # skip run if we analyzed this run already
131 if os.path.exists(outFileName) or os.path.exists(outFileName.replace('masks/', 'masks/imported/')):
132 logging.debug('Output file exists for run ='+str(run)+', skipping')
133 continue
134
135 # open file and check if it is a good physics run
136 file = TFile(fileName)
137 if not file:
138 logging.error(fileName + ' ... cannot open')
139 continue
140 file.ReadAll()
141 nev = checkIfRunUsable(file)
142 if nev == -1:
143 file.Close()
144 continue
145
146 masked = makeChannelMasks(file, outFileName)
147 logging.info(fileName + 'run =' + str(run) + 'events =' + str(nev) + 'masked channels =' + str(masked))
148 if masked > 3000:
149 logging.critical("Result looks completely wrong, over 3000 masked channels. Inspect run "+str(run))
150 file.Close()