CMS 3D CMS Logo

PhaseITreeProducer.py
Go to the documentation of this file.
1 #!/usr/bin/env python
2 
3 import sys
4 from ROOT import *
5 from array import array
6 from copy import deepcopy
7 
8 gROOT.SetBatch() # don't pop up canvases
9 
10 #get data files
11 
12 def getFileInPath(rfile):
13  import os
14  for dir in os.environ['CMSSW_SEARCH_PATH'].split(":"):
15  if os.path.exists(os.path.join(dir,rfile)): return os.path.join(dir,rfile)
16  return None
17 
18 # Default values
19 inputFileName = "DQM.root"
20 outputFileName = "DQMTree.root"
21 #detIDsFileName = "DATA/detids.dat"
22 detIDsFileName = getFileInPath('DQM/SiStripMonitorClient/data/detids.dat')
23 
25 
26  ############################################################################
27 
28  def __TraverseDirTree(self, dir):
29 
30  currPath = (dir.GetPath().split(":/"))[1]
31 
32  for obj in dir.GetListOfKeys():
33  if not obj.IsFolder():
34  if obj.ReadObjectAny(TClass.GetClass("TH2")):
35  th2 = deepcopy(obj.ReadObj())
36  name = th2.GetName()
37  if 6 < th2.GetNbinsX() < 10 and name.find("per") != -1 and name.find("Lumisection") == -1: #take only module lvl plots
38  print(''.join([dir.GetPath(), '/', name]))
39 
40  # fix when there are plots starting with the same strings in different directories
41  prefix = ""
42  for i in self.dirs:
43  if currPath.startswith(i):
44  prefix = self.dirsAliases[i]
45  break
46  # print(currPath, prefix)
47  th2.SetName(prefix + th2.GetName())
48  self.listOfNumHistograms.append(th2)
49  else:
50  self.__TraverseDirTree(obj.ReadObj())
51 
52  def __GetPartStr(self, isXlowerThanZero, isYlowerThanZero):
53  if isXlowerThanZero and isYlowerThanZero:
54  return "mO"
55  if isXlowerThanZero and isYlowerThanZero == False:
56  return "mI"
57  if isXlowerThanZero == False and isYlowerThanZero:
58  return "pO"
59  if isXlowerThanZero == False and isYlowerThanZero == False:
60  return "pI"
61 
62  def __GetBarrelSector(self, layer, signedLadder, signedModule): #adapted from PixelBarrelName
63  theLadder = abs(signedLadder)
64  theModule = abs(signedModule)
65 
66  sector = 0
67 
68  if layer == 1:
69 
70  if theLadder == 1:
71  if theModule >= 2:
72  return 1
73  else:
74  return 2
75  if theLadder == 2:
76  if theModule >= 3:
77  return 2
78  else:
79  return 3
80  if theLadder == 3:
81  if theModule >= 4:
82  return 3
83  else:
84  return 4
85  if theLadder == 4:
86  if theModule >= 2:
87  return 5
88  else:
89  return 6
90  if theLadder == 5:
91  if theModule >= 3:
92  return 6
93  else:
94  return 7
95  if theLadder == 6:
96  if theModule >= 4:
97  return 7
98  else:
99  return 8
100  # here is used simplified form of assignment, see source file for reference
101  elif layer == 2:
102  i = theLadder // 5
103  sector = i * 3
104  shortLadder = theLadder - 5 * i
105  for i in range(0, shortLadder, 2):
106  sector = sector + 1
107  return sector
108  elif layer == 3:
109  sector = 1
110  for i in range(2, theLadder, 3):
111  if (i + 1) % 3 == 0:
112  sector = sector + 1
113  return sector
114  elif layer == 4:
115  sector = (theLadder + 3) // 4
116  return sector
117 
118  def __BuildOnlineBarrelName(self, signedModule, signedLadder, layer): #in Phase1 it is assumed that there are only full modules
119  thePart = self.__GetPartStr(signedModule < 0, signedLadder < 0)
120  theSector = str(self.__GetBarrelSector(layer, signedLadder, signedModule))
121  return "BPix_B" + thePart + "_SEC" + theSector + "_LYR" + str(layer) + "_LDR" + str(abs(signedLadder)) + "F_MOD" + str(abs(signedModule))
122 
123  def __BuildOnlineDiskName(self, signedDisk, signedBlade, panel, ring):
124  thePart = self.__GetPartStr(signedDisk < 0, signedBlade < 0)
125  return "FPix_B" + thePart + "_D" + str(abs(signedDisk)) + "_BLD" + str(abs(signedBlade)) + "_PNL" + str(panel) + "_RNG" + str(ring)
126 
127  def __GroupHistograms(self):
128  currentGroupName = ""
129  groupOfHists = []
131 
132  ##### GROUP ALL LAYERS/RINGS HAVING SIMILAR INFORMATION
133  for obj in self.listOfNumHistograms:
134  objName = obj.GetName()
135  objNameSplit = objName.split("_")
136  objNameCollected = ''.join(objNameSplit[0:-1])
137  if objNameCollected != currentGroupName:
138  if len(groupOfHists):
139  self.groupedHistograms.append(groupOfHists)
140  groupOfHists = []
141 
142  currentGroupName = objNameCollected
143  groupOfHists.append(obj)
144  self.groupedHistograms.append(groupOfHists) #the last group
145 
146  def __CreateDummyStructAsStr(self, dicData):
147  str = "struct MyStruct{"
148  str = str + "Int_t key;"
149  leafStr = "key/I"
150  for k in dicData:
151  str = str + "Float_t " + k + ";"
152  leafStr = leafStr + ":" + k + "/F"
153  str = str + "};"
154  return str, leafStr
155 
156  ############################################################################
157 
158  def __init__(self, inputDQMName, outputDQMName, modDicName):
159  self.inputFileName = inputDQMName
160  self.outputFileName = outputDQMName
161  self.detIDsFileName = modDicName
162 
163  index = self.inputFileName.find('R000')
164  runNumber = self.inputFileName[index+4:index+10]
165 
166  self.dirs = ['DQMData/Run %s/PixelPhase1/Run summary/Phase1_MechanicalView' % (runNumber),
167  'DQMData/Run %s/PixelPhase1/Run summary/Tracks' % (runNumber)]
168  self.dirsAliases = {self.dirs[0]:"", self.dirs[1]: "T"}
169 
170  self.inputFile = TFile(self.inputFileName)
172  self.availableNames = []
173 
174  self.maxLadderToLayer = {6:1, 14:2, 22:3, 32:4}
175  self.maxBladeToRing = {11:1, 17:2}
176 
177  self.internalData = {}
178 
179  if self.inputFile.IsOpen():
180  print("%s opened successfully!" % (self.inputFileName))
181  #Get all neeeded histograms
182  for dir in self.dirs:
183  self.__TraverseDirTree(self.inputFile.Get(dir))
184  print("Histograms to read %d" % (len(self.listOfNumHistograms)))
185 
186  self.detDict = {}
187 
188  with open(self.detIDsFileName, "r") as detIDs: # create dictionary online -> rawid
189  for entry in detIDs:
190  items = entry.replace("\n", " ").split(" ")
191  self.detDict.update({items[1] : int(items[0])})
192  # init internal data structure
193  self.internalData.update({int(items[0]) : {}})
194 
195  self.__GroupHistograms()
196 
197  else:
198  print("Unable to open file %s" % (self.inputFileName))
199 
200  def ReadHistograms(self):
201  if self.inputFile.IsOpen():
202  for group in self.groupedHistograms:
203  # name = ''.join(group[0].GetName().split("_")[0:-1])
204  name = ''.join(group[0].GetName().split("_per_")[0])
205  self.availableNames.append(name)
206  # print(name)
207  for obj in group:
208  nbinsX = obj.GetNbinsX()
209  nbinsY = obj.GetNbinsY()
210 
211  if nbinsX == 9: # BARREL
212  maxX = nbinsX // 2
213  maxY = nbinsY // 2
214 
215  for x in range(-maxX, maxX + 1):
216  if x == 0:
217  continue
218  for y in range(-maxY, maxY + 1, 1):
219  if y == 0:
220  continue
221  onlineName = self.__BuildOnlineBarrelName(x, y, self.maxLadderToLayer[maxY])
222  self.internalData[self.detDict[onlineName]].update({name : obj.GetBinContent(x + maxX + 1, y + maxY + 1)})
223 
224  elif nbinsX == 7: # FORWARD
225  maxX = nbinsX // 2
226  maxY = nbinsY // 4
227 
228  for x in range(-maxX, maxX + 1):
229  if x == 0:
230  continue
231  for y in range(-maxY, maxY + 1):
232  if int(y) == 0:
233  continue
234  for panel in range(1, 3):
235  onlineName = self.__BuildOnlineDiskName(x, y, panel, self.maxBladeToRing[maxY])
236  self.internalData[self.detDict[onlineName]].update({name : obj.GetBinContent(x + maxX + 1, (y + maxY) * 2 + panel)})
237  else:
238  print("Unrecognized plot")
239  else:
240  print("Histograms saved to internal data structure")
241  def CreateTree(self): # too much complication lvl, use CreateTree2
242  if len(self.internalData):
243  # <- TOTAL WORKAROUND: START -> #
244  # SEE: http://wlav.web.cern.ch/wlav/pyroot/tpytree.html
245 
246  leafStr = ""
247  for key in self.internalData:
248  s, leafStr = self.__CreateDummyStructAsStr(self.internalData[key])
249  print(s)
250  print(leafStr)
251  gROOT.ProcessLine(s)
252  break #all modules are assumed to have the same set of measured parameters
253 
254  from ROOT import MyStruct
255  ms = MyStruct()
256  # <- TOTAL WORKAROUND: END -> #
257 
258  self.outputFile = TFile(self.outputFileName, "recreate")
259  tree = TTree("tree", "readData")
260  tree.Branch("b", ms, leafStr)
261 
262  tree.SetBranchAddress("b", tree)
263 
264  for key in self.internalData:
265  setattr(ms, "key", key)
266  for d in self.internalData[key]:
267  setattr(ms, d, (self.internalData[key])[d])
268  tree.Fill()
269 
270  # break
271  tree.Write()
272  self.outputFile.Close()
273 
274  print("Data saved as TTree object")
275 
276  def CreateTree2(self): # use for TkCommissioner
277  if len(self.internalData):
278 
279  self.outputFile = TFile(self.outputFileName, "recreate")
280  tree = TTree("tree", "readData")
281 
282  key = array("i", [0])
283  tree.Branch("detid", key, "detid/i")
284 
285  dat = {}
286  for k in self.internalData:
287  for d in self.internalData[k]:
288  dat.update({d : array("f", [0])})
289  tree.Branch(d, dat[d], d + "/F")
290  break
291 
292  for k in self.internalData:
293  key[0] = k
294  for d in self.internalData[k]:
295  (dat[d])[0] = (self.internalData[k])[d]
296  tree.Fill()
297 
298  tree.Write()
299  self.outputFile.Close()
300 
301  print("Data saved as TTree object")
302 
303  def DumpData(self):
304  for key in self.internalData:
305  print("#"*20)
306  print(key)
307  module = self.internalData[key]
308  for d in module:
309  print(d, module[d])
310 
311  for i in self.availableNames:
312  print(i)
313  print(len(self.availableNames))
314 
315  def __del__(self):
316  if self.inputFile.IsOpen():
317  self.inputFile.Close()
318 
319 #--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--
320 for i in range(1, len(sys.argv), 1):
321  if i == 1:
322  inputFileName = sys.argv[i]
323  elif i == 2:
324 # detIDsFileName = sys.argv[i]
325 # elif i == 3:
326  outputFileName = sys.argv[i]
327 
328 #readerObj = ModuleLvlValuesReader(inputFileName, outputFileName)
329 readerObj = ModuleLvlValuesReader(inputFileName, outputFileName, detIDsFileName)
330 readerObj.ReadHistograms()
331 readerObj.CreateTree2()
332 # readerObj.DumpData()
S & print(S &os, JobReport::InputFile const &f)
Definition: JobReport.cc:65
def __init__(self, inputDQMName, outputDQMName, modDicName)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
def __GetBarrelSector(self, layer, signedLadder, signedModule)
static std::string join(char **cmd)
Definition: RemoteFile.cc:18
def __GetPartStr(self, isXlowerThanZero, isYlowerThanZero)
#define update(a, b)
def __BuildOnlineDiskName(self, signedDisk, signedBlade, panel, ring)
double split
Definition: MVATrainer.cc:139
def __BuildOnlineBarrelName(self, signedModule, signedLadder, layer)