Belle II Software development
SVDHotStripsCalibrationsAlgorithm.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#include <svd/calibration/SVDHotStripsCalibrationsAlgorithm.h>
9
10
11#include <svd/calibration/SVDHotStripsCalibrations.h>
12
13#include <TH1F.h>
14#include <framework/logging/Logger.h>
15#include <iostream>
16#include <TString.h>
17
18using namespace Belle2;
19
21 CalibrationAlgorithm("SVDOccupancyCalibrationsCollector")
22{
23 setDescription("SVDHotStripsCalibrations calibration algorithm");
24 m_id = str;
25 m_base = -1;
26 m_relOccPrec = 5.;
27}
28
30{
31
32 int isHotStrip = 0;
33 int vecHS[768];
34 double stripOccAfterAbsCut[768];
35
36 for (int i = 0; i < 768; i++) { vecHS[i] = 0; stripOccAfterAbsCut[i] = 0;}
37
38 // float occCal = 1.;
39 float occThr = 0.2; //first absolute cut for strips with occupancy >20%
40
41 // auto HSBitmap = new Belle2::SVDCalibrationsBitmap();
42 auto payload = new Belle2::SVDHotStripsCalibrations::t_payload(isHotStrip, m_id);
43
44 auto tree = getObjectPtr<TTree>("HTreeOccupancyCalib");
45
46 TH1F* hocc = new TH1F("", "", 768, 0, 768);
47
48 if (!tree) {
49 B2WARNING("No tree object.");
50 } else if (!tree->GetEntries()) {
51 B2WARNING("No data in the tree.");
52 }
53
54 int layer = 0;
55 int ladder = 0;
56 int sensor = 0;
57 int side = 0;
58
59 tree->SetBranchAddress("hist", &hocc);
60 tree->SetBranchAddress("layer", &layer);
61 tree->SetBranchAddress("ladder", &ladder);
62 tree->SetBranchAddress("sensor", &sensor);
63 tree->SetBranchAddress("view", &side);
64
65
66 for (int i = 0; i < tree->GetEntries(); i++) {
67 tree->GetEntry(i);
68 int nstrips = 768;
69 if (!side && layer != 3) nstrips = 512;
70
71 for (int iterStrip = 0; iterStrip < nstrips; iterStrip++) {
72 float occCal = hocc->GetBinContent(iterStrip + 1);
73 // B2INFO("Occupancy for" << layer << "." << ladder << "." << sensor << "." << side << ", strip:" << iterStrip << ": " << occCal);
74
75 if (occCal > occThr) {
76 vecHS[iterStrip] = 1;
77 stripOccAfterAbsCut[iterStrip] = 0;
78 } else stripOccAfterAbsCut[iterStrip] = occCal;
79 }
80
81 //iterative procedure
82 // B2INFO("Starting iterative procedure for hot strips finding");
83 bool moreHS = true;
84
85 while (moreHS && theHSFinder(stripOccAfterAbsCut, vecHS, nstrips)) {
86 moreHS = theHSFinder(stripOccAfterAbsCut, vecHS, nstrips);
87 }
88
89 for (int l = 0; l < nstrips; l++) {
90 isHotStrip = (int) vecHS[l];
91
92 payload->set(layer, ladder, sensor, bool(side), l, isHotStrip);
93 }
94 }
95
96 saveCalibration(payload, "SVDHotStripsCalibrations");
97
98 // probably not needed - would trigger re-doing the collection
99 //if ( ... too large corrections ... ) return c_Iterate;
100 return c_OK;
101}
102
103
104bool SVDHotStripsCalibrationsAlgorithm::theHSFinder(double* stripOccAfterAbsCut, int* hsflag, int nstrips)
105{
106 bool found = false;
107
108 if (m_base == -1)
109 m_base = nstrips;
110
111 int N = nstrips / m_base;
112
113 for (int sector = 0; sector < N; sector++) {
114
115 int nafter = 0;
116 double sensorOccAverage = 0;
117
118 for (int l = sector * m_base; l < sector * m_base + m_base; l++) {
119 sensorOccAverage = sensorOccAverage + stripOccAfterAbsCut[l];
120 if (stripOccAfterAbsCut[l] > 0) nafter++;
121 }
122 sensorOccAverage = sensorOccAverage / nafter;
123
124 B2DEBUG(1, "Average occupancy: " << sensorOccAverage);
125
126 for (int l = sector * m_base; l < sector * m_base + m_base; l++) {
127
128 // flag additional HS by comparing each strip occupancy with the sensor-based average occupancy
129
130 if (stripOccAfterAbsCut[l] > sensorOccAverage * m_relOccPrec) {
131 hsflag[l] = 1;
132 found = true;
133 stripOccAfterAbsCut[l] = 0;
134 }
135 // else hsflag[l]=0;
136 }
137 }
138
139 return found;
140}
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 successfully =0 in Python.
std::string m_id
Parameter given to set the UniqueID of the payload.
float m_relOccPrec
parameter to tue the finder algorithm
bool theHSFinder(double *stripOccAfterAbsCut, int *hsflag, int nstrips)
returns true if the strip is hot
SVDHotStripsCalibrationsAlgorithm(const std::string &str)
Constructor set the prefix to SVDHotStripsCalibrationsCollector.
virtual EResult calibrate() override
Run algo on data.
SVDCalibrationsBase< SVDCalibrationsBitmap > t_payload
typedef of the SVDHotStripsCalibrations payload for all SVD strips
Abstract base class for different kinds of events.