Belle II Software development
PXDHitRateCounter Class Reference

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

#include <PXDHitRateCounter.h>

Inheritance diagram for PXDHitRateCounter:
HitRateBase

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.
 
virtual void clear () override
 Clear time-stamp buffer to prepare for 'accumulate'.
 
virtual void accumulate (unsigned timeStamp) override
 Accumulate hits.
 
virtual void normalize (unsigned timeStamp) override
 Normalize accumulated hits (e.g.
 

Private Member Functions

void setActivePixels ()
 Sets fractions of active channels.
 

Static Private Member Functions

static const PXD::SensorInfogetInfo (VxdID sensorID)
 Get PXD::SensorInfo.
 

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 segments
 
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.

Constructor & Destructor Documentation

◆ PXDHitRateCounter()

PXDHitRateCounter ( )
inline

Constructor.

Definition at line 59 of file PXDHitRateCounter.h.

60 {}

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 }
TString getInfo(const TObject *obj)
Get object info HTML (e.g.
Definition ObjectInfo.cc:55

◆ clear()

void clear ( )
overridevirtual

Clear time-stamp buffer to prepare for 'accumulate'.

Implements HitRateBase.

Definition at line 56 of file PXDHitRateCounter.cc.

57 {
58 m_buffer.clear();
59 }

◆ getInfo()

const PXD::SensorInfo & getInfo ( VxdID sensorID)
inlinestaticprivate

Get PXD::SensorInfo.

Definition at line 121 of file PXDHitRateCounter.h.

122 {
123 return dynamic_cast<const PXD::SensorInfo&>(VXD::GeoCache::getInstance().getSensorInfo(sensorID));
124 }

◆ 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.

28 {
29 //Pointer to GeoTools instance
30 auto gTools = VXD::GeoCache::getInstance().getGeoTools();
31 if (gTools->getNumberOfPXDLayers() == 0) {
32 B2FATAL("Missing geometry for PXD.");
33 }
34
35 m_nPXDSensors = gTools->getNumberOfPXDSensors();
36
37 // register collection(s) as optional, your detector might be excluded in DAQ
38 m_digits.isOptional();
39
40 // register collection(s) as optional, your detector might be excluded in DAQ
41 m_clusters.isOptional();
42
43 string leaflist =
44 "meanOccupancies[40]/F:maxOccupancies[40]/F:doseRates[40]/F:softPhotonFluxes[40]/F:hardPhotonFluxes[40]/F:chargedFluxes[40]/F:segmentDoseRates[240]/F:segmentSoftPhotonFluxes[240]/F:segmentHardPhotonFluxes[240]/F:segmentChargedFluxes[240]/F:averageRate/F:numEvents/I:valid/O";
45
46 // set branch address
47 tree->Branch("pxd", &m_rates, leaflist.c_str());
48
49 // check parameters
50 if (m_integrationTime <= 0) B2FATAL("invalid integration time window for PXD: " << m_integrationTime);
51
52 // set fractions of active channels
53 setActivePixels();
54 }

◆ 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.

132 {
133 // copy buffer element
134 m_rates = m_buffer[timeStamp];
135
136 if (not m_rates.valid) return;
137
138 if (m_rates.numEvents == 0) return;
139
140 // Average number of PXDDigits per 1Hz
141 m_rates.averageRate /= m_rates.numEvents;
142
143 // Total integration time of PXD per 1Hz
144 double currentComponentTime = m_rates.numEvents * m_integrationTime;
145
146 //Pointer to GeoTools instance
147 auto gTools = VXD::GeoCache::getInstance().getGeoTools();
148
149 // Compute sensor level averages per 1Hz
150 for (int index = 0; index < m_nPXDSensors; index++) {
151 VxdID sensorID = gTools->getSensorIDFromPXDIndex(index);
152 auto info = getInfo(sensorID);
153 double currentSensorMass = m_activeAreas[index] * info.getThickness() * c_densitySi;
154 double currentSensorArea = m_activeAreas[index];
155
156 m_rates.meanOccupancies[index] /= m_rates.numEvents;
157 if (currentSensorArea > 0) {
158 m_rates.doseRates[index] *= (1.0 / (currentComponentTime / Unit::s)) * (1000 / currentSensorMass);
159 m_rates.softPhotonFluxes[index] *= (1.0 / currentSensorArea) * (1.0 / (currentComponentTime / Unit::s));
160 m_rates.hardPhotonFluxes[index] *= (1.0 / currentSensorArea) * (1.0 / (currentComponentTime / Unit::s));
161 m_rates.chargedFluxes[index] *= (1.0 / currentSensorArea) * (1.0 / (currentComponentTime / Unit::s));
162
163 // Loop over v segments
164 for (int vBin = 0; vBin < 6; ++vBin) {
165 double currentSegmentMass = m_segmentActiveAreas[vBin + index * 6] * info.getThickness() * c_densitySi;
166 double currentSegmentArea = m_segmentActiveAreas[vBin + index * 6];
167 if (currentSegmentArea > 0) {
168 m_rates.segmentDoseRates[vBin + index * 6] *= (1.0 / (currentComponentTime / Unit::s)) * (1000 / currentSegmentMass);
169 m_rates.segmentSoftPhotonFluxes[vBin + index * 6] *= (1.0 / currentSegmentArea) * (1.0 / (currentComponentTime / Unit::s));
170 m_rates.segmentHardPhotonFluxes[vBin + index * 6] *= (1.0 / currentSegmentArea) * (1.0 / (currentComponentTime / Unit::s));
171 m_rates.segmentChargedFluxes[vBin + index * 6] *= (1.0 / currentSegmentArea) * (1.0 /
172 (currentComponentTime / Unit::s));
173 }
174 }
175 }
176 }
177 }

◆ setActivePixels()

void setActivePixels ( )
private

Sets fractions of active channels.

Definition at line 179 of file PXDHitRateCounter.cc.

180 {
181 //Pointer to GeoTools instance
182 auto gTools = VXD::GeoCache::getInstance().getGeoTools();
183
184 // Initialize active areas and active fractions for all sensors
185 for (int index = 0; index < m_nPXDSensors; index++) {
186 VxdID sensorID = gTools->getSensorIDFromPXDIndex(index);
187 auto info = getInfo(sensorID);
188
189 // Compute nominal number of pixel per sensor
190 m_activePixels[index] = info.getUCells() * info.getVCells();
191 // Compute nominal area of sensor
192 m_activeAreas[index] = info.getWidth() * info.getLength();
193
194 // Loop over v segments
195 for (int vBin = 0; vBin < 6; ++vBin) {
196 // Compute nominal number of pixel per segment
197 // Sgements have same number of rows, but due to pitch change a different length
198 m_segmentActivePixels[vBin + index * 6] = info.getUCells() * info.getVCells() / 6;
199 // Compute position of v segment
200 double v = info.getVCellPosition(vBin * info.getVCells() / 6);
201 // Compute nominal area of segment
202 m_segmentActiveAreas[vBin + index * 6] = info.getWidth() * info.getVPitch(v) * info.getVCells() / 6;
203 }
204
205 if (m_maskDeadPixels) {
206 for (int ui = 0; ui < info.getUCells(); ++ui) {
207 for (int vi = 0; vi < info.getVCells(); ++vi) {
208 if (PXD::PXDPixelMasker::getInstance().pixelDead(sensorID, ui, vi)
209 || !PXD::PXDPixelMasker::getInstance().pixelOK(sensorID, ui, vi)) {
210 m_activePixels[index] -= 1;
211 m_activeAreas[index] -= info.getVPitch(info.getVCellPosition(vi)) * info.getUPitch();
212 int vBin = PXD::PXDGainCalibrator::getInstance().getBinV(sensorID, vi, 6);
213 m_segmentActivePixels[vBin + index * 6] -= 1;
214 m_segmentActiveAreas[vBin + index * 6] -= info.getVPitch(info.getVCellPosition(vi)) * info.getUPitch();
215 }
216 }
217 }
218 }
219 }
220 }

Member Data Documentation

◆ c_densitySi

const double c_densitySi = 2.3290 * Unit::g_cm3
private

Density of crystalline Silicon.

Definition at line 102 of file PXDHitRateCounter.h.

◆ m_activeAreas

double m_activeAreas[40] = {0}
private

area of active pixels in sensor

Definition at line 116 of file PXDHitRateCounter.h.

116{0};

◆ m_activePixels

double m_activePixels[40] = {0}
private

number of active pixels in sensor

Definition at line 115 of file PXDHitRateCounter.h.

115{0};

◆ m_buffer

std::map<unsigned, TreeStruct> m_buffer
private

average rates in time stamps

Definition at line 108 of file PXDHitRateCounter.h.

◆ m_clusters

StoreArray<PXDCluster> m_clusters
private

collection of clusters

Definition at line 112 of file PXDHitRateCounter.h.

◆ m_digits

StoreArray<PXDDigit> m_digits
private

collection of digits

Definition at line 111 of file PXDHitRateCounter.h.

◆ m_integrationTime

double m_integrationTime = 20 * Unit::us
private

Integration time of PXD in ns.

Definition at line 98 of file PXDHitRateCounter.h.

◆ m_maskDeadPixels

bool m_maskDeadPixels = true
private

Correct bg rates by taking into account masked pixels.

Definition at line 99 of file PXDHitRateCounter.h.

◆ m_nPXDSensors

int m_nPXDSensors = 0
private

number of sensors

Definition at line 100 of file PXDHitRateCounter.h.

◆ m_rates

TreeStruct m_rates
private

tree variables

Definition at line 105 of file PXDHitRateCounter.h.

◆ m_segmentActiveAreas

double m_segmentActiveAreas[240] = {0}
private

area of active pixels in v segments

Definition at line 118 of file PXDHitRateCounter.h.

118{0};

◆ m_segmentActivePixels

double m_segmentActivePixels[240] = {0}
private

number of active pixels in v segments

Definition at line 117 of file PXDHitRateCounter.h.

117{0};

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