Belle II Software  release-05-01-25
EclPainterCommon.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2015 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Milkail Remnev, Dmitry Matvienko *
7  * *
8  * This software is provided "as is" without any warranty. *
9  ***************************************************************************/
10 //This module
11 #include <ecl/modules/eclDisplay/EclPainterCommon.h>
12 
13 //Root
14 #include <TH1.h>
15 
16 using namespace Belle2;
17 
19  EclPainter(data)
20 {
21  m_type = type;
22 
23  initHisto();
24  m_hist->GetXaxis()->CenterTitle();
25  m_hist->GetXaxis()->SetTitleOffset(1.1);
26  m_hist->GetYaxis()->SetTitleOffset(1.1);
27 }
28 
30 {
31  delete m_hist;
32 }
33 
35 {
36  switch (m_type) {
37  case ENERGY:
38  return 0;
39  case ENERGY_SUM:
40  return 0;
41  case TIME:
42  return -2048;
43  }
44 
45  return 0;
46 }
47 
49 {
50  switch (m_type) {
51  case ENERGY:
52  return 150;
53  case ENERGY_SUM:
54  return 2500;
55  case TIME:
56  return 2048;
57  }
58 
59  return 0;
60 }
61 
63 {
64 
65  char obj_name[255];
66  getNewRootObjectName(obj_name, 255);
67  m_hist = new TH1F(obj_name, "title", getMaxX() / 10,
68  getMinX(), getMaxX());
69 }
70 
72 {
73  EclPainter::getInformation(px, py, panel);
74 }
75 
77 {
78  EclData* data = getData();
79 
80  setTitles();
81  m_hist->Reset();
82 
83  switch (getType()) {
84  case ENERGY:
85  data->fillEnergyHistogram(m_hist, getMinX(), getMaxX(), getDisplayedSubsystem());
86  break;
87  case ENERGY_SUM:
88  data->fillEnergySumHistogram(m_hist, getMinX(), getMaxX(), getDisplayedSubsystem());
89  break;
90  case TIME:
91  data->fillTimeHistogram(m_hist, getMinX(), getMaxX(), getDisplayedSubsystem());
92  break;
93  }
94 
95  m_hist->Draw();
96 }
97 
98 
100 {
101  const char* name[3] = {
102  "Energy per channel",
103  "Energy sum per event",
104  "Time"
105  };
106  const char* xname[3] = {
107  "Energy (MeV)",
108  "Energy (MeV)",
109  "Time (ns)"
110  };
111 
112  int type = (int)getType();
113 
114  TString title = TString(name[type]) + " (" +
115  getSubsystemTitle(getDisplayedSubsystem()) + ")";
116  m_hist->SetTitle(title);
117 
118  m_hist->SetXTitle(xname[type]);
119  m_hist->SetYTitle("");
120 }
121 
123 {
124  return m_type;
125 }
Belle2::EclPainterCommon::getMinX
int getMinX()
Return m_x_min.
Definition: EclPainterCommon.cc:34
Belle2::EclPainterCommon::ENERGY
@ ENERGY
Energy per channel distribution.
Definition: EclPainterCommon.h:33
Belle2::EclPainterCommon::Draw
void Draw() override
Redraw the canvas.
Definition: EclPainterCommon.cc:76
Belle2::EclPainterCommon::Type
Type
Subtype of histogram to draw.
Definition: EclPainterCommon.h:32
Belle2::EclPainterCommon::m_type
Type m_type
Display subtypes of this class.
Definition: EclPainterCommon.h:49
Belle2::EclPainter::getData
EclData * getData()
Return currently displayed EclData.
Definition: EclPainter.h:48
Belle2::EclPainterCommon::getMaxX
int getMaxX()
Return m_x_max.
Definition: EclPainterCommon.cc:48
Belle2::EclPainterCommon::TIME
@ TIME
Time distribution.
Definition: EclPainterCommon.h:35
Belle2::EclPainter::getDisplayedSubsystem
EclData::EclSubsystem getDisplayedSubsystem()
Get currently displayed ECL subsystem.
Definition: EclPainter.cc:52
Belle2::EclPainterCommon::setTitles
void setTitles()
Update titles of the histogram.
Definition: EclPainterCommon.cc:99
Belle2::EclPainter::getInformation
virtual void getInformation(int px, int py, MultilineWidget *panel)
Sets the information to be displayed in the provided MultilineWidget.
Definition: EclPainter.cc:73
Belle2::EclPainterCommon::~EclPainterCommon
virtual ~EclPainterCommon()
Destructor for EclPainter subclass.
Definition: EclPainterCommon.cc:29
Belle2::EclPainter
Painter for EclData, parent class, created with EclPainterFactory.
Definition: EclPainter.h:31
Belle2::EclPainterCommon::ENERGY_SUM
@ ENERGY_SUM
Energy per event distribution.
Definition: EclPainterCommon.h:34
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::EclPainter::getNewRootObjectName
void getNewRootObjectName(char *buf, int size)
Make unique name for next root object.
Definition: EclPainter.cc:89
Belle2::MultilineWidget
Widget which contains the dynamic amount of TGLabel objects.
Definition: MultilineWidget.h:32
Belle2::EclPainterCommon::getType
Type getType()
Return subtype of ECLPainterCommon.
Definition: EclPainterCommon.cc:122
Belle2::EclPainterCommon::m_hist
TH1F * m_hist
Histogram for energy distribution.
Definition: EclPainterCommon.h:51
Belle2::EclPainterCommon::EclPainterCommon
EclPainterCommon(EclData *data, Type type)
Constructor for EclPainter subclass.
Definition: EclPainterCommon.cc:18
Belle2::EclData
This class contains data for ECLSimHit's and provides several relevant conversion functions for bette...
Definition: EclData.h:41
Belle2::EclPainterCommon::getInformation
virtual void getInformation(int px, int py, MultilineWidget *panel) override
Sets the information to be displayed in the provided MultilineWidget.
Definition: EclPainterCommon.cc:71
Belle2::EclPainterCommon::initHisto
void initHisto()
Initialize histogram.
Definition: EclPainterCommon.cc:62