Belle II Software  release-05-02-19
TOPPulseHeightAlgorithm.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2019 - Belle II Collaboration *
4  * *
5  * *
6  * Author: The Belle II Collaboration *
7  * Contributors: Marko Staric *
8  * *
9  * This software is provided "as is" without any warranty. *
10  **************************************************************************/
11 
12 #include <top/calibration/TOPPulseHeightAlgorithm.h>
13 #include "TROOT.h"
14 #include <TDirectory.h>
15 #include <TF1.h>
16 #include <TMath.h>
17 #include <string>
18 #include <algorithm>
19 
20 using namespace std;
21 
22 namespace Belle2 {
27  namespace TOP {
28 
29  TOPPulseHeightAlgorithm::TOPPulseHeightAlgorithm():
30  CalibrationAlgorithm("TOPPulseHeightCollector")
31  {
32  setDescription("Calibration algorithm for pulse-height and threshold efficiency calibration");
33  }
34 
36  {
37  gROOT->SetBatch();
38 
39  // construct file name, open output root file and book output tree
40 
41  const auto& expRun = getRunList();
42  string expNo = to_string(expRun[0].first);
43  while (expNo.length() < 4) expNo.insert(0, "0");
44  string runNo = to_string(expRun[0].second);
45  while (runNo.length() < 5) runNo.insert(0, "0");
46  string outputFileName = "pulseHeight-e" + expNo + "-r" + runNo + ".root";
47  m_file = TFile::Open(outputFileName.c_str(), "recreate");
48 
49  m_tree = new TTree("tree", "pulse-height calibration results");
50  m_tree->Branch<int>("slot", &m_slot);
51  m_tree->Branch<unsigned>("channel", &m_channel);
52  m_tree->Branch<int>("numPhot", &m_numPhot);
53  m_tree->Branch<float>("x0", &m_x0);
54  m_tree->Branch<float>("p1", &m_p1);
55  m_tree->Branch<float>("p2", &m_p2);
56  m_tree->Branch<float>("x0err", &m_x0err);
57  m_tree->Branch<float>("p1err", &m_p1err);
58  m_tree->Branch<float>("p2err", &m_p2err);
59  m_tree->Branch<float>("effi", &m_effi);
60  m_tree->Branch<float>("meanHist", &m_meanHist);
61  m_tree->Branch<float>("meanFunc", &m_meanFunc);
62  m_tree->Branch<float>("chi2", &m_chi2);
63  m_tree->Branch<int>("ndf", &m_ndf);
64  m_tree->Branch<int>("fitStatus", &m_fitStatus);
65  m_tree->Branch<bool>("good", &m_good);
66 
67  // write-out control histograms
68 
69  auto h1a = getObjectPtr<TH1F>("time");
70  if (h1a) h1a->Write();
71  auto h1b = getObjectPtr<TH1F>("time_sel");
72  if (h1b) h1b->Write();
73  auto h2a = getObjectPtr<TH2F>("ph_vs_width");
74  if (h2a) h2a->Write();
75  auto h2b = getObjectPtr<TH2F>("ph_vs_width_sel");
76  if (h2b) h2b->Write();
77 
78  // create payloads for storing results
79 
82 
83  // fit histograms
84 
85  TDirectory* oldDir = gDirectory;
86  int nsucc = 0;
87  for (int slot = 1; slot <= 16; slot++) {
88  m_slot = slot;
89  string name = "ph_slot_" + to_string(slot);
90  auto h = getObjectPtr<TH2F>(name);
91  if (not h) continue;
92 
93  name = "slot_" + to_string(slot);
94  oldDir->mkdir(name.c_str())->cd();
95  h->Write();
96 
97  int n = fitChannels(h);
98  B2INFO("slot " << slot << ": " << n << "/512 channels successfully fitted");
99  nsucc += n;
100  }
101 
102  // write the results and close the file
103 
104  m_file->Write();
105  m_file->Close();
106 
107  // check the results and return if fraction of well fitted too small
108 
109  if (nsucc / 8192.0 < m_minCalibrated) {
110  delete m_pulseHeight;
111  delete m_thresholdEffi;
112  return c_NotEnoughData;
113  }
114 
115  // otherwise import payloads to DB
116 
119 
120  return c_OK;
121  }
122 
123 
124  int TOPPulseHeightAlgorithm::fitChannels(std::shared_ptr<TH2F> h2d)
125  {
126  int n = 0;
127 
128  for (int ch = 0; ch < h2d->GetNbinsX(); ch++) {
129  m_channel = ch;
130  string chan = to_string(ch);
131  while (chan.length() < 3) chan.insert(0, "0");
132  string name = "chan_" + chan;
133  auto* h = h2d->ProjectionY(name.c_str(), ch + 1, ch + 1);
134  string title = "slot " + to_string(m_slot) + " channel " + to_string(ch);
135  h->SetTitle(title.c_str());
136 
137  auto status = fitPulseHeight(h);
138 
139  if (status == 0 and m_good) {
141  m_thresholdEffi->setThrEff(m_slot, m_channel, std::min(m_effi, 1.0f), m_xmin);
142  n++;
143  }
144  }
145 
146  return n;
147  }
148 
149 
151  {
152  m_numPhot = h->GetSumOfWeights();
153  if (m_numPhot < std::max(m_minEntries, 5)) return -1;
154 
155  double xmin = h->GetXaxis()->GetXmin();
156  double xmax = h->GetXaxis()->GetXmax();
157 
158  auto* func = new TF1("func", "[0]*x/[1]*exp(-pow(x/[1],[2]))", xmin, xmax);
159 
160  double p2 = 1.0;
161  m_meanHist = h->GetMean();
162  double x0 = m_meanHist / 3; // for p2 = 1 only
163  func->SetParameter(0, m_numPhot / x0); // for p2 = 1 only
164  func->SetParameter(1, x0);
165  func->SetParameter(2, p2);
166 
167  int status = h->Fit(func, "LRSQ", "", m_xmin, xmax);
168 
169  m_x0 = func->GetParameter(1);
170  m_p1 = 1.0;
171  m_p2 = func->GetParameter(2);
172  m_x0err = func->GetParError(1);
173  m_p1err = 0;
174  m_p2err = func->GetParError(2);
175  m_chi2 = func->GetChisquare();
176  m_effi = getEfficiency(h, func);
177  m_meanFunc = TMath::Gamma((m_p1 + 2) / m_p2) / TMath::Gamma((m_p1 + 1) / m_p2) * m_x0;
178  m_fitStatus = status;
179  m_good = (m_chi2 / m_ndf < 10 and m_p2 > 0.5 and m_p2 < 5);
180  m_tree->Fill();
181 
182  return status;
183  }
184 
185 
186  double TOPPulseHeightAlgorithm::getEfficiency(TH1D* h, TF1* func)
187  {
188  double fbelow = 0;
189  double fabove = 0;
190  double count = 0;
191  m_ndf = -func->GetNumberFreeParameters();
192  for (int i = 0; i < h->GetNbinsX(); i++) {
193  double x = h->GetBinCenter(i + 1);
194  if (x > m_xmin) {
195  fabove += func->Eval(x);
196  count += h->GetBinContent(i);
197  if (h->GetBinContent(i) > 0) m_ndf++;
198  } else {
199  fbelow += func->Eval(x);
200  }
201  }
202  double integral = (fbelow + fabove) * count / fabove;
203  return h->GetSumOfWeights() / integral;
204  }
205 
206 
207  } // end namespace TOP
209 } // end namespace Belle2
Belle2::TOP::TOPPulseHeightAlgorithm::calibrate
virtual EResult calibrate() final
algorithm implementation
Definition: TOPPulseHeightAlgorithm.cc:35
Belle2::TOP::TOPPulseHeightAlgorithm::m_x0err
float m_x0err
error on x0 [ADC counts]
Definition: TOPPulseHeightAlgorithm.h:109
Belle2::TOP::TOPPulseHeightAlgorithm::fitPulseHeight
int fitPulseHeight(TH1D *h)
Fit pulse-height distribution of a single channel with P(x) = (x/x0)^p1 * exp(-(x/x0)^p2) and p1 = 1.
Definition: TOPPulseHeightAlgorithm.cc:150
Belle2::TOPCalChannelPulseHeight
Pulse height parameterizations for all 512 channels of 16 modules.
Definition: TOPCalChannelPulseHeight.h:36
Belle2::TOPCalChannelThresholdEff::setThrEff
void setThrEff(int moduleID, unsigned channel, float ThrEff, short offlineThreshold)
Sets the threshold efficiency and correspolding threshold for a single channel and switches status to...
Definition: TOPCalChannelThresholdEff.h:64
Belle2::TOP::TOPPulseHeightAlgorithm::m_x0
float m_x0
fitted distribution parameter x0 [ADC counts]
Definition: TOPPulseHeightAlgorithm.h:106
Belle2::TOP::TOPPulseHeightAlgorithm::m_fitStatus
int m_fitStatus
fit status
Definition: TOPPulseHeightAlgorithm.h:117
Belle2::TOP::TOPPulseHeightAlgorithm::m_pulseHeight
TOPCalChannelPulseHeight * m_pulseHeight
parametrized PH distributions
Definition: TOPPulseHeightAlgorithm.h:125
Belle2::CalibrationAlgorithm::saveCalibration
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
Definition: CalibrationAlgorithm.cc:290
Belle2::TOP::TOPPulseHeightAlgorithm::m_xmin
double m_xmin
fit range lower limit [ADC counts]
Definition: TOPPulseHeightAlgorithm.h:99
Belle2::TOP::TOPPulseHeightAlgorithm::m_thresholdEffi
TOPCalChannelThresholdEff * m_thresholdEffi
threshold efficiencies
Definition: TOPPulseHeightAlgorithm.h:126
Belle2::TOP::TOPPulseHeightAlgorithm::m_p1err
float m_p1err
error on p1
Definition: TOPPulseHeightAlgorithm.h:110
Belle2::TOP::TOPPulseHeightAlgorithm::m_p2
float m_p2
fitted distribution parameter p2
Definition: TOPPulseHeightAlgorithm.h:108
Belle2::TOP::TOPPulseHeightAlgorithm::m_chi2
float m_chi2
chi^2
Definition: TOPPulseHeightAlgorithm.h:115
Belle2::TOP::TOPPulseHeightAlgorithm::m_minEntries
int m_minEntries
minimal number of histogram entries to perform fit
Definition: TOPPulseHeightAlgorithm.h:98
Belle2::CalibrationAlgorithm::c_OK
@ c_OK
Finished successfuly =0 in Python.
Definition: CalibrationAlgorithm.h:51
Belle2::CalibrationAlgorithm::setDescription
void setDescription(const std::string &description)
Set algorithm description (in constructor)
Definition: CalibrationAlgorithm.h:331
Belle2::TOP::TOPPulseHeightAlgorithm::m_slot
int m_slot
module ID
Definition: TOPPulseHeightAlgorithm.h:103
Belle2::TOP::TOPPulseHeightAlgorithm::m_p2err
float m_p2err
error on p2
Definition: TOPPulseHeightAlgorithm.h:111
Belle2::TOP::TOPPulseHeightAlgorithm::m_p1
float m_p1
fitted distribution parameter p1
Definition: TOPPulseHeightAlgorithm.h:107
Belle2::TOP::TOPPulseHeightAlgorithm::m_meanHist
float m_meanHist
histogram mean value
Definition: TOPPulseHeightAlgorithm.h:113
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::TOP::TOPPulseHeightAlgorithm::m_ndf
int m_ndf
NDF.
Definition: TOPPulseHeightAlgorithm.h:116
Belle2::TOPCalChannelThresholdEff
Class to store the threshold efficiency (i.e.
Definition: TOPCalChannelThresholdEff.h:38
Belle2::TOP::TOPPulseHeightAlgorithm::m_channel
unsigned m_channel
channel ID
Definition: TOPPulseHeightAlgorithm.h:104
Belle2::TOP::TOPPulseHeightAlgorithm::m_minCalibrated
double m_minCalibrated
min.
Definition: TOPPulseHeightAlgorithm.h:100
Belle2::CalibrationAlgorithm::EResult
EResult
The result of calibration.
Definition: CalibrationAlgorithm.h:50
Belle2::TOP::TOPPulseHeightAlgorithm::m_tree
TTree * m_tree
output ntuple
Definition: TOPPulseHeightAlgorithm.h:122
Belle2::CalibrationAlgorithm::getRunList
const std::vector< Calibration::ExpRun > & getRunList() const
Get the list of runs for which calibration is called.
Definition: CalibrationAlgorithm.h:276
Belle2::TOP::TOPPulseHeightAlgorithm::m_meanFunc
float m_meanFunc
fitted distribution mean value
Definition: TOPPulseHeightAlgorithm.h:114
Belle2::TOP::TOPPulseHeightAlgorithm::m_effi
float m_effi
threshold efficiency
Definition: TOPPulseHeightAlgorithm.h:112
Belle2::TOPCalChannelPulseHeight::setParameters
void setParameters(int moduleID, unsigned channel, double x0, double p1, double p2)
Sets calibration for a single channel and switches status to calibrated.
Definition: TOPCalChannelPulseHeight.h:62
Belle2::CalibrationAlgorithm::c_NotEnoughData
@ c_NotEnoughData
Needs more data =2 in Python.
Definition: CalibrationAlgorithm.h:53
Belle2::TOP::TOPPulseHeightAlgorithm::m_numPhot
int m_numPhot
number of photons (histogram entries)
Definition: TOPPulseHeightAlgorithm.h:105
Belle2::CalibrationAlgorithm
Base class for calibration algorithms.
Definition: CalibrationAlgorithm.h:47
Belle2::TOP::TOPPulseHeightAlgorithm::getEfficiency
double getEfficiency(TH1D *h, TF1 *func)
Calculate and return threshold efficiency.
Definition: TOPPulseHeightAlgorithm.cc:186
Belle2::TOP::TOPPulseHeightAlgorithm::m_file
TFile * m_file
root file
Definition: TOPPulseHeightAlgorithm.h:121
Belle2::TOP::TOPPulseHeightAlgorithm::fitChannels
int fitChannels(std::shared_ptr< TH2F > h)
Fit pulse-height distributions of a single module.
Definition: TOPPulseHeightAlgorithm.cc:124
Belle2::TOP::TOPPulseHeightAlgorithm::m_good
bool m_good
on true fit is good
Definition: TOPPulseHeightAlgorithm.h:118