Belle II Software development
PXDHitRateCounter.cc
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// Own header.
10#include <background/modules/BeamBkgHitRateMonitor/PXDHitRateCounter.h>
11#include <pxd/reconstruction/PXDGainCalibrator.h>
12#include <pxd/reconstruction/PXDPixelMasker.h>
13#include <string>
14
15// framework aux
16#include <framework/logging/Logger.h>
17
18using namespace std;
19
20namespace Belle2 {
25 namespace Background {
26
28 {
29 //Pointer to GeoTools instance
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
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
54 }
55
57 {
58 m_buffer.clear();
59 }
60
61 void PXDHitRateCounter::accumulate(unsigned timeStamp)
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
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 }
130
131 void PXDHitRateCounter::normalize(unsigned timeStamp)
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
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
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 segements
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 }
178
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 segements
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 }
221
222 } // Background namespace
224} // Belle2 namespace
225
double m_segmentActivePixels[240]
number of active pixels in v segements
std::map< unsigned, TreeStruct > m_buffer
average rates in time stamps
bool m_maskDeadPixels
Correct bg rates by taking into account masked pixels.
double m_integrationTime
Integration time of PXD in ns.
virtual void initialize(TTree *tree) override
Class initializer: set branch addresses and other staf.
const double c_densitySi
Density of crystalline Silicon.
StoreArray< PXDCluster > m_clusters
collection of clusters
double m_activePixels[40]
number of active pixels in sensor
virtual void accumulate(unsigned timeStamp) override
Accumulate hits.
const PXD::SensorInfo & getInfo(VxdID sensorID) const
Get PXD::SensorInfo.
StoreArray< PXDDigit > m_digits
collection of digits
void setActivePixels()
Sets fractions of active channels.
double m_activeAreas[40]
area of active pixels in sensor
double m_segmentActiveAreas[240]
area of active pixels in v segments
virtual void normalize(unsigned timeStamp) override
Normalize accumulated hits (e.g.
virtual void clear() override
Clear time-stamp buffer to prepare for 'accumulate'.
The PXD Cluster class This class stores all information about reconstructed PXD clusters The position...
Definition: PXDCluster.h:30
The PXD digit class.
Definition: PXDDigit.h:27
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.
static PXDPixelMasker & getInstance()
Main (and only) way to access the PXDPixelMasker.
bool isOptional(const std::string &name="")
Tell the DataStore about an optional input.
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 const double s
[second]
Definition: Unit.h:95
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:142
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
Abstract base class for different kinds of events.
STL namespace.
float doseRates[40]
mean dose rate from PXDDigits per sensor [Gy/s]
float meanOccupancies[40]
mean hit occupancy from PXDDigits per sensor [Hits/Channel]
float softPhotonFluxes[40]
mean soft photon flux per sensor (Single pixel cluster <10keV) [clusters/cm2/s]
float chargedFluxes[40]
mean charged particle flux per sensor (Multi pixel cluster >10keV) [clusters/cm2/s]
float segmentChargedFluxes[240]
mean charged particle flux per v segment of sensor (Multi pixel cluster >10keV) [clusters/cm2/s]
float segmentHardPhotonFluxes[240]
mean hard photon flux per v segment of sensor (Single pixel cluster >10keV) [clusters/cm2/s]
float hardPhotonFluxes[40]
mean hard photon flux per sensor (Single pixel cluster >10keV) [clusters/cm2/s]
float averageRate
total detector average hit rate
float segmentSoftPhotonFluxes[240]
mean soft photon flux per v segment of sensor (Single pixel cluster <10keV) [clusters/cm2/s]
float segmentDoseRates[240]
mean dose rate from PXDDigits per v segment of sensor [Gy/s]