Belle II Software development
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 <pxd/geometry/SensorInfo.h>
15#include <vxd/geometry/GeoCache.h>
16#include <string>
17#include <vector>
18#include <map>
19
20
21namespace Belle2 {
26
27 namespace PXD {
28
37
38 public:
39
41 struct SensorData {
45 double m_expo;
47 double m_dose;
54
56 std::vector<double> m_regionExpoMap;
58 std::vector<double> m_regionDoseMap;
60 std::vector<double> m_regionSoftPhotonFluxMap;
62 std::vector<double> m_regionChargedParticleFluxMap;
64 std::vector<double> m_regionHardPhotonFluxMap;
65 };
66
69
71 void initialize() override final;
73 void beginRun() override final;
75 void event() override final;
77 void terminate() override final;
78
79 private:
80
81 // General
82 const double c_densitySi = 2.3290 * Unit::g_cm3;
83
88 inline const PXD::SensorInfo& getInfo(VxdID sensorID) const;
90 inline double getSensorArea(VxdID sensorID) const;
92 inline int getRegionID(int uBin, int vBin) const;
94 inline double getRegionArea(VxdID sensorID, int vBin) const;
95
96 // Output directory
98
99 // StoreArrays
107
116 };
117
118 inline int PXDMCBgTupleProducerModule::getRegionID(int uBin, int vBin) const
119 {
120 return uBin * m_nBinsV + vBin;
121 }
122
124 {
125 return dynamic_cast<const PXD::SensorInfo&>(VXD::GeoCache::getInstance().getSensorInfo(sensorID));
126 }
127
129 {
130 const PXD::SensorInfo& info = getInfo(sensorID);
131 return info.getWidth() * info.getLength();
132 }
133
134 inline double PXDMCBgTupleProducerModule::getRegionArea(VxdID sensorID, int vBin) const
135 {
136 const PXD::SensorInfo& info = getInfo(sensorID);
137 int vi = vBin * info.getVCells() / m_nBinsV;
138 double length = info.getVPitch(info.getVCellPosition(vi)) * info.getVCells() / m_nBinsV;
139 return length * info.getWidth() / m_nBinsU;
140 }
141
142 } // namespace PXD
144} // namespace Belle2
145
Module()
Constructor.
Definition Module.cc:30
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.
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
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a reference to the SensorInfo of a given SensorID.
Definition GeoCache.cc:67
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition GeoCache.cc:214
Class to uniquely identify a any structure of the PXD and SVD.
Definition VxdID.h:32
STL class.
STL class.
Namespace to encapsulate code needed for simulation and reconstrucion of the PXD.
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)