Belle II Software  release-08-01-10
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 33 of file PXDHitRateCounter.h.

Member Function Documentation

◆ accumulate()

void accumulate ( unsigned  timeStamp)
overridevirtual

Accumulate hits.

Parameters
timeStamptime stamp

Implements HitRateBase.

Definition at line 61 of file PXDHitRateCounter.cc.

62  {
63  // check if data are available
64  if ((not m_digits.isValid()) or (not m_clusters.isValid())) return;
65 
66  // get buffer element
67  auto& rates = m_buffer[timeStamp];
68 
69  // increment event counter
70  rates.numEvents++;
71 
72  // accumulate hits
73  rates.averageRate += m_digits.getEntries();
74 
75  // Use the same indexing of sensors as for PXDDQM
76  auto gTools = VXD::GeoCache::getInstance().getGeoTools();
77 
78  // accumulate dose rate and event occupancy
79  float occupancies[40] = {0};
80  for (const PXDDigit& storeDigit : m_digits) {
81  VxdID sensorID = storeDigit.getSensorID();
82  int index = gTools->getPXDSensorIndex(storeDigit.getSensorID());
83  int vBin = PXD::PXDGainCalibrator::getInstance().getBinV(sensorID, storeDigit.getVCellID(), 6);
84  double ADUToEnergy = PXD::PXDGainCalibrator::getInstance().getADUToEnergy(sensorID, storeDigit.getUCellID(),
85  storeDigit.getVCellID());
86  double hitEnergy = storeDigit.getCharge() * ADUToEnergy;
87  rates.doseRates[index] += (hitEnergy / Unit::J);
88  occupancies[index] += 1.0;
89  rates.segmentDoseRates[vBin + index * 6] += (hitEnergy / Unit::J);
90  }
91 
92  for (int index = 0; index < m_nPXDSensors; index++) {
93  if (m_activePixels[index] > 0) {
94  occupancies[index] /= m_activePixels[index];
95  }
96  rates.meanOccupancies[index] += occupancies[index];
97  if (rates.maxOccupancies[index] < occupancies[index]) {
98  rates.maxOccupancies[index] = occupancies[index];
99  }
100  }
101 
102  // accumulate fluxes
103  for (const PXDCluster& cluster : m_clusters) {
104  // Update if we have a new sensor
105  VxdID sensorID = cluster.getSensorID();
106  int index = gTools->getPXDSensorIndex(sensorID);
107  auto info = getInfo(sensorID);
108 
109  auto cluster_uID = info.getUCellID(cluster.getU());
110  auto cluster_vID = info.getVCellID(cluster.getV());
111  int vBin = PXD::PXDGainCalibrator::getInstance().getBinV(sensorID, cluster_vID, 6);
112  double ADUToEnergy = PXD::PXDGainCalibrator::getInstance().getADUToEnergy(sensorID, cluster_uID, cluster_vID);
113  double clusterEnergy = cluster.getCharge() * ADUToEnergy;
114 
115  if (cluster.getSize() == 1 && clusterEnergy < 10000 * Unit::eV && clusterEnergy > 6000 * Unit::eV) {
116  rates.softPhotonFluxes[index] += 1.0;
117  rates.segmentSoftPhotonFluxes[vBin + index * 6] += 1.0;
118  } else if (cluster.getSize() == 1 && clusterEnergy > 10000 * Unit::eV) {
119  rates.hardPhotonFluxes[index] += 1.0;
120  rates.segmentHardPhotonFluxes[vBin + index * 6] += 1.0;
121  } else if (cluster.getSize() > 1 && clusterEnergy > 10000 * Unit::eV) {
122  rates.chargedFluxes[index] += 1.0;
123  rates.segmentChargedFluxes[vBin + index * 6] += 1.0;
124  }
125  }
126 
127  // set flag to true to indicate the rates are valid
128  rates.valid = true;
129  }
std::map< unsigned, TreeStruct > m_buffer
average rates in time stamps
StoreArray< PXDCluster > m_clusters
collection of clusters
double m_activePixels[40]
number of active pixels in sensor
const PXD::SensorInfo & getInfo(VxdID sensorID) const
Get PXD::SensorInfo.
StoreArray< PXDDigit > m_digits
collection of digits
unsigned short getBinV(VxdID id, unsigned int vid) const
Get gain correction bin along sensor v side.
float getADUToEnergy(VxdID id, unsigned int uid, unsigned int vid) const
Get conversion factor from ADU to energy.
static PXDGainCalibrator & getInstance()
Main (and only) way to access the PXDGainCalibrator.
bool isValid() const
Check wether the array was registered.
Definition: StoreArray.h:288
static const double eV
[electronvolt]
Definition: Unit.h:112
static const double J
[joule]
Definition: Unit.h:116
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:214
const GeoTools * getGeoTools()
Return a raw pointer to a GeoTools object.
Definition: GeoCache.h:147

◆ initialize()

void initialize ( TTree *  tree)
overridevirtual

Class initializer: set branch addresses and other staf.

Parameters
treea valid TTree pointer

Implements HitRateBase.

Definition at line 27 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 131 of file PXDHitRateCounter.cc.


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