Belle II Software  release-06-01-15
PXDMCBgTupleProducerModule.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/PXDMCBgTupleProducerModule.h>
12 #include <framework/datastore/DataStore.h>
13 #include <framework/datastore/StoreObjPtr.h>
14 #include <framework/datastore/StoreArray.h>
15 #include <framework/dataobjects/BackgroundMetaData.h>
16 #include <pxd/reconstruction/PXDGainCalibrator.h>
17 #include <pxd/reconstruction/PXDPixelMasker.h>
18 
19 #include <pxd/dataobjects/PXDDigit.h>
20 #include <pxd/dataobjects/PXDCluster.h>
21 #include <boost/format.hpp>
22 
23 #include <TFile.h>
24 #include <TTree.h>
25 
26 using namespace std;
27 using boost::format;
28 using namespace Belle2;
29 using namespace Belle2::PXD;
30 
31 
32 
33 //-----------------------------------------------------------------
34 // Register the Module
35 //-----------------------------------------------------------------
36 REG_MODULE(PXDMCBgTupleProducer)
37 
38 
39 //-----------------------------------------------------------------
40 // Implementation
41 //-----------------------------------------------------------------
42 
44  , m_hasPXDData(false), m_componentTime(0)
45 {
46  //Set module properties
47  setDescription("PXD background tuple producer module");
48  addParam("integrationTime", m_integrationTime, "PXD integration time in micro seconds", double(20));
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  addParam("overrideComponentTime", m_overrideComponentTime, "User specified component time in micro seconds", double(0.0));
54 }
55 
56 
57 
58 void PXDMCBgTupleProducerModule::initialize()
59 {
60  //Register collections
61  StoreArray<PXDCluster> storeClusters(m_storeClustersName);
62  StoreArray<PXDDigit> storeDigits(m_storeDigitsName);
63  StoreObjPtr<BackgroundMetaData> storeBgMetaData(m_storeBgMetaDataName, DataStore::c_Persistent);
64 
65  //Store names to speed up creation later
66  m_storeDigitsName = storeDigits.getName();
67  m_storeBgMetaDataName = storeBgMetaData.getName();
68 
69  // PXD integration time
70  m_integrationTime *= Unit::us;
71 
72  m_overrideComponentTime *= Unit::us;
73 
74  // So far, we did not see PXD data
75  m_hasPXDData = false;
76 
77  //Pointer to GeoTools instance
78  auto gTools = VXD::GeoCache::getInstance().getGeoTools();
79  if (gTools->getNumberOfPXDLayers() == 0) {
80  B2WARNING("Missing geometry for PXD, PXD-masking is skiped.");
81  }
82  int nPXDSensors = gTools->getNumberOfPXDSensors();
83 
84  // Initialize m_sensorData with empty sensorData for all sensors
85  for (int i = 0; i < nPXDSensors; i++) {
86  VxdID sensorID = gTools->getSensorIDFromPXDIndex(i);
87  m_sensorData[sensorID] = SensorData();
88  // Initialize counters for subdivisions per sensor
89  m_sensorData[sensorID].m_regionExpoMap = vector<double>(m_nBinsU * m_nBinsV, 0);
90  m_sensorData[sensorID].m_regionDoseMap = vector<double>(m_nBinsU * m_nBinsV, 0);
91  m_sensorData[sensorID].m_regionSoftPhotonFluxMap = vector<double>(m_nBinsU * m_nBinsV, 0);
92  m_sensorData[sensorID].m_regionChargedParticleFluxMap = vector<double>(m_nBinsU * m_nBinsV, 0);
93  m_sensorData[sensorID].m_regionHardPhotonFluxMap = vector<double>(m_nBinsU * m_nBinsV, 0);
94  }
95 }
96 
97 void PXDMCBgTupleProducerModule::beginRun()
98 {
99  // Compute the sensitive area for all PXD sensors
100  for (auto const& pair2 : m_sensorData) {
101  auto const& sensorID = pair2.first;
102  auto info = getInfo(sensorID);
103 
104  // Compute nominal number of pixel per sensor
105  m_sensitivePixelMap[sensorID] = info.getUCells() * info.getVCells();
106  // Compute nominal area per sensor
107  m_sensitiveAreaMap[sensorID] = getSensorArea(sensorID);
108 
109  for (int uBin = 0; uBin < m_nBinsU; ++uBin) {
110  for (int vBin = 0; vBin < m_nBinsV; ++vBin) {
111  std::pair<VxdID, int> key(sensorID, getRegionID(uBin, vBin));
112  // Compute nominal number of pixel per sensor subregion
113  m_regionSensitivePixelMap[key] = info.getUCells() * info.getVCells() / m_nBinsU / m_nBinsV;
114  // Compute nominal area per sensor subregion
115  m_regionSensitiveAreaMap[key] = getRegionArea(sensorID, vBin);
116  }
117  }
118 
119  if (m_maskDeadPixels) {
120  for (int ui = 0; ui < info.getUCells(); ++ui) {
121  for (int vi = 0; vi < info.getVCells(); ++vi) {
122  if (PXD::PXDPixelMasker::getInstance().pixelDead(sensorID, ui, vi)
123  || !PXD::PXDPixelMasker::getInstance().pixelOK(sensorID, ui, vi)) {
124  m_sensitivePixelMap[sensorID] -= 1;
125  m_sensitiveAreaMap[sensorID] -= info.getVPitch(info.getVCellPosition(vi)) * info.getUPitch();
126  int uBin = PXDGainCalibrator::getInstance().getBinU(sensorID, ui, vi, m_nBinsU);
127  int vBin = PXDGainCalibrator::getInstance().getBinV(sensorID, vi, m_nBinsV);
128  std::pair<VxdID, int> key(sensorID, getRegionID(uBin, vBin));
129  m_regionSensitivePixelMap[key] -= 1;
130  m_regionSensitiveAreaMap[key] -= info.getVPitch(info.getVCellPosition(vi)) * info.getUPitch();
131  }
132  }
133  }
134  }
135 
136  if (m_sensitivePixelMap[sensorID] == 0) {
137  B2WARNING("All pixels from Sensor=" << sensorID << " are masked.");
138  }
139 
140  for (int uBin = 0; uBin < m_nBinsU; ++uBin) {
141  for (int vBin = 0; vBin < m_nBinsV; ++vBin) {
142  std::pair<VxdID, int> key(sensorID, getRegionID(uBin, vBin));
143  if (m_regionSensitivePixelMap[key] == 0) {
144  B2WARNING("All pixels from subregion uBin=" << uBin << " vBin=" << vBin << " on Sensor=" << sensorID << " are masked.");
145  }
146  }
147  }
148 
149  }
150 }
151 
152 void PXDMCBgTupleProducerModule::event()
153 {
154  //Register collections
155  const StoreArray<PXDCluster> storeClusters(m_storeClustersName);
156  const StoreArray<PXDDigit> storeDigits(m_storeDigitsName);
157  const StoreObjPtr<BackgroundMetaData> storeBgMetaData(m_storeBgMetaDataName, DataStore::c_Persistent);
158 
159  // Set the real time
160  m_componentTime = storeBgMetaData->getRealTime();
161 
162  // Empty map for computing event wise occupancy
163  std::map<VxdID, double> occupancyMap;
164 
165  // Check if there is PXD data
166  if (storeDigits.getEntries() > 0) {
167  m_hasPXDData = true;
168  }
169 
170  for (const PXDDigit& storeDigit : storeDigits) {
171  VxdID sensorID = storeDigit.getSensorID();
172  double ADUToEnergy = PXDGainCalibrator::getInstance().getADUToEnergy(sensorID, storeDigit.getUCellID(), storeDigit.getVCellID());
173  double hitEnergy = storeDigit.getCharge() * ADUToEnergy;
174 
175  if (m_sensitivePixelMap[sensorID] != 0) {
176  occupancyMap[sensorID] += 1.0 / m_sensitivePixelMap[sensorID];
177  }
178  m_sensorData[sensorID].m_dose += (hitEnergy / Unit::J);
179  m_sensorData[sensorID].m_expo += hitEnergy;
180 
181  int uBin = PXDGainCalibrator::getInstance().getBinU(sensorID, storeDigit.getUCellID(), storeDigit.getVCellID(), m_nBinsU);
182  int vBin = PXDGainCalibrator::getInstance().getBinV(sensorID, storeDigit.getVCellID(), m_nBinsV);
183  int regionID = getRegionID(uBin, vBin);
184  m_sensorData[sensorID].m_regionDoseMap[regionID] += (hitEnergy / Unit::J);
185  m_sensorData[sensorID].m_regionExpoMap[regionID] += hitEnergy;
186  }
187 
188  for (auto& pair : m_sensorData) {
189  auto& sensorID = pair.first;
190  auto& bgdata = pair.second;
191 
192  // Check if there is actually data for this sensor
193  if (occupancyMap.find(sensorID) != occupancyMap.end()) {
194  bgdata.m_meanOccupancy += occupancyMap[sensorID];
195  }
196  }
197 
198  for (const PXDCluster& cluster : storeClusters) {
199  // Update if we have a new sensor
200  VxdID sensorID = cluster.getSensorID();
201  auto info = getInfo(sensorID);
202 
203  auto cluster_uID = info.getUCellID(cluster.getU());
204  auto cluster_vID = info.getVCellID(cluster.getV());
205  int uBin = PXDGainCalibrator::getInstance().getBinU(sensorID, cluster_uID, cluster_vID, m_nBinsU);
206  int vBin = PXDGainCalibrator::getInstance().getBinV(sensorID, cluster_vID, m_nBinsV);
207  int regionID = getRegionID(uBin, vBin);
208  double ADUToEnergy = PXDGainCalibrator::getInstance().getADUToEnergy(sensorID, cluster_uID, cluster_vID);
209  double clusterEnergy = cluster.getCharge() * ADUToEnergy;
210 
211  if (cluster.getSize() == 1 && clusterEnergy < 10000 * Unit::eV && clusterEnergy > 6000 * Unit::eV) {
212  m_sensorData[sensorID].m_softPhotonFlux += 1.0;
213  m_sensorData[sensorID].m_regionSoftPhotonFluxMap[regionID] += 1.0;
214  } else if (cluster.getSize() == 1 && clusterEnergy > 10000 * Unit::eV) {
215  m_sensorData[sensorID].m_hardPhotonFlux += 1.0;
216  m_sensorData[sensorID].m_regionHardPhotonFluxMap[regionID] += 1.0;
217  } else if (cluster.getSize() > 1 && clusterEnergy > 10000 * Unit::eV) {
218  m_sensorData[sensorID].m_chargedParticleFlux += 1.0;
219  m_sensorData[sensorID].m_regionChargedParticleFluxMap[regionID] += 1.0;
220  }
221  }
222 }
223 
224 void PXDMCBgTupleProducerModule::terminate()
225 {
226  // Create beast tuple
227  if (m_hasPXDData) {
228  TFile* rfile = new TFile(m_outputFileName.c_str(), "RECREATE");
229  TTree* treeBEAST = new TTree("tout", "BEAST data tree");
230 
231  double currentComponentTime = m_componentTime;
232  if (m_overrideComponentTime > 0.0) currentComponentTime = m_overrideComponentTime;
233 
234  B2RESULT("Total real time is " << currentComponentTime / Unit::us << " microseconds.");
235  B2RESULT("This is equivalent to " << currentComponentTime / m_integrationTime << " random triggered events.");
236 
237  for (auto& pair : m_sensorData) {
238  auto& sensorID = pair.first;
239  auto& bgdata = pair.second;
240  const PXD::SensorInfo& info = getInfo(sensorID);
241  double currentSensorMass = m_sensitiveAreaMap[sensorID] * info.getThickness() * c_densitySi;
242  double currentSensorArea = m_sensitiveAreaMap[sensorID];
243 
244  // Finalize computation of rates
245  m_sensorData[sensorID].m_meanOccupancy = bgdata.m_meanOccupancy * (m_integrationTime / currentComponentTime);
246 
247  if (currentSensorArea > 0) {
248  m_sensorData[sensorID].m_dose *= (1.0 / (currentComponentTime / Unit::s)) * (1000 / currentSensorMass);
249  m_sensorData[sensorID].m_expo *= (1.0 / currentSensorArea) * (1.0 / (currentComponentTime / Unit::s));
250  m_sensorData[sensorID].m_softPhotonFlux *= (1.0 / currentSensorArea) * (1.0 / (currentComponentTime / Unit::s));
251  m_sensorData[sensorID].m_hardPhotonFlux *= (1.0 / currentSensorArea) * (1.0 / (currentComponentTime / Unit::s));
252  m_sensorData[sensorID].m_chargedParticleFlux *= (1.0 / currentSensorArea) * (1.0 / (currentComponentTime / Unit::s));
253 
254  for (int regionID = 0; regionID < m_nBinsU * m_nBinsV; ++regionID) {
255  std::pair<VxdID, int> key(sensorID, regionID);
256  double currentRegionMass = m_regionSensitiveAreaMap[key] * info.getThickness() * c_densitySi;
257  double currentRegionArea = m_regionSensitiveAreaMap[key];
258  if (currentRegionArea > 0) {
259  m_sensorData[sensorID].m_regionDoseMap[regionID] *= (1.0 / currentComponentTime) * (1000 / currentRegionMass);
260  m_sensorData[sensorID].m_regionExpoMap[regionID] *= (1.0 / currentRegionArea) * (1.0 / (currentComponentTime / Unit::s));
261  m_sensorData[sensorID].m_regionSoftPhotonFluxMap[regionID] *= (1.0 / currentRegionArea) * (1.0 / (currentComponentTime / Unit::s));
262  m_sensorData[sensorID].m_regionHardPhotonFluxMap[regionID] *= (1.0 / currentRegionArea) * (1.0 / (currentComponentTime / Unit::s));
263  m_sensorData[sensorID].m_regionChargedParticleFluxMap[regionID] *= (1.0 / currentRegionArea) * (1.0 /
264  (currentComponentTime / Unit::s));
265  }
266  }
267  }
268 
269  // Prepare output tree
270  string sensorDescr = sensorID;
271  treeBEAST->Branch(str(format("pxd_%1%_meanOccupancy") % sensorDescr).c_str(), &(bgdata.m_meanOccupancy));
272  treeBEAST->Branch(str(format("pxd_%1%_exposition") % sensorDescr).c_str(), &(bgdata.m_expo));
273  treeBEAST->Branch(str(format("pxd_%1%_dose") % sensorDescr).c_str(), &(bgdata.m_dose));
274  treeBEAST->Branch(str(format("pxd_%1%_softPhotonFlux") % sensorDescr).c_str(), &(bgdata.m_softPhotonFlux));
275  treeBEAST->Branch(str(format("pxd_%1%_hardPhotonFlux") % sensorDescr).c_str(), &(bgdata.m_hardPhotonFlux));
276  treeBEAST->Branch(str(format("pxd_%1%_chargedParticleFlux") % sensorDescr).c_str(),
277  &(bgdata.m_chargedParticleFlux));
278 
279  for (int uBin = 0; uBin < m_nBinsU; ++uBin) {
280  for (int vBin = 0; vBin < m_nBinsV; ++vBin) {
281  int regionID = getRegionID(uBin, vBin);
282  treeBEAST->Branch(str(format("pxd_%1%_region_%2%_%3%_exposition") % sensorDescr % uBin % vBin).c_str(),
283  &(bgdata.m_regionExpoMap[regionID]));
284  treeBEAST->Branch(str(format("pxd_%1%_region_%2%_%3%_dose") % sensorDescr % uBin % vBin).c_str(),
285  &(bgdata.m_regionDoseMap[regionID]));
286  treeBEAST->Branch(str(format("pxd_%1%_region_%2%_%3%_softPhotonFlux") % sensorDescr % uBin % vBin).c_str(),
287  &(bgdata.m_regionSoftPhotonFluxMap[regionID]));
288  treeBEAST->Branch(str(format("pxd_%1%_region_%2%_%3%_hardPhotonFlux") % sensorDescr % uBin % vBin).c_str(),
289  &(bgdata.m_regionHardPhotonFluxMap[regionID]));
290  treeBEAST->Branch(str(format("pxd_%1%_region_%2%_%3%_chargedParticleFlux") % sensorDescr % uBin % vBin).c_str(),
291  &(bgdata.m_regionChargedParticleFluxMap[regionID]));
292  }
293  }
294  }
295 
296  // Dump variables into tree
297  treeBEAST->Fill();
298  // Write output tuple
299  rfile->cd();
300  treeBEAST->Write();
301  rfile->Close();
302  }
303 }
Base class for Modules.
Definition: Module.h:72
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
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.
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:95
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
#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.