Belle II Software release-09-00-00
PXDMCBgTupleProducerModule.h
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#pragma once
10
11#include <framework/core/Module.h>
12#include <framework/gearbox/Unit.h>
13#include <vxd/dataobjects/VxdID.h>
14#include <vxd/geometry/SensorInfoBase.h>
15#include <pxd/geometry/SensorInfo.h>
16#include <vxd/geometry/GeoCache.h>
17#include <string>
18#include <vector>
19#include <map>
20
21
22namespace Belle2 {
28 namespace PXD {
29
41
42 public:
43
45 struct SensorData {
49 double m_expo;
51 double m_dose;
58
60 std::vector<double> m_regionExpoMap;
62 std::vector<double> m_regionDoseMap;
64 std::vector<double> m_regionSoftPhotonFluxMap;
66 std::vector<double> m_regionChargedParticleFluxMap;
68 std::vector<double> m_regionHardPhotonFluxMap;
69 };
70
73
75 void initialize() override final;
77 void beginRun() override final;
79 void event() override final;
81 void terminate() override final;
82
83 private:
84
85 // General
86 const double c_densitySi = 2.3290 * Unit::g_cm3;
92 inline const PXD::SensorInfo& getInfo(VxdID sensorID) const;
94 inline double getSensorArea(VxdID sensorID) const;
96 inline int getRegionID(int uBin, int vBin) const;
98 inline double getRegionArea(VxdID sensorID, int vBin) const;
99
100 // Output directory
103 // StoreArrays
115 std::map<std::pair<VxdID, int>, int> m_regionSensitivePixelMap;
116 std::map<std::pair<VxdID, int>, double> m_regionSensitiveAreaMap;
120 };
121
122 inline int PXDMCBgTupleProducerModule::getRegionID(int uBin, int vBin) const
123 {
124 return uBin * m_nBinsV + vBin;
125 }
126
128 {
129 return dynamic_cast<const PXD::SensorInfo&>(VXD::GeoCache::get(sensorID));
130 }
131
133 {
134 const PXD::SensorInfo& info = getInfo(sensorID);
135 return info.getWidth() * info.getLength();
136 }
137
138 inline double PXDMCBgTupleProducerModule::getRegionArea(VxdID sensorID, int vBin) const
139 {
140 const PXD::SensorInfo& info = getInfo(sensorID);
141 int vi = vBin * info.getVCells() / m_nBinsV;
142 double length = info.getVPitch(info.getVCellPosition(vi)) * info.getVCells() / m_nBinsV;
143 return length * info.getWidth() / m_nBinsU;
144 }
145
146 } // namespace PXD
148} // namespace Belle2
149
Base class for Modules.
Definition: Module.h:72
std::map< std::pair< VxdID, int >, int > m_regionSensitivePixelMap
Struct to hold region-wise number of sensitive pixels.
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::string m_storeBgMetaDataName
Name of the persistent BackgroundMetaDta object.
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.
double m_overrideComponentTime
Time of current component given by user.
bool m_hasPXDData
Flag to indicate there was at least one PXDDigit in the run.
void beginRun() override final
Start-of-run initializations.
double m_componentTime
Time of current component.
std::string m_outputFileName
output tuple file name
Specific implementation of SensorInfo for PXD Sensors which provides additional pixel specific inform...
Definition: SensorInfo.h:23
The Unit class.
Definition: Unit.h:40
static const SensorInfoBase & get(Belle2::VxdID id)
Return a reference to the SensorInfo of a given SensorID.
Definition: GeoCache.h:139
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.
double m_expo
Exposition (energy deposited per cm2 and second)
double m_meanOccupancy
Average occupancy of all events inside one second block.
std::vector< double > m_regionChargedParticleFluxMap
Charged particle flux (selected clusters per cm and second) for sensor regions.
double m_hardPhotonFlux
Hard photon flux (selected clusters per cm and second)
double m_softPhotonFlux
Soft photon flux (selected clusters per cm and second)
std::vector< double > m_regionDoseMap
Dose (Gy per second) for sensor regions.
std::vector< double > m_regionHardPhotonFluxMap
Hard photon flux (selected clusters per cm and second) for sensor regions.
std::vector< double > m_regionExpoMap
Expositions (energy deposited per cm2 and second) for sensor regions.
std::vector< double > m_regionSoftPhotonFluxMap
Soft photon flux (selected clusters per cm and second) for sensor regions.
double m_chargedParticleFlux
Charged particle flux (selected clusters per cm and second)