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>
18 #include <pxd/dataobjects/PXDDigit.h>
19 #include <pxd/dataobjects/PXDCluster.h>
20 #include <boost/format.hpp>
43 , m_nPXDSensors(0), m_hasPXDData(false)
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));
56 void PXDBgTupleProducerModule::initialize()
66 m_storeDigitsName = storeDigits.
getName();
69 m_integrationTime *= Unit::us;
72 m_timePeriod *= Unit::s;
78 auto gTools = VXD::GeoCache::getInstance().getGeoTools();
79 if (gTools->getNumberOfPXDLayers() == 0) {
80 B2WARNING(
"Missing geometry for PXD, PXD-masking is skiped.");
82 m_nPXDSensors = gTools->getNumberOfPXDSensors();
85 for (
int i = 0; i < m_nPXDSensors; i++) {
86 VxdID sensorID = gTools->getSensorIDFromPXDIndex(i);
89 m_sensorData[sensorID].m_minOccupancy = 1.0;
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);
99 void PXDBgTupleProducerModule::beginRun()
102 for (
auto const& pair2 : m_sensorData) {
103 auto const& sensorID = pair2.first;
104 auto info = getInfo(sensorID);
107 m_sensitivePixelMap[sensorID] = info.getUCells() * info.getVCells();
109 m_sensitiveAreaMap[sensorID] = getSensorArea(sensorID);
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));
115 m_regionSensitivePixelMap[key] = info.getUCells() * info.getVCells() / m_nBinsU / m_nBinsV;
117 m_regionSensitiveAreaMap[key] = getRegionArea(sensorID, vBin);
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));
131 m_regionSensitivePixelMap[key] -= 1;
132 m_regionSensitiveAreaMap[key] -= info.getVPitch(info.getVCellPosition(vi)) * info.getUPitch();
138 if (m_sensitivePixelMap[sensorID] == 0) {
139 B2WARNING(
"All pixels from Sensor=" << sensorID <<
" are masked.");
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.");
154 void PXDBgTupleProducerModule::event()
164 unsigned long long int ts = eventMetaDataPtr->getTime() / m_timePeriod;
167 auto iter = m_buffer.find(ts);
168 if (iter == m_buffer.end()) {
169 m_buffer[ts] = m_sensorData;
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;
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;
190 if (m_sensitivePixelMap[sensorID] != 0) {
191 occupancyMap[sensorID] += 1.0 / m_sensitivePixelMap[sensorID];
193 m_buffer[ts][sensorID].m_dose += (hitEnergy / Unit::J);
194 m_buffer[ts][sensorID].m_expo += hitEnergy;
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;
203 for (
auto& pair : m_buffer[ts]) {
204 auto& sensorID = pair.first;
205 auto& bgdata = pair.second;
206 bgdata.m_run = eventMetaDataPtr->getRun();
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];
214 if (occupancyMap[sensorID] < bgdata.m_minOccupancy) {
215 bgdata.m_minOccupancy = occupancyMap[sensorID];
220 for (
const PXDCluster& cluster : storeClusters) {
222 VxdID sensorID = cluster.getSensorID();
223 auto info = getInfo(sensorID);
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;
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;
246 void PXDBgTupleProducerModule::terminate()
250 TFile* rfile =
new TFile(m_outputFileName.c_str(),
"RECREATE");
251 TTree* treeBEAST =
new TTree(
"tout",
"BEAST data tree");
254 treeBEAST->Branch(
"ts", &(ts));
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));
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]));
290 for (
auto const& pair1 : m_buffer) {
291 auto const& timestamp = pair1.first;
292 auto const& sensors = pair1.second;
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;
301 double currentSensorMass = m_sensitiveAreaMap[sensorID] * info.getThickness() * c_densitySi;
302 double currentSensorArea = m_sensitiveAreaMap[sensorID];
303 m_sensorData[sensorID] = bgdata;
305 m_sensorData[sensorID].m_meanOccupancy = bgdata.m_meanOccupancy / bgdata.m_nEvents;
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));
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));
The PXD Cluster class This class stores all information about reconstructed PXD clusters The position...
PXD Background Tuple Producer.
Specific implementation of SensorInfo for PXD Sensors which provides additional pixel specific inform...
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
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.
Type-safe access to single objects in the data store.
Class to uniquely identify a any structure of the PXD and SVD.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Namespace to encapsulate code needed for simulation and reconstrucion of the PXD.
Abstract base class for different kinds of events.
Struct to hold data of an PXD sensor.