Belle II Software  release-06-00-14
EclData.cc
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 //This module
9 #include <ecl/modules/eclDisplay/EclData.h>
10 
11 //ROOT
12 #include <TLeaf.h>
13 #include <TFile.h>
14 
15 //Framework
16 #include <framework/logging/Logger.h>
17 
18 //ECL
19 #include <ecl/dataobjects/ECLCalDigit.h>
20 
21 using namespace Belle2;
22 
23 const int EclData::ring_start_id[70] = {
24  // forward (0-12)
25  1, 49, 97, 161, 225, 289, 385, 481, 577, 673,
26  769, 865, 1009,
27  // barrel (13-58)
28  1153, 1297, 1441, 1585, 1729, 1873, 2017, 2161, 2305, 2449,
29  2593, 2737, 2881, 3025, 3169, 3313, 3457, 3601, 3745, 3889,
30  4033, 4177, 4321, 4465, 4609, 4753, 4897, 5041, 5185, 5329,
31  5473, 5617, 5761, 5905, 6049, 6193, 6337, 6481, 6625, 6769,
32  6913, 7057, 7201, 7345, 7489, 7633,
33  // forward (59-68)
34  7777, 7921, 8065, 8161, 8257, 8353, 8449, 8545, 8609, 8673,
35  // last_crystal+1
37 };
38 
40 {
41  m_tree = new TTree("tree", "tree");
42  m_tree->Branch("ch", &m_branch_ch, "ch/I");
43  m_tree->Branch("energy", &m_branch_energy, "energy/D");
44  m_tree->Branch("time", &m_branch_time, "time/D");
45  m_tree->Branch("evtn", &m_branch_evtn, "evtn/I");
46 
47  initVariables();
49 }
50 
52 {
53  initVariables();
55 
56  cloneFrom(data);
57 }
58 
60 {
61  if (this != &other) { // self-assignment check
62  initVariables();
64 
65  cloneFrom(other);
66  }
67  return *this;
68 }
69 
70 EclData::~EclData()
71 {
72  delete m_event_counts;
73  delete m_energy_sums;
74  delete m_tree;
75 }
76 
77 void EclData::cloneFrom(const EclData& other)
78 {
79  m_tree = other.m_tree->CloneTree();
80 
81  m_tree = new TTree("tree", "tree");
82  m_tree->Branch("ch", &m_branch_ch, "ch/I");
83  m_tree->Branch("energy", &m_branch_energy, "energy/I");
84  m_tree->Branch("time", &m_branch_time, "time/I");
85  m_tree->Branch("evtn", &m_branch_evtn, "evtn/I");
86 
88  m_time_max = other.m_time_max;
89  // m_excluded_ch is not copied.
90  m_excluded_ch.clear();
91 
94 
95  for (int i = 0; i < getCrystalCount() + 1; i++) {
96  m_event_counts[i] = other.m_event_counts[i];
97  m_energy_sums[i] = other.m_energy_sums[i];
98  }
99 
100  // TODO: Vectors surely can be copied in more simple way.
101  for (unsigned int i = 0; i < other.m_event_entry.size(); i++) {
102  m_event_entry.push_back(other.m_event_entry[i]);
103  }
104 
107 
110 
113 }
114 
116 {
117  // Setting
119 
120  m_event_counts = new int[getCrystalCount() + 1];
121  m_energy_sums = new float[getCrystalCount() + 1];
122 
123  m_event_count_max = 0;
124  m_energy_sums_max = 0;
125  m_energy_total = 0;
126 
127  m_last_event_id = -1;
128 
129  m_time_max = 0;
130 
131  m_excluded_ch.clear();
132 }
133 
135 {
136  m_en_range_min = 0;
137  m_en_range_max = -1;
138 
139  m_time_range_min = -2048;
141 
142  m_ev_range_min = 0;
144 }
145 
147 {
148  return 8736;
149 }
150 
152 {
153  return m_tree;
154 }
155 
157 {
158  return m_event_counts;
159 }
160 // Alias for getEventCounts()
162 {
163  return getEventCounts();
164 }
165 
167 {
168  return m_event_count_max;
169 }
170 
172 {
173  return m_energy_sums;
174 }
175 // Alias for getEnergySums()
177 {
178  return getEnergySums();
179 }
180 
182 {
183  return m_energy_sums_max;
184 }
185 
187 {
188  return m_energy_total;
189 }
190 
192 {
193  switch (subsys) {
194  case ALL:
195  return true;
196  case BARR:
197  return crystal >= 1153 && crystal <= 7776;
198  case FORW:
199  return crystal < 1153;
200  case BACKW:
201  return crystal > 7776;
202  default:
203  return false;
204  }
205 }
206 
208 {
209  return m_time_range_min;
210 }
212 {
213  return m_time_range_max;
214 }
215 
216 void EclData::setTimeRange(int time_min, int time_max, bool do_update)
217 {
218  if (m_time_range_min == time_min && m_time_range_max == time_max)
219  return;
220 
221  m_time_range_min = time_min;
222  m_time_range_max = time_max;
223  if (do_update)
224  update();
225 }
226 
228 {
229  return m_ev_range_min;
230 }
232 {
233  return m_ev_range_max;
234 }
235 void EclData::setEventRange(int ev_min, int ev_max, bool do_update)
236 {
237  if (ev_min == m_ev_range_min && ev_max == m_ev_range_max)
238  return;
239 
240  m_ev_range_min = ev_min;
241  m_ev_range_max = ev_max;
242  if (do_update)
243  update();
244 }
245 
246 void EclData::setEnergyThreshold(int en_min, int en_max, bool do_update)
247 {
248  if (en_min == m_en_range_min && en_max == m_en_range_max)
249  return;
250 
251  m_en_range_min = en_min;
252  m_en_range_max = en_max;
253  if (do_update)
254  update();
255 }
256 
258 {
259  return m_time_range_min;
260 }
262 {
263  return m_time_max;
264 }
265 
267 {
268  return m_last_event_id;
269 }
270 
271 int EclData::getChannel(int phi_id, int theta_id)
272 {
273  return ring_start_id[theta_id] + phi_id;
274 }
275 
276 int EclData::getPhiId(int _ch)
277 {
278  for (int i = 0; i < 69; i++) {
279  if (_ch < ring_start_id[i + 1])
280  return _ch - ring_start_id[i];
281  }
282 
283  return -1;
284 }
285 
287 {
288  for (int i = 0; i < 69; i++) {
289  if (_ch < ring_start_id[i + 1])
290  return i;
291  }
292 
293  return -1;
294 }
295 
296 void EclData::excludeChannel(int _ch, bool do_update)
297 {
298  m_excluded_ch.insert(_ch);
299  if (do_update)
300  update();
301 }
302 
303 void EclData::includeChannel(int _ch, bool do_update)
304 {
305  m_excluded_ch.erase(_ch);
306  if (do_update)
307  update();
308 }
309 
310 void EclData::loadRootFile(const char* path)
311 {
312  // TODO: Move this process to some new reset() method.
313  m_excluded_ch.clear();
314  m_event_entry.clear();
315  m_tree->Reset();
316  m_last_event_id = -1;
317  m_en_range_max = -1;
318  m_time_max = 0;
319 
320  TFile* file = new TFile(path, "READ");
321  TTree* tree = 0;
322  file->GetObject("tree", tree);
323  long nentries = tree->GetEntries();
324 
325  tree->SetBranchAddress("ECLCalDigits.m_CellId", &m_branch_ch);
326  tree->SetBranchAddress("ECLCalDigits.m_Energy", &m_branch_energy);
327  tree->SetBranchAddress("ECLCalDigits.m_Time" , &m_branch_time);
328 
329  TLeaf* leafCellId;
330  TLeaf* leafEnergy;
331  TLeaf* leafTimeFit;
332 
333  leafCellId = tree->GetLeaf("ECLCalDigits.m_CellId");
334  leafEnergy = tree->GetLeaf("ECLCalDigits.m_Energy");
335  leafTimeFit = tree->GetLeaf("ECLCalDigits.m_Time");
336 
337  for (long i = 0; i < nentries; i++) {
338  tree->GetEntry(i);
339 
340  const int len = leafCellId->GetLen();
341 
342  if (len > 0) {
343  m_last_event_id++;
345 
346  m_event_entry.push_back(m_tree->GetEntries());
347  }
348 
349  for (int j = 0; j < len; j++) {
350  m_branch_ch = leafCellId ->GetTypedValue<int>(j);
351  m_branch_energy = leafEnergy ->GetTypedValue<double>(j) * 1e3;
352  m_branch_time = leafTimeFit->GetTypedValue<double>(j);
353 
355 
356  B2DEBUG(10, "Added event ECLCalDigits for crystal " << m_branch_ch);
357 
358  // TODO: Add energy cut.
359  m_tree->Fill();
360  }
361  }
362 
363  update(true);
364 }
365 
366 void EclData::update(bool reset_event_ranges)
367 {
368  int start, end;
369 
370  if (reset_event_ranges) initEventRanges();
371 
372  if (m_ev_range_min < 0)
373  return;
375  return;
376 
377  start = m_event_entry[m_ev_range_min];
379  end = m_event_entry[m_ev_range_max + 1];
380  else
381  end = m_tree->GetEntries();
382 
383  for (int i = 1; i <= getCrystalCount(); i++) {
384  m_event_counts[i] = 0;
385  m_energy_sums[i] = 0;
386  }
387  m_event_count_max = 0;
388  m_energy_sums_max = 0;
389  m_energy_total = 0;
390 
391  for (int i = start; i < end; i++) {
392  m_tree->GetEntry(i);
393 
394  if (m_excluded_ch.count(m_branch_ch) > 0) continue;
395 
396  if (m_branch_ch >= 1 && m_branch_ch <= getCrystalCount()) {
397  // Check if current time belongs to time_range, set by user.
398  if (m_time_range_max >= 0)
399  if (m_branch_time < m_time_range_min || m_branch_time > m_time_range_max)
400  continue;
401 
402  if (m_en_range_max >= 0)
403  if (m_branch_energy < m_en_range_min || m_branch_energy > m_en_range_max)
404  continue;
405 
409 
412 
415  }
416  }
417 
418  if (m_time_range_max < 0)
420  if (m_ev_range_max < 0)
422 
423  B2DEBUG(250, end - start << " events handled.");
424 }
425 
426 int EclData::addEvent(ECLCalDigit* event, int _evtn)
427 {
428  if (event->getEnergy() <= 0 || event->getCellId() <= 0) {
429  return -1;
430  }
431 
432  m_branch_ch = event->getCellId();
433  m_branch_energy = event->getEnergy() * 1e3; // To MeV
434  m_branch_time = event->getTime();
435  m_branch_evtn = _evtn;
436 
438  if (m_last_event_id < _evtn) {
439  m_last_event_id = _evtn;
440  m_event_entry.push_back(m_tree->GetEntries());
441  }
442 
443  B2DEBUG(200, "Added event #" << _evtn << ", ch:" << m_branch_ch
444  << ", energy:" << m_branch_energy << ", time:" << m_branch_time);
445 
446  m_tree->Fill();
447 
448  return 0;
449 }
450 
451 void EclData::fillEnergyHistogram(TH1F* hist, int energy_min, int energy_max, EclSubsystem subsys)
452 {
453  int start, end;
454 
455  if (m_ev_range_min < 0)
456  return;
458  return;
459 
460  start = m_event_entry[m_ev_range_min];
462  end = m_event_entry[m_ev_range_max + 1];
463  else
464  end = m_tree->GetEntries();
465 
466  for (int i = start; i < end; i++) {
467  m_tree->GetEntry(i);
468  if (isCrystalInSubsystem(m_branch_ch, subsys)) {
469  if (m_branch_energy >= energy_min && m_branch_energy <= energy_max)
470  hist->Fill(m_branch_energy);
471  }
472  }
473 }
474 
475 void EclData::fillEnergySumHistogram(TH1F* hist, int energy_min, int energy_max, EclData::EclSubsystem subsys)
476 {
477  int start, end;
478 
479  if (m_ev_range_min < 0)
480  return;
482  return;
483 
484  start = m_event_entry[m_ev_range_min];
486  end = m_event_entry[m_ev_range_max + 1];
487  else
488  end = m_tree->GetEntries();
489 
490  m_tree->GetEntry(start);
491  int cur_evtn = m_branch_evtn;
492  int energy_sum = 0;
493 
494  for (int i = start; i < end; i++) {
495  m_tree->GetEntry(i);
496 
497  // After reading all of the data from event cur_evtn,
498  // save the sum in histogram.
499  if (m_branch_evtn > cur_evtn) {
500  if (energy_sum > 0) hist->Fill(energy_sum);
501  energy_sum = 0;
502  cur_evtn = m_branch_evtn;
503  }
504 
505  if (isCrystalInSubsystem(m_branch_ch, subsys)) {
506  if (m_branch_energy >= energy_min && m_branch_energy <= energy_max)
507  energy_sum += m_branch_energy;
508  }
509  }
510 
511  // Add last selected event.
512  if (energy_sum > 0) hist->Fill(energy_sum);
513 }
514 
515 void EclData::fillTimeHistogram(TH1F* hist, int time_min, int time_max, EclData::EclSubsystem subsys)
516 {
517  int start, end;
518 
519  if (m_ev_range_min < 0)
520  return;
522  return;
523 
524  start = m_event_entry[m_ev_range_min];
526  end = m_event_entry[m_ev_range_max + 1];
527  else
528  end = m_tree->GetEntries();
529 
530  for (int i = start; i < end; i++) {
531  m_tree->GetEntry(i);
532 
533  if (m_en_range_max >= 0)
534  if (m_branch_energy < m_en_range_min || m_branch_energy > m_en_range_max)
535  continue;
536 
537  if (isCrystalInSubsystem(m_branch_ch, subsys)) {
538  if (m_branch_time >= time_min && m_branch_time <= time_max)
539  hist->Fill(m_branch_time);
540  }
541  }
542 }
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:266
int addEvent(ECLCalDigit *event, int evtn)
Add ECLDigit event to inner TTree (m_tree).
Definition: EclData.cc:426
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:176
float * getEnergySums()
Get array of total energy for each channel in the specified time and event range.
Definition: EclData.cc:171
int m_time_range_min
Time range (min) for display.
Definition: EclData.h:62
int getEventCountsMax()
Definition: EclData.cc:166
int getThetaId(int ch)
ECL CellId -> theta_id.
Definition: EclData.cc:286
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:156
int * getEventCountsPerCrystal()
Alias for GetEventCounts()
Definition: EclData.cc:161
void initVariables()
Initialization of arrays.
Definition: EclData.cc:115
int getPhiId(int ch)
ECL CellId -> phi_id.
Definition: EclData.cc:276
int getChannel(int phi_id, int theta_id)
Converts (phi_id, theta_id) pair to ECL CellId.
Definition: EclData.cc:271
int getTimeMax()
Definition: EclData.cc:261
void setEventRange(int ev_min, int ev_max, bool do_update=true)
Sets event range to (ev_min, ev_max)
Definition: EclData.cc:235
void update(bool reset_event_ranges=false)
Update time_min, time_max, event_counts and energy_sums.
Definition: EclData.cc:366
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:77
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:515
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:134
void setTimeRange(int time_min, int time_max, bool do_update=true)
Display only events in the specified time range.
Definition: EclData.cc:216
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:186
void includeChannel(int ch, bool do_update=false)
Includes specific channel in the count of events and energy.
Definition: EclData.cc:303
int getTimeRangeMax()
Return max time in time range.
Definition: EclData.cc:211
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:246
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:475
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:181
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:451
static int getCrystalCount()
Get number of crystals in ECL.
Definition: EclData.cc:146
int getTimeRangeMin()
Return min time in time range.
Definition: EclData.cc:207
int m_branch_ch
Tree channel field.
Definition: EclData.h:36
int getEventRangeMax()
Return max event number to display.
Definition: EclData.cc:231
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:191
EclData & operator=(const EclData &other)
Assignment operator: utilizes copy constructor.
Definition: EclData.cc:59
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:227
EclData()
Default constructor.
Definition: EclData.cc:39
void excludeChannel(int ch, bool do_update=false)
Excludes specific channel from the count of events and energy.
Definition: EclData.cc:296
void loadRootFile(const char *path)
Load root file containing ECLCalDigit data from the specified path.
Definition: EclData.cc:310
int getTimeMin()
Definition: EclData.cc:257
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:151
Abstract base class for different kinds of events.