Belle II Software  release-05-01-25
EclPainter1D.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 
11 //This module
12 #include <ecl/modules/eclDisplay/EclPainter1D.h>
13 
14 //Root
15 #include <TPad.h>
16 #include <TH1.h>
17 
18 //ECL
19 #include <ecl/modules/eclDisplay/geometry.h>
20 #include <ecl/utility/ECLChannelMapper.h>
21 #include <ecl/modules/eclDisplay/MultilineWidget.h>
22 
23 using namespace Belle2;
24 using namespace ECLDisplayUtility;
25 
27  EclPainter(data)
28 {
29  m_type = type;
30 
31  initHisto();
32  m_hist->GetXaxis()->CenterTitle();
33  m_hist->GetXaxis()->SetTitleOffset(1.1);
34  m_hist->GetYaxis()->SetTitleOffset(1.1);
35 
36  m_shaper = -1;
37  m_crate = -1;
38 }
39 
41 {
42  delete m_hist;
43 }
44 
46 {
47  switch (m_type) {
48  case CHANNEL:
49  return channel;
50  case SHAPER:
51  return 12 * getMapper()->getCrateID(channel) +
52  getMapper()->getShaperPosition(channel) - 12;
53  case CRATE:
54  return getMapper()->getCrateID(channel);
55  case PHI:
56  return getData()->getPhiId(channel);
57  case THETA:
58  return getData()->getThetaId(channel);
59  }
60 
61  return 0;
62 }
63 
65 {
66  switch (m_type) {
67  case CHANNEL:
68  return getData()->getCrystalCount();
69  case SHAPER:
70  return 52 * 12;
71  case CRATE:
72  return 52;
73  case PHI:
74  return 144;
75  case THETA:
76  return 69;
77  }
78 
79  return 1;
80 }
81 
83 {
84  char obj_name[255];
85  getNewRootObjectName(obj_name, 255);
86 
87  m_hist = new TH1F(obj_name, "title", getMaxX() + 1, 0, getMaxX() + 1);
88 }
89 
91 {
92  // TODO: These extra labels might not be necessary.
93  const char* name[3][5] = {
94  {
95  "Events per channel", "Events per shaper", "Events per collector",
96  "Events per phi_id", "Events per theta_id"
97  },
98 
99  {
100  "Energy per channel (MeV)", "Energy per shaper (MeV)", "Energy per collector (MeV)"
101  "Energy per phi_id (MeV)", "Energy per theta_id (MeV)"
102  },
103 
104  {
105  "Time per channel (ns)", "Time per shaper (ns)", "Time per collector (ns)",
106  "Time per phi_id (ns)", "Time per theta_id (ns)"
107  }
108  };
109  const char* xname[3] = {
110  "Channel id", "Shaper id", "Collector id"
111  };
112  const char* yname[3] = {
113  "Events", "Energy", "Time"
114  };
115 
116  TString title = TString(name[GetMode()][(int)m_type]) + " (" +
118  const char* xtitle = xname[(int)m_type];
119  const char* ytitle = yname[GetMode()];
120 
121  m_hist->SetTitle(title);
122  m_hist->SetXTitle(xtitle);
123  m_hist->SetYTitle(ytitle);
124 }
125 
126 void EclPainter1D::getInformation(int px, int py, MultilineWidget* panel)
127 {
128  EclPainter::getInformation(px, py, panel);
129 
130  char info[255];
131 
132  Float_t upx = gPad->AbsPixeltoX(px);
133  Float_t x = gPad->PadtoX(upx);
134  int binx = m_hist->GetXaxis()->FindBin(x) - 1;
135 
136  if (m_type == CHANNEL) {
137  sprintf(info, "channel_id = %d (%d)", binx,
138  getMapper()->getShaperChannel(binx));
139  panel->setLine(1, info);
140  sprintf(info, "shaper_id = %d", getMapper()->getShaperPosition(binx));
141  panel->setLine(2, info);
142  sprintf(info, "crate_id = %d", getMapper()->getCrateID(binx));
143  panel->setLine(3, info);
144  }
145  if (m_type == SHAPER) {
146  sprintf(info, "shaper_id = %d (%d)", binx, (binx - 1) % 12 + 1);
147  panel->setLine(1, info);
148  sprintf(info, "crate_id = %d", (binx - 1) / 12 + 1);
149  panel->setLine(2, info);
150  }
151  if (m_type == CRATE) {
152  sprintf(info, "crate_id = %d", binx);
153  panel->setLine(1, info);
154  }
155  if (m_type == PHI) {
156  sprintf(info, "phi_id = %d", binx);
157  panel->setLine(1, info);
158  }
159  if (m_type == THETA) {
160  sprintf(info, "theta_id = %d", binx);
161  panel->setLine(1, info);
162  }
163 }
164 
166 {
167  return m_type;
168 }
169 
171 {
172  Float_t upx = gPad->AbsPixeltoX(px);
173  Float_t x = gPad->PadtoX(upx);
174  int binx = m_hist->GetXaxis()->FindBin(x) - 1;
175 
176  Float_t upy = gPad->AbsPixeltoY(py);
177  Float_t y = gPad->PadtoY(upy);
178 
179  if (y < 0) return NULL;
180 
181  if (m_type == SHAPER) {
182  EclPainter1D* ret = new EclPainter1D(getData(), CHANNEL);
183  ret->setMapper(getMapper());
184  ret->setDisplayedSubsystem(getDisplayedSubsystem());
185  //ret->setXRange(binx * 16, binx * 16 + 15);
186  ret->setShaper((binx - 1) / 12 + 1, (binx - 1) % 12 + 1);
187  return ret;
188  }
189  if (m_type == CRATE) {
190  EclPainter1D* ret = new EclPainter1D(getData(), SHAPER);
191  ret->setMapper(getMapper());
192  ret->setDisplayedSubsystem(getDisplayedSubsystem());
193  ret->setXRange((binx - 1) * 12 + 1, (binx - 1) * 12 + 12);
194  //ret->setCrate(binx);
195  return ret;
196  }
197 
198  return NULL;
199 }
200 
201 void EclPainter1D::setXRange(int xmin, int xmax)
202 {
203  m_hist->GetXaxis()->SetRange(xmin + 1, xmax + 1);
204 }
205 
206 void EclPainter1D::setShaper(int crate, int shaper)
207 {
208  m_crate = crate;
209  m_shaper = shaper;
210 }
211 
212 void EclPainter1D::setCrate(int crate)
213 {
214  m_crate = crate;
215 }
216 
218 {
219  setTitles();
220 
221  EclData* data = getData();
222 
223  const int* ev_counts = data->getEventCounts();
224  const float* energy_sums = data->getEnergySums();
225 
226  m_hist->Reset();
227  for (int i = 1; i <= getData()->getCrystalCount(); i++) {
228  if (!data->isCrystalInSubsystem(i, getDisplayedSubsystem())) continue;
229 
230  // If filter is set, display only specified channels.
231  if (m_type == CHANNEL && m_crate > 0 && m_shaper > 0) {
232  if (m_crate != getMapper()->getCrateID(i) ||
234  continue;
235  }
236  }
237  if (m_type == SHAPER && m_crate > 0) {
238  if (m_crate != getMapper()->getCrateID(i)) {
239  continue;
240  }
241  }
242  //
243  int id = channelToSegId(i);
244  if (GetMode())
245  m_hist->Fill(id, energy_sums[i]);
246  else
247  m_hist->Fill(id, ev_counts[i]);
248  }
249 
250  m_hist->Draw("HIST");
251 }
Belle2::EclPainter1D::CHANNEL
@ CHANNEL
Events/energy per channel.
Definition: EclPainter1D.h:30
Belle2::EclPainter1D::SHAPER
@ SHAPER
Events/energy per ShaperDSP.
Definition: EclPainter1D.h:31
Belle2::MultilineWidget::setLine
void setLine(int line_id, const char *text)
Set content of the specified line to 'text'.
Definition: MultilineWidget.cc:68
Belle2::EclPainter1D::getMaxX
int getMaxX()
Returns number of X bins.
Definition: EclPainter1D.cc:64
Belle2::EclPainter1D::setTitles
void setTitles()
Update titles of the histogram.
Definition: EclPainter1D.cc:90
Belle2::EclPainter1D::~EclPainter1D
virtual ~EclPainter1D()
Destructor for EclPainter subclass.
Definition: EclPainter1D.cc:40
Belle2::EclPainter1D::CRATE
@ CRATE
Events/energy per crate/ECLCollector.
Definition: EclPainter1D.h:32
Belle2::EclPainter1D::initHisto
void initHisto()
Initialize histogram.
Definition: EclPainter1D.cc:82
Belle2::ECL::ECLChannelMapper::getCrateID
int getCrateID(int iCOPPERNode, int iFINESSE)
get crate number by given COPPER node number and FINESSE number
Definition: ECLChannelMapper.cc:211
Belle2::ECL::ECLChannelMapper::getShaperPosition
int getShaperPosition(int cellID)
get position of the shaper in the crate by given CellId
Definition: ECLChannelMapper.cc:289
Belle2::EclPainter1D::m_hist
TH1F * m_hist
Displayed histogram.
Definition: EclPainter1D.h:50
Belle2::EclPainter::getData
EclData * getData()
Return currently displayed EclData.
Definition: EclPainter.h:48
Belle2::EclPainter1D::Draw
void Draw() override
Redraw the canvas.
Definition: EclPainter1D.cc:217
Belle2::EclPainter::getDisplayedSubsystem
EclData::EclSubsystem getDisplayedSubsystem()
Get currently displayed ECL subsystem.
Definition: EclPainter.cc:52
Belle2::EclPainter1D::m_shaper
int m_shaper
ID of currently selected shaper.
Definition: EclPainter1D.h:53
Belle2::EclPainter1D::setXRange
void setXRange(int xmin, int xmax) override
Set XRange for histogram.
Definition: EclPainter1D.cc:201
Belle2::EclData::getPhiId
int getPhiId(int ch)
ECL CellId -> phi_id.
Definition: EclData.cc:278
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::EclPainter1D::getType
Type getType()
Return subtype of ECLPainter1D.
Definition: EclPainter1D.cc:165
Belle2::EclPainter
Painter for EclData, parent class, created with EclPainterFactory.
Definition: EclPainter.h:31
Belle2::EclPainter1D::m_type
Type m_type
Display subtypes of this class.
Definition: EclPainter1D.h:48
Belle2::EclPainter1D::handleClick
virtual EclPainter * handleClick(int px, int py) override
Creates sub-histogram for crates and shapers.
Definition: EclPainter1D.cc:170
Belle2::EclPainter1D::getInformation
virtual void getInformation(int px, int py, MultilineWidget *panel) override
Sets the information to be displayed in the provided MultilineWidget.
Definition: EclPainter1D.cc:126
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::EclData::getCrystalCount
static int getCrystalCount()
Get number of crystals in ECL.
Definition: EclData.cc:148
Belle2::EclPainter::getSubsystemTitle
TString getSubsystemTitle(EclData::EclSubsystem subsys)
Return title of ECL subsystem to use in painter.
Definition: EclPainter.cc:57
Belle2::EclPainter::getMapper
ECL::ECLChannelMapper * getMapper()
Return currently set ECLChannelMapper.
Definition: EclPainter.cc:42
Belle2::EclData::getThetaId
int getThetaId(int ch)
ECL CellId -> theta_id.
Definition: EclData.cc:288
Belle2::EclPainter1D::Type
Type
Subtype of histogram to draw.
Definition: EclPainter1D.h:29
Belle2::EclData
This class contains data for ECLSimHit's and provides several relevant conversion functions for bette...
Definition: EclData.h:41
Belle2::EclPainter1D::channelToSegId
int channelToSegId(int channel)
Convert channel id to X bin number.
Definition: EclPainter1D.cc:45
Belle2::EclPainter1D::m_crate
int m_crate
ID of currently selected crate.
Definition: EclPainter1D.h:55
Belle2::EclPainter1D::EclPainter1D
EclPainter1D(EclData *data, Type type)
Constructor for EclPainter subclass.
Definition: EclPainter1D.cc:26
Belle2::EclPainter1D
Painter for EclData, 1D histograms.
Definition: EclPainter1D.h:26
Belle2::EclPainter1D::setCrate
void setCrate(int crate)
Show data only from specific crate.
Definition: EclPainter1D.cc:212
Belle2::EclPainter1D::setShaper
void setShaper(int crate, int shaper)
Show data only from specific shaper.
Definition: EclPainter1D.cc:206