Belle II Software development
SVDHitRateCounter.h
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#pragma once
10
11#include <background/modules/BeamBkgHitRateMonitor/HitRateBase.h>
12#include <framework/datastore/StoreArray.h>
13#include <framework/database/DBObjPtr.h>
14#include <framework/dbobjects/HardwareClockSettings.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
26namespace Belle2 {
31 namespace Background {
32
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
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
208 // other
212 int m_layerActiveStrips[4] = {0};
213 int m_layerLadderActiveStrips[4][16] = {{0}};
214 int m_layerSensorActiveStrips[4][5] = {{0}};
218 int m_layerActiveStripsU[4] = {0};
219 int m_layerLadderActiveStripsU[4][16] = {{0}};
220 int m_layerSensorActiveStripsU[4][5] = {{0}};
224 int m_layerActiveStripsV[4] = {0};
225 int m_layerLadderActiveStripsV[4][16] = {{0}};
226 int m_layerSensorActiveStripsV[4][5] = {{0}};
229 double m_thrCharge = 0;
233 double m_massKg = 0;
234 double m_layerMassKg[4] = {0};
235 double m_layerLadderMassKg[4][16] = {{0}};
236 double m_layerSensorMassKg[4][5] = {{0}};
237 // No need for m_l3LadderSensorMassKg, massOfSensor() can be used.
238
245 static constexpr double c_SVDSamplingClockFrequency = 0.03175;
246 };
247
248 } // Background namespace
250} // Belle2 namespace
Abstract base class for monitoring beam background hit rates All the monitoring classes must inherit ...
Definition: HitRateBase.h:24
Class for monitoring beam background hit rates of SVD.
int m_layerLadderActiveStripsV[4][16]
number of active V-strips in each layer, ladder
std::map< unsigned, TreeStruct > m_buffer
average strip occupancies in time stamps
void normalizeRates(TreeStruct &rates, bool isU=false, bool isV=false)
Normalize TreeStruct.
int m_activeStripsV
number of active V-strips
int m_activeStrips
number of active strips
double m_layerLadderMassKg[4][16]
Active mass of each ladder of each layer, in Kg.
double m_layerMassKg[4]
Active mass of each layer in Kg.
TreeStruct m_rates
tree variables for fired strips
std::map< unsigned, TreeStruct > m_buffer_energyV
Average deposited energy (V-side) per event in time stamps.
TreeStruct m_rates_highE
tree variables for high-energy clusters
std::string m_svdShaperDigitsName
name of the input SVDShaperDigits collection
int m_layerSensorActiveStripsV[4][5]
number of active V-strips in each layer, sensor position
int m_l3LadderSensorActiveStripsU[7][2]
number of active U-strips in each sensor in Layer 3
StoreObjPtr< SVDEventInfo > m_eventInfo
For number of APV samples taken (3 or 6).
TreeStruct m_rates_energyU
Tree variables for deposited charge per unit time, then converted to dose rate (U-side)
std::map< unsigned, TreeStruct > m_buffer_clustersV
average cluster occupancies (V-side) in time stamps
std::map< unsigned, TreeStruct > m_buffer_lowE
average cluster occupancies (low energy) in time stamps
TreeStruct m_ratesU
tree variables for fired U-strips
bool m_ignoreHotStripsPayload
count hot strips as active
int m_nLadders[4]
number of ladders on each layer
int m_l3LadderSensorActiveStrips[7][2]
number of active strips in each sensor in Layer 3
virtual void initialize(TTree *tree) override
Class initializer: set branch addresses and other staf.
int nStripsOnLayerSide(int layer, bool isU)
Return number of strips on a sensor.
int m_layerSensorActiveStrips[4][5]
number of active strips in each layer, sensor position
SVDFADCMaskedStrips m_FADCMaskedStrips
payload for strips masked on FADC level
bool m_ignoreMaskedStripsPayload
SVD: count FAD-masked strips as active.
static constexpr double c_SVDSamplingClockFrequency
SVD Sampling Clock frequency (approximated) in GHz (standard unit of frequency in basf2).
StoreArray< SVDCluster > m_clusters
collection of clusters
virtual void accumulate(unsigned timeStamp) override
Accumulate hits.
double m_layerSensorMassKg[4][5]
Active mass of each ladder/sensor position, in Kg.
std::map< unsigned, TreeStruct > m_buffer_highE
average cluster occupancies (high energy) in time stamps
double massOfSensor(int layer, int ladder, int sensor)
Returns the (active) mass of the given sensor in Kg.
double m_thrCharge
cut on cluster energy in electrons
SVDHotStripsCalibrations m_HotStripsCalib
payload for hot strips
int m_l3LadderSensorActiveStripsV[7][2]
number of active V-strips in each sensor in Layer 3
TreeStruct m_rates_energyV
Tree variables for deposited charge per unit time, then converted to dose rate (V-side)
OptionalDBObjPtr< HardwareClockSettings > m_clockSettings
hardware clock settings
double m_massKg
Active mass of the whole SVD in Kg.
int m_nSensors[4]
number of sensors on a ladder on each layer
int m_layerLadderActiveStrips[4][16]
number of active strips in each layer, ladder
StoreArray< SVDShaperDigit > m_digits
collection of digits
TreeStruct m_rates_lowE
tree variables for low-energy clusters
std::map< unsigned, TreeStruct > m_buffer_clustersU
average cluster occupancies (U-side) in time stamps
void normalizeEnergyRates(TreeStruct &rates)
Normalize a TreeStruct that stores charge/energy, not hits.
int m_activeStripsU
number of active U-strips
int m_layerSensorActiveStripsU[4][5]
number of active U-strips in each layer, sensor position
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...
std::map< unsigned, TreeStruct > m_bufferU
average U-strip occupancies in time stamps
int m_layerActiveStripsV[4]
number of active V-strips in each layer
SVDHitRateCounter(const std::string &svdShaperDigitsName, double thrCharge, bool ignoreHotStripsPayload=false, bool ignoreMaskedStripsPayload=false)
Constructor.
StoreObjPtr< SoftwareTriggerResult > m_resultStoreObjectPointer
trigger decision
int m_layerLadderActiveStripsU[4][16]
number of active U-strips in each layer, ladder
TreeStruct m_clustersV
tree variables for V-side clusters
TreeStruct m_clustersU
tree variables for U-side clusters
std::map< unsigned, TreeStruct > m_buffer_energyU
Average deposited energy (U-side) per event in time stamps.
std::map< unsigned, TreeStruct > m_bufferV
average V-strip occupancies in time stamps
virtual void normalize(unsigned timeStamp) override
Normalize accumulated hits (e.g.
virtual void clear() override
Clear time-stamp buffer to prepare for 'accumulate'.
int m_layerActiveStrips[4]
number of active strips in each layer
int m_layerActiveStripsU[4]
number of active U-strips in each layer
TreeStruct m_ratesV
tree variables for fired V-strips
Optional DBObjPtr: This class behaves the same as the DBObjPtr except that it will not raise errors w...
Definition: DBObjPtr.h:48
This class defines the dbobject and the method to access strips which are masked at FADC level.
This class defines the wrapper to retrieve the the list of the hot strips flgged offline.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
Abstract base class for different kinds of events.
float layerAverageRates[4]
layer average occupancy
float layerLadderAverageRates[4][16]
[#layer][#ladder]
float layerSensorAverageRates[4][5]
[#layer][#sensor]
float averageRate
total SVD average occupancy
void normalize()
normalize accumulated hits to single event
float l3LadderSensorAverageRates[7][2]
Layer 3 sensors [#ladder][#sensor].