Belle II Software development
ECLHitRateCounter.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 <ecl/dataobjects/ECLDigit.h>
14#include <ecl/dataobjects/ECLDsp.h>
15#include <ecl/geometry/ECLGeometryPar.h>
16#include <TTree.h>
17
18namespace Belle2 {
23 namespace Background {
24
29
30 public:
31
35 struct TreeStruct {
36
37 float averageDspBkgRate[16] = {0};
38 int numEvents = 0;
39 bool validDspRate = false;
44 void normalize()
45 {
46 if (numEvents == 0) return;
47
48 for (int i = 0; i < 16; i++) {
50 }
51 }
52
53
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
87 // class parameters: to be set via constructor or setters
88
89 //functions
102 void segmentECL();
103
108 int findECLSegment(int cellid)
109 {
110 return m_segmentMap.find(cellid)->second;
111 }
112
118
119 // tree structure
122 // buffer
123 std::map<unsigned, TreeStruct> m_buffer;
125 // collections
129 std::vector<float> m_ADCtoEnergy;
130 std::vector<float> m_waveformNoise;
132 // other
134 std::map<int, int> m_segmentMap;
135 int m_crystalsInSegment[16] = {0};
136 };
137 }
139}
140
Class for monitoring beam background hit rates of ECL.
std::map< int, int > m_segmentMap
map with keys containing ECL CellID and values containing segment number
int m_crystalsInSegment[16]
array cotaining the number of crystals in given segment
std::map< unsigned, TreeStruct > m_buffer
average rates in time stamps
void findElectronicsNoise()
Find the electronics noise correction for each cellID Reads a file with a histogram containing electr...
StoreArray< ECLDsp > m_dsps
collection of ECL waveforms
virtual void initialize(TTree *tree) override
Class initializer: set branch addresses and other staf.
Belle2::ECL::ECLGeometryPar * m_geometry
pointer to ECLGeometryPar
std::vector< float > m_waveformNoise
vector used to store ECL electronic noise constants foe each crystal
virtual void accumulate(unsigned timeStamp) override
Accumulate hits.
StoreArray< ECLDigit > m_digits
collection of digits
int findECLSegment(int cellid)
Find the correcsponding ECL segment based on the cellID.
void segmentECL()
Performs ECL segmentation; Done once per run; Populates a map which connects each ECL crystal with a ...
virtual void normalize(unsigned timeStamp) override
Normalize accumulated hits (e.g.
virtual void clear() override
Clear time-stamp buffer to prepare for 'accumulate'.
std::vector< float > m_ADCtoEnergy
vector used to store ECL calibration constants for each crystal
Abstract base class for monitoring beam background hit rates All the monitoring classes must inherit ...
Definition: HitRateBase.h:24
The Class for ECL Geometry Parameters.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
Abstract base class for different kinds of events.
bool validDspRate
status for rates calculated from waveforms, true if waveforms for all crystals are recorded
float averageDspBkgRate[16]
average background rate per crystal in given segment, calculated using ECL waveforms [hits/second]
void normalize()
normalize accumulated rates based on ECL waveforms