Belle II Software  release-08-01-10
TOPModuleT0LLAlgorithm.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/calibration/TOPModuleT0LLAlgorithm.h>
10 #include <top/dbobjects/TOPCalModuleT0.h>
11 #include <top/utilities/Chi2MinimumFinder1D.h>
12 
13 #include <math.h>
14 #include <TFile.h>
15 #include <TH1F.h>
16 #include <TH1D.h>
17 #include <TH2F.h>
18 
19 #include <string>
20 #include <vector>
21 
22 using namespace std;
23 
24 namespace Belle2 {
29  namespace TOP {
30 
31  TOPModuleT0LLAlgorithm::TOPModuleT0LLAlgorithm():
32  CalibrationAlgorithm("TOPModuleT0LLCollector")
33  {
34  setDescription("Calibration algorithm for method LL");
35  }
36 
38  {
39  // get basic histogram
40 
41  auto h1 = getObjectPtr<TH2F>("tracks_per_slot");
42  if (not h1) {
43  B2ERROR("TOPModuleT0LLAlgorithm: histogram 'tracks_per_slot' not found");
44  return c_Failure;
45  }
46  unsigned numModules = h1->GetNbinsX();
47  unsigned numSets = h1->GetNbinsY();
48 
49  // construct file name, open output root file
50 
51  const auto& expRun = getRunList();
52  string expNo = to_string(expRun[0].first);
53  while (expNo.length() < 4) expNo.insert(0, "0");
54  string runNo = to_string(expRun[0].second);
55  while (runNo.length() < 5) runNo.insert(0, "0");
56  string outputFileName = "moduleT0_final-e" + expNo + "-r" + runNo + ".root";
57  auto* file = TFile::Open(outputFileName.c_str(), "recreate");
58 
59  // write-out input histograms
60 
61  h1->Write();
62  auto h2 = getObjectPtr<TH1F>("numHits");
63  if (h2) h2->Write();
64  auto h3 = getObjectPtr<TH2F>("timeHits");
65  if (h3) h3->Write();
66  auto h4 = getObjectPtr<TH1F>("offset");
67  if (h4) h4->Write();
68 
69  // create histogram of pulls which is needed for error rescaling
70 
71  auto h_pulls = new TH1F("pulls", "Pulls of statistically independent results",
72  200, -15.0, 15.0);
73  h_pulls->SetXTitle("pulls");
74 
75  // create vectors for keeping the results
76 
77  std::vector<double> T0(numModules, 0);
78  std::vector<double> errT0(numModules, 0);
79 
80  // find minimum for each slot
81 
82  bool ok = true;
83  for (unsigned slot = 1; slot <= numModules; slot++) {
84 
85  Chi2MinimumFinder1D slotFinder;
86  std::vector<double> pos, err;
87 
88  // loop over statistically independent subsamples, add them
89  for (unsigned set = 0; set < numSets; set++) {
90 
91  if (h1->GetBinContent(slot, set + 1) == 0) continue;
92 
93  string name = "chi2_set" + to_string(set) + "_slot" + to_string(slot);
94  auto h = getObjectPtr<TH1D>(name);
95  if (not h) {
96  B2ERROR("TOPModuleT0LLAlgorithm: histogram '" << name << "' not found");
97  continue;
98  }
99  Chi2MinimumFinder1D setFinder(h);
100  if (slotFinder.getNbins() == 0) {
101  slotFinder = setFinder;
102  } else {
103  slotFinder.add(setFinder);
104  }
105  const auto& minimum = setFinder.getMinimum();
106  if (not minimum.valid) continue;
107  pos.push_back(minimum.position);
108  err.push_back(minimum.error);
109  }
110 
111  // fill pulls
112  for (unsigned i = 0; i < pos.size(); i++) {
113  for (unsigned j = i + 1; j < pos.size(); j++) {
114  double pull = (pos[i] - pos[j]) / sqrt(err[i] * err[i] + err[j] * err[j]);
115  h_pulls->Fill(pull);
116  }
117  }
118 
119  // determine minimum and store it
120  const auto& minimum = slotFinder.getMinimum();
121  if (minimum.valid) {
122  T0[slot - 1] = minimum.position;
123  errT0[slot - 1] = minimum.error;
124  } else {
125  B2ERROR("TOPModuleT0Algorithm: no minimum found for slot " << slot);
126  ok = false;
127  }
128 
129  // write-out histogram of a chi2 scan
130  string name = "chi2_slot_" + to_string(slot);
131  string title = "chi2 scan, slot " + to_string(slot) + "; t0 [ns]; -2 logL";
132  auto h = slotFinder.getHistogram(name, title);
133  h.Write();
134  }
135 
136  // rescale the errors
137 
138  if (h_pulls->GetEntries() > 1) {
139  double scaleError = h_pulls->GetRMS();
140  for (auto& err : errT0) err *= scaleError;
141  }
142 
143  // store the results in a histogram
144 
145  auto* h_moduleT0 = new TH1F("moduleT0", "Module T0", 16, 0.5, 16.5);
146  h_moduleT0->SetXTitle("slot number");
147  h_moduleT0->SetYTitle("module T0 [ns]");
148  for (unsigned i = 0; i < T0.size(); i++) {
149  h_moduleT0->SetBinContent(i + 1, T0[i]);
150  h_moduleT0->SetBinError(i + 1, errT0[i]);
151  }
152 
153  // write the results and close the file
154 
155  file->Write();
156  file->Close();
157 
158  // check the results and return if not good enough
159 
160  if (not ok) return c_Failure;
161  for (auto err : errT0) {
162  if (err > m_minError) return c_NotEnoughData;
163  }
164 
165  // otherwise create and import payload to DB
166 
167  auto* moduleT0 = new TOPCalModuleT0();
168  for (unsigned i = 0; i < T0.size(); i++) {
169  moduleT0->setT0(i + 1, T0[i], errT0[i]);
170  }
171  moduleT0->suppressAverage();
172  saveCalibration(moduleT0);
173 
174  return c_OK;
175  }
176 
177  } // end namespace TOP
179 } // end namespace Belle2
Base class for calibration algorithms.
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
void setDescription(const std::string &description)
Set algorithm description (in constructor)
EResult
The result of calibration.
@ c_OK
Finished successfuly =0 in Python.
@ c_NotEnoughData
Needs more data =2 in Python.
@ c_Failure
Failed =3 in Python.
const std::vector< Calibration::ExpRun > & getRunList() const
Get the list of runs for which calibration is called.
Module T0 calibration constants for all 16 modules.
Minimum finder using tabulated chi^2 values in one dimension.
TH1F getHistogram(std::string name, std::string title) const
Returns chi^2 packed into 1D histogram.
void add(unsigned i, double chi2)
Add chi^2 value to bin i.
const Minimum & getMinimum()
Returns parabolic minimum.
int getNbins() const
Returns number of bins.
virtual EResult calibrate() final
algorithm implementation
double m_minError
minimal moduleT0 uncertainty [ns] to declare c_OK
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
Abstract base class for different kinds of events.
double position
position of the minimum