Belle II Software  release-05-01-25
PXDHitRateCounter Class Reference

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

#include <PXDHitRateCounter.h>

Inheritance diagram for PXDHitRateCounter:
Collaboration diagram for PXDHitRateCounter:

Classes

struct  TreeStruct
 tree structure More...
 

Public Member Functions

 PXDHitRateCounter ()
 Constructor.
 
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...
 

Private Member Functions

const PXD::SensorInfogetInfo (VxdID sensorID) const
 Get PXD::SensorInfo.
 
void setActivePixels ()
 Sets fractions of active channels.
 

Private Attributes

double m_integrationTime = 20 * Unit::us
 Integration time of PXD in ns.
 
bool m_maskDeadPixels = true
 Correct bg rates by taking into account masked pixels.
 
int m_nPXDSensors = 0
 number of sensors
 
const double c_densitySi = 2.3290 * Unit::g_cm3
 Density of crystalline Silicon.
 
TreeStruct m_rates
 tree variables
 
std::map< unsigned, TreeStructm_buffer
 average rates in time stamps
 
StoreArray< PXDDigitm_digits
 collection of digits
 
StoreArray< PXDClusterm_clusters
 collection of clusters
 
double m_activePixels [40] = {0}
 number of active pixels in sensor
 
double m_activeAreas [40] = {0}
 area of active pixels in sensor
 
double m_segmentActivePixels [240] = {0}
 number of active pixels in v segements
 
double m_segmentActiveAreas [240] = {0}
 area of active pixels in v segments
 

Detailed Description

Class for monitoring beam background hit rates of PXD.

Definition at line 43 of file PXDHitRateCounter.h.

Member Function Documentation

◆ accumulate()

void accumulate ( unsigned  timeStamp)
overridevirtual

Accumulate hits.

Parameters
timeStamptime stamp

Implements HitRateBase.

Definition at line 63 of file PXDHitRateCounter.cc.

64  {
65  // check if data are available
66  if ((not m_digits.isValid()) or (not m_clusters.isValid())) return;
67 
68  // get buffer element
69  auto& rates = m_buffer[timeStamp];
70 
71  // increment event counter
72  rates.numEvents++;
73 
74  // accumulate hits
75  rates.averageRate += m_digits.getEntries();
76 
77  // Use the same indexing of sensors as for PXDDQM
78  auto gTools = VXD::GeoCache::getInstance().getGeoTools();
79 
80  // accumulate dose rate and event occupancy
81  float occupancies[40] = {0};
82  for (const PXDDigit& storeDigit : m_digits) {
83  VxdID sensorID = storeDigit.getSensorID();
84  int index = gTools->getPXDSensorIndex(storeDigit.getSensorID());
85  int vBin = PXD::PXDGainCalibrator::getInstance().getBinV(sensorID, storeDigit.getVCellID(), 6);
86  double ADUToEnergy = PXD::PXDGainCalibrator::getInstance().getADUToEnergy(sensorID, storeDigit.getUCellID(),
87  storeDigit.getVCellID());
88  double hitEnergy = storeDigit.getCharge() * ADUToEnergy;
89  rates.doseRates[index] += (hitEnergy / Unit::J);
90  occupancies[index] += 1.0;
91  rates.segmentDoseRates[vBin + index * 6] += (hitEnergy / Unit::J);
92  }
93 
94  for (int index = 0; index < m_nPXDSensors; index++) {
95  if (m_activePixels[index] > 0) {
96  occupancies[index] /= m_activePixels[index];
97  }
98  rates.meanOccupancies[index] += occupancies[index];
99  if (rates.maxOccupancies[index] < occupancies[index]) {
100  rates.maxOccupancies[index] = occupancies[index];
101  }
102  }
103 
104  // accumulate fluxes
105  for (const PXDCluster& cluster : m_clusters) {
106  // Update if we have a new sensor
107  VxdID sensorID = cluster.getSensorID();
108  int index = gTools->getPXDSensorIndex(sensorID);
109  auto info = getInfo(sensorID);
110 
111  auto cluster_uID = info.getUCellID(cluster.getU());
112  auto cluster_vID = info.getVCellID(cluster.getV());
113  int vBin = PXD::PXDGainCalibrator::getInstance().getBinV(sensorID, cluster_vID, 6);
114  double ADUToEnergy = PXD::PXDGainCalibrator::getInstance().getADUToEnergy(sensorID, cluster_uID, cluster_vID);
115  double clusterEnergy = cluster.getCharge() * ADUToEnergy;
116 
117  if (cluster.getSize() == 1 && clusterEnergy < 10000 * Unit::eV && clusterEnergy > 6000 * Unit::eV) {
118  rates.softPhotonFluxes[index] += 1.0;
119  rates.segmentSoftPhotonFluxes[vBin + index * 6] += 1.0;
120  } else if (cluster.getSize() == 1 && clusterEnergy > 10000 * Unit::eV) {
121  rates.hardPhotonFluxes[index] += 1.0;
122  rates.segmentHardPhotonFluxes[vBin + index * 6] += 1.0;
123  } else if (cluster.getSize() > 1 && clusterEnergy > 10000 * Unit::eV) {
124  rates.chargedFluxes[index] += 1.0;
125  rates.segmentChargedFluxes[vBin + index * 6] += 1.0;
126  }
127  }
128 
129  // set flag to true to indicate the rates are valid
130  rates.valid = true;
131  }

◆ initialize()

void initialize ( TTree *  tree)
overridevirtual

Class initializer: set branch addresses and other staf.

Parameters
treea valid TTree pointer

Implements HitRateBase.

Definition at line 29 of file PXDHitRateCounter.cc.

◆ normalize()

void normalize ( unsigned  timeStamp)
overridevirtual

Normalize accumulated hits (e.g.

transform to rates)

Parameters
timeStamptime stamp

Implements HitRateBase.

Definition at line 133 of file PXDHitRateCounter.cc.


The documentation for this class was generated from the following files:
Belle2::Background::PXDHitRateCounter::getInfo
const PXD::SensorInfo & getInfo(VxdID sensorID) const
Get PXD::SensorInfo.
Definition: PXDHitRateCounter.h:131
Belle2::PXD::PXDGainCalibrator::getBinV
unsigned short getBinV(VxdID id, unsigned int vid) const
Get gain correction bin along sensor v side.
Definition: PXDGainCalibrator.cc:73
Belle2::Background::PXDHitRateCounter::m_digits
StoreArray< PXDDigit > m_digits
collection of digits
Definition: PXDHitRateCounter.h:121
Belle2::Unit::J
static const double J
[joule]
Definition: Unit.h:126
Belle2::Background::PXDHitRateCounter::m_clusters
StoreArray< PXDCluster > m_clusters
collection of clusters
Definition: PXDHitRateCounter.h:122
Belle2::Unit::eV
static const double eV
[electronvolt]
Definition: Unit.h:122
Belle2::PXD::PXDGainCalibrator::getADUToEnergy
float getADUToEnergy(VxdID id, unsigned int uid, unsigned int vid) const
Get conversion factor from ADU to energy.
Definition: PXDGainCalibrator.cc:45
Belle2::VXD::GeoCache::getInstance
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:215
Belle2::VXD::GeoCache::getGeoTools
const GeoTools * getGeoTools()
Return a raw pointer to a GeoTools object.
Definition: GeoCache.h:149
Belle2::StoreArray::isValid
bool isValid() const
Check wether the array was registered.
Definition: StoreArray.h:298
Belle2::Background::PXDHitRateCounter::m_buffer
std::map< unsigned, TreeStruct > m_buffer
average rates in time stamps
Definition: PXDHitRateCounter.h:118
Belle2::PXD::PXDGainCalibrator::getInstance
static PXDGainCalibrator & getInstance()
Main (and only) way to access the PXDGainCalibrator.
Definition: PXDGainCalibrator.cc:29
Belle2::Background::PXDHitRateCounter::m_activePixels
double m_activePixels[40]
number of active pixels in sensor
Definition: PXDHitRateCounter.h:125
Belle2::Background::PXDHitRateCounter::m_nPXDSensors
int m_nPXDSensors
number of sensors
Definition: PXDHitRateCounter.h:110