Belle II Software development
TOPCalPhotonYields.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9#include <top/dbobjects/TOPCalPhotonYields.h>
10#include <framework/logging/Logger.h>
11
12using namespace std;
13
14namespace Belle2 {
20 void TOPCalPhotonYields::set(int slot, const TH1F* photonYields, const TH1F* backgroundYields, const TH1F* alphaRatio,
21 const TH1F* activePixels, const TH2F* pulseHeights, const TH1F* muonZ)
22 {
23 string slotName = (slot < 10) ? "_0" + to_string(slot) : "_" + to_string(slot);
24
25 m_photonYields.push_back(TH2F(("photonYields" + slotName).c_str(),
26 ("Photon yields for slot " + to_string(slot) + "; pixel column; pixel row").c_str(),
27 c_numCols, 0.5, c_numCols + 0.5, c_numRows, 0.5, c_numRows + 0.5));
28 copyContent(photonYields, m_photonYields.back());
29
30 m_backgroundYields.push_back(TH2F(("backgroundYields" + slotName).c_str(),
31 ("Background yields for slot " + to_string(slot) + "; pixel column; pixel row").c_str(),
32 c_numCols, 0.5, c_numCols + 0.5, c_numRows, 0.5, c_numRows + 0.5));
33 copyContent(backgroundYields, m_backgroundYields.back());
34
35 m_alphaRatio.push_back(TH2F(("alphaRatio" + slotName).c_str(),
36 ("Equalized alpha ratio for slot " + to_string(slot) + "; pixel column; pixel row").c_str(),
37 c_numCols, 0.5, c_numCols + 0.5, c_numRows, 0.5, c_numRows + 0.5));
38 copyContent(alphaRatio, m_alphaRatio.back());
39
40 m_activePixels.push_back(TH2F(("activePixels" + slotName).c_str(),
41 ("Active pixels for slot " + to_string(slot) + "; pixel column; pixel row").c_str(),
42 c_numCols, 0.5, c_numCols + 0.5, c_numRows, 0.5, c_numRows + 0.5));
43 copyContent(activePixels, m_activePixels.back());
44 m_activePixels.back().Scale(1 / muonZ->GetEntries());
45
46 m_pulseHeights.push_back(*pulseHeights);
47 m_muonZ.push_back(*muonZ);
48 }
49
50
51 void TOPCalPhotonYields::copyContent(const TH1F* input, TH2F& output)
52 {
53 for (int bin = 1; bin <= input->GetNbinsX(); bin++) {
54 int row = (bin - 1) / 64 + 1;
55 int col = (bin - 1) % 64 + 1;
56 output.SetBinContent(col, row, input->GetBinContent(bin));
57 output.SetBinError(col, row, input->GetBinError(bin));
58 }
59 }
60
61
62 const TH2F* TOPCalPhotonYields::getPhotonYields(int slot) const
63 {
64 unsigned index = slot - 1;
65 if (index < m_photonYields.size()) return &m_photonYields[index];
66 return 0;
67 }
68
69
70 const TH2F* TOPCalPhotonYields::getBackgroundYields(int slot) const
71 {
72 unsigned index = slot - 1;
73 if (index < m_backgroundYields.size()) return &m_backgroundYields[index];
74 return 0;
75 }
76
77
78 const TH2F* TOPCalPhotonYields::getAlphaRatio(int slot) const
79 {
80 unsigned index = slot - 1;
81 if (index < m_alphaRatio.size()) return &m_alphaRatio[index];
82 return 0;
83 }
84
85
86 const TH2F* TOPCalPhotonYields::getActivePixels(int slot) const
87 {
88 unsigned index = slot - 1;
89 if (index < m_activePixels.size()) return &m_activePixels[index];
90 return 0;
91 }
92
93
94 const TH2F* TOPCalPhotonYields::getPulseHeights(int slot) const
95 {
96 unsigned index = slot - 1;
97 if (index < m_pulseHeights.size()) return &m_pulseHeights[index];
98 return 0;
99 }
100
101
102 const TH1F* TOPCalPhotonYields::getMuonZ(int slot) const
103 {
104 unsigned index = slot - 1;
105 if (index < m_muonZ.size()) return &m_muonZ[index];
106 return 0;
107 }
108
109
111}
std::vector< TH2F > m_activePixels
active pixels (index = slot - 1)
std::vector< TH2F > m_pulseHeights
pixel pulse-heights (index = slot - 1)
std::vector< TH1F > m_muonZ
local z distribution of tracks (index = slot - 1)
std::vector< TH2F > m_alphaRatio
equalized alpha ratio per pixel (index = slot - 1)
std::vector< TH2F > m_photonYields
photon yields per pixel (index = slot - 1)
std::vector< TH2F > m_backgroundYields
background yields per pixel (index = slot - 1)
@ c_numCols
number of pixel columns
@ c_numRows
number of pixel rows
void set(int slot, const TH1F *photonYields, const TH1F *backgroundYields, const TH1F *alphaRatio, const TH1F *activePixels, const TH2F *pulseHeights, const TH1F *muonZ)
Sets the data of a given slot.
const TH2F * getAlphaRatio(int slot) const
Returns a 2D histogram of equalized pixel alpha ratio.
const TH1F * getMuonZ(int slot) const
Returns z distribution of tracks used to determine pixel yields.
const TH2F * getActivePixels(int slot) const
Returns a 2D histogram of active pixels.
const TH2F * getPulseHeights(int slot) const
Returns a 2D histogram of pixel pulse-heights.
const TH2F * getBackgroundYields(int slot) const
Returns a 2D histogram of background pixel yields.
void copyContent(const TH1F *input, TH2F &output)
Copy content of 1D histogram into 2D histogram.
const TH2F * getPhotonYields(int slot) const
Returns a 2D histogram of photon pixel yields.
Abstract base class for different kinds of events.
STL namespace.