Belle II Software  release-05-01-25
SVDHitRateCounter.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2019 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Marko Staric, Hikaru Tanigawa *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 // Own include
12 #include <background/modules/BeamBkgHitRateMonitor/SVDHitRateCounter.h>
13 
14 // framework aux
15 #include <framework/logging/Logger.h>
16 #include <hlt/softwaretrigger/core/FinalTriggerDecisionCalculator.h>
17 
18 using namespace std;
19 
20 namespace Belle2 {
25  namespace Background {
26 
27  void SVDHitRateCounter::initialize(TTree* tree)
28  {
29  // register collection(s) as optional, your detector might be excluded in DAQ
30  m_digits.isOptional(m_svdShaperDigitsName);
31  m_clusters.isOptional();
32 
33  B2DEBUG(10, "SVDHitRateCounter: initialize()");
34  // set branch address
35  tree->Branch("svd", &m_rates,
36  "layerAverageRates[4]/F:layerLadderAverageRates[4][16]/F:layerSensorAverageRates[4][5]:averageRate/F:l3LadderSensorAverageRates[7][2]/F:numEvents/I:valid/O");
37  tree->Branch("svd_highE", &m_rates_highE,
38  "layerAverageRates[4]/F:layerLadderAverageRates[4][16]/F:layerSensorAverageRates[4][5]:averageRate/F:l3LadderSensorAverageRates[7][2]/F:numEvents/I:valid/O");
39  tree->Branch("svd_lowE", &m_rates_lowE,
40  "layerAverageRates[4]/F:layerLadderAverageRates[4][16]/F:layerSensorAverageRates[4][5]:averageRate/F:l3LadderSensorAverageRates[7][2]/F:numEvents/I:valid/O");
41 
42 
43  // count active strips
44  for (int layer = 0; layer < m_nLayers; layer++) {
45  for (int ladder = 0; ladder < m_nLadders[layer]; ladder++) {
46  for (int sensor = 0; sensor < m_nSensors[layer]; sensor++) {
47  VxdID sensorID(layer + 3, ladder + 1, sensor + 1);
48  for (bool isU : {true, false}) {
49  int nStrips = nStripsOnLayerSide(layer, isU);
50  for (int strip = 0; strip < nStrips; strip++) {
51  if (!m_HotStripsCalib.isHot(sensorID, isU, strip) && !m_FADCMaskedStrips.isMasked(sensorID, isU, strip)) {
52  m_activeStrips ++;
53  m_layerActiveStrips[layer] ++;
54  m_layerLadderActiveStrips[layer][ladder] ++;
55  m_layerSensorActiveStrips[layer][sensor] ++;
56  }
57  }
58  }
59  }
60  }
61  }
62  int layer = 0;
63  for (int ladder = 0; ladder < m_nLadders[layer]; ladder++) {
64  for (int sensor = 0; sensor < m_nSensors[layer]; sensor++) {
65  VxdID sensorID(layer + 3, ladder + 1, sensor + 1);
66  for (bool isU : {true, false}) {
67  int nStrips = nStripsOnLayerSide(layer, isU);
68  for (int strip = 0; strip < nStrips; strip++) {
69  if (!m_HotStripsCalib.isHot(sensorID, isU, strip) && !m_FADCMaskedStrips.isMasked(sensorID, isU, strip)) {
70  m_l3LadderSensorActiveStrips[ladder][sensor] ++;
71  }
72  }
73  }
74  }
75  }
76 
77  }
78 
79  void SVDHitRateCounter::clear()
80  {
81  m_buffer.clear();
82  }
83 
84  void SVDHitRateCounter::accumulate(unsigned timeStamp)
85  {
86  B2DEBUG(10, "SVDHitRateCounter: accumulate()");
87 
88  // check if the event has passed HLT filter
89  if (m_resultStoreObjectPointer.isValid()) {
90  const bool eventAccepted = SoftwareTrigger::FinalTriggerDecisionCalculator::getFinalTriggerDecision(*m_resultStoreObjectPointer);
91  if (!eventAccepted) return;
92  }
93 
94  // check if data are available
95  if (m_digits.isValid()) {
96 
97  // get buffer element
98  auto& rates = m_buffer[timeStamp];
99 
100  // increment event counter
101  rates.numEvents++;
102 
103  // accumulate hits
104  for (const auto& digit : m_digits) {
105  // select digits to count (usualy only good ones)
106  VxdID sensorID = digit.getSensorID();
107  int layer = sensorID.getLayerNumber() - 3;
108  int ladder = sensorID.getLadderNumber() - 1;
109  int sensor = sensorID.getSensorNumber() - 1;
110  rates.layerAverageRates[layer] ++;
111  rates.layerLadderAverageRates[layer][ladder] ++;
112  rates.layerSensorAverageRates[layer][sensor] ++;
113  rates.averageRate ++;
114  if (layer == 0)
115  rates.l3LadderSensorAverageRates[ladder][sensor] ++;
116  }
117 
118  // set flag to true to indicate the rates are valid
119  rates.valid = true;
120  }
121 
122  // check if data are available
123  if (m_clusters.isValid()) {
124 
125  // get buffer element
126  auto& rates_highE = m_buffer_highE[timeStamp];
127  auto& rates_lowE = m_buffer_lowE[timeStamp];
128 
129  // increment event counter
130  rates_highE.numEvents++;
131  rates_lowE.numEvents++;
132 
133  // accumulate clusters
134  for (const auto& cluster : m_clusters) {
135  VxdID sensorID = cluster.getSensorID();
136  int layer = sensorID.getLayerNumber() - 3;
137  int ladder = sensorID.getLadderNumber() - 1;
138  int sensor = sensorID.getSensorNumber() - 1;
139  if (cluster.getCharge() > m_thrCharge) {
140  rates_highE.layerAverageRates[layer] ++;
141  rates_highE.layerLadderAverageRates[layer][ladder] ++;
142  rates_highE.layerSensorAverageRates[layer][sensor] ++;
143  rates_highE.averageRate ++;
144  if (layer == 0)
145  rates_highE.l3LadderSensorAverageRates[ladder][sensor] ++;
146  } else {
147  rates_lowE.layerAverageRates[layer] ++;
148  rates_lowE.layerLadderAverageRates[layer][ladder] ++;
149  rates_lowE.layerSensorAverageRates[layer][sensor] ++;
150  rates_lowE.averageRate ++;
151  if (layer == 0)
152  rates_lowE.l3LadderSensorAverageRates[ladder][sensor] ++;
153  }
154  }
155 
156  // set flag to true to indicate the rates are valid
157  rates_highE.valid = true;
158  rates_lowE.valid = true;
159  }
160 
161  }
162 
163  void SVDHitRateCounter::normalize(unsigned timeStamp)
164  {
165  B2DEBUG(10, "SVDHitRateCounter: normalize()");
166  // copy buffer element
167  m_rates = m_buffer[timeStamp];
168  m_rates_highE = m_buffer_highE[timeStamp];
169  m_rates_lowE = m_buffer_lowE[timeStamp];
170 
171  SVDHitRateCounter::normalize_rates(m_rates);
172  SVDHitRateCounter::normalize_rates(m_rates_highE);
173  SVDHitRateCounter::normalize_rates(m_rates_lowE);
174  }
175 
176  void SVDHitRateCounter::normalize_rates(TreeStruct& rates)
177  {
178  if (not rates.valid) return;
179 
180  // normalize -> nHits on each segment in single event
181  rates.normalize();
182 
183  // convert to occupancy [%]
184  rates.averageRate /= m_activeStrips / 100;
185  for (int layer = 0; layer < m_nLayers; layer++) {
186  rates.layerAverageRates[layer] /= m_layerActiveStrips[layer] / 100;
187  for (int ladder = 0; ladder < m_nLadders[layer]; ladder++) {
188  rates.layerLadderAverageRates[layer][ladder] /= m_layerLadderActiveStrips[layer][ladder] / 100;
189  }
190  for (int sensor = 0; sensor < m_nSensors[layer]; sensor++) {
191  rates.layerSensorAverageRates[layer][sensor] /= m_layerSensorActiveStrips[layer][sensor] / 100;
192  }
193  }
194  int layer = 0;
195  for (int ladder = 0; ladder < m_nLadders[layer]; ladder++) {
196  for (int sensor = 0; sensor < m_nSensors[layer]; sensor++) {
197  rates.l3LadderSensorAverageRates[ladder][sensor] /= m_l3LadderSensorActiveStrips[ladder][sensor] / 100;
198  }
199  }
200  }
201 
202  } // background namespace
204 } // Belle2 namespace
205 
Belle2::VxdID
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:43
Belle2::Background::SVDHitRateCounter::TreeStruct
tree structure
Definition: SVDHitRateCounter.h:50
Belle2::Background::SVDHitRateCounter::TreeStruct::l3LadderSensorAverageRates
float l3LadderSensorAverageRates[7][2]
Layer 3 sensors [#ladder][#sensor].
Definition: SVDHitRateCounter.h:55
Belle2::VxdID::getLadderNumber
baseType getLadderNumber() const
Get the ladder id.
Definition: VxdID.h:108
Belle2::Background::SVDHitRateCounter::TreeStruct::normalize
void normalize()
normalize accumulated hits to single event
Definition: SVDHitRateCounter.h:62
Belle2::Background::SVDHitRateCounter::TreeStruct::averageRate
float averageRate
total SVD average occupancy
Definition: SVDHitRateCounter.h:54
Belle2::Background::SVDHitRateCounter::TreeStruct::layerLadderAverageRates
float layerLadderAverageRates[4][16]
[#layer][#ladder]
Definition: SVDHitRateCounter.h:52
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::Background::SVDHitRateCounter::TreeStruct::layerSensorAverageRates
float layerSensorAverageRates[4][5]
[#layer][#sensor]
Definition: SVDHitRateCounter.h:53
Belle2::VxdID::getSensorNumber
baseType getSensorNumber() const
Get the sensor id.
Definition: VxdID.h:110
Belle2::Background::SVDHitRateCounter::TreeStruct::valid
bool valid
status: true = rates valid
Definition: SVDHitRateCounter.h:57
Belle2::VxdID::getLayerNumber
baseType getLayerNumber() const
Get the layer id.
Definition: VxdID.h:106
Belle2::Background::SVDHitRateCounter::TreeStruct::layerAverageRates
float layerAverageRates[4]
layer average occupancy
Definition: SVDHitRateCounter.h:51