Belle II Software  release-08-01-10
SVDHitRateCounter Class Reference

Class for monitoring beam background hit rates of SVD. More...

#include <SVDHitRateCounter.h>

Inheritance diagram for SVDHitRateCounter:
Collaboration diagram for SVDHitRateCounter:

Classes

struct  TreeStruct
 tree structure More...
 

Public Member Functions

 SVDHitRateCounter (const std::string &svdShaperDigitsName, double thrCharge, bool ignoreHotStripsPayload=false, bool ignoreMaskedStripsPayload=false)
 Constructor. More...
 
virtual void initialize (TTree *tree) override
 Class initializer: set branch addresses and other staf. More...
 
virtual void clear () override
 Clear time-stamp buffer to prepare for 'accumulate'.
 
virtual void accumulate (unsigned timeStamp) override
 Accumulate hits. More...
 
virtual void normalize (unsigned timeStamp) override
 Normalize accumulated hits (e.g. More...
 
void normalizeRates (TreeStruct &rates, bool isU=false, bool isV=false)
 Normalize TreeStruct. More...
 
void normalizeEnergyRates (TreeStruct &rates)
 Normalize a TreeStruct that stores charge/energy, not hits. More...
 
int nStripsOnLayerSide (int layer, bool isU)
 Return number of strips on a sensor. More...
 
double massOfSensor (int layer, int ladder, int sensor)
 Returns the (active) mass of the given sensor in Kg. More...
 
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 and ignoreMaskedStrips settings. More...
 

Private Attributes

int m_nLayers = 4
 number of layers
 
int m_nLadders [4] = {7, 10, 12, 16}
 number of ladders on each layer
 
int m_nSensors [4] = {2, 3, 4, 5}
 number of sensors on a ladder on each layer
 
TreeStruct m_rates
 tree variables for fired strips
 
TreeStruct m_ratesU
 tree variables for fired U-strips
 
TreeStruct m_ratesV
 tree variables for fired V-strips
 
TreeStruct m_rates_highE
 tree variables for high-energy clusters
 
TreeStruct m_rates_lowE
 tree variables for low-energy clusters
 
TreeStruct m_clustersU
 tree variables for U-side clusters
 
TreeStruct m_clustersV
 tree variables for V-side clusters
 
TreeStruct m_rates_energyU
 Tree variables for deposited charge per unit time, then converted to dose rate (U-side)
 
TreeStruct m_rates_energyV
 Tree variables for deposited charge per unit time, then converted to dose rate (V-side)
 
std::map< unsigned, TreeStructm_buffer
 average strip occupancies in time stamps
 
std::map< unsigned, TreeStructm_bufferU
 average U-strip occupancies in time stamps
 
std::map< unsigned, TreeStructm_bufferV
 average V-strip occupancies in time stamps
 
std::map< unsigned, TreeStructm_buffer_highE
 average cluster occupancies (high energy) in time stamps
 
std::map< unsigned, TreeStructm_buffer_lowE
 average cluster occupancies (low energy) in time stamps
 
std::map< unsigned, TreeStructm_buffer_clustersU
 average cluster occupancies (U-side) in time stamps
 
std::map< unsigned, TreeStructm_buffer_clustersV
 average cluster occupancies (V-side) in time stamps
 
std::map< unsigned, TreeStructm_buffer_energyU
 Average deposited energy (U-side) per event in time stamps.
 
std::map< unsigned, TreeStructm_buffer_energyV
 Average deposited energy (V-side) per event in time stamps.
 
StoreArray< SVDShaperDigitm_digits
 collection of digits
 
StoreArray< SVDClusterm_clusters
 collection of clusters
 
StoreObjPtr< SoftwareTriggerResultm_resultStoreObjectPointer
 trigger decision
 
StoreObjPtr< SVDEventInfom_eventInfo
 For number of APV samples taken (3 or 6).
 
SVDHotStripsCalibrations m_HotStripsCalib
 payload for hot strips
 
SVDFADCMaskedStrips m_FADCMaskedStrips
 payload for strips masked on FADC level
 
OptionalDBObjPtr< HardwareClockSettingsm_clockSettings
 hardware clock settings
 
std::string m_svdShaperDigitsName
 name of the input SVDShaperDigits collection
 
int m_activeStrips = 0
 number of active strips
 
int m_layerActiveStrips [4] = {0}
 number of active strips in each layer
 
int m_layerLadderActiveStrips [4][16] = {{0}}
 number of active strips in each layer, ladder
 
int m_layerSensorActiveStrips [4][5] = {{0}}
 number of active strips in each layer, sensor position
 
int m_l3LadderSensorActiveStrips [7][2] = {{0}}
 number of active strips in each sensor in Layer 3
 
int m_activeStripsU = 0
 number of active U-strips
 
int m_layerActiveStripsU [4] = {0}
 number of active U-strips in each layer
 
int m_layerLadderActiveStripsU [4][16] = {{0}}
 number of active U-strips in each layer, ladder
 
int m_layerSensorActiveStripsU [4][5] = {{0}}
 number of active U-strips in each layer, sensor position
 
int m_l3LadderSensorActiveStripsU [7][2] = {{0}}
 number of active U-strips in each sensor in Layer 3
 
int m_activeStripsV = 0
 number of active V-strips
 
int m_layerActiveStripsV [4] = {0}
 number of active V-strips in each layer
 
int m_layerLadderActiveStripsV [4][16] = {{0}}
 number of active V-strips in each layer, ladder
 
int m_layerSensorActiveStripsV [4][5] = {{0}}
 number of active V-strips in each layer, sensor position
 
int m_l3LadderSensorActiveStripsV [7][2] = {{0}}
 number of active V-strips in each sensor in Layer 3
 
double m_thrCharge = 0
 cut on cluster energy in electrons
 
bool m_ignoreHotStripsPayload
 count hot strips as active
 
bool m_ignoreMaskedStripsPayload
 SVD: count FAD-masked strips as active.
 
double m_massKg = 0
 Active mass of the whole SVD in Kg.
 
double m_layerMassKg [4] = {0}
 Active mass of each layer in Kg.
 
double m_layerLadderMassKg [4][16] = {{0}}
 Active mass of each ladder of each layer, in Kg.
 
double m_layerSensorMassKg [4][5] = {{0}}
 Active mass of each ladder/sensor position, in Kg.
 

Static Private Attributes

static constexpr double c_SVDSamplingClockFrequency = 0.03175
 SVD Sampling Clock frequency (approximated) in GHz (standard unit of frequency in basf2). More...
 

Detailed Description

Class for monitoring beam background hit rates of SVD.

Definition at line 36 of file SVDHitRateCounter.h.

Constructor & Destructor Documentation

◆ SVDHitRateCounter()

SVDHitRateCounter ( const std::string &  svdShaperDigitsName,
double  thrCharge,
bool  ignoreHotStripsPayload = false,
bool  ignoreMaskedStripsPayload = false 
)
inline

Constructor.

Parameters
svdShaperDigitsNamename of the input SVDShaperDigits collection
thrChargecut on cluster energy in electrons
ignoreHotStripsPayloadCount also hot strips as active
ignoreMaskedStripsPayloadCount also FADC-masked strips as active

Definition at line 85 of file SVDHitRateCounter.h.

87  :
88  m_svdShaperDigitsName(svdShaperDigitsName), m_thrCharge(thrCharge),
89  m_ignoreHotStripsPayload(ignoreHotStripsPayload),
90  m_ignoreMaskedStripsPayload(ignoreMaskedStripsPayload)
91  {}
std::string m_svdShaperDigitsName
name of the input SVDShaperDigits collection
bool m_ignoreHotStripsPayload
count hot strips as active
bool m_ignoreMaskedStripsPayload
SVD: count FAD-masked strips as active.
double m_thrCharge
cut on cluster energy in electrons

Member Function Documentation

◆ accumulate()

void accumulate ( unsigned  timeStamp)
overridevirtual

Accumulate hits.

Parameters
timeStamptime stamp

Implements HitRateBase.

Definition at line 144 of file SVDHitRateCounter.cc.

145  {
146  B2DEBUG(10, "SVDHitRateCounter: accumulate()");
147 
148  // check if the event has passed HLT filter
149  if (m_resultStoreObjectPointer.isValid()) {
151  if (!eventAccepted) return;
152  }
153 
154  // check if data are available
155  if (m_digits.isValid()) {
156 
157  // get buffer element
158  auto& rates = m_buffer[timeStamp];
159  auto& ratesU = m_bufferU[timeStamp];
160  auto& ratesV = m_bufferV[timeStamp];
161 
162  // increment event counter
163  rates.numEvents++;
164  ratesU.numEvents++;
165  ratesV.numEvents++;
166 
167  // accumulate hits
168  for (const auto& digit : m_digits) {
169  // select digits to count (usualy only good ones)
170  VxdID sensorID = digit.getSensorID();
171  int layer = sensorID.getLayerNumber() - 3;
172  int ladder = sensorID.getLadderNumber() - 1;
173  int sensor = sensorID.getSensorNumber() - 1;
174  rates.layerAverageRates[layer] ++;
175  rates.layerLadderAverageRates[layer][ladder] ++;
176  rates.layerSensorAverageRates[layer][sensor] ++;
177  rates.averageRate ++;
178  if (layer == 0)
179  rates.l3LadderSensorAverageRates[ladder][sensor] ++;
180 
181  if (digit.isUStrip()) {
182  ratesU.layerAverageRates[layer]++;
183  ratesU.layerLadderAverageRates[layer][ladder]++;
184  ratesU.layerSensorAverageRates[layer][sensor]++;
185  ratesU.averageRate++;
186  if (layer == 0)
187  ratesU.l3LadderSensorAverageRates[ladder][sensor]++;
188  } else {
189  ratesV.layerAverageRates[layer]++;
190  ratesV.layerLadderAverageRates[layer][ladder]++;
191  ratesV.layerSensorAverageRates[layer][sensor]++;
192  ratesV.averageRate++;
193  if (layer == 0)
194  ratesV.l3LadderSensorAverageRates[ladder][sensor]++;
195  }
196  }
197 
198  // set flag to true to indicate the rates are valid
199  rates.valid = true;
200  ratesU.valid = true;
201  ratesV.valid = true;
202  }
203 
204  // check if data are available
205  if (m_clusters.isValid()) {
206 
207  // get buffer element
208  auto& rates_highE = m_buffer_highE[timeStamp];
209  auto& rates_lowE = m_buffer_lowE[timeStamp];
210  auto& clustersU = m_buffer_clustersU[timeStamp];
211  auto& clustersV = m_buffer_clustersV[timeStamp];
212  auto& rates_energyU = m_buffer_energyU[timeStamp];
213  auto& rates_energyV = m_buffer_energyV[timeStamp];
214 
215  // increment event counter
216  rates_highE.numEvents++;
217  rates_lowE.numEvents++;
218  clustersU.numEvents++;
219  clustersV.numEvents++;
220  rates_energyU.numEvents++;
221  rates_energyV.numEvents++;
222 
223  // accumulate clusters
224  for (const auto& cluster : m_clusters) {
225  VxdID sensorID = cluster.getSensorID();
226  int layer = sensorID.getLayerNumber() - 3;
227  int ladder = sensorID.getLadderNumber() - 1;
228  int sensor = sensorID.getSensorNumber() - 1;
229  if (cluster.getCharge() > m_thrCharge) {
230  rates_highE.layerAverageRates[layer] ++;
231  rates_highE.layerLadderAverageRates[layer][ladder] ++;
232  rates_highE.layerSensorAverageRates[layer][sensor] ++;
233  rates_highE.averageRate ++;
234  if (layer == 0)
235  rates_highE.l3LadderSensorAverageRates[ladder][sensor] ++;
236  } else {
237  rates_lowE.layerAverageRates[layer] ++;
238  rates_lowE.layerLadderAverageRates[layer][ladder] ++;
239  rates_lowE.layerSensorAverageRates[layer][sensor] ++;
240  rates_lowE.averageRate ++;
241  if (layer == 0)
242  rates_lowE.l3LadderSensorAverageRates[ladder][sensor] ++;
243  }
244 
245  // Compute integration time for obtaining the dose rate
246  int nSamp = 6; // Fallback value
247  double svdSamplingClock = c_SVDSamplingClockFrequency; // Fallback value
248  // Replace number of samples with real value if available
249  if (m_eventInfo.isValid())
250  nSamp = m_eventInfo->getNSamples();
251  else
252  B2WARNING("SVDEventInfo not available: assuming 6 samples were taken.");
253  // Avoid division by zero if nSamp = 1. Valid values are only 6 and 3,
254  // so this line is only meant to prevent a crash in case of corrupted data.
255  if (nSamp < 2) nSamp = 2;
256  // Replace sampling frequency with real value if available
257  if (m_clockSettings.isValid())
258  svdSamplingClock = m_clockSettings->getClockFrequency(Const::SVD, "sampling");
259  else
260  B2WARNING("HardwareClockSettings not available: using approximated SVD sampling clock frequency.");
261  double integrationTimeSeconds = (nSamp - 1) / svdSamplingClock / Unit::s;
262  double chargePerUnitTime = cluster.getCharge() / integrationTimeSeconds;
263  if (cluster.isUCluster()) {
264  rates_energyU.layerAverageRates[layer] += chargePerUnitTime;
265  rates_energyU.layerLadderAverageRates[layer][ladder] += chargePerUnitTime;
266  rates_energyU.layerSensorAverageRates[layer][sensor] += chargePerUnitTime;
267  rates_energyU.averageRate += chargePerUnitTime;
268  if (layer == 0)
269  rates_energyU.l3LadderSensorAverageRates[ladder][sensor] += chargePerUnitTime;
270  clustersU.layerAverageRates[layer]++;
271  clustersU.layerLadderAverageRates[layer][ladder]++;
272  clustersU.layerSensorAverageRates[layer][sensor]++;
273  clustersU.averageRate++;
274  if (layer == 0)
275  clustersU.l3LadderSensorAverageRates[ladder][sensor]++;
276  } else {
277  rates_energyV.layerAverageRates[layer] += chargePerUnitTime;
278  rates_energyV.layerLadderAverageRates[layer][ladder] += chargePerUnitTime;
279  rates_energyV.layerSensorAverageRates[layer][sensor] += chargePerUnitTime;
280  rates_energyV.averageRate += chargePerUnitTime;
281  if (layer == 0)
282  rates_energyV.l3LadderSensorAverageRates[ladder][sensor] += chargePerUnitTime;
283  clustersV.layerAverageRates[layer]++;
284  clustersV.layerLadderAverageRates[layer][ladder]++;
285  clustersV.layerSensorAverageRates[layer][sensor]++;
286  clustersV.averageRate++;
287  if (layer == 0)
288  clustersV.l3LadderSensorAverageRates[ladder][sensor]++;
289  }
290  }
291 
292  // set flag to true to indicate the rates are valid
293  rates_highE.valid = true;
294  rates_lowE.valid = true;
295  clustersU.valid = true;
296  clustersV.valid = true;
297  rates_energyU.valid = true;
298  rates_energyV.valid = true;
299  }
300 
301  }
std::map< unsigned, TreeStruct > m_buffer
average strip occupancies in time stamps
std::map< unsigned, TreeStruct > m_buffer_energyV
Average deposited energy (V-side) per event in time stamps.
StoreObjPtr< SVDEventInfo > m_eventInfo
For number of APV samples taken (3 or 6).
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
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
std::map< unsigned, TreeStruct > m_buffer_highE
average cluster occupancies (high energy) in time stamps
OptionalDBObjPtr< HardwareClockSettings > m_clockSettings
hardware clock settings
StoreArray< SVDShaperDigit > m_digits
collection of digits
std::map< unsigned, TreeStruct > m_buffer_clustersU
average cluster occupancies (U-side) in time stamps
std::map< unsigned, TreeStruct > m_bufferU
average U-strip occupancies in time stamps
StoreObjPtr< SoftwareTriggerResult > m_resultStoreObjectPointer
trigger decision
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
static bool getFinalTriggerDecision(const SoftwareTriggerResult &result, bool forgetTotalResult=false)
Calculate the final cut decision using all "total_results" of all sub triggers in the software trigge...
bool isValid() const
Check wether the array was registered.
Definition: StoreArray.h:288
static const double s
[second]
Definition: Unit.h:95

◆ initialize()

void initialize ( TTree *  tree)
overridevirtual

Class initializer: set branch addresses and other staf.

Parameters
treea valid TTree pointer

Implements HitRateBase.

Definition at line 30 of file SVDHitRateCounter.cc.

◆ 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 and ignoreMaskedStrips settings.

Parameters
sensorIDThe VxdID of the sensor
isUThe side of the sensor
stripThe strip index

Definition at line 398 of file SVDHitRateCounter.cc.

◆ massOfSensor()

double massOfSensor ( int  layer,
int  ladder,
int  sensor 
)

Returns the (active) mass of the given sensor in Kg.

Parameters
layertha layer number (starting from 0, not 3)
ladderthe ladder number (starting from 0, not 1)
sensorthe sensor number (starting from 0, not 1)

Definition at line 388 of file SVDHitRateCounter.cc.

◆ normalize()

void normalize ( unsigned  timeStamp)
overridevirtual

Normalize accumulated hits (e.g.

transform to rates)

Parameters
timeStamptime stamp

Implements HitRateBase.

Definition at line 303 of file SVDHitRateCounter.cc.

◆ normalizeEnergyRates()

void normalizeEnergyRates ( TreeStruct rates)

Normalize a TreeStruct that stores charge/energy, not hits.

Parameters
ratesTreeStruct to be normalized.

Assumes the TreeStruct contains the total accumulated cluster charge (in e-), and uses Const::ehEnergy and the mass of the sensors to convert to the average dose rate per event in mrad/s.

Assumes the integration time is 155 ns (6 samples).

Definition at line 364 of file SVDHitRateCounter.cc.

◆ normalizeRates()

void normalizeRates ( TreeStruct rates,
bool  isU = false,
bool  isV = false 
)

Normalize TreeStruct.

Parameters
ratesTreeStruct to be normalized
isUWhether the TreeStruct only contains U-side hits
isVWhether the TreeStruct only contains V-side hits

Definition at line 328 of file SVDHitRateCounter.cc.

◆ nStripsOnLayerSide()

int nStripsOnLayerSide ( int  layer,
bool  isU 
)
inline

Return number of strips on a sensor.

Parameters
layerlayer number of the sensor (starting from 0)
isUtrue if the sensor is U side, false if V.

Definition at line 142 of file SVDHitRateCounter.h.

Member Data Documentation

◆ c_SVDSamplingClockFrequency

constexpr double c_SVDSamplingClockFrequency = 0.03175
staticconstexprprivate

SVD Sampling Clock frequency (approximated) in GHz (standard unit of frequency in basf2).

This is used as a fallback value when HardwareClockSettings' payload is not available.

Definition at line 245 of file SVDHitRateCounter.h.


The documentation for this class was generated from the following files: