Belle II Software development
EclData.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/* ROOT headers. */
12#include <TH1F.h>
13#include <TTree.h>
14
15/* C++ headers. */
16#include <set>
17#include <vector>
18
19namespace Belle2 {
25 class ECLCalDigit;
26
31 class EclData {
32 private:
34 TTree* m_tree;
43
54
57
66
75
80
82 std::set<int> m_excluded_ch;
84 std::vector<int> m_event_entry;
85
89 static const int ring_start_id[70];
90
91 public:
98 enum EclSubsystem {ALL, BARR, FORW, BACKW};
99
103 EclData();
104
108 EclData(const EclData& data);
112 EclData& operator=(const EclData& other);
113
114
115 ~EclData();
116
117 private:
121 void cloneFrom(const EclData& other);
122
126 void initVariables();
130 void initEventRanges();
131
132 public:
136 static int getCrystalCount();
137
141 TTree* getTree();
142
146 int* getEventCounts();
154 int getEventCountsMax();
155
160 float* getEnergySums();
164 float* getEnergySumPerCrystal();
168 float getEnergySumsMax();
172 float getEnergyTotal();
173
179 bool isCrystalInSubsystem(int crystal, EclSubsystem subsys);
180
184 int getTimeRangeMin();
188 int getTimeRangeMax();
193 void setTimeRange(int time_min, int time_max, bool do_update = true);
194
198 int getEventRangeMin();
202 int getEventRangeMax();
210 void setEventRange(int ev_min, int ev_max, bool do_update = true);
211
219 void setEnergyThreshold(int en_min, int en_max, bool do_update = true);
220
224 int getTimeMin();
228 int getTimeMax();
229
233 int getLastEventId();
234
238 int getChannel(int phi_id, int theta_id);
242 int getPhiId(int ch);
246 int getThetaId(int ch);
247
251 void excludeChannel(int ch, bool do_update = false);
255 void includeChannel(int ch, bool do_update = false);
256
260 void loadRootFile(const char* path);
261
265 void update(bool reset_event_ranges = false);
266
273 int addEvent(ECLCalDigit* event, int evtn);
278 void fillEnergyHistogram(TH1F* hist, int energy_min, int energy_max, EclSubsystem subsys);
283 void fillEnergySumHistogram(TH1F* hist, int energy_min, int energy_max, EclSubsystem subsys);
288 void fillTimeHistogram(TH1F* hist, int time_min, int time_max, EclSubsystem subsys);
289 };
291}
Class to store calibrated ECLDigits: ECLCalDigits.
Definition: ECLCalDigit.h:23
This class contains data for ECLSimHit's and provides several relevant conversion functions for bette...
Definition: EclData.h:31
int getLastEventId()
Definition: EclData.cc:268
int addEvent(ECLCalDigit *event, int evtn)
Add ECLDigit event to inner TTree (m_tree).
Definition: EclData.cc:428
int m_event_count_max
Max value in event_counts array.
Definition: EclData.h:47
TTree * m_tree
Tree with loaded events.
Definition: EclData.h:34
float m_energy_total
Total energy for last displayed range of events.
Definition: EclData.h:53
float * getEnergySumPerCrystal()
Alias for GetEnergySums()
Definition: EclData.cc:178
float * getEnergySums()
Get array of total energy for each channel in the specified time and event range.
Definition: EclData.cc:173
int m_time_range_min
Time range (min) for display.
Definition: EclData.h:62
int getEventCountsMax()
Definition: EclData.cc:168
int getThetaId(int ch)
ECL CellId -> theta_id.
Definition: EclData.cc:288
int m_en_range_max
Upper boundary of energy threshold.
Definition: EclData.h:79
int m_branch_evtn
Tree event number branch.
Definition: EclData.h:42
int * getEventCounts()
Returns array of event counts per crystal[getCrystalsCount()].
Definition: EclData.cc:158
int * getEventCountsPerCrystal()
Alias for GetEventCounts()
Definition: EclData.cc:163
void initVariables()
Initialization of arrays.
Definition: EclData.cc:117
int getPhiId(int ch)
ECL CellId -> phi_id.
Definition: EclData.cc:278
int getChannel(int phi_id, int theta_id)
Converts (phi_id, theta_id) pair to ECL CellId.
Definition: EclData.cc:273
int getTimeMax()
Definition: EclData.cc:263
void setEventRange(int ev_min, int ev_max, bool do_update=true)
Sets event range to (ev_min, ev_max)
Definition: EclData.cc:237
void update(bool reset_event_ranges=false)
Update time_min, time_max, event_counts and energy_sums.
Definition: EclData.cc:368
static const int ring_start_id[70]
First crystal id in the beginning of i-th ECL ring.
Definition: EclData.h:89
void cloneFrom(const EclData &other)
Clone attributes from other EclData.
Definition: EclData.cc:79
void fillTimeHistogram(TH1F *hist, int time_min, int time_max, EclSubsystem subsys)
Fill time histogram for the specified EclSubsystem (Barrel, forward endcap, backward endcap,...
Definition: EclData.cc:517
int * m_event_counts
Number of events for each crystal.
Definition: EclData.h:45
int m_ev_range_min
Events from ev_range_min will be counted in energy_sums and event_counts.
Definition: EclData.h:70
void initEventRanges()
Set initial values for time and event range.
Definition: EclData.cc:136
void setTimeRange(int time_min, int time_max, bool do_update=true)
Display only events in the specified time range.
Definition: EclData.cc:218
std::vector< int > m_event_entry
This vector holds the position of each entry which starts an event.
Definition: EclData.h:84
float getEnergyTotal()
Definition: EclData.cc:188
void includeChannel(int ch, bool do_update=false)
Includes specific channel in the count of events and energy.
Definition: EclData.cc:305
int getTimeRangeMax()
Return max time in time range.
Definition: EclData.cc:213
float * m_energy_sums
Sum of energies of every event captured by crystal (MeV).
Definition: EclData.h:49
int m_time_max
Maximum possible time.
Definition: EclData.h:59
void setEnergyThreshold(int en_min, int en_max, bool do_update=true)
Sets energy range to (en_min, en_max).
Definition: EclData.cc:248
void fillEnergySumHistogram(TH1F *hist, int energy_min, int energy_max, EclSubsystem subsys)
Fill energy per event histogram for the specified EclSubsystem (Barrel, forward endcap,...
Definition: EclData.cc:477
int m_last_event_id
Id of the event with max recorded event number (evtn).
Definition: EclData.h:56
double m_branch_energy
Tree energy branch.
Definition: EclData.h:38
float getEnergySumsMax()
Definition: EclData.cc:183
void fillEnergyHistogram(TH1F *hist, int energy_min, int energy_max, EclSubsystem subsys)
Fill energy per channel histogram for the specified EclSubsystem (Barrel, forward endcap,...
Definition: EclData.cc:453
static int getCrystalCount()
Get number of crystals in ECL.
Definition: EclData.cc:148
int getTimeRangeMin()
Return min time in time range.
Definition: EclData.cc:209
int m_branch_ch
Tree channel field.
Definition: EclData.h:36
int getEventRangeMax()
Return max event number to display.
Definition: EclData.cc:233
std::set< int > m_excluded_ch
Channels which are filtered out from the count of events and energy.
Definition: EclData.h:82
double m_branch_time
Tree time branch.
Definition: EclData.h:40
bool isCrystalInSubsystem(int crystal, EclSubsystem subsys)
Definition: EclData.cc:193
EclData & operator=(const EclData &other)
Assignment operator: utilizes copy constructor.
Definition: EclData.cc:61
float m_energy_sums_max
Max value in m_energy_sums array.
Definition: EclData.h:51
int m_time_range_max
Time range (max) for display.
Definition: EclData.h:65
int getEventRangeMin()
Return min event number to display.
Definition: EclData.cc:229
EclData()
Default constructor.
Definition: EclData.cc:41
void excludeChannel(int ch, bool do_update=false)
Excludes specific channel from the count of events and energy.
Definition: EclData.cc:298
void loadRootFile(const char *path)
Load root file containing ECLCalDigit data from the specified path.
Definition: EclData.cc:312
int getTimeMin()
Definition: EclData.cc:259
int m_en_range_min
Lower boundary of energy threshold.
Definition: EclData.h:77
EclSubsystem
Subsystems of ECL: ALL all subsystems BARR barrel only FORW forward endcap only BACKW backward endcap...
Definition: EclData.h:98
int m_ev_range_max
Events up to ev_range_max will be counted in energy_sums and event_counts.
Definition: EclData.h:74
TTree * getTree()
Returns data contained in EclDisplay.
Definition: EclData.cc:153
Abstract base class for different kinds of events.