Belle II Software development
CDCHitRateCounter.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 <cdc/dataobjects/CDCHit.h>
14#include <TTree.h>
15#include <map>
16
17
18namespace Belle2 {
23 namespace Background {
24
29
30 private:
31 static const int f_nSuperLayer = 9;
32 static const int f_nLayer = 56;
33 static const int f_nPhiDivision = 8;
35 public:
36
40 struct TreeStruct {
41 float averageRate = 0;
43 float layerHitRate[f_nLayer] = {0};
53 int numEvents = 0;
54 bool valid = false;
60 {
61 for (int i = 0 ; i < f_nLayer ; ++i) {
62 for (int j = 0 ; j < f_nPhiDivision ; ++j) {
63 layerPhiHitRate[i][j] = 0;
64 nActiveWireInLayerPhi[i][j] = 0;
65 }
66 }
67 }
68
72 void normalize()
73 {
74 if (numEvents == 0) return;
75 averageRate /= (numEvents * timeWindowForNormalCell * 0.982536 * 1e-6);
76
77 for (int iSL = 0 ; iSL < f_nSuperLayer ; ++iSL) {
78 if (iSL == 0)
79 superLayerHitRate[iSL] /= (numEvents * timeWindowForSmallCell * 0.982536 * 1e-6);
80 else
81 superLayerHitRate[iSL] /= (numEvents * timeWindowForNormalCell * 0.982536 * 1e-6);
82 }
83
84 for (int iL = 0 ; iL < f_nLayer ; ++iL) {
85 if (iL <= 7) {
86 layerHitRate[iL] /= (numEvents * timeWindowForSmallCell * 0.982536 * 1e-6);
87 for (int iPhi = 0 ; iPhi < f_nPhiDivision ; ++iPhi)
88 layerPhiHitRate[iL][iPhi] /= (numEvents * timeWindowForSmallCell * 0.982536 * 1e-6);
89 } else {
90 layerHitRate[iL] /= (numEvents * timeWindowForNormalCell * 0.982536 * 1e-6);
91 for (int iPhi = 0 ; iPhi < f_nPhiDivision ; ++iPhi)
92 layerPhiHitRate[iL][iPhi] /= (numEvents * timeWindowForNormalCell * 0.982536 * 1e-6);
93 }
94 }
95 }
96 };
97
108 CDCHitRateCounter(const int timeWindowLowerEdge_smallCell,
109 const int timeWindowUpperEdge_smallCell,
110 const int timeWindowLowerEdge_normalCell,
111 const int timeWindowUpperEdge_normalCell,
112 const bool enableBadWireTreatment,
113 const bool enableBackgroundHitFilter,
114 const bool enableMarkBackgroundHit):
115 m_timeWindowLowerEdge_smallCell(timeWindowLowerEdge_smallCell),
116 m_timeWindowUpperEdge_smallCell(timeWindowUpperEdge_smallCell),
117 m_timeWindowLowerEdge_normalCell(timeWindowLowerEdge_normalCell),
118 m_timeWindowUpperEdge_normalCell(timeWindowUpperEdge_normalCell),
119 m_enableBadWireTreatment(enableBadWireTreatment),
120 m_enableBackgroundHitFilter(enableBackgroundHitFilter),
121 m_enableMarkBackgroundHit(enableMarkBackgroundHit)
122 {
125 B2FATAL("invalid seting of CDC time window");
126 }
127 }
132 virtual void initialize(TTree* tree) override;
133
137 virtual void clear() override;
138
143 virtual void accumulate(unsigned timeStamp) override;
144
149 virtual void normalize(unsigned timeStamp) override;
150
151 private:
152
153 // class parameters: to be set via constructor or setters
154
155 // tree structure
158 // buffer
159 std::map<unsigned, TreeStruct> m_buffer;
161 // collections
175 bool isInTimeWindow(const int SL, const short tdc)
176 {
177 if (SL == 0) {
179 } else {
181 }
182 }
183
184 // DB payloads
185
186 // other
189 const bool
193 0;
204
210 void countActiveWires();
211
212
216 unsigned short getIPhiBin(unsigned short iSL, unsigned short iWireInLayer);
217
218 };
219
220 } // Background namespace
222} // Belle2 namespace
Class for monitoring beam background hit rates of CDC.
static const int f_nLayer
the number of layers
std::map< unsigned, TreeStruct > m_buffer
average rates in time stamps
int m_nActiveWireInLayer[f_nLayer]
the number of wires used in this hit-rate calculation in each layer
void countActiveWires_countAll()
set m_nActiveWireInTotal, m_nActiveWireInLayer[] and m_nActiveWireInSuperLayer[].
const bool m_enableBackgroundHitFilter
flag to enable the CDC background hit (crosstakl, noise) filter.
bool isInTimeWindow(const int SL, const short tdc)
return true if the hit is in the given time window
virtual void initialize(TTree *tree) override
Class initializer: set branch addresses and other staf.
const int m_timeWindowLowerEdge_normalCell
lower edge of the time window for normal cells [tdc count = 0.982536 ns]
int m_nActiveWireInTotal
the number of wires used in this hit-rate calculation in the whole CDC
StoreArray< CDCHit > m_digits
collection of digits
virtual void accumulate(unsigned timeStamp) override
Accumulate hits.
unsigned short getIPhiBin(unsigned short iSL, unsigned short iWireInLayer)
get the bin ID of the division.
int m_nActiveWireInLayerPhi[f_nLayer][f_nPhiDivision]
the number of wires used in this hit-rate calculation in each phi bin in each layer
CDCHitRateCounter(const int timeWindowLowerEdge_smallCell, const int timeWindowUpperEdge_smallCell, const int timeWindowLowerEdge_normalCell, const int timeWindowUpperEdge_normalCell, const bool enableBadWireTreatment, const bool enableBackgroundHitFilter, const bool enableMarkBackgroundHit)
Constructor.
const int m_timeWindowUpperEdge_smallCell
upper edge of the time window for small cells [tdc count = 0.982536 ns]
const bool m_enableBadWireTreatment
flag to enable the bad wire treatment.
const int m_timeWindowUpperEdge_normalCell
upper edge of the time window for normal cells [tdc count = 0.982536 ns]
static const int f_nSuperLayer
the number of super layers
const bool m_enableMarkBackgroundHit
flag to enable to mark background flag on CDCHit (set 0x100 bit for CDCHit::m_status).
static const int f_nPhiDivision
the number of division in phi
void countActiveWires()
set m_nActiveWireInTotal, m_nActiveWireInLayer[] and m_nActiveWireInSuperLayer[].
const int m_timeWindowLowerEdge_smallCell
lower edge of the time window for small cells [tdc count = 0.982536 ns]
int m_nActiveWireInSuperLayer[f_nSuperLayer]
the number of wires used in this hit-rate calculation in each suler layer
virtual void normalize(unsigned timeStamp) override
Normalize accumulated hits (e.g.
virtual void clear() override
Clear time-stamp buffer to prepare for 'accumulate'.
Abstract base class for monitoring beam background hit rates All the monitoring classes must inherit ...
Definition: HitRateBase.h:24
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
Abstract base class for different kinds of events.
float superLayerHitRate[f_nSuperLayer]
SuperLayer average hit rate in kHz.
float layerPhiHitRate[f_nLayer][f_nPhiDivision]
Layer (in the phi bin) average hit rate in kHz.
int nActiveWireInLayerPhi[f_nLayer][f_nPhiDivision]
number of wires used in this analysis in each phi bin in each layer
float layerHitRate[f_nLayer]
Layer average hit rate in kHz.
int nActiveWireInTotal
number of wires used in this analysis in the whole CDC
int nActiveWireInSuperLayer[f_nSuperLayer]
number of wires used in this analysis in each super layer
int nActiveWireInLayer[f_nLayer]
number of wires used in this analysis in each layer
float averageRate
total detector average hit rate in KHz
void normalize()
normalize accumulated hits to hit rate in kHz
int timeWindowForSmallCell
time window for the small cells in 2*508.887 MHz clock ( 1 clock = 0.982536 ns)
int timeWindowForNormalCell
time window for the normal cells in 2*508.887 MHz clock ( 1 clock = 0.982536 ns)