Belle II Software development
PXDHitRateCounter.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 <background/modules/BeamBkgHitRateMonitor/HitRateBase.h>
12#include <framework/datastore/StoreArray.h>
13#include <pxd/dataobjects/PXDDigit.h>
14#include <pxd/dataobjects/PXDCluster.h>
15#include <vxd/dataobjects/VxdID.h>
16#include <pxd/geometry/SensorInfo.h>
17#include <vxd/geometry/GeoCache.h>
18#include <TTree.h>
19#include <map>
20
21#include <framework/gearbox/Unit.h>
22
23namespace Belle2 {
28 namespace Background {
29
34
35 public:
36
40 struct TreeStruct {
41 float meanOccupancies[40] = {0};
42 float maxOccupancies[40] = {0};
43 float doseRates[40] = {0};
44 float softPhotonFluxes[40] = {0};
45 float hardPhotonFluxes[40] = {0};
46 float chargedFluxes[40] = {0};
47 float segmentDoseRates[240] = {0};
48 float segmentSoftPhotonFluxes[240] = {0};
49 float segmentHardPhotonFluxes[240] = {0};
50 float segmentChargedFluxes[240] = {0};
51 float averageRate = 0;
52 int numEvents = 0;
53 bool valid = false;
54 };
55
60 {}
61
66 virtual void initialize(TTree* tree) override;
67
71 virtual void clear() override;
72
77 virtual void accumulate(unsigned timeStamp) override;
78
83 virtual void normalize(unsigned timeStamp) override;
84
85 private:
86
90 inline const PXD::SensorInfo& getInfo(VxdID sensorID) const;
91
95 void setActivePixels();
96
97 // class parameters: to be set via constructor or setters
99 bool m_maskDeadPixels = true;
102 const double c_densitySi = 2.3290 * Unit::g_cm3;
104 // tree structure
107 // buffer
108 std::map<unsigned, TreeStruct> m_buffer;
110 // collections
114 // other
115 double m_activePixels[40] = {0};
116 double m_activeAreas[40] = {0};
117 double m_segmentActivePixels[240] = {0};
118 double m_segmentActiveAreas[240] = {0};
119 };
120
122 {
123 return dynamic_cast<const PXD::SensorInfo&>(VXD::GeoCache::getInstance().getSensorInfo(sensorID));
124 }
125
126 } // Background namespace
128} // Belle2 namespace
Abstract base class for monitoring beam background hit rates All the monitoring classes must inherit ...
Definition: HitRateBase.h:24
Class for monitoring beam background hit rates of PXD.
double m_segmentActivePixels[240]
number of active pixels in v segements
std::map< unsigned, TreeStruct > m_buffer
average rates in time stamps
bool m_maskDeadPixels
Correct bg rates by taking into account masked pixels.
double m_integrationTime
Integration time of PXD in ns.
virtual void initialize(TTree *tree) override
Class initializer: set branch addresses and other staf.
const double c_densitySi
Density of crystalline Silicon.
StoreArray< PXDCluster > m_clusters
collection of clusters
double m_activePixels[40]
number of active pixels in sensor
virtual void accumulate(unsigned timeStamp) override
Accumulate hits.
const PXD::SensorInfo & getInfo(VxdID sensorID) const
Get PXD::SensorInfo.
StoreArray< PXDDigit > m_digits
collection of digits
void setActivePixels()
Sets fractions of active channels.
double m_activeAreas[40]
area of active pixels in sensor
double m_segmentActiveAreas[240]
area of active pixels in v segments
virtual void normalize(unsigned timeStamp) override
Normalize accumulated hits (e.g.
virtual void clear() override
Clear time-stamp buffer to prepare for 'accumulate'.
Specific implementation of SensorInfo for PXD Sensors which provides additional pixel specific inform...
Definition: SensorInfo.h:23
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
static const double us
[microsecond]
Definition: Unit.h:97
static const double g_cm3
Practical units with the value set at 1.
Definition: Unit.h:60
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne 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:33
Abstract base class for different kinds of events.
float doseRates[40]
mean dose rate from PXDDigits per sensor [Gy/s]
float maxOccupancies[40]
max hit occupancy from PXDDigits per sensor [Hits/Channel]
float meanOccupancies[40]
mean hit occupancy from PXDDigits per sensor [Hits/Channel]
float softPhotonFluxes[40]
mean soft photon flux per sensor (Single pixel cluster <10keV) [clusters/cm2/s]
float chargedFluxes[40]
mean charged particle flux per sensor (Multi pixel cluster >10keV) [clusters/cm2/s]
float segmentChargedFluxes[240]
mean charged particle flux per v segment of sensor (Multi pixel cluster >10keV) [clusters/cm2/s]
float segmentHardPhotonFluxes[240]
mean hard photon flux per v segment of sensor (Single pixel cluster >10keV) [clusters/cm2/s]
float hardPhotonFluxes[40]
mean hard photon flux per sensor (Single pixel cluster >10keV) [clusters/cm2/s]
float averageRate
total detector average hit rate
float segmentSoftPhotonFluxes[240]
mean soft photon flux per v segment of sensor (Single pixel cluster <10keV) [clusters/cm2/s]
float segmentDoseRates[240]
mean dose rate from PXDDigits per v segment of sensor [Gy/s]