Belle II Software  release-08-01-10
SVDCrossTalkCalibrationsAlgorithm.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 <svd/calibration/SVDCrossTalkCalibrationsAlgorithm.h>
10 #include <svd/calibration/SVDCrossTalkStripsCalibrations.h>
11 
12 #include <TH1F.h>
13 
14 #include <iostream>
15 
16 using namespace std;
17 using namespace Belle2;
18 
19 SVDCrossTalkCalibrationsAlgorithm::SVDCrossTalkCalibrationsAlgorithm(const std::string& str) :
20  CalibrationAlgorithm("SVDCrossTalkCalibrationsCollector")
21 {
22  setDescription("Calibration algorithm for SVDCrossTalkCalibrations payload");
23  m_id = str;
24 }
25 
27 {
28  int isCrossTalkCal = 0;
29  auto payload = new Belle2::SVDCrossTalkStripsCalibrations::t_payload(isCrossTalkCal, m_id);
30 
31  auto tree = getObjectPtr<TTree>("HTreeCrossTalkCalib");
32 
33  TH1F* hist = new TH1F("", "", 768, 0, 768);
34 
35  if (!tree) {
36  B2WARNING("No tree");
37  } else if (!tree->GetEntries()) {
38  B2WARNING("Empty tree");
39  }
40 
41  int layer = 0;
42  int ladder = 0;
43  int sensor = 0;
44  int side = 0;
45 
46  tree->SetBranchAddress("hist", &hist);
47  tree->SetBranchAddress("layer", &layer);
48  tree->SetBranchAddress("ladder", &ladder);
49  tree->SetBranchAddress("sensor", &sensor);
50  tree->SetBranchAddress("side", &side);
51 
52 
53 
54  for (int i = 0; i < tree->GetEntries(); i++) {
55  tree->GetEntry(i);
56 
57  int nstrips = 768;
58  if (side == 0 && layer != 3) nstrips = 512;
59 
60  for (int strip = 0; strip < nstrips; strip++) {
61 
62  isCrossTalkCal = hist->GetBinContent(strip + 1);
63 
64  if (layer == 4 && ladder == 1 && sensor == 2 && side == 1
65  && hist->GetEntries() < m_minEntries) { //Check that we have enough events populating the calibration
66  cout << "Not enough Data: " << hist->GetEntries() << " entries found" << endl;
67  return c_NotEnoughData;
68  }
69 
70  payload->set(layer, ladder, sensor, side, strip, isCrossTalkCal);
71  }
72  }
73 
74 
75  saveCalibration(payload, "SVDCrossTalkStripsCalibrations");
76 
77  return c_OK;
78 
79 }
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.
int m_minEntries
Minimum number of required entries for collector histogram L4.1.2 u-side.
virtual EResult calibrate() override
Run algo on data.
SVDCalibrationsBase< SVDCalibrationsBitmap > t_payload
typedef of the SVDCrossTalkStripsCalibrations payload for all SVD strips
Abstract base class for different kinds of events.