Belle II Software  release-05-01-25
CDCHitRateCounter.h
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2019 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Marko Staric *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #pragma once
12 
13 #include <background/modules/BeamBkgHitRateMonitor/HitRateBase.h>
14 #include <framework/datastore/StoreArray.h>
15 #include <cdc/dataobjects/CDCHit.h>
16 #include <TTree.h>
17 #include <map>
18 
19 
20 namespace Belle2 {
25  namespace Background {
26 
30  class CDCHitRateCounter: public HitRateBase {
31 
32  private:
33  static const int f_nLayer = 56;
34  static const int f_nSuperLayer = 9;
36  public:
37 
41  struct TreeStruct {
42  float layerHitRate[f_nLayer] = {0};
43  float superLayerHitRate[f_nSuperLayer] = {0};
44  float averageRate = 0;
45  int timeWindowForSmallCell = -1;
46  int timeWindowForNormalCell = -1;
48  int nActiveWireInLayer[f_nLayer] = {0};
51  int numEvents = 0;
52  bool valid = false;
57  void normalize()
58  {
59  if (numEvents == 0) return;
61 
62  for (int iSL = 0 ; iSL < f_nSuperLayer ; ++iSL) {
63  if (iSL == 0)
65  else
67  }
68 
69  for (int iL = 0 ; iL < f_nLayer ; ++iL) {
70  if (iL <= 7)
72  else
74  }
75  }
76  };
77 
88  CDCHitRateCounter(const int timeWindowLowerEdge_smallCell,
89  const int timeWindowUpperEdge_smallCell,
90  const int timeWindowLowerEdge_normalCell,
91  const int timeWindowUpperEdge_normalCell,
92  const bool enableBadWireTreatment,
93  const bool enableBackgroundHitFilter,
94  const bool enableMarkBackgroundHit):
95  m_timeWindowLowerEdge_smallCell(timeWindowLowerEdge_smallCell),
96  m_timeWindowUpperEdge_smallCell(timeWindowUpperEdge_smallCell),
97  m_timeWindowLowerEdge_normalCell(timeWindowLowerEdge_normalCell),
98  m_timeWindowUpperEdge_normalCell(timeWindowUpperEdge_normalCell),
99  m_enableBadWireTreatment(enableBadWireTreatment),
100  m_enableBackgroundHitFilter(enableBackgroundHitFilter),
101  m_enableMarkBackgroundHit(enableMarkBackgroundHit)
102  {
105  B2FATAL("invalid seting of CDC time window");
106  }
107  }
112  virtual void initialize(TTree* tree) override;
113 
117  virtual void clear() override;
118 
123  virtual void accumulate(unsigned timeStamp) override;
124 
129  virtual void normalize(unsigned timeStamp) override;
130 
131  private:
132 
133  // class parameters: to be set via constructor or setters
134 
135  // tree structure
136  TreeStruct m_rates;
138  // buffer
139  std::map<unsigned, TreeStruct> m_buffer;
141  // collections
155  bool isInTimeWindow(const int SL, const short tdc)
156  {
157  if (SL == 0) {
159  } else {
161  }
162  }
163 
164  // DB payloads
165 
166  // other
167  const bool m_enableBadWireTreatment;
168  const bool m_enableBackgroundHitFilter;
169  const bool
172  int m_nActiveWireInTotal = 0;
174  int m_nActiveWireInLayer[f_nLayer] = {0};
182 
188  void countActiveWires();
189 
190 
191  };
192 
193  } // Background namespace
195 } // Belle2 namespace
Belle2::Background::CDCHitRateCounter::m_timeWindowUpperEdge_smallCell
const int m_timeWindowUpperEdge_smallCell
upper edge of the time window for small cells [tdc count = ns]
Definition: CDCHitRateCounter.h:154
Belle2::Background::CDCHitRateCounter::m_digits
StoreArray< CDCHit > m_digits
collection of digits
Definition: CDCHitRateCounter.h:150
Belle2::Background::CDCHitRateCounter::TreeStruct::numEvents
int numEvents
number of events accumulated
Definition: CDCHitRateCounter.h:59
Belle2::Background::CDCHitRateCounter::accumulate
virtual void accumulate(unsigned timeStamp) override
Accumulate hits.
Definition: CDCHitRateCounter.cc:71
Belle2::Background::CDCHitRateCounter::f_nSuperLayer
static const int f_nSuperLayer
the number of super layers
Definition: CDCHitRateCounter.h:42
Belle2::Background::CDCHitRateCounter::clear
virtual void clear() override
Clear time-stamp buffer to prepare for 'accumulate'.
Definition: CDCHitRateCounter.cc:66
Belle2::Background::CDCHitRateCounter::m_nActiveWireInLayer
int m_nActiveWireInLayer[f_nLayer]
the number of wires used in this hit-rate calculation in each layer
Definition: CDCHitRateCounter.h:182
Belle2::Background::CDCHitRateCounter::TreeStruct::timeWindowForSmallCell
int timeWindowForSmallCell
time window for the small cells in ns
Definition: CDCHitRateCounter.h:53
Belle2::Background::CDCHitRateCounter::m_enableBadWireTreatment
const bool m_enableBadWireTreatment
flag to enable the bad wire treatment.
Definition: CDCHitRateCounter.h:175
Belle2::Background::CDCHitRateCounter::TreeStruct::nActiveWireInTotal
int nActiveWireInTotal
number of wires used in this analysis in the whole CDC
Definition: CDCHitRateCounter.h:58
Belle2::Background::CDCHitRateCounter::m_timeWindowUpperEdge_normalCell
const int m_timeWindowUpperEdge_normalCell
upper edge of the time window for normal cells [tdc count = ns]
Definition: CDCHitRateCounter.h:156
Belle2::Background::CDCHitRateCounter::m_timeWindowLowerEdge_normalCell
const int m_timeWindowLowerEdge_normalCell
lower edge of the time window for normal cells [tdc count = ns]
Definition: CDCHitRateCounter.h:155
Belle2::Background::CDCHitRateCounter::m_buffer
std::map< unsigned, TreeStruct > m_buffer
average rates in time stamps
Definition: CDCHitRateCounter.h:147
Belle2::Background::CDCHitRateCounter::initialize
virtual void initialize(TTree *tree) override
Class initializer: set branch addresses and other staf.
Definition: CDCHitRateCounter.cc:38
Belle2::Background::CDCHitRateCounter::m_enableMarkBackgroundHit
const bool m_enableMarkBackgroundHit
flag to enable to mark background flag on CDCHit (set 0x100 bit for CDCHit::m_status).
Definition: CDCHitRateCounter.h:178
Belle2::Background::CDCHitRateCounter::TreeStruct::timeWindowForNormalCell
int timeWindowForNormalCell
time window for the normal cells in ns
Definition: CDCHitRateCounter.h:54
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::Background::CDCHitRateCounter::TreeStruct::normalize
void normalize()
normalize accumulated hits to hit rate in kHz
Definition: CDCHitRateCounter.h:65
Belle2::Background::CDCHitRateCounter::TreeStruct::superLayerHitRate
float superLayerHitRate[f_nSuperLayer]
SuperLayer average hit rate in kHz.
Definition: CDCHitRateCounter.h:51
Belle2::Background::CDCHitRateCounter::TreeStruct
tree structure
Definition: CDCHitRateCounter.h:49
Belle2::Background::CDCHitRateCounter::TreeStruct::valid
bool valid
status: true = rates valid
Definition: CDCHitRateCounter.h:60
Belle2::Background::CDCHitRateCounter::TreeStruct::averageRate
float averageRate
total detector average hit rate in KHz
Definition: CDCHitRateCounter.h:52
Belle2::Background::CDCHitRateCounter::TreeStruct::nActiveWireInSuperLayer
int nActiveWireInSuperLayer[f_nSuperLayer]
number of wires used in this analysis in each super layer
Definition: CDCHitRateCounter.h:57
Belle2::Background::CDCHitRateCounter::countActiveWires_countAll
void countActiveWires_countAll()
set m_nActiveWireInTotal, m_nActiveWireInLayer[] and m_nActiveWireInSuperLayer[].
Definition: CDCHitRateCounter.cc:186
Belle2::Background::CDCHitRateCounter::m_rates
TreeStruct m_rates
tree variables
Definition: CDCHitRateCounter.h:144
Belle2::Background::CDCHitRateCounter::isInTimeWindow
bool isInTimeWindow(const int SL, const short tdc)
return true if the hit is in the given time window
Definition: CDCHitRateCounter.h:163
Belle2::Background::CDCHitRateCounter::CDCHitRateCounter
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.
Definition: CDCHitRateCounter.h:96
Belle2::Background::CDCHitRateCounter::m_enableBackgroundHitFilter
const bool m_enableBackgroundHitFilter
flag to enable the CDC background hit (crosstakl, noise) filter.
Definition: CDCHitRateCounter.h:176
Belle2::Background::CDCHitRateCounter::countActiveWires
void countActiveWires()
set m_nActiveWireInTotal, m_nActiveWireInLayer[] and m_nActiveWireInSuperLayer[].
Definition: CDCHitRateCounter.cc:212
Belle2::Background::CDCHitRateCounter::m_nActiveWireInTotal
int m_nActiveWireInTotal
the number of wires used in this hit-rate calculation in the whole CDC
Definition: CDCHitRateCounter.h:180
Belle2::Background::CDCHitRateCounter::m_timeWindowLowerEdge_smallCell
const int m_timeWindowLowerEdge_smallCell
lower edge of the time window for small cells [tdc count = ns]
Definition: CDCHitRateCounter.h:153
Belle2::StoreArray
Accessor to arrays stored in the data store.
Definition: ECLMatchingPerformanceExpertModule.h:33
Belle2::Background::CDCHitRateCounter::TreeStruct::layerHitRate
float layerHitRate[f_nLayer]
Layer average hit rate in kHz.
Definition: CDCHitRateCounter.h:50
Belle2::Background::CDCHitRateCounter::TreeStruct::nActiveWireInLayer
int nActiveWireInLayer[f_nLayer]
number of wires used in this analysis in each layer
Definition: CDCHitRateCounter.h:56
Belle2::Background::CDCHitRateCounter::normalize
virtual void normalize(unsigned timeStamp) override
Normalize accumulated hits (e.g.
Definition: CDCHitRateCounter.cc:144
Belle2::Background::CDCHitRateCounter::m_nActiveWireInSuperLayer
int m_nActiveWireInSuperLayer[f_nSuperLayer]
the number of wires used in this hit-rate calculation in each suler layer
Definition: CDCHitRateCounter.h:181
Belle2::Background::CDCHitRateCounter::f_nLayer
static const int f_nLayer
the number of layers
Definition: CDCHitRateCounter.h:41