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

Public Member Functions

def initialize (self)
 

Detailed Description

 determination of average quantum and collection efficiency from TOPPmtQEs 

Definition at line 24 of file averageQE.py.

Member Function Documentation

◆ initialize()

def initialize (   self)
 run the procedure, store the results in a root file and print in xml format 

Definition at line 27 of file averageQE.py.

27 def initialize(self):
28 ''' run the procedure, store the results in a root file and print in xml format '''
29
30 dbarray = Belle2.PyDBArray('TOPPmtQEs')
31 if dbarray.getEntries() == 0:
32 b2.B2ERROR("No entries in TOPPmtQEs")
33 return
34
35 geo = Belle2.PyDBObj('TOPGeometry')
36 if not geo:
37 b2.B2ERROR("No TOPGeometry")
38 return
39
40 Bfield_on = True # use collection efficiencies for 1.5 T
41
42 lambdaStep = dbarray[0].getLambdaStep()
43 lambdaFirst = dbarray[0].getLambdaFirst()
44 lambdaLast = dbarray[0].getLambdaLast()
45 n = int((lambdaLast - lambdaFirst) / lambdaStep) + 1
46 print(lambdaFirst, lambdaLast, lambdaStep, n)
47 if Bfield_on:
48 Bfield = '(1.5 T)'
49 else:
50 Bfield = '(0 T)'
51
52 pde_2d = TH2F("pde_2d", "quantum times collection efficiencies " + Bfield,
53 n, lambdaFirst - lambdaStep / 2, lambdaLast + lambdaStep / 2,
54 300, 0.0, 1.0)
55 pde_prof = TProfile("pde_prof", "quantum times collection efficiencies " + Bfield,
56 n, lambdaFirst - lambdaStep / 2, lambdaLast + lambdaStep / 2,
57 0.0, 1.0)
58 pde_nom = TH1F("pde_nom", "nominal quantum times collection efficiency (as found in DB)",
59 n, lambdaFirst - lambdaStep / 2, lambdaLast + lambdaStep / 2)
60 ce_1d = TH1F("ce_1d", "collection efficiencies " + Bfield, 100, 0.0, 1.0)
61
62 ce = 0
63 for pmtQE in dbarray:
64 ce_1d.Fill(pmtQE.getCE(Bfield_on))
65 ce += pmtQE.getCE(Bfield_on)
66 for pmtPixel in range(1, 17):
67 for i in range(n):
68 lam = lambdaFirst + lambdaStep * i
69 effi = pmtQE.getEfficiency(pmtPixel, lam, Bfield_on)
70 pde_2d.Fill(lam, effi)
71 pde_prof.Fill(lam, effi)
72
73 nominalQE = geo.getNominalQE()
74 for i in range(n):
75 lam = lambdaFirst + lambdaStep * i
76 effi = nominalQE.getEfficiency(lam)
77 pde_nom.SetBinContent(i+1, effi)
78
79 ce /= dbarray.getEntries()
80 print('average CE =', ce)
81 file = TFile('averageQE.root', 'recreate')
82 pde_2d.Write()
83 pde_prof.Write()
84 pde_nom.Write()
85 ce_1d.Write()
86 file.Close()
87
88 print()
89 if Bfield_on:
90 print('<ColEffi descr="average over all PMTs for 1.5 T">' + str(round(ce, 4)) + '</ColEffi>')
91 else:
92 print('<ColEffi descr="average over all PMTs for 0 T">' + str(round(ce, 4)) + '</ColEffi>')
93 print('<LambdaFirst unit="nm">' + str(lambdaFirst) + '</LambdaFirst>')
94 print('<LambdaStep unit="nm">' + str(lambdaStep) + '</LambdaStep>')
95 for i in range(n):
96 qe = pde_prof.GetBinContent(i+1) / ce
97 print('<Qeffi>' + str(round(qe, 3)) + '</Qeffi>')
98 print()
99
100
101# Central database
Class to access a DB Array from Python.
Definition: PyDBArray.h:43
Class to access a DBObjPtr from Python.
Definition: PyDBObj.h:50

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