Belle II Software  release-05-02-19
ECLHitRateCounter.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, Andrea Fodor *
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 <ecl/dataobjects/ECLDigit.h>
16 #include <ecl/dataobjects/ECLDsp.h>
17 #include <ecl/geometry/ECLGeometryPar.h>
18 #include <TTree.h>
19 
20 namespace Belle2 {
25  namespace Background {
26 
30  class ECLHitRateCounter: public HitRateBase {
31 
32  public:
33 
37  struct TreeStruct {
38 
39  float averageDspBkgRate[16] = {0};
40  int numEvents = 0;
41  bool validDspRate = false;
46  void normalize()
47  {
48  if (numEvents == 0) return;
49 
50  for (int i = 0; i < 16; i++) {
52  }
53  }
54 
55 
56  };
57 
62  {}
63 
68  virtual void initialize(TTree* tree) override;
69 
73  virtual void clear() override;
74 
79  virtual void accumulate(unsigned timeStamp) override;
80 
85  virtual void normalize(unsigned timeStamp) override;
86 
87  private:
88 
89  // class parameters: to be set via constructor or setters
90 
91  //functions
104  void segmentECL();
105 
110  int findECLSegment(int cellid)
111  {
112  return m_segmentMap.find(cellid)->second;
113  }
114 
119  void findElectronicsNoise();
120 
121  // tree structure
122  TreeStruct m_rates;
124  // buffer
125  std::map<unsigned, TreeStruct> m_buffer;
127  // collections
131  std::vector<float> m_ADCtoEnergy;
132  std::vector<float> m_waveformNoise;
134  // other
136  std::map<int, int> m_segmentMap;
137  int m_crystalsInSegment[16] = {0};
138  };
139  }
141 }
142 
Belle2::Background::ECLHitRateCounter::m_geometry
Belle2::ECL::ECLGeometryPar * m_geometry
pointer to ECLGeometryPar
Definition: ECLHitRateCounter.h:143
Belle2::Background::ECLHitRateCounter::TreeStruct::normalize
void normalize()
normalize accumulated rates based on ECL waveforms
Definition: ECLHitRateCounter.h:54
Belle2::Background::ECLHitRateCounter::initialize
virtual void initialize(TTree *tree) override
Class initializer: set branch addresses and other staf.
Definition: ECLHitRateCounter.cc:30
Belle2::Background::ECLHitRateCounter::normalize
virtual void normalize(unsigned timeStamp) override
Normalize accumulated hits (e.g.
Definition: ECLHitRateCounter.cc:104
Belle2::Background::ECLHitRateCounter::ECLHitRateCounter
ECLHitRateCounter()
Constructor.
Definition: ECLHitRateCounter.h:69
Belle2::Background::ECLHitRateCounter::m_segmentMap
std::map< int, int > m_segmentMap
map with keys containing ECL CellID and values containing segment number
Definition: ECLHitRateCounter.h:144
Belle2::Background::ECLHitRateCounter::m_waveformNoise
std::vector< float > m_waveformNoise
vector used to store ECL electronic noise constants foe each crystal
Definition: ECLHitRateCounter.h:140
Belle2::Background::ECLHitRateCounter::TreeStruct::validDspRate
bool validDspRate
status for rates calculated from waveforms, true if waveforms for all crystals are recorded
Definition: ECLHitRateCounter.h:49
Belle2::Background::ECLHitRateCounter::m_dsps
StoreArray< ECLDsp > m_dsps
collection of ECL waveforms
Definition: ECLHitRateCounter.h:137
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::Background::ECLHitRateCounter::m_rates
TreeStruct m_rates
tree variables
Definition: ECLHitRateCounter.h:130
Belle2::Background::ECLHitRateCounter::m_buffer
std::map< unsigned, TreeStruct > m_buffer
average rates in time stamps
Definition: ECLHitRateCounter.h:133
Belle2::Background::ECLHitRateCounter::m_ADCtoEnergy
std::vector< float > m_ADCtoEnergy
vector used to store ECL calibration constants for each crystal
Definition: ECLHitRateCounter.h:139
Belle2::Background::ECLHitRateCounter::TreeStruct::averageDspBkgRate
float averageDspBkgRate[16]
average background rate per crystal in given segment, calculated using ECL waveforms [hits/second]
Definition: ECLHitRateCounter.h:47
Belle2::Background::ECLHitRateCounter::findElectronicsNoise
void findElectronicsNoise()
Find the electronics noise correction for each cellID Reads a file with a histogram containing electr...
Belle2::Background::ECLHitRateCounter::m_digits
StoreArray< ECLDigit > m_digits
collection of digits
Definition: ECLHitRateCounter.h:136
Belle2::Background::ECLHitRateCounter::TreeStruct::numEvents
int numEvents
number of valid events
Definition: ECLHitRateCounter.h:48
Belle2::ECL::ECLGeometryPar
The Class for ECL Geometry Parameters.
Definition: ECLGeometryPar.h:45
Belle2::StoreArray
Accessor to arrays stored in the data store.
Definition: ECLMatchingPerformanceExpertModule.h:33
Belle2::Background::ECLHitRateCounter::findECLSegment
int findECLSegment(int cellid)
Find the correcsponding ECL segment based on the cellID.
Definition: ECLHitRateCounter.h:118
Belle2::Background::ECLHitRateCounter::clear
virtual void clear() override
Clear time-stamp buffer to prepare for 'accumulate'.
Definition: ECLHitRateCounter.cc:49
Belle2::Background::ECLHitRateCounter::m_crystalsInSegment
int m_crystalsInSegment[16]
array cotaining the number of crystals in given segment
Definition: ECLHitRateCounter.h:145
Belle2::Background::ECLHitRateCounter::accumulate
virtual void accumulate(unsigned timeStamp) override
Accumulate hits.
Definition: ECLHitRateCounter.cc:54
Belle2::Background::ECLHitRateCounter::segmentECL
void segmentECL()
Performs ECL segmentation; Done once per run; Populates a map which connects each ECL crystal with a ...
Definition: ECLHitRateCounter.cc:121