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