Belle II Software  release-05-02-19
SVDHitRateCounter.h
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, Ludovico Massaccesi *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #pragma once
12 
13 #include <background/modules/BeamBkgHitRateMonitor/HitRateBase.h>
14 #include <framework/datastore/StoreArray.h>
15 #include <mdst/dataobjects/SoftwareTriggerResult.h>
16 #include <svd/dataobjects/SVDShaperDigit.h>
17 #include <svd/dataobjects/SVDCluster.h>
18 #include <svd/dataobjects/SVDEventInfo.h>
19 #include <svd/calibration/SVDHotStripsCalibrations.h>
20 #include <svd/calibration/SVDFADCMaskedStrips.h>
21 #include <TTree.h>
22 #include <map>
23 #include <string>
24 
25 
26 namespace Belle2 {
31  namespace Background {
32 
36  class SVDHitRateCounter: public HitRateBase {
37 
38  public:
39 
43  struct TreeStruct {
44  float layerAverageRates[4] = {0};
45  float layerLadderAverageRates[4][16] = {0};
46  float layerSensorAverageRates[4][5] = {0};
47  float averageRate = 0;
48  float l3LadderSensorAverageRates[7][2] = {0};
49  int numEvents = 0;
50  bool valid = false;
55  void normalize()
56  {
57  if (numEvents == 0) return;
59  for (auto& rate : layerAverageRates) rate /= numEvents;
60  for (auto& row : layerLadderAverageRates) {
61  for (auto& rate : row) {
62  rate /= numEvents;
63  }
64  }
65  for (auto& row : layerSensorAverageRates) {
66  for (auto& rate : row) {
67  rate /= numEvents;
68  }
69  }
70  for (auto& row : l3LadderSensorAverageRates) {
71  for (auto& rate : row) {
72  rate /= numEvents;
73  }
74  }
75  }
76  };
77 
85  SVDHitRateCounter(const std::string& svdShaperDigitsName, double thrCharge,
86  bool ignoreHotStripsPayload = false,
87  bool ignoreMaskedStripsPayload = false):
88  m_svdShaperDigitsName(svdShaperDigitsName), m_thrCharge(thrCharge),
89  m_ignoreHotStripsPayload(ignoreHotStripsPayload),
90  m_ignoreMaskedStripsPayload(ignoreMaskedStripsPayload)
91  {}
92 
97  virtual void initialize(TTree* tree) override;
98 
102  virtual void clear() override;
103 
108  virtual void accumulate(unsigned timeStamp) override;
109 
114  virtual void normalize(unsigned timeStamp) override;
115 
122  void normalizeRates(TreeStruct& rates, bool isU = false, bool isV = false);
123 
135  void normalizeEnergyRates(TreeStruct& rates);
136 
142  int nStripsOnLayerSide(int layer, bool isU)
143  {
144  if (!isU && layer > 0) return 512; // V side on Layer 4,5,6
145  else return 768;
146  }
147 
154  double massOfSensor(int layer, int ladder, int sensor);
155 
164  bool isStripActive(const VxdID& sensorID, const bool& isU, const unsigned short& strip);
165 
166  private:
167 
168  // class parameters: to be set via constructor or setters
169  int m_nLayers = 4;
170  int m_nLadders[4] = {7, 10, 12, 16};
171  int m_nSensors[4] = {2, 3, 4, 5};
173  // tree structure
174  TreeStruct m_rates;
175  TreeStruct m_ratesU;
176  TreeStruct m_ratesV;
184  // buffer
185  std::map<unsigned, TreeStruct> m_buffer;
186  std::map<unsigned, TreeStruct> m_bufferU;
187  std::map<unsigned, TreeStruct> m_bufferV;
188  std::map<unsigned, TreeStruct> m_buffer_highE;
189  std::map<unsigned, TreeStruct> m_buffer_lowE;
190  std::map<unsigned, TreeStruct> m_buffer_clustersU;
191  std::map<unsigned, TreeStruct> m_buffer_clustersV;
192  std::map<unsigned, TreeStruct> m_buffer_energyU;
193  std::map<unsigned, TreeStruct> m_buffer_energyV;
195  // collections
199  // store object
203  // DB payloads
207  // other
208  std::string m_svdShaperDigitsName;
210  int m_activeStrips = 0;
211  int m_layerActiveStrips[4] = {0};
212  int m_layerLadderActiveStrips[4][16] = {0};
213  int m_layerSensorActiveStrips[4][5] = {0};
214  int m_l3LadderSensorActiveStrips[7][2] = {0};
216  int m_activeStripsU = 0;
217  int m_layerActiveStripsU[4] = {0};
218  int m_layerLadderActiveStripsU[4][16] = {0};
219  int m_layerSensorActiveStripsU[4][5] = {0};
222  int m_activeStripsV = 0;
223  int m_layerActiveStripsV[4] = {0};
224  int m_layerLadderActiveStripsV[4][16] = {0};
225  int m_layerSensorActiveStripsV[4][5] = {0};
228  double m_thrCharge = 0;
232  double m_massKg = 0;
233  double m_layerMassKg[4] = {0};
234  double m_layerLadderMassKg[4][16] = {0};
235  double m_layerSensorMassKg[4][5] = {0};
236  // No need for m_l3LadderSensorMassKg, massOfSensor() can be used.
237 
245  static constexpr double c_SVDSamplingClockFrequency = 0.03175;
246  };
247 
248  } // Background namespace
250 } // Belle2 namespace
Belle2::Background::SVDHitRateCounter::m_buffer_energyV
std::map< unsigned, TreeStruct > m_buffer_energyV
Average deposited energy (V-side) per event in time stamps.
Definition: SVDHitRateCounter.h:201
Belle2::Background::SVDHitRateCounter::nStripsOnLayerSide
int nStripsOnLayerSide(int layer, bool isU)
Return number of strips on a sensor.
Definition: SVDHitRateCounter.h:150
Belle2::SVDFADCMaskedStrips
This class defines the dbobject and the method to access strips which are masked at FADC level.
Definition: SVDFADCMaskedStrips.h:44
Belle2::Background::SVDHitRateCounter::m_resultStoreObjectPointer
StoreObjPtr< SoftwareTriggerResult > m_resultStoreObjectPointer
trigger decision
Definition: SVDHitRateCounter.h:208
Belle2::Background::SVDHitRateCounter::m_nLayers
int m_nLayers
number of layers
Definition: SVDHitRateCounter.h:177
Belle2::Background::SVDHitRateCounter::TreeStruct
tree structure
Definition: SVDHitRateCounter.h:51
Belle2::Background::SVDHitRateCounter::m_layerLadderActiveStripsU
int m_layerLadderActiveStripsU[4][16]
number of active U-strips in each layer, ladder
Definition: SVDHitRateCounter.h:226
Belle2::Background::SVDHitRateCounter::m_layerSensorActiveStripsU
int m_layerSensorActiveStripsU[4][5]
number of active U-strips in each layer, sensor position
Definition: SVDHitRateCounter.h:227
Belle2::Background::SVDHitRateCounter::m_layerLadderActiveStripsV
int m_layerLadderActiveStripsV[4][16]
number of active V-strips in each layer, ladder
Definition: SVDHitRateCounter.h:232
Belle2::Background::SVDHitRateCounter::c_SVDSamplingClockFrequency
static constexpr double c_SVDSamplingClockFrequency
SVD Sampling Clock frequency (approximated) in GHz (standard unit of frequency in basf2).
Definition: SVDHitRateCounter.h:253
Belle2::Background::SVDHitRateCounter::m_eventInfo
StoreObjPtr< SVDEventInfo > m_eventInfo
For number of APV samples taken (3 or 6).
Definition: SVDHitRateCounter.h:209
Belle2::Background::SVDHitRateCounter::TreeStruct::numEvents
int numEvents
number of events accumulated
Definition: SVDHitRateCounter.h:57
Belle2::Background::SVDHitRateCounter::m_ratesV
TreeStruct m_ratesV
tree variables for fired V-strips
Definition: SVDHitRateCounter.h:184
Belle2::Background::SVDHitRateCounter::m_thrCharge
double m_thrCharge
cut on cluster energy in electrons
Definition: SVDHitRateCounter.h:236
Belle2::Background::SVDHitRateCounter::m_layerLadderMassKg
double m_layerLadderMassKg[4][16]
Active mass of each ladder of each layer, in Kg.
Definition: SVDHitRateCounter.h:242
Belle2::Background::SVDHitRateCounter::TreeStruct::l3LadderSensorAverageRates
float l3LadderSensorAverageRates[7][2]
Layer 3 sensors [#ladder][#sensor].
Definition: SVDHitRateCounter.h:56
Belle2::Background::SVDHitRateCounter::m_l3LadderSensorActiveStripsV
int m_l3LadderSensorActiveStripsV[7][2]
number of active V-strips in each sensor in Layer 3
Definition: SVDHitRateCounter.h:234
Belle2::Background::SVDHitRateCounter::m_nSensors
int m_nSensors[4]
number of sensors on a ladder on each layer
Definition: SVDHitRateCounter.h:179
Belle2::Background::SVDHitRateCounter::m_clusters
StoreArray< SVDCluster > m_clusters
collection of clusters
Definition: SVDHitRateCounter.h:205
Belle2::Background::SVDHitRateCounter::normalize
virtual void normalize(unsigned timeStamp) override
Normalize accumulated hits (e.g.
Definition: SVDHitRateCounter.cc:288
Belle2::Background::SVDHitRateCounter::m_layerActiveStripsU
int m_layerActiveStripsU[4]
number of active U-strips in each layer
Definition: SVDHitRateCounter.h:225
Belle2::Background::SVDHitRateCounter::m_buffer_clustersV
std::map< unsigned, TreeStruct > m_buffer_clustersV
average cluster occupancies (V-side) in time stamps
Definition: SVDHitRateCounter.h:199
Belle2::Background::SVDHitRateCounter::m_l3LadderSensorActiveStripsU
int m_l3LadderSensorActiveStripsU[7][2]
number of active U-strips in each sensor in Layer 3
Definition: SVDHitRateCounter.h:228
Belle2::Background::SVDHitRateCounter::TreeStruct::normalize
void normalize()
normalize accumulated hits to single event
Definition: SVDHitRateCounter.h:63
Belle2::Background::SVDHitRateCounter::normalizeRates
void normalizeRates(TreeStruct &rates, bool isU=false, bool isV=false)
Normalize TreeStruct.
Definition: SVDHitRateCounter.cc:313
Belle2::Background::SVDHitRateCounter::m_ignoreMaskedStripsPayload
bool m_ignoreMaskedStripsPayload
SVD: count FAD-masked strips as active.
Definition: SVDHitRateCounter.h:238
Belle2::Background::SVDHitRateCounter::m_layerActiveStripsV
int m_layerActiveStripsV[4]
number of active V-strips in each layer
Definition: SVDHitRateCounter.h:231
Belle2::Background::SVDHitRateCounter::m_rates_highE
TreeStruct m_rates_highE
tree variables for high-energy clusters
Definition: SVDHitRateCounter.h:185
Belle2::Background::SVDHitRateCounter::TreeStruct::averageRate
float averageRate
total SVD average occupancy
Definition: SVDHitRateCounter.h:55
Belle2::Background::SVDHitRateCounter::m_digits
StoreArray< SVDShaperDigit > m_digits
collection of digits
Definition: SVDHitRateCounter.h:204
Belle2::Background::SVDHitRateCounter::m_rates_energyU
TreeStruct m_rates_energyU
Tree variables for deposited charge per unit time, then converted to dose rate (U-side)
Definition: SVDHitRateCounter.h:189
Belle2::Background::SVDHitRateCounter::m_buffer_highE
std::map< unsigned, TreeStruct > m_buffer_highE
average cluster occupancies (high energy) in time stamps
Definition: SVDHitRateCounter.h:196
Belle2::Background::SVDHitRateCounter::isStripActive
bool isStripActive(const VxdID &sensorID, const bool &isU, const unsigned short &strip)
Returns wether a strips is active (neither hot nor masked), taking into account the ignoreHotStrips a...
Definition: SVDHitRateCounter.cc:383
Belle2::Background::SVDHitRateCounter::m_buffer_clustersU
std::map< unsigned, TreeStruct > m_buffer_clustersU
average cluster occupancies (U-side) in time stamps
Definition: SVDHitRateCounter.h:198
Belle2::Background::SVDHitRateCounter::m_activeStripsV
int m_activeStripsV
number of active V-strips
Definition: SVDHitRateCounter.h:230
Belle2::SVDHotStripsCalibrations
This class defines the wrapper to retrieve the the list of the hot strips flgged offline.
Definition: SVDHotStripsCalibrations.h:44
Belle2::Background::SVDHitRateCounter::m_layerSensorMassKg
double m_layerSensorMassKg[4][5]
Active mass of each ladder/sensor position, in Kg.
Definition: SVDHitRateCounter.h:243
Belle2::Background::SVDHitRateCounter::TreeStruct::layerLadderAverageRates
float layerLadderAverageRates[4][16]
[#layer][#ladder]
Definition: SVDHitRateCounter.h:53
Belle2::Background::SVDHitRateCounter::m_l3LadderSensorActiveStrips
int m_l3LadderSensorActiveStrips[7][2]
number of active strips in each sensor in Layer 3
Definition: SVDHitRateCounter.h:222
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::StoreObjPtr
Type-safe access to single objects in the data store.
Definition: ParticleList.h:33
Belle2::Background::SVDHitRateCounter::m_ignoreHotStripsPayload
bool m_ignoreHotStripsPayload
count hot strips as active
Definition: SVDHitRateCounter.h:237
Belle2::Background::SVDHitRateCounter::m_activeStrips
int m_activeStrips
number of active strips
Definition: SVDHitRateCounter.h:218
Belle2::Background::SVDHitRateCounter::clear
virtual void clear() override
Clear time-stamp buffer to prepare for 'accumulate'.
Definition: SVDHitRateCounter.cc:139
Belle2::Background::SVDHitRateCounter::m_layerSensorActiveStripsV
int m_layerSensorActiveStripsV[4][5]
number of active V-strips in each layer, sensor position
Definition: SVDHitRateCounter.h:233
Belle2::Background::SVDHitRateCounter::normalizeEnergyRates
void normalizeEnergyRates(TreeStruct &rates)
Normalize a TreeStruct that stores charge/energy, not hits.
Definition: SVDHitRateCounter.cc:349
Belle2::Background::SVDHitRateCounter::TreeStruct::layerSensorAverageRates
float layerSensorAverageRates[4][5]
[#layer][#sensor]
Definition: SVDHitRateCounter.h:54
Belle2::Background::SVDHitRateCounter::m_svdShaperDigitsName
std::string m_svdShaperDigitsName
name of the input SVDShaperDigits collection
Definition: SVDHitRateCounter.h:216
Belle2::Background::SVDHitRateCounter::m_buffer
std::map< unsigned, TreeStruct > m_buffer
average strip occupancies in time stamps
Definition: SVDHitRateCounter.h:193
Belle2::Background::SVDHitRateCounter::m_rates
TreeStruct m_rates
tree variables for fired strips
Definition: SVDHitRateCounter.h:182
Belle2::Background::SVDHitRateCounter::m_ratesU
TreeStruct m_ratesU
tree variables for fired U-strips
Definition: SVDHitRateCounter.h:183
Belle2::Background::SVDHitRateCounter::massOfSensor
double massOfSensor(int layer, int ladder, int sensor)
Returns the (active) mass of the given sensor in Kg.
Definition: SVDHitRateCounter.cc:373
Belle2::Background::SVDHitRateCounter::m_layerSensorActiveStrips
int m_layerSensorActiveStrips[4][5]
number of active strips in each layer, sensor position
Definition: SVDHitRateCounter.h:221
Belle2::Background::SVDHitRateCounter::m_layerMassKg
double m_layerMassKg[4]
Active mass of each layer in Kg.
Definition: SVDHitRateCounter.h:241
Belle2::Background::SVDHitRateCounter::m_clustersV
TreeStruct m_clustersV
tree variables for V-side clusters
Definition: SVDHitRateCounter.h:188
Belle2::Background::SVDHitRateCounter::TreeStruct::valid
bool valid
status: true = rates valid
Definition: SVDHitRateCounter.h:58
Belle2::Background::SVDHitRateCounter::m_massKg
double m_massKg
Active mass of the whole SVD in Kg.
Definition: SVDHitRateCounter.h:240
Belle2::Background::SVDHitRateCounter::m_clustersU
TreeStruct m_clustersU
tree variables for U-side clusters
Definition: SVDHitRateCounter.h:187
Belle2::Background::SVDHitRateCounter::m_layerLadderActiveStrips
int m_layerLadderActiveStrips[4][16]
number of active strips in each layer, ladder
Definition: SVDHitRateCounter.h:220
Belle2::Background::SVDHitRateCounter::m_HotStripsCalib
SVDHotStripsCalibrations m_HotStripsCalib
payload for hot strips
Definition: SVDHitRateCounter.h:212
Belle2::Background::SVDHitRateCounter::m_activeStripsU
int m_activeStripsU
number of active U-strips
Definition: SVDHitRateCounter.h:224
Belle2::Background::SVDHitRateCounter::accumulate
virtual void accumulate(unsigned timeStamp) override
Accumulate hits.
Definition: SVDHitRateCounter.cc:144
Belle2::Background::SVDHitRateCounter::m_bufferV
std::map< unsigned, TreeStruct > m_bufferV
average V-strip occupancies in time stamps
Definition: SVDHitRateCounter.h:195
Belle2::Background::SVDHitRateCounter::m_buffer_energyU
std::map< unsigned, TreeStruct > m_buffer_energyU
Average deposited energy (U-side) per event in time stamps.
Definition: SVDHitRateCounter.h:200
Belle2::Background::SVDHitRateCounter::m_FADCMaskedStrips
SVDFADCMaskedStrips m_FADCMaskedStrips
payload for strips masked on FADC level
Definition: SVDHitRateCounter.h:213
Belle2::StoreArray
Accessor to arrays stored in the data store.
Definition: ECLMatchingPerformanceExpertModule.h:33
Belle2::Background::SVDHitRateCounter::m_nLadders
int m_nLadders[4]
number of ladders on each layer
Definition: SVDHitRateCounter.h:178
Belle2::Background::SVDHitRateCounter::m_rates_energyV
TreeStruct m_rates_energyV
Tree variables for deposited charge per unit time, then converted to dose rate (V-side)
Definition: SVDHitRateCounter.h:190
Belle2::Background::SVDHitRateCounter::m_buffer_lowE
std::map< unsigned, TreeStruct > m_buffer_lowE
average cluster occupancies (low energy) in time stamps
Definition: SVDHitRateCounter.h:197
Belle2::Background::SVDHitRateCounter::m_rates_lowE
TreeStruct m_rates_lowE
tree variables for low-energy clusters
Definition: SVDHitRateCounter.h:186
Belle2::Background::SVDHitRateCounter::m_layerActiveStrips
int m_layerActiveStrips[4]
number of active strips in each layer
Definition: SVDHitRateCounter.h:219
Belle2::Background::SVDHitRateCounter::initialize
virtual void initialize(TTree *tree) override
Class initializer: set branch addresses and other staf.
Definition: SVDHitRateCounter.cc:32
Belle2::Background::SVDHitRateCounter::TreeStruct::layerAverageRates
float layerAverageRates[4]
layer average occupancy
Definition: SVDHitRateCounter.h:52
Belle2::Background::SVDHitRateCounter::m_bufferU
std::map< unsigned, TreeStruct > m_bufferU
average U-strip occupancies in time stamps
Definition: SVDHitRateCounter.h:194
Belle2::Background::SVDHitRateCounter::SVDHitRateCounter
SVDHitRateCounter(const std::string &svdShaperDigitsName, double thrCharge, bool ignoreHotStripsPayload=false, bool ignoreMaskedStripsPayload=false)
Constructor.
Definition: SVDHitRateCounter.h:93