Belle II Software  release-06-01-15
validatePhase2Cosmics.py
1 #!/usr/bin/env python3
2 # -*- coding: utf-8 -*-
3 
4 
11 
12 # *****************************************************************************
13 
14 # title : 3_ValidationCosmics.py
15 # description : Validation of cosmic tracks in phase II
16 
17 # *****************************************************************************
18 
19 import math
20 import basf2 as b2
21 
22 import ROOT
23 from ROOT import Belle2
24 from ROOT import gROOT, AddressOf
25 import reconstruction as reco
26 
27 gROOT.ProcessLine('struct TrackData {\
28  float chi2;\
29  float ndf;\
30  float chiSquaredOverNdf;\
31  float momentum;\
32 };')
33 
34 gROOT.ProcessLine('struct DEDXData {\
35  float acdCount;\
36  float dE;\
37  float path;\
38  float momentumCDC;\
39  int numberOfHits;\
40  float distanceKLM;\
41 };')
42 
43 from ROOT import TrackData, DEDXData # noqa
44 
45 
46 class CosmicAnalysis(b2.Module):
47 
48  """A module to analyse residuals in overlaps of ladders."""
49 
50  def __init__(self):
51  """Initialize the module"""
52 
53  super(CosmicAnalysis, self).__init__()
54 
55 
56  self.rootfilerootfile = ROOT.TFile('cosmicAnalysis.root', 'recreate')
57 
58  self.tree_tracktree_track = ROOT.TTree('track', '')
59 
60  self.tree_DEDXtree_DEDX = ROOT.TTree('dedx', '')
61 
62  ROOT.gStyle.Reset()
63  ROOT.gStyle.SetCanvasColor(0)
64  ROOT.gStyle.SetStatBorderSize(1)
65  ROOT.gStyle.SetStatColor(0)
66  ROOT.gStyle.SetTitleColor(1)
67  ROOT.gStyle.SetCanvasColor(0)
68  ROOT.gStyle.SetPadColor(0)
69  ROOT.gStyle.SetPadBorderMode(0)
70  ROOT.gStyle.SetCanvasBorderMode(0)
71  ROOT.gStyle.SetFrameBorderMode(0)
72  ROOT.gStyle.SetOptStat(0)
73 
74 
75  self.trackDatatrackData = TrackData()
76  # Declare tree branches
77  for key in TrackData.__dict__:
78  if '__' not in key:
79  formstring = '/F'
80  if isinstance(self.trackDatatrackData.__getattribute__(key), int):
81  formstring = '/I'
82  self.tree_tracktree_track.Branch(key, AddressOf(self.trackDatatrackData, key), key + formstring)
83 
84  self.dedxDatadedxData = DEDXData()
85  # Declare tree branches
86  for key in DEDXData.__dict__:
87  if '__' not in key:
88  formstring = '/F'
89  if isinstance(self.dedxDatadedxData.__getattribute__(key), int):
90  formstring = '/I'
91  self.tree_DEDXtree_DEDX.Branch(key, AddressOf(self.dedxDatadedxData, key), key + formstring)
92 
93  self.TotalNumberOfHitsTotalNumberOfHits = ROOT.TH1F('TotalNumberOfHits', '', 6, 0.5, 6.5)
94 
95  self.HitsVsLayerHitsVsLayer = ROOT.TH2F('HitsVsLayer', '', 6, 0.5, 6.5, 6, 0.5, 6.5)
96 
97  self.HitsVsSensorHitsVsSensor = ROOT.TH2F('HitsVsSensor', '', 6, 0.5, 6.5, 5, 0.5, 5.5)
98 
99  self.LayerVsSensorLayerVsSensor = ROOT.TH2F('LayerVsSensor', '', 6, 0.5, 6.5, 5, 0.5, 5.5)
100 
101  self.PXDClusterSizePXDClusterSize = ROOT.TH1F('PXDClusterSize', '', 20, 0.5, 20.5)
102 
103  self.Chi2Chi2 = ROOT.TH1F('Chi2', '', 300, 0.0, 500)
104 
105  self.NDFNDF = ROOT.TH1F('NDF', '', 200, 0.0, 200)
106 
107  self.Chi2OverNDFChi2OverNDF = ROOT.TH1F('Chi2OverNDF', '', 300, 0.0, 5)
108 
109  self.MomentumMomentum = ROOT.TH1F('Momentum', '', 500, 0.0, 1000)
110 
111  self.ADCCountOverNumberOfHitsInCDCADCCountOverNumberOfHitsInCDC = ROOT.TH1F('ADCCountOverNumberOfHitsInCDC', '', 200, 0.0, 300)
112 
113  self.ADCCountOverNumberOfHitsInCDCVsMomentumADCCountOverNumberOfHitsInCDCVsMomentum = ROOT.TProfile(
114  'ADCCountOverNumberOfHitsInCDCVsMomentum', '', 200, 0.0, 300, 0.0, 1000)
115 
116  self.MomentumVsADCCountOverNumberOfHitsInCDCMomentumVsADCCountOverNumberOfHitsInCDC = ROOT.TProfile(
117  'MomentumVsADCCountOverNumberOfHitsInCDC', '', 100, 0.0, 100, 0.0, 1000)
118 
119  def beginRun(self):
120  """Do nothing"""
121 
122  def event(self):
123  """ Fill histograms """
124 
125  # Study dEdx (CDC) as prediction of momentum
126  cdcDedxTracks = Belle2.PyStoreArray('CDCDedxTracks')
127  nCDCDedxTracks = cdcDedxTracks.getEntries()
128  for DEDXTrack_index in range(nCDCDedxTracks):
129  dedxTrack = cdcDedxTracks[DEDXTrack_index]
130  self.dedxData.momentumCDC = dedxTrack.getMomentum()
131  self.dedxData.acdCount = 0
132  self.dedxData.dE = 0
133  self.dedxData.path = 0
134  self.dedxData.numberOfHits = dedxTrack.size()
135  KLMClusters = Belle2.PyStoreArray('KLMClusters')
136  nKLMClusters = KLMClusters.getEntries()
137  if nKLMClusters == 2:
138  for nHit in range(dedxTrack.size()):
139  self.dedxData.acdCount = self.dedxData.acdCount + dedxTrack.getADCCount(nHit)
140  self.dedxData.dE = self.dedxData.dE + dedxTrack.getDE(nHit)
141  self.dedxData.path = self.dedxData.path + dedxTrack.getPath(nHit)
142  self.dedxData.distanceKLM = math.sqrt(
143  math.pow(
144  KLMClusters[0].getClusterPosition().Z() -
145  KLMClusters[1].getClusterPosition().Z(),
146  2))
147 
148  self.ADCCountOverNumberOfHitsInCDC.Fill(self.dedxData.acdCount / dedxTrack.size())
149  self.MomentumVsADCCountOverNumberOfHitsInCDC.Fill(
150  self.dedxData.momentumCDC, self.dedxData.acdCount / dedxTrack.size())
151  self.ADCCountOverNumberOfHitsInCDCVsMomentum.Fill(
152  self.dedxData.acdCount / dedxTrack.size(), self.dedxData.momentumCDC)
153 
154  self.rootfile.cd()
155  self.tree_DEDX.Fill()
156 
157  # We are using RecoTracks to finding overlaps.
158  RecoTracks = Belle2.PyStoreArray('RecoTracks')
159  nRecoTracks = RecoTracks.getEntries()
161  for track_index in range(nRecoTracks):
162  track = RecoTracks[track_index]
163 
164  if track.wasFitSuccessful():
165 
166  TrackFitResults = Belle2.PyStoreArray('TrackFitResults')
167  nTrackFitResults = TrackFitResults.getEntries()
168  if nTrackFitResults == 1:
169  self.trackData.momentum = TrackFitResults[0].getMomentum().Mag()
170  self.Momentum.Fill(TrackFitResults[0].getMomentum().Mag())
171 
172  if track.hasPXDHits() or track.hasSVDHits():
173 
174  if track.getNumberOfSVDHits() % 2 == 0:
175  totalNumberOfHits = track.getNumberOfPXDHits() + track.getNumberOfSVDHits() / 2
176  else:
177  totalNumberOfHits = track.getNumberOfPXDHits() + (track.getNumberOfSVDHits() - 1) / 2
178 
179  self.TotalNumberOfHits.Fill(totalNumberOfHits)
180  # print('Total number of hits:', totalNumberOfHits)
181 
182  self.Chi2.Fill(track.getTrackFitStatus().getChi2())
183  self.NDF.Fill(track.getTrackFitStatus().getNdf())
184  self.Chi2OverNDF.Fill(track.getTrackFitStatus().getChi2() / track.getTrackFitStatus().getNdf())
185 
186  self.trackData.chi2 = track.getTrackFitStatus().getChi2()
187  self.trackData.ndf = track.getTrackFitStatus().getNdf()
188  self.trackData.chiSquaredOverNdf = track.getTrackFitStatus().getChi2() / track.getTrackFitStatus().getNdf()
189  # print('Chi2/NDF:', self.trackData.chiSquaredOverNdf, 'Chi2:', self.trackData.chi2)
190  self.rootfile.cd()
191  self.tree_track.Fill()
192 
193  # Check overlaps in PXD ladders & Check track in RecoTracks
194  if track.hasPXDHits():
195 
196  # Print number of PXD hits
197  # EventMetaData = Belle2.PyStoreObj('EventMetaData')
198  # event = EventMetaData.getEvent()
199  # nPxdHits = track.getNumberOfPXDHits()
200  # print('Event', event, 'has PXD', nPxdHits, 'hit(s):')
201 
202  # First loop over PXD Hits
203  for n in range(0, len(track.getPXDHitList())):
204  pxdHit = track.getPXDHitList()[n]
205  sensorID = Belle2.VxdID(pxdHit.getSensorID())
206  info = geoCache.get(sensorID)
207  layer = sensorID.getLayerNumber()
208  sensor = sensorID.getSensorNumber()
209  # print('Hit information: #hit:', n, 'sensor information:', layer, ladder, sensor)
210  self.HitsVsLayer.Fill(totalNumberOfHits, layer)
211  self.HitsVsSensor.Fill(totalNumberOfHits, sensor)
212  self.LayerVsSensor.Fill(layer, sensor)
213  self.PXDClusterSize.Fill(track.getPXDHitList()[n].getSize())
214 
215  # Check overlaps in SVD ladders
216  if track.hasSVDHits():
217 
218  # Print number of SVD hits
219  # EventMetaData = Belle2.PyStoreObj('EventMetaData')
220  # event = EventMetaData.getEvent()
221  # nSvdHits = track.getNumberOfSVDHits()
222  # print('Event', event, 'has SVD', nSvdHits, 'hit(s):')
223 
224  # First loop over SVD Hits
225  for n in range(0, len(track.getSVDHitList())):
226  svdHit = track.getSVDHitList()[n]
227  sensorID = Belle2.VxdID(svdHit.getSensorID())
228  info = geoCache.get(sensorID)
229  layer = sensorID.getLayerNumber()
230  sensor = sensorID.getSensorNumber()
231 
232  if svdHit.isUCluster():
233  # print('Hit information: #hit:', n, 'sensor information:', layer, ladder, sensor, 'isU')
234  for lst in range(0, len(track.getSVDHitList())):
235  svdHitA = track.getSVDHitList()[lst]
236  sensorIDA = Belle2.VxdID(svdHitA.getSensorID())
237  if info == geoCache.get(sensorIDA):
238  if svdHitA.isUCluster() == 0:
239  # print('Hit information: #hit:', n, 'sensor information:', layer, ladder, sensor)
240  self.HitsVsLayer.Fill(totalNumberOfHits, layer)
241  self.HitsVsSensor.Fill(totalNumberOfHits, sensor)
242  self.LayerVsSensor.Fill(layer, sensor)
243  # else:
244  # print('Hit information: #hit:', n, 'sensor information:', layer, ladder, sensor, 'isNotU')
245  else:
246  print('Fit was not successful')
247 
248  def terminate(self):
249  """Close & write output files"""
250  self.rootfilerootfile.cd()
251 
252  labelHits = ('', 'One', 'Two', 'Three', 'Four', 'Five', 'Six')
253  labelLayers = ('', 'First', 'Second', 'Third', 'Fourth', 'Fifth', 'Sixth')
254  labelSensors = ('', 'First', 'Second', 'Third', 'Fourth', 'Fifth')
255 
256  for i in range(1, 7):
257  self.TotalNumberOfHitsTotalNumberOfHits.GetXaxis().SetBinLabel(i, labelHits[i])
258  self.TotalNumberOfHitsTotalNumberOfHits.GetXaxis().SetLabelSize(0.05)
259  self.TotalNumberOfHitsTotalNumberOfHits.GetYaxis().SetTitle("Number of tracks in VXD")
260  self.TotalNumberOfHitsTotalNumberOfHits.GetXaxis().SetTitle("Number of hits in VXD")
261  self.TotalNumberOfHitsTotalNumberOfHits.GetYaxis().CenterTitle()
262  self.TotalNumberOfHitsTotalNumberOfHits.GetXaxis().CenterTitle()
263  self.TotalNumberOfHitsTotalNumberOfHits.GetYaxis().SetTitleOffset(1.3)
264  self.TotalNumberOfHitsTotalNumberOfHits.GetXaxis().SetTitleOffset(1.3)
265  self.TotalNumberOfHitsTotalNumberOfHits.SetFillStyle(3365)
266  self.TotalNumberOfHitsTotalNumberOfHits.SetFillColor(9)
267  self.TotalNumberOfHitsTotalNumberOfHits.SetLineColor(9)
268  self.TotalNumberOfHitsTotalNumberOfHits.Draw()
269  # self.TotalNumberOfHits.Write()
270 
271  for i in range(1, 7):
272  self.HitsVsLayerHitsVsLayer.GetXaxis().SetBinLabel(i, labelHits[i])
273  self.HitsVsLayerHitsVsLayer.GetYaxis().SetBinLabel(i, labelLayers[i])
274  self.HitsVsLayerHitsVsLayer.GetXaxis().SetLabelSize(0.05)
275  self.HitsVsLayerHitsVsLayer.GetYaxis().SetLabelSize(0.05)
276  self.HitsVsLayerHitsVsLayer.GetXaxis().SetTitle("Number hits in VXD")
277  self.HitsVsLayerHitsVsLayer.GetYaxis().SetTitle("Layer")
278  self.HitsVsLayerHitsVsLayer.GetYaxis().CenterTitle()
279  self.HitsVsLayerHitsVsLayer.GetXaxis().CenterTitle()
280  self.HitsVsLayerHitsVsLayer.GetYaxis().SetTitleOffset(1.6)
281  self.HitsVsLayerHitsVsLayer.GetXaxis().SetTitleOffset(1.3)
282  ROOT.gPad.SetMargin(0.11, 0.11, 0.1, 0.1)
283  self.HitsVsLayerHitsVsLayer.Draw()
284  # self.HitsVsLayer.Write()
285 
286  for i in range(1, 7):
287  self.HitsVsSensorHitsVsSensor.GetXaxis().SetBinLabel(i, labelHits[i])
288  for i in range(1, 6):
289  self.HitsVsSensorHitsVsSensor.GetYaxis().SetBinLabel(i, labelSensors[i])
290  self.HitsVsSensorHitsVsSensor.GetXaxis().SetLabelSize(0.05)
291  self.HitsVsSensorHitsVsSensor.GetYaxis().SetLabelSize(0.05)
292  self.HitsVsSensorHitsVsSensor.GetXaxis().SetTitle("Number hits in VXD")
293  self.HitsVsSensorHitsVsSensor.GetYaxis().SetTitle("Sensor")
294  self.HitsVsSensorHitsVsSensor.GetYaxis().CenterTitle()
295  self.HitsVsSensorHitsVsSensor.GetXaxis().CenterTitle()
296  self.HitsVsSensorHitsVsSensor.GetYaxis().SetTitleOffset(1.6)
297  self.HitsVsSensorHitsVsSensor.GetXaxis().SetTitleOffset(1.3)
298  ROOT.gPad.SetMargin(0.11, 0.11, 0.1, 0.1)
299  self.HitsVsSensorHitsVsSensor.Draw()
300  # self.HitsVsSensor.Write()
301 
302  for i in range(1, 7):
303  self.LayerVsSensorLayerVsSensor.GetXaxis().SetBinLabel(i, labelLayers[i])
304  for i in range(1, 6):
305  self.LayerVsSensorLayerVsSensor.GetYaxis().SetBinLabel(i, labelSensors[i])
306  self.LayerVsSensorLayerVsSensor.GetXaxis().SetLabelSize(0.05)
307  self.LayerVsSensorLayerVsSensor.GetYaxis().SetLabelSize(0.05)
308  self.LayerVsSensorLayerVsSensor.GetXaxis().SetTitle("Layer")
309  self.LayerVsSensorLayerVsSensor.GetYaxis().SetTitle("Sensor")
310  self.LayerVsSensorLayerVsSensor.GetYaxis().CenterTitle()
311  self.LayerVsSensorLayerVsSensor.GetXaxis().CenterTitle()
312  self.LayerVsSensorLayerVsSensor.GetYaxis().SetTitleOffset(1.6)
313  self.LayerVsSensorLayerVsSensor.GetXaxis().SetTitleOffset(1.3)
314  ROOT.gPad.SetMargin(0.11, 0.11, 0.1, 0.1)
315  self.LayerVsSensorLayerVsSensor.Draw()
316  # self.LayerVsSensor.Write()
317 
318  self.Chi2Chi2.GetXaxis().SetTitle("#Chi^{2}")
319  self.Chi2Chi2.GetYaxis().SetTitle("Number of tracks")
320  self.Chi2Chi2.GetYaxis().CenterTitle()
321  self.Chi2Chi2.GetXaxis().CenterTitle()
322  self.Chi2Chi2.GetYaxis().SetTitleOffset(1.3)
323  self.Chi2Chi2.GetXaxis().SetTitleOffset(1.3)
324  self.Chi2Chi2.SetFillStyle(3365)
325  self.Chi2Chi2.SetFillColor(9)
326  self.Chi2Chi2.SetLineColor(9)
327  self.Chi2Chi2.Draw()
328  # self.Chi2.Write()
329 
330  self.NDFNDF.GetXaxis().SetTitle("Degrees of freedom")
331  self.NDFNDF.GetYaxis().SetTitle("Number of tracks")
332  self.NDFNDF.GetYaxis().CenterTitle()
333  self.NDFNDF.GetXaxis().CenterTitle()
334  self.NDFNDF.GetYaxis().SetTitleOffset(1.3)
335  self.NDFNDF.GetXaxis().SetTitleOffset(1.3)
336  self.NDFNDF.SetFillStyle(3365)
337  self.NDFNDF.SetFillColor(9)
338  self.NDFNDF.SetLineColor(9)
339  self.NDFNDF.Draw()
340  # self.NDF.Write()
341 
342  self.Chi2OverNDFChi2OverNDF.SetTitle("#Chi^{2}/Degrees of freedom")
343  self.Chi2OverNDFChi2OverNDF.GetYaxis().SetTitle("Number of tracks")
344  self.Chi2OverNDFChi2OverNDF.GetYaxis().CenterTitle()
345  self.Chi2OverNDFChi2OverNDF.GetXaxis().CenterTitle()
346  self.Chi2OverNDFChi2OverNDF.GetYaxis().SetTitleOffset(1.3)
347  self.Chi2OverNDFChi2OverNDF.GetXaxis().SetTitleOffset(1.3)
348  self.Chi2OverNDFChi2OverNDF.SetFillStyle(3365)
349  self.Chi2OverNDFChi2OverNDF.SetFillColor(9)
350  self.Chi2OverNDFChi2OverNDF.SetLineColor(9)
351  self.Chi2OverNDFChi2OverNDF.Draw()
352  # self.Chi2OverNDF.Write()
353 
354  self.MomentumMomentum.GetXaxis().SetTitle("Momentum [GeVc^{-1}]")
355  self.MomentumMomentum.GetYaxis().SetTitle("Number of tracks")
356  self.MomentumMomentum.GetYaxis().CenterTitle()
357  self.MomentumMomentum.GetXaxis().CenterTitle()
358  self.MomentumMomentum.GetYaxis().SetTitleOffset(1.3)
359  self.MomentumMomentum.GetXaxis().SetTitleOffset(1.3)
360  self.MomentumMomentum.SetFillStyle(3365)
361  self.MomentumMomentum.SetFillColor(9)
362  self.MomentumMomentum.SetLineColor(9)
363  self.MomentumMomentum.Draw()
364  # self.Momentum.Write()
365 
366  self.ADCCountOverNumberOfHitsInCDCADCCountOverNumberOfHitsInCDC.GetXaxis().SetTitle("ADC Count [ADC] / number of hit")
367  self.ADCCountOverNumberOfHitsInCDCADCCountOverNumberOfHitsInCDC.GetYaxis().SetTitle("Number of tracks")
368  self.ADCCountOverNumberOfHitsInCDCADCCountOverNumberOfHitsInCDC.GetYaxis().CenterTitle()
369  self.ADCCountOverNumberOfHitsInCDCADCCountOverNumberOfHitsInCDC.GetXaxis().CenterTitle()
370  self.ADCCountOverNumberOfHitsInCDCADCCountOverNumberOfHitsInCDC.GetYaxis().SetTitleOffset(1.3)
371  self.ADCCountOverNumberOfHitsInCDCADCCountOverNumberOfHitsInCDC.GetXaxis().SetTitleOffset(1.3)
372  self.ADCCountOverNumberOfHitsInCDCADCCountOverNumberOfHitsInCDC.SetFillStyle(3365)
373  self.ADCCountOverNumberOfHitsInCDCADCCountOverNumberOfHitsInCDC.SetFillColor(9)
374  self.ADCCountOverNumberOfHitsInCDCADCCountOverNumberOfHitsInCDC.SetLineColor(9)
375  self.ADCCountOverNumberOfHitsInCDCADCCountOverNumberOfHitsInCDC.Draw()
376  # self.ADCCountOverNumberOfHitsInCDC.Write()
377 
378  function = ROOT.TF1("function", "[0]+[1]*log(x)")
379  function.SetParameters(0.0, 0.0)
380  self.MomentumVsADCCountOverNumberOfHitsInCDCMomentumVsADCCountOverNumberOfHitsInCDC.GetYaxis().SetTitle("ADC Count [ADC] / number of hit")
381  self.MomentumVsADCCountOverNumberOfHitsInCDCMomentumVsADCCountOverNumberOfHitsInCDC.GetXaxis().SetTitle("Momentum [GeVc^{-1}]")
382  self.MomentumVsADCCountOverNumberOfHitsInCDCMomentumVsADCCountOverNumberOfHitsInCDC.GetYaxis().CenterTitle()
383  self.MomentumVsADCCountOverNumberOfHitsInCDCMomentumVsADCCountOverNumberOfHitsInCDC.GetXaxis().CenterTitle()
384  self.MomentumVsADCCountOverNumberOfHitsInCDCMomentumVsADCCountOverNumberOfHitsInCDC.GetYaxis().SetTitleOffset(1.3)
385  self.MomentumVsADCCountOverNumberOfHitsInCDCMomentumVsADCCountOverNumberOfHitsInCDC.GetXaxis().SetTitleOffset(1.3)
386  self.MomentumVsADCCountOverNumberOfHitsInCDCMomentumVsADCCountOverNumberOfHitsInCDC.Fit(function)
387  self.MomentumVsADCCountOverNumberOfHitsInCDCMomentumVsADCCountOverNumberOfHitsInCDC.Draw()
388  # self.MomentumVsADCCountOverNumberOfHitsInCDC.Write()
389 
390  inverse_function = ROOT.TF1("inverse_function", "exp((-[0]+x)/[1])", 0, 115)
391  inverse_function.SetParameter(0, function.GetParameter(0))
392  inverse_function.SetParameter(1, function.GetParameter(1))
393  constant = ROOT.TF1("constant", "[0]", 115, 300)
394  constant.SetParameter(0, 20)
395  self.ADCCountOverNumberOfHitsInCDCVsMomentumADCCountOverNumberOfHitsInCDCVsMomentum.Fit(inverse_function, "R+")
396  self.ADCCountOverNumberOfHitsInCDCVsMomentumADCCountOverNumberOfHitsInCDCVsMomentum.Fit(constant, "R+")
397  self.ADCCountOverNumberOfHitsInCDCVsMomentumADCCountOverNumberOfHitsInCDCVsMomentum.GetXaxis().SetTitle("ADC Count [ADC] / number of hit")
398  self.ADCCountOverNumberOfHitsInCDCVsMomentumADCCountOverNumberOfHitsInCDCVsMomentum.GetYaxis().SetTitle("Momentum [GeVc^{-1}]")
399  self.ADCCountOverNumberOfHitsInCDCVsMomentumADCCountOverNumberOfHitsInCDCVsMomentum.GetYaxis().CenterTitle()
400  self.ADCCountOverNumberOfHitsInCDCVsMomentumADCCountOverNumberOfHitsInCDCVsMomentum.GetXaxis().CenterTitle()
401  self.ADCCountOverNumberOfHitsInCDCVsMomentumADCCountOverNumberOfHitsInCDCVsMomentum.GetYaxis().SetTitleOffset(1.3)
402  self.ADCCountOverNumberOfHitsInCDCVsMomentumADCCountOverNumberOfHitsInCDCVsMomentum.GetXaxis().SetTitleOffset(1.3)
403  self.ADCCountOverNumberOfHitsInCDCVsMomentumADCCountOverNumberOfHitsInCDCVsMomentum.Draw()
404  # self.ADCCountOverNumberOfHitsInCDCVsMomentum.Write()
405 
406  print('Parameters of fitted function y = exp((-p0+x)/p1) + p2, where x ~ ADCCounts/NumberOfCDCHits, y ~ momentum')
407  print('p0 =', inverse_function.GetParameter(0), 'p1 =', inverse_function.GetParameter(1), 'p2 =', constant.GetParameter(0))
408 
409  self.PXDClusterSizePXDClusterSize.GetXaxis().SetTitle("Size of PXD Clusters")
410  self.PXDClusterSizePXDClusterSize.GetYaxis().SetTitle("Number of tracks")
411  self.PXDClusterSizePXDClusterSize.GetYaxis().CenterTitle()
412  self.PXDClusterSizePXDClusterSize.GetXaxis().CenterTitle()
413  self.PXDClusterSizePXDClusterSize.GetYaxis().SetTitleOffset(1.3)
414  self.PXDClusterSizePXDClusterSize.GetXaxis().SetTitleOffset(1.3)
415  self.PXDClusterSizePXDClusterSize.SetFillStyle(3365)
416  self.PXDClusterSizePXDClusterSize.SetFillColor(9)
417  self.PXDClusterSizePXDClusterSize.SetLineColor(9)
418  self.PXDClusterSizePXDClusterSize.Draw()
419 
420  self.rootfilerootfile.Write()
421  self.rootfilerootfile.Close()
422 
423 
424 main = b2.create_path()
425 
426 main.add_module('RootInput')
427 main.add_module('Gearbox')
428 main.add_module('Geometry')
429 
430 reco.add_cosmics_reconstruction(main, pruneTracks=False, merge_tracks=True)
431 main.add_module('DAFRecoFitter', pdgCodesToUseForFitting=[13])
432 main.add_module("CDCDedxPID")
433 
434 CosmicAnalysis = CosmicAnalysis()
435 main.add_module(CosmicAnalysis)
436 
437 progress = b2.register_module('ProgressBar')
438 main.add_module(progress)
439 
440 b2.process(main)
441 
442 # Print call statistics
443 print(b2.statistics)
a (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:56
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:213
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
Momentum
Histogram with track momentum.
HitsVsSensor
Histogram with number of hits vs.
TotalNumberOfHits
Histogram with total number of hits on track.
HitsVsLayer
Histogram with number of hits vs.
Chi2OverNDF
Histogram with track Chi2/NDF.
ADCCountOverNumberOfHitsInCDC
Histogram with ADC count vs.
LayerVsSensor
Histogram with hitted layer vs.
PXDClusterSize
Histogram with PXD cluster size.
ADCCountOverNumberOfHitsInCDCVsMomentum
Profile with ADC count vs.
MomentumVsADCCountOverNumberOfHitsInCDC
Profile with ADC count vs.