Belle II Software development
PXDBgTupleProducerModule.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
10
11#include <pxd/modules/pxdBackground/PXDBgTupleProducerModule.h>
12#include <framework/datastore/StoreObjPtr.h>
13#include <framework/datastore/StoreArray.h>
14#include <framework/dataobjects/EventMetaData.h>
15#include <pxd/reconstruction/PXDGainCalibrator.h>
16#include <pxd/reconstruction/PXDPixelMasker.h>
17
18#include <pxd/dataobjects/PXDDigit.h>
19#include <pxd/dataobjects/PXDCluster.h>
20#include <boost/format.hpp>
21
22#include <TFile.h>
23#include <TTree.h>
24
25using namespace std;
26using boost::format;
27using namespace Belle2;
28using namespace Belle2::PXD;
29
30
31
32//-----------------------------------------------------------------
33// Register the Module
34//-----------------------------------------------------------------
35REG_MODULE(PXDBgTupleProducer);
36
37
38//-----------------------------------------------------------------
39// Implementation
40//-----------------------------------------------------------------
41
43 , m_nPXDSensors(0), m_hasPXDData(false)
44{
45 //Set module properties
46 setDescription("PXD background tuple producer module");
47 addParam("integrationTime", m_integrationTime, "PXD integration time in micro seconds", double(20));
48 addParam("timePeriod", m_timePeriod, "Period for background time series in seconds.", double(1));
49 addParam("outputFileName", m_outputFileName, "Output file name", string("beast_tuple.root"));
50 addParam("maskDeadPixels", m_maskDeadPixels, "Correct bg rates by known dead pixels", bool(true));
51 addParam("nBinsU", m_nBinsU, "Number of regions per sensor along u side", int(1));
52 addParam("nBinsV", m_nBinsV, "Number of regions per sensor along v side", int(6));
53}
54
55
57{
58 //Register collections
61
62 //Make sure the EventMetaData already exists.
64
65 //Store names to speed up creation later
66 m_storeDigitsName = storeDigits.getName();
67
68 // PXD integration time
70
71 // Period for time series
73
74 // So far, we did not see PXD data
75 m_hasPXDData = false;
76
77 //Pointer to GeoTools instance
79 if (gTools->getNumberOfPXDLayers() == 0) {
80 B2WARNING("Missing geometry for PXD, PXD-masking is skiped.");
81 }
82 m_nPXDSensors = gTools->getNumberOfPXDSensors();
83
84 // Initialize m_sensorData with empty sensorData for all sensors
85 for (int i = 0; i < m_nPXDSensors; i++) {
86 VxdID sensorID = gTools->getSensorIDFromPXDIndex(i);
87 m_sensorData[sensorID] = SensorData();
88 // Start value for minOccupancy should be one not zero
89 m_sensorData[sensorID].m_minOccupancy = 1.0;
90 // Initialize counters for subdivisions per sensor
91 m_sensorData[sensorID].m_regionExpoMap = vector<double>(m_nBinsU * m_nBinsV, 0);
92 m_sensorData[sensorID].m_regionDoseMap = vector<double>(m_nBinsU * m_nBinsV, 0);
93 m_sensorData[sensorID].m_regionSoftPhotonFluxMap = vector<double>(m_nBinsU * m_nBinsV, 0);
94 m_sensorData[sensorID].m_regionChargedParticleFluxMap = vector<double>(m_nBinsU * m_nBinsV, 0);
95 m_sensorData[sensorID].m_regionHardPhotonFluxMap = vector<double>(m_nBinsU * m_nBinsV, 0);
96 }
97}
98
100{
101 // Compute the sensitive area for all PXD sensors
102 for (auto const& pair2 : m_sensorData) {
103 auto const& sensorID = pair2.first;
104 auto info = getInfo(sensorID);
105
106 // Compute nominal number of pixel per sensor
107 m_sensitivePixelMap[sensorID] = info.getUCells() * info.getVCells();
108 // Compute nominal area per sensor
109 m_sensitiveAreaMap[sensorID] = getSensorArea(sensorID);
110
111 for (int uBin = 0; uBin < m_nBinsU; ++uBin) {
112 for (int vBin = 0; vBin < m_nBinsV; ++vBin) {
113 std::pair<VxdID, int> key(sensorID, getRegionID(uBin, vBin));
114 // Compute nominal number of pixel per sensor subregion
115 m_regionSensitivePixelMap[key] = info.getUCells() * info.getVCells() / m_nBinsU / m_nBinsV;
116 // Compute nominal area per sensor subregion
117 m_regionSensitiveAreaMap[key] = getRegionArea(sensorID, vBin);
118 }
119 }
120
121 if (m_maskDeadPixels) {
122 for (int ui = 0; ui < info.getUCells(); ++ui) {
123 for (int vi = 0; vi < info.getVCells(); ++vi) {
124 if (PXD::PXDPixelMasker::getInstance().pixelDead(sensorID, ui, vi)
125 || !PXD::PXDPixelMasker::getInstance().pixelOK(sensorID, ui, vi)) {
126 m_sensitivePixelMap[sensorID] -= 1;
127 m_sensitiveAreaMap[sensorID] -= info.getVPitch(info.getVCellPosition(vi)) * info.getUPitch();
128 int uBin = PXDGainCalibrator::getInstance().getBinU(sensorID, ui, vi, m_nBinsU);
129 int vBin = PXDGainCalibrator::getInstance().getBinV(sensorID, vi, m_nBinsV);
130 std::pair<VxdID, int> key(sensorID, getRegionID(uBin, vBin));
132 m_regionSensitiveAreaMap[key] -= info.getVPitch(info.getVCellPosition(vi)) * info.getUPitch();
133 }
134 }
135 }
136 }
137
138 if (m_sensitivePixelMap[sensorID] == 0) {
139 B2WARNING("All pixels from Sensor=" << sensorID << " are masked.");
140 }
141
142 for (int uBin = 0; uBin < m_nBinsU; ++uBin) {
143 for (int vBin = 0; vBin < m_nBinsV; ++vBin) {
144 std::pair<VxdID, int> key(sensorID, getRegionID(uBin, vBin));
145 if (m_regionSensitivePixelMap[key] == 0) {
146 B2WARNING("All pixels from subregion uBin=" << uBin << " vBin=" << vBin << " on Sensor=" << sensorID << " are masked.");
147 }
148 }
149 }
150
151 }
152}
153
155{
156 //Register collections
157 const StoreArray<PXDCluster> storeClusters(m_storeClustersName);
158 const StoreArray<PXDDigit> storeDigits(m_storeDigitsName);
159
160 //Get the event meta data
161 StoreObjPtr<EventMetaData> eventMetaDataPtr;
162
163 // Compute the curent one second timestamp
164 unsigned long long int ts = eventMetaDataPtr->getTime() / m_timePeriod;
165
166 // If needed, add a new one second block to buffer
167 auto iter = m_buffer.find(ts);
168 if (iter == m_buffer.end()) {
170 }
171
172 // Empty map for computing event wise occupancy
173 std::map<VxdID, double> occupancyMap;
174 auto gTools = VXD::GeoCache::getInstance().getGeoTools();
175 for (int i = 0; i < m_nPXDSensors; i++) {
176 VxdID sensorID = gTools->getSensorIDFromPXDIndex(i);
177 occupancyMap[sensorID] = 0.0;
178 }
179
180 // Check if there is PXD data
181 if (storeDigits.getEntries() > 0) {
182 m_hasPXDData = true;
183 }
184
185 for (const PXDDigit& storeDigit : storeDigits) {
186 VxdID sensorID = storeDigit.getSensorID();
187 double ADUToEnergy = PXDGainCalibrator::getInstance().getADUToEnergy(sensorID, storeDigit.getUCellID(), storeDigit.getVCellID());
188 double hitEnergy = storeDigit.getCharge() * ADUToEnergy;
189
190 if (m_sensitivePixelMap[sensorID] != 0) {
191 occupancyMap[sensorID] += 1.0 / m_sensitivePixelMap[sensorID];
192 }
193 m_buffer[ts][sensorID].m_dose += (hitEnergy / Unit::J);
194 m_buffer[ts][sensorID].m_expo += hitEnergy;
195
196 int uBin = PXDGainCalibrator::getInstance().getBinU(sensorID, storeDigit.getUCellID(), storeDigit.getVCellID(), m_nBinsU);
197 int vBin = PXDGainCalibrator::getInstance().getBinV(sensorID, storeDigit.getVCellID(), m_nBinsV);
198 int regionID = getRegionID(uBin, vBin);
199 m_buffer[ts][sensorID].m_regionDoseMap[regionID] += (hitEnergy / Unit::J);
200 m_buffer[ts][sensorID].m_regionExpoMap[regionID] += hitEnergy;
201 }
202
203 for (auto& pair : m_buffer[ts]) {
204 auto& sensorID = pair.first;
205 auto& bgdata = pair.second;
206 bgdata.m_run = eventMetaDataPtr->getRun();
207 // Check if there is actually data for this sensor
208 if (occupancyMap.find(sensorID) != occupancyMap.end()) {
209 bgdata.m_nEvents += 1;
210 bgdata.m_meanOccupancy += occupancyMap[sensorID];
211 if (occupancyMap[sensorID] > bgdata.m_maxOccupancy) {
212 bgdata.m_maxOccupancy = occupancyMap[sensorID];
213 }
214 if (occupancyMap[sensorID] < bgdata.m_minOccupancy) {
215 bgdata.m_minOccupancy = occupancyMap[sensorID];
216 }
217 }
218 }
219
220 for (const PXDCluster& cluster : storeClusters) {
221 // Update if we have a new sensor
222 VxdID sensorID = cluster.getSensorID();
223 auto info = getInfo(sensorID);
224
225 auto cluster_uID = info.getUCellID(cluster.getU());
226 auto cluster_vID = info.getVCellID(cluster.getV());
227 int uBin = PXDGainCalibrator::getInstance().getBinU(sensorID, cluster_uID, cluster_vID, m_nBinsU);
228 int vBin = PXDGainCalibrator::getInstance().getBinV(sensorID, cluster_vID, m_nBinsV);
229 int regionID = getRegionID(uBin, vBin);
230 double ADUToEnergy = PXDGainCalibrator::getInstance().getADUToEnergy(sensorID, cluster_uID, cluster_vID);
231 double clusterEnergy = cluster.getCharge() * ADUToEnergy;
232
233 if (cluster.getSize() == 1 && clusterEnergy < 10000 * Unit::eV && clusterEnergy > 6000 * Unit::eV) {
234 m_buffer[ts][sensorID].m_softPhotonFlux += 1.0;
235 m_buffer[ts][sensorID].m_regionSoftPhotonFluxMap[regionID] += 1.0;
236 } else if (cluster.getSize() == 1 && clusterEnergy > 10000 * Unit::eV) {
237 m_buffer[ts][sensorID].m_hardPhotonFlux += 1.0;
238 m_buffer[ts][sensorID].m_regionHardPhotonFluxMap[regionID] += 1.0;
239 } else if (cluster.getSize() > 1 && clusterEnergy > 10000 * Unit::eV) {
240 m_buffer[ts][sensorID].m_chargedParticleFlux += 1.0;
241 m_buffer[ts][sensorID].m_regionChargedParticleFluxMap[regionID] += 1.0;
242 }
243 }
244}
245
247{
248 // Create beast tuple
249 if (m_hasPXDData) {
250 TFile* rfile = new TFile(m_outputFileName.c_str(), "RECREATE");
251 TTree* treeBEAST = new TTree("tout", "BEAST data tree");
252
253 unsigned int ts = 0;
254 treeBEAST->Branch("ts", &(ts));
255
256 for (auto& pair : m_sensorData) {
257 auto& sensorID = pair.first;
258 auto& bgdata = pair.second;
259 string sensorDescr = sensorID;
260 treeBEAST->Branch(str(format("pxd_%1%_run") % sensorDescr).c_str(), &(bgdata.m_run));
261 treeBEAST->Branch(str(format("pxd_%1%_nEvents") % sensorDescr).c_str(), &(bgdata.m_nEvents));
262 treeBEAST->Branch(str(format("pxd_%1%_minOccupancy") % sensorDescr).c_str(), &(bgdata.m_minOccupancy));
263 treeBEAST->Branch(str(format("pxd_%1%_maxOccupancy") % sensorDescr).c_str(), &(bgdata.m_maxOccupancy));
264 treeBEAST->Branch(str(format("pxd_%1%_meanOccupancy") % sensorDescr).c_str(), &(bgdata.m_meanOccupancy));
265 treeBEAST->Branch(str(format("pxd_%1%_exposition") % sensorDescr).c_str(), &(bgdata.m_expo));
266 treeBEAST->Branch(str(format("pxd_%1%_dose") % sensorDescr).c_str(), &(bgdata.m_dose));
267 treeBEAST->Branch(str(format("pxd_%1%_softPhotonFlux") % sensorDescr).c_str(), &(bgdata.m_softPhotonFlux));
268 treeBEAST->Branch(str(format("pxd_%1%_hardPhotonFlux") % sensorDescr).c_str(), &(bgdata.m_hardPhotonFlux));
269 treeBEAST->Branch(str(format("pxd_%1%_chargedParticleFlux") % sensorDescr).c_str(),
270 &(bgdata.m_chargedParticleFlux));
271
272 for (int uBin = 0; uBin < m_nBinsU; ++uBin) {
273 for (int vBin = 0; vBin < m_nBinsV; ++vBin) {
274 int regionID = getRegionID(uBin, vBin);
275 treeBEAST->Branch(str(format("pxd_%1%_region_%2%_%3%_exposition") % sensorDescr % uBin % vBin).c_str(),
276 &(bgdata.m_regionExpoMap[regionID]));
277 treeBEAST->Branch(str(format("pxd_%1%_region_%2%_%3%_dose") % sensorDescr % uBin % vBin).c_str(),
278 &(bgdata.m_regionDoseMap[regionID]));
279 treeBEAST->Branch(str(format("pxd_%1%_region_%2%_%3%_softPhotonFlux") % sensorDescr % uBin % vBin).c_str(),
280 &(bgdata.m_regionSoftPhotonFluxMap[regionID]));
281 treeBEAST->Branch(str(format("pxd_%1%_region_%2%_%3%_hardPhotonFlux") % sensorDescr % uBin % vBin).c_str(),
282 &(bgdata.m_regionHardPhotonFluxMap[regionID]));
283 treeBEAST->Branch(str(format("pxd_%1%_region_%2%_%3%_chargedParticleFlux") % sensorDescr % uBin % vBin).c_str(),
284 &(bgdata.m_regionChargedParticleFluxMap[regionID]));
285 }
286 }
287 }
288
289 // Write timestamp and background rates into TTree
290 for (auto const& pair1 : m_buffer) {
291 auto const& timestamp = pair1.first;
292 auto const& sensors = pair1.second;
293
294 // Set variables for dumping into tree
295 ts = timestamp;
296 for (auto const& pair2 : sensors) {
297 auto const& sensorID = pair2.first;
298 auto const& bgdata = pair2.second;
299 double currentComponentTime = bgdata.m_nEvents * m_integrationTime;
300 const PXD::SensorInfo& info = getInfo(sensorID);
301 double currentSensorMass = m_sensitiveAreaMap[sensorID] * info.getThickness() * c_densitySi;
302 double currentSensorArea = m_sensitiveAreaMap[sensorID];
303 m_sensorData[sensorID] = bgdata;
304 // Some bg rates are still in wrong units. We have to fix this now.
305 m_sensorData[sensorID].m_meanOccupancy = bgdata.m_meanOccupancy / bgdata.m_nEvents;
306
307 if (currentSensorArea > 0) {
308 m_sensorData[sensorID].m_dose *= (1.0 / (currentComponentTime / Unit::s)) * (1000 / currentSensorMass);
309 m_sensorData[sensorID].m_expo *= (1.0 / currentSensorArea) * (1.0 / (currentComponentTime / Unit::s));
310 m_sensorData[sensorID].m_softPhotonFlux *= (1.0 / currentSensorArea) * (1.0 / (currentComponentTime / Unit::s));
311 m_sensorData[sensorID].m_hardPhotonFlux *= (1.0 / currentSensorArea) * (1.0 / (currentComponentTime / Unit::s));
312 m_sensorData[sensorID].m_chargedParticleFlux *= (1.0 / currentSensorArea) * (1.0 / (currentComponentTime / Unit::s));
313
314 for (int regionID = 0; regionID < m_nBinsU * m_nBinsV; ++regionID) {
315 std::pair<VxdID, int> key(sensorID, regionID);
316 double currentRegionMass = m_regionSensitiveAreaMap[key] * info.getThickness() * c_densitySi;
317 double currentRegionArea = m_regionSensitiveAreaMap[key];
318 if (currentRegionArea > 0) {
319 m_sensorData[sensorID].m_regionDoseMap[regionID] *= (1.0 / currentComponentTime) * (1000 / currentRegionMass);
320 m_sensorData[sensorID].m_regionExpoMap[regionID] *= (1.0 / currentRegionArea) * (1.0 / (currentComponentTime / Unit::s));
321 m_sensorData[sensorID].m_regionSoftPhotonFluxMap[regionID] *= (1.0 / currentRegionArea) * (1.0 / (currentComponentTime / Unit::s));
322 m_sensorData[sensorID].m_regionHardPhotonFluxMap[regionID] *= (1.0 / currentRegionArea) * (1.0 / (currentComponentTime / Unit::s));
323 m_sensorData[sensorID].m_regionChargedParticleFluxMap[regionID] *= (1.0 / currentRegionArea) * (1.0 /
324 (currentComponentTime / Unit::s));
325 }
326 }
327 }
328 }
329 // Dump variables into tree
330 treeBEAST->Fill();
331 }
332
333 // Write output tuple
334 rfile->cd();
335 treeBEAST->Write();
336 rfile->Close();
337 }
338}
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
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
std::map< std::pair< VxdID, int >, int > m_regionSensitivePixelMap
Struct to hold region-wise number of sensitive pixels.
double m_timePeriod
Period for background time series.
void initialize() override final
Initialize module.
bool m_maskDeadPixels
Correct bg rates by taking into account masked pixels.
std::map< VxdID, SensorData > m_sensorData
Struct to hold sensor-wise background data.
std::map< VxdID, double > m_sensitiveAreaMap
Struct to hold sensor-wise sensitive area.
double m_integrationTime
Integration time of PXD.
int m_nBinsV
Number of regions per sensor along v side.
const double c_densitySi
Density of crystalline Silicon.
void terminate() override final
Final summary and cleanup.
int getRegionID(int uBin, int vBin) const
Get region id from region uBin and vBin.
const PXD::SensorInfo & getInfo(VxdID sensorID) const
This is a shortcut to getting PXD::SensorInfo from the GeoCache.
void event() override final
Event processing.
int m_nBinsU
Number of regions per sensor along u side.
std::map< std::pair< VxdID, int >, double > m_regionSensitiveAreaMap
Struct to hold region-wise sensitive area.
std::map< VxdID, int > m_sensitivePixelMap
Struct to hold sensor-wise number of sensitive pixels.
double getRegionArea(VxdID sensorID, int vBin) const
Return area of the region with the given sensor ID and region vBin.
std::map< unsigned long long int, std::map< VxdID, SensorData > > m_buffer
Struct to hold sensor-wise background data.
double getSensorArea(VxdID sensorID) const
Return area of the sensor with the given sensor ID.
std::string m_storeDigitsName
PXDDigits StoreArray name.
std::string m_storeClustersName
PXDClusters StoreArray name.
bool m_hasPXDData
Flag to indicate there was at least one PXDDigit in the run.
void beginRun() override final
Start-of-run initializations.
std::string m_outputFileName
output tuple file name
int m_nPXDSensors
Total number of PXD sensors.
unsigned short getBinV(VxdID id, unsigned int vid) const
Get gain correction bin along sensor v side.
unsigned short getBinU(VxdID id, unsigned int uid, unsigned int vid) const
Get gain correction bin along sensor u 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.
Specific implementation of SensorInfo for PXD Sensors which provides additional pixel specific inform...
Definition: SensorInfo.h:23
const std::string & getName() const
Return name under which the object is saved in the DataStore.
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:216
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
static const double us
[microsecond]
Definition: Unit.h:97
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
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Namespace to encapsulate code needed for simulation and reconstrucion of the PXD.
Abstract base class for different kinds of events.
STL namespace.