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