Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TEveCaloLegoGL.cxx
Go to the documentation of this file.
1// @(#)root/eve:$Id$
2// Author: Alja Mrak-Tadel 2007
3
4/*************************************************************************
5 * Copyright (C) 1995-2007, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12#include "TAxis.h"
13#include "THLimitsFinder.h"
14
15#include "TGLViewer.h"
16#include "TGLIncludes.h"
17#include "TGLPhysicalShape.h"
18#include "TGLRnrCtx.h"
19#include "TGLSelectRecord.h"
20#include "TGLScene.h"
21#include "TGLCamera.h"
22#include "TGLUtil.h"
23#include "TColor.h"
24#include "TROOT.h"
25
26
27#include "TEveCaloLegoGL.h"
28#include "TEveCalo.h"
29#include "TEveManager.h"
30#include "TEveRGBAPalette.h"
31
32#include <algorithm>
33
34/** \class TEveCaloLegoGL
35\ingroup TEve
36OpenGL renderer class for TEveCaloLego.
37*/
38
40
41////////////////////////////////////////////////////////////////////////////////
42/// Constructor.
43
45 TGLObject(),
46
47 fGridColor(-1),
48 fFontColor(-1),
49
50 fEtaAxis(nullptr),
51 fPhiAxis(nullptr),
52 fZAxis(nullptr),
53 fM(nullptr),
54 fDLCacheOK(kFALSE),
55 fMaxVal(0),
56 fValToPixel(0),
57 fCurrentPixelsPerBin(0),
58 fCells3D(kTRUE),
59 fBinStep(-1)
60{
62
63 fEtaAxis = new TAxis();
64 fPhiAxis = new TAxis();
65 fZAxis = new TAxis();
66
68}
69
70////////////////////////////////////////////////////////////////////////////////
71/// Destructor.
72
74{
76
77 delete fEtaAxis;
78 delete fPhiAxis;
79 delete fZAxis;
80}
81
82////////////////////////////////////////////////////////////////////////////////
83/// Set model object.
84
86{
87 fM = SetModelDynCast<TEveCaloLego>(obj);
88 return kTRUE;
89}
90
91////////////////////////////////////////////////////////////////////////////////
92/// Set bounding box.
93
95{
97}
98
99////////////////////////////////////////////////////////////////////////////////
100/// Drop all display-list definitions.
101
103{
105 for (SliceDLMap_i i = fDLMap.begin(); i != fDLMap.end(); ++i)
106 i->second = 0;
107
109}
110
111////////////////////////////////////////////////////////////////////////////////
112/// Unregister all display-lists.
113
115{
116 // all lego cells
118 if (! fDLMap.empty()) {
119 for (SliceDLMap_i i = fDLMap.begin(); i != fDLMap.end(); ++i) {
120 if (i->second) {
121 PurgeDLRange(i->second, 1);
122 i->second = 0;
123 }
124 }
125 }
127}
128
129////////////////////////////////////////////////////////////////////////////////
130/// Draw an axis-aligned box using quads.
131/// ~~~ {.cpp}
132/// z
133/// |
134/// |
135/// |________y
136/// / 6-------7
137/// / /| /|
138/// x 5-------4 |
139/// | 2-----|-3
140/// |/ |/
141/// 1-------0
142/// ~~~
143
145 Float_t xw, Float_t yw, Float_t h) const
146{
147 Float_t x2 = x1 + xw;
148 Float_t y2 = y1 + yw;
149 Float_t z2 = z1 + h;
150
151 if (x1 < fM->GetEtaMin()) x1 = fM->GetEtaMin();
152 if (x2 > fM->GetEtaMax()) x2 = fM->GetEtaMax();
153
154 if (y1 < fM->GetPhiMin()) y1 = fM->GetPhiMin();
155 if (y2 > fM->GetPhiMax()) y2 = fM->GetPhiMax();
156
157 glBegin(GL_QUADS);
158 {
159 // bottom 0123
160 glNormal3f(0, 0, -1);
161 glVertex3f(x2, y2, z1);
162 glVertex3f(x2, y1, z1);
163 glVertex3f(x1, y1, z1);
164 glVertex3f(x1, y2, z1);
165 // top 4765
166 glNormal3f(0, 0, 1);
167 glVertex3f(x2, y2, z2);
168 glVertex3f(x1, y2, z2);
169 glVertex3f(x1, y1, z2);
170 glVertex3f(x2, y1, z2);
171
172 // back 0451
173 glNormal3f(1, 0, 0);
174 glVertex3f(x2, y2, z1);
175 glVertex3f(x2, y2, z2);
176 glVertex3f(x2, y1, z2);
177 glVertex3f(x2, y1, z1);
178 // front 3267
179 glNormal3f(-1, 0, 0);
180 glVertex3f(x1, y2, z1);
181 glVertex3f(x1, y1, z1);
182 glVertex3f(x1, y1, z2);
183 glVertex3f(x1, y2, z2);
184
185 // left 0374
186 glNormal3f(0, 1, 0);
187 glVertex3f(x2, y2, z1);
188 glVertex3f(x1, y2, z1);
189 glVertex3f(x1, y2, z2);
190 glVertex3f(x2, y2, z2);
191 // right 1562
192 glNormal3f(0, -1, 0);
193 glVertex3f(x2, y1, z1);
194 glVertex3f(x2, y1, z2);
195 glVertex3f(x1, y1, z2);
196 glVertex3f(x1, y1, z1);
197 }
198 glEnd();
199}
200
201////////////////////////////////////////////////////////////////////////////////
202/// Create display-list that draws histogram bars for non-rebinned data.
203/// It is used for filled and outline passes.
204
206{
208 Int_t prevTower = 0;
209 Float_t offset = 0;
210
211 // ids in eta phi rng
212 Int_t nSlices = fM->fData->GetNSlices();
213 for (Int_t s = 0; s < nSlices; ++s)
214 {
215 if (dlMap.empty() || dlMap[s] == 0)
216 dlMap[s] = glGenLists(1);
217
218 glNewList(dlMap[s], GL_COMPILE);
219
220 for (UInt_t i = 0; i < cellList.size(); ++i)
221 {
222 if (cellList[i].fSlice > s) continue;
223 if (cellList[i].fTower != prevTower) {
224 offset = 0;
225 prevTower = cellList[i].fTower;
226 }
227
228 fM->fData->GetCellData(cellList[i], cellData);
229 if (s == cellList[i].fSlice)
230 {
231 if (selection) glLoadName(i);
232
233 WrapTwoPi(cellData.fPhiMin, cellData.fPhiMax);
234 MakeQuad(cellData.EtaMin(), cellData.PhiMin(), offset,
235 cellData.EtaDelta(), cellData.PhiDelta(), cellData.Value(fM->fPlotEt));
236 }
237 offset += cellData.Value(fM->fPlotEt);
238 }
239 glEndList();
240 }
241}
242
243////////////////////////////////////////////////////////////////////////////////
244/// Create display-list that draws histogram bars for rebinned data.
245/// It is used for filled and outline passes.
246
248{
249 Int_t nSlices = fM->fData->GetNSlices();
250 Float_t *vals;
252 Float_t y0, y1;
253
254 for (Int_t s = 0; s < nSlices; ++s)
255 {
256 if (dlMap.empty() || dlMap[s] == 0)
257 dlMap[s] = glGenLists(1);
258
259 glNewList(dlMap[s], GL_COMPILE);
260
261 for (Int_t i = 1; i <= fEtaAxis->GetNbins(); ++i)
262 {
263 for (Int_t j = 1; j <= fPhiAxis->GetNbins(); ++j)
264 {
265 const Int_t bin = (i)+(j)*(fEtaAxis->GetNbins()+2);
266
267 if (rebinData.fBinData[bin] !=-1)
268 {
269 vals = rebinData.GetSliceVals(bin);
270 offset =0;
271 for (Int_t t = 0; t < s; ++t)
272 offset += vals[t];
273
274 y0 = fPhiAxis->GetBinLowEdge(j);
276 WrapTwoPi(y0, y1);
277
278 if (selection) glLoadName(bin);
279
281 fEtaAxis->GetBinWidth(i), y1-y0, vals[s]);
282 }
283 }
284 }
285 glEndList();
286 }
287}
288
289////////////////////////////////////////////////////////////////////////////////
290/// Set the axis 3D title position.
291
293{
294 const GLdouble *pm = rnrCtx.RefCamera().RefLastNoPickProjM().CArr();
295 GLdouble mm[16];
296 GLint vp[4];
297 glGetDoublev(GL_MODELVIEW_MATRIX, mm);
298 glGetIntegerv(GL_VIEWPORT, vp);
299 GLdouble projX[4], projY[4], projZ[4];
300
301 GLdouble cornerX[4];
302 GLdouble cornerY[4];
303 cornerX[0] = x0; cornerY[0] = y0;
304 cornerX[1] = x1; cornerY[1] = y0;
305 cornerX[2] = x1; cornerY[2] = y1;
306 cornerX[3] = x0; cornerY[3] = y1;
307
308 gluProject(cornerX[0], cornerY[0], 0, mm, pm, vp, &projX[0], &projY[0], &projZ[0]);
309 gluProject(cornerX[1], cornerY[1], 0, mm, pm, vp, &projX[1], &projY[1], &projZ[1]);
310 gluProject(cornerX[2], cornerY[2], 0, mm, pm, vp, &projX[2], &projY[2], &projZ[2]);
311 gluProject(cornerX[3], cornerY[3], 0, mm, pm, vp, &projX[3], &projY[3], &projZ[3]);
312
313
314 // Z axis location (left most corner)
315 //
316 Int_t idxLeft = 0;
317 Float_t xt = projX[0];
318 for (Int_t i = 1; i < 4; ++i) {
319 if (projX[i] < xt) {
320 xt = projX[i];
321 idxLeft = i;
322 }
323 }
324 fZAxisTitlePos.Set(cornerX[idxLeft], cornerY[idxLeft], 1.05 * fMaxVal);
325
326
327 // XY axis location (closest to eye) first
328 //
329 Float_t zt = 1.f;
330 Float_t zMin = 0.f;
331 Int_t idxFront = 0;
332 for (Int_t i = 0; i < 4; ++i) {
333 if (projZ[i] < zt) {
334 zt = projZ[i];
335 idxFront = i;
336 }
337 if (projZ[i] > zMin) zMin = projZ[i];
338 }
339
340
341 Int_t xyIdx = idxFront;
342 if (zMin - zt < 1e-2) xyIdx = 0; // avoid flipping in front view
343
344
345 switch (xyIdx) {
346 case 0:
348 fXAxisTitlePos.fY = y0;
349 fYAxisTitlePos.fX = x0;
351 break;
352 case 1:
353 fXAxisTitlePos.fX = x0;
354 fXAxisTitlePos.fY = y0;
357 break;
358 case 2:
359 fXAxisTitlePos.fX = x0;
362 fYAxisTitlePos.fY = y0;
363 break;
364 case 3:
367 fYAxisTitlePos.fX = x0;
368 fYAxisTitlePos.fY = y0;
369 break;
370 }
371
372 // move title 5% over the axis length
373 Float_t off = 0.05;
374 Float_t tOffX = (x1-x0) * off; if (fYAxisTitlePos.fX > x0) tOffX = -tOffX;
375 Float_t tOffY = (y1-y0) * off; if (fXAxisTitlePos.fY > y0) tOffY = -tOffY;
376 fXAxisTitlePos.fX += tOffX;
377 fYAxisTitlePos.fY += tOffY;
378
379
380 // frame box
381 //
382 if (fM->fBoxMode)
383 {
384 // get corner closest to eye excluding left corner
385 Double_t zm = 1.f;
386 Int_t idxDepthT = 0;
387 for (Int_t i = 0; i < 4; ++i)
388 {
389 if (projZ[i] < zm && projZ[i] >= zt && i != idxFront )
390 {
391 zm = projZ[i];
392 idxDepthT = i;
393 }
394 }
395 if (idxFront == idxLeft) idxFront =idxDepthT;
396
397 switch (idxFront)
398 {
399 case 0:
400 fBackPlaneXConst[0].Set(x1, y0, 0); fBackPlaneXConst[1].Set(x1, y1, 0);
401 fBackPlaneYConst[0].Set(x0, y1, 0); fBackPlaneYConst[1].Set(x1, y1, 0);
402 break;
403 case 1:
404 fBackPlaneXConst[0].Set(x0, y0, 0); fBackPlaneXConst[1].Set(x0, y1, 0);
405 fBackPlaneYConst[0].Set(x0, y1, 0); fBackPlaneYConst[1].Set(x1, y1, 0);
406 break;
407 case 2:
408 fBackPlaneXConst[0].Set(x0, y0, 0); fBackPlaneXConst[1].Set(x0, y1, 0);
409 fBackPlaneYConst[0].Set(x0, y0, 0); fBackPlaneYConst[1].Set(x1, y0, 0);
410 break;
411 case 3:
412 fBackPlaneXConst[0].Set(x1, y0, 0); fBackPlaneXConst[1].Set(x1, y1, 0);
413 fBackPlaneYConst[0].Set(x0, y0, 0); fBackPlaneYConst[1].Set(x1, y0, 0);
414 break;
415 }
416 }
417}
418
419////////////////////////////////////////////////////////////////////////////////
420/// Draw z-axis and z-box at the appropriate grid corner-point including
421/// tick-marks and labels.
422
424{
425 // set font size first depending on size of projected axis
426
427 TGLMatrix mm;
428 GLdouble pm[16];
429 glGetDoublev(GL_MODELVIEW_MATRIX, mm.Arr());
430 glGetDoublev(GL_PROJECTION_MATRIX, pm);
431 Int_t* vp = rnrCtx.RefCamera().RefViewport().CArr();
432 GLdouble dn[3];
433 GLdouble up[3];
434 gluProject(fXAxisTitlePos.fX, fXAxisTitlePos.fY, fXAxisTitlePos.fZ, mm.Arr(), pm, vp, &up[0], &up[1], &up[2]);
435 gluProject(fYAxisTitlePos.fX, fYAxisTitlePos.fY, fYAxisTitlePos.fZ, mm.Arr(), pm, vp, &dn[0], &dn[1], &dn[2]);
436 Float_t len = TMath::Sqrt((up[0] - dn[0]) * (up[0] - dn[0])
437 + (up[1] - dn[1]) * (up[1] - dn[1])
438 + (up[2] - dn[2]) * (up[2] - dn[2]));
439 len = TMath::Min(len, rnrCtx.RefCamera().RefViewport().Diagonal()*0.7f);
440 len /= TMath::Sqrt2();
441
443 fAxisPainter.RefTMOff(0) = rnrCtx.RefCamera().ViewportDeltaToWorld(worldRef, -len, 0, &mm);
446
447 Float_t tickLength = TMath::Max(fM->GetData()->GetEtaBins()->GetTickLength(), 0.02f);
448 Float_t labelOffset = TMath::Max(fM->GetData()->GetEtaBins()->GetLabelOffset(), 0.02f);
449
450 // Z axis
451 //
452 if (fM->fData->Empty() == kFALSE)
453 {
454 Int_t ondiv;
455 Double_t omin=0, omax=0, bw1;
456 THLimitsFinder::Optimize(0, fMaxVal, fM->fNZSteps, omin, omax, ondiv, bw1);
459 // check z axis title does not overalp with label
461 fZAxisTitlePos.fZ = omax + zto.Z();
462
463
467 fZAxis->SetNdivisions(fM->fNZSteps*100 + 10);
469 fZAxis->SetTitle(fM->GetPlotEt() ? "Et[GeV]" : "E[GeV]");
470
472 fAxisPainter.RefDir().Set(0., 0., 1.);
474 glPushMatrix();
475 glTranslatef(fZAxisTitlePos.fX, fZAxisTitlePos.fY, 0);
476
477 // tickmark vector = 10 pixels left
479 fZAxis->SetLabelOffset(labelOffset);
480 fZAxis->SetTickLength(tickLength);
482 glPopMatrix();
483
484 // draw box frame
485 //
486 if (fM->fBoxMode) {
487
488 glPushAttrib(GL_ENABLE_BIT | GL_LINE_BIT);
489
490 // box verticals
492 glBegin(GL_LINES);
494
495 glVertex3f(fBackPlaneXConst[0].fX ,fBackPlaneXConst[0].fY ,0);
496 glVertex3f(fBackPlaneXConst[0].fX ,fBackPlaneXConst[0].fY ,fMaxVal);
497 glVertex3f(fBackPlaneXConst[1].fX ,fBackPlaneXConst[1].fY ,0);
498 glVertex3f(fBackPlaneXConst[1].fX ,fBackPlaneXConst[1].fY ,fMaxVal);
499
500
501 glVertex3f(fBackPlaneYConst[0].fX ,fBackPlaneYConst[0].fY ,0);
502 glVertex3f(fBackPlaneYConst[0].fX ,fBackPlaneYConst[0].fY ,fMaxVal);
503 glVertex3f(fBackPlaneYConst[1].fX ,fBackPlaneYConst[1].fY ,0);
504 glVertex3f(fBackPlaneYConst[1].fX ,fBackPlaneYConst[1].fY ,fMaxVal);
505
506 // box top
507 glVertex3f(fBackPlaneXConst[0].fX ,fBackPlaneXConst[0].fY ,fMaxVal);
508 glVertex3f(fBackPlaneXConst[1].fX ,fBackPlaneXConst[1].fY ,fMaxVal);
509 glVertex3f(fBackPlaneYConst[0].fX ,fBackPlaneYConst[0].fY ,fMaxVal);
510 glVertex3f(fBackPlaneYConst[1].fX ,fBackPlaneYConst[1].fY ,fMaxVal);
511
512 glEnd();
513
514 // box horizontals stippled
515 glEnable(GL_LINE_STIPPLE);
516 glLineStipple(1, 0x5555);
517 glBegin(GL_LINES);
518 Float_t hz = bw1;
519 for (Int_t i = 1; i <= ondiv; ++i, hz += bw1) {
520 glVertex3f(fBackPlaneXConst[0].fX ,fBackPlaneXConst[0].fY ,hz);
521 glVertex3f(fBackPlaneXConst[1].fX ,fBackPlaneXConst[1].fY ,hz);
522 glVertex3f(fBackPlaneYConst[0].fX ,fBackPlaneYConst[0].fY ,hz);
523 glVertex3f(fBackPlaneYConst[1].fX ,fBackPlaneYConst[1].fY ,hz);
524 }
525 glEnd();
526
527 glPopAttrib();
528 }
529 }
530
531 // XY Axis
532 //
533
534 Float_t yOff = fM->GetPhiRng();
535 if (fXAxisTitlePos.fY < fM->GetPhiMax()) yOff = -yOff;
536
537 Float_t xOff = fM->GetEtaRng();
538 if (fYAxisTitlePos.fX < fM->GetEtaMax()) xOff = -xOff;
539
540 TAxis ax;
545 ax.SetLabelOffset(labelOffset);
546 ax.SetTickLength(tickLength);
548 fAxisPainter.RefTMOff(1).Set(0, 0, -fMaxVal);
550
551 // eta
552 glPushMatrix();
553 fAxisPainter.RefDir().Set(1, 0, 0);
554 fAxisPainter.RefTMOff(0).Set(0, yOff, 0);
555 glTranslatef(0, fXAxisTitlePos.fY, 0);
556
558 ax.SetLimits(fM->GetEtaMin(), fM->GetEtaMax());
561 fAxisPainter.PaintAxis(rnrCtx, &ax);
562 glPopMatrix();
563
564 // phi
565 fAxisPainter.RefDir().Set(0, 1, 0);
566 fAxisPainter.RefTMOff(0).Set(xOff, 0, 0);
568 ax.SetLimits(fM->GetPhiMin(), fM->GetPhiMax());
570 glPushMatrix();
571 glTranslatef(fYAxisTitlePos.fX, 0, 0);
573 fAxisPainter.PaintAxis(rnrCtx, &ax);
574 glPopMatrix();
575
576} // DrawAxis3D
577
578////////////////////////////////////////////////////////////////////////////////
579/// Get scale for matrix.
580
582{
583 Double_t em, eM, pm, pM;
584 fM->fData->GetEtaLimits(em, eM);
585 fM->fData->GetPhiLimits(pm, pM);
586 Double_t unit = ((eM - em) < (pM - pm)) ? (eM - em) : (pM - pm);
587 sx = (eM - em) / (fM->GetEtaRng() * unit);
588 sy = (pM - pm) / (fM->GetPhiRng() * unit);
589
590 sz = 1;
591 if (fM->fScaleAbs)
592 {
593 sz = fM->GetMaxTowerH() / fM->fMaxValAbs;
594 }
595 else if (!fM->fData->Empty())
596 {
597 sz = fM->GetMaxTowerH() / fMaxVal;
598 }
599}
600
601////////////////////////////////////////////////////////////////////////////////
602/// Draw XY axis.
603
605{
606 if (fM->GetData()->Empty())
608
609 TGLCamera& cam = rnrCtx.RefCamera();
610
611 TAxis ax;
619
620 // set fonts
622
623 // get projected length of diagonal to determine
624 TGLMatrix mm;
625 GLdouble pm[16];
626 GLint vp[4];
627 glGetDoublev(GL_MODELVIEW_MATRIX, mm.Arr());
628 glGetDoublev(GL_PROJECTION_MATRIX, pm);
629 glGetIntegerv(GL_VIEWPORT, vp);
630
631 GLdouble dn[3];
632 GLdouble up[3];
633 gluProject(fM->GetEtaMin(), fM->GetPhiMin(), 0, mm.Arr(), pm, vp, &dn[0], &dn[1], &dn[2]);
634 gluProject(fM->GetEtaMax(), fM->GetPhiMax(), 0, mm.Arr(), pm, vp, &up[0], &up[1], &up[2]);
635 Double_t len = TMath::Sqrt((up[0] - dn[0]) * (up[0] - dn[0])
636 + (up[1] - dn[1]) * (up[1] - dn[1])
637 + (up[2] - dn[2]) * (up[2] - dn[2]));
638
639 // lock upper limit to of relative font size relative to viewport diagonal
640 Double_t vpLimit = cam.RefViewport().Diagonal()*0.5/TMath::Sqrt2();
641 len = TMath::Min(len, vpLimit);
642
643 // eta
647 ax.SetLimits(fM->GetEtaMin(), fM->GetEtaMax());
650 fAxisPainter.RefDir().Set(1, 0, 0);
651
653 fAxisPainter.RefTMOff(0).Set(0, -TMath::Min(fM->GetPhiRng(), tmOffFrustX), 0);
655
656 glPushMatrix();
657 glTranslatef(0, fM->GetPhiMin(), 0);
658 fAxisPainter.PaintAxis(rnrCtx, &ax);
659 glPopMatrix();
660
661 // phi
663 ax.SetLimits(fM->GetPhiMin(), fM->GetPhiMax());
666 fAxisPainter.RefDir().Set(0, 1, 0);
668 fAxisPainter.RefTMOff(0).Set(-TMath::Min(fM->GetEtaRng(), tmOffFrustY), 0, 0);
670
671 glPushMatrix();
672 glTranslatef(fM->GetEtaMin(), 0, 0);
673 fAxisPainter.PaintAxis(rnrCtx, &ax);
674 glPopMatrix();
675
677}
678
679////////////////////////////////////////////////////////////////////////////////
680/// Calculate view-dependent grid density.
681
683{
684 TGLCamera &camera = rnrCtx.RefCamera();
689 Float_t frustD = TMath::Hypot(r-l, t-b);
690
691 GLint vp[4]; glGetIntegerv(GL_VIEWPORT, vp);
692 Float_t viewportD = TMath::Sqrt((vp[1] - vp[0]) * (vp[1] - vp[0]) + (vp[3] - vp[1]) * (vp[3] - vp[1]));
693 Float_t deltaToViewport = viewportD/frustD;
694
695 // average bin width
696 GLdouble em, eM, pm, pM;
697 fM->GetData()->GetEtaLimits(pm, pM);
698 fM->GetData()->GetPhiLimits(em, eM);
703
704 Float_t averageBinWidth = TMath::Hypot(eM - em, pM - pm)/TMath::Sqrt((i0 - i1) * (i0 - i1) + (j0 - j1) * (j0 - j1));
705 Float_t ppb = deltaToViewport*averageBinWidth;
706
707 Int_t ngroup = 1;
708 if (fM->fAutoRebin && fM->fPixelsPerBin > ppb)
709 {
710 // limit rebin realtive to number of axis bins
711 Int_t maxGroup = TMath::Min(fM->fData->GetEtaBins()->GetNbins(), fM->fData->GetPhiBins()->GetNbins())/4;
712 if (maxGroup > 1) {
713 ngroup = TMath::Nint(fM->fPixelsPerBin*0.5/ppb); // symetrical rebin factor 2
714 if (ngroup > maxGroup) ngroup = maxGroup;
715 }
716 }
718
719 return ngroup;
720}
721
722////////////////////////////////////////////////////////////////////////////////
723/// Rebin eta, phi axis.
724
725void TEveCaloLegoGL::RebinAxis(TAxis *orig, TAxis *curr) const
726{
727 Double_t center = 0.5 * (orig->GetXmin() + orig->GetXmax());
728 Int_t idx0 = orig->FindBin(center);
729 Double_t bc = orig->GetBinCenter(idx0);
730 if (bc > center) --idx0;
731
732 Int_t nbR = TMath::FloorNint(idx0/fBinStep) + TMath::FloorNint((orig->GetNbins() - idx0)/fBinStep);
733 Int_t off = idx0 - TMath::FloorNint(idx0/fBinStep)*fBinStep;
734 std::vector<Double_t> bins(nbR + 1);
735 for (Int_t i = 0; i <= nbR; ++i)
736 {
737 bins[i] = orig->GetBinUpEdge(off + i*fBinStep);
738 }
739 curr->Set(nbR, &bins[0]);
740}
741
742////////////////////////////////////////////////////////////////////////////////
743/// Draw basic histogram components: x-y grid
744
746{
747 Float_t eta0 = fM->fEtaMin;
748 Float_t eta1 = fM->fEtaMax;
749 Float_t phi0 = fM->GetPhiMin();
750 Float_t phi1 = fM->GetPhiMax();
751
752 // XY grid
753 //
756 glBegin(GL_LINES);
757 glVertex2f(eta0, phi0);
758 glVertex2f(eta0, phi1);
759 glVertex2f(eta1, phi0);
760 glVertex2f(eta1, phi1);
761
762 glVertex2f(eta0, phi0);
763 glVertex2f(eta1, phi0);
764 glVertex2f(eta0, phi1);
765 glVertex2f(eta1, phi1);
766
767 // eta grid
768 Float_t val;
769 Int_t neb = fEtaAxis->GetNbins();
770 for (Int_t i = 0; i<= neb; i++)
771 {
772 val = fEtaAxis->GetBinUpEdge(i);
773 if (val > eta0 && val < eta1 )
774 {
775 glVertex2f(val, phi0);
776 glVertex2f(val, phi1);
777 }
778 }
779
780 // phi grid
781 Int_t npb = fPhiAxis->GetNbins();
782 for (Int_t i = 1; i <= npb; i++) {
783 val = fPhiAxis->GetBinUpEdge(i);
784 if (val > phi0 && val < phi1)
785 {
786 glVertex2f(eta0, val);
787 glVertex2f(eta1, val);
788 }
789 }
790
791 glEnd();
792
793 // XYZ axes
794 //
795 glPushAttrib(GL_ENABLE_BIT | GL_LINE_BIT | GL_POLYGON_BIT);
797 if (fCells3D)
798 {
799 SetAxis3DTitlePos(rnrCtx, eta0, eta1, phi0, phi1);
800 DrawAxis3D(rnrCtx);
801 }
802 else
803 {
804 DrawAxis2D(rnrCtx);
805 }
806 glPopAttrib();
807}
808
809////////////////////////////////////////////////////////////////////////////////
810/// Render the calo lego-plot with OpenGL.
811
813{
814 // quads
815 {
816 for (SliceDLMap_i i = fDLMap.begin(); i != fDLMap.end(); ++i) {
818 glLoadName(i->first);
819 glPushName(0);
820 glCallList(i->second);
821 glPopName();
822 }
823 }
824 // outlines
825 {
826 if (rnrCtx.SceneStyle() == TGLRnrCtx::kFill) {
827 glPolygonMode(GL_FRONT_AND_BACK, GL_LINE);
828 glDisable(GL_POLYGON_OFFSET_FILL);
830 for (SliceDLMap_i i = fDLMap.begin(); i != fDLMap.end(); ++i)
831 glCallList(i->second);
832 }
833 }
834}
835
836////////////////////////////////////////////////////////////////////////////////
837/// Prepare cells 2D data non-rebinned for drawing.
838
840{
841 Int_t max_energy_slice, cellID=0;
842 Float_t sum, max_energy;
843
844 TEveCaloData::vCellId_t::iterator currentCell = cellList.begin();
845 TEveCaloData::vCellId_t::iterator nextCell = currentCell;
846 ++nextCell;
847
848 while (true)
849 {
850 TEveCaloData::CellData_t currentCellData;
851 TEveCaloData::CellData_t nextCellData;
852
853 fM->fData->GetCellData(*currentCell, currentCellData);
854 sum = max_energy = currentCellData.Value(fM->fPlotEt);
855 max_energy_slice = currentCell->fSlice;
856 while (nextCell != cellList.end() && currentCell->fTower == nextCell->fTower)
857 {
858 fM->fData->GetCellData(*nextCell, nextCellData);
859 Float_t energy = nextCellData.Value(fM->fPlotEt);
860 sum += energy;
861 if (energy > max_energy)
862 {
863 max_energy = energy;
864 max_energy_slice = nextCell->fSlice;
865 }
866 ++nextCell;
867 ++cellID;
868 }
869
870 WrapTwoPi(currentCellData.fPhiMin, currentCellData.fPhiMax);
871 cells2D.push_back(Cell2D_t(cellID, sum, max_energy_slice));
872 cells2D.back().SetGeom(currentCellData.fEtaMin, currentCellData.fEtaMax,
873 currentCellData.fPhiMin, currentCellData.fPhiMax);
874
875 if (nextCell == cellList.end())
876 break;
877
878 currentCell = nextCell;
879 ++nextCell;
880 ++cellID;
881 }
882}
883
884////////////////////////////////////////////////////////////////////////////////
885/// Prepare cells 2D rebinned data for drawing.
886
888{
889 const Int_t nEta = fEtaAxis->GetNbins();
890 const Int_t nPhi = fPhiAxis->GetNbins();
891 std::vector<Float_t> vec;
892 vec.assign((nEta + 2)*(nPhi + 2), 0.f);
893 std::vector<Float_t> max_e;
894 std::vector<Int_t> max_e_slice;
895 max_e.assign((nEta + 2) * (nPhi + 2), 0.f);
896 max_e_slice.assign((nEta + 2) * (nPhi + 2), -1);
897
898 for (UInt_t bin = 0; bin < rebinData.fBinData.size(); ++bin) {
899 Float_t ssum = 0;
900 if (rebinData.fBinData[bin] != -1) {
901 Float_t *val = rebinData.GetSliceVals(bin);
902 for (Int_t s = 0; s < rebinData.fNSlices; ++s) {
903 ssum += val[s];
904 if (val[s] > max_e[bin]) {
905 max_e[bin] = val[s];
906 max_e_slice[bin] = s;
907 }
908 }
909 }
910 vec[bin] = ssum;
911 }
912
913 // smallest threshold
914 Float_t threshold = fM->GetDataSliceThreshold(0);
915 for (Int_t s = 1; s < fM->fData->GetNSlices(); ++s) {
916 if (threshold > fM->GetDataSliceThreshold(s))
917 threshold = fM->GetDataSliceThreshold(s);
918 }
919
920 // write cells
921 for (Int_t i = 1; i <= fEtaAxis->GetNbins(); ++i) {
922 for (Int_t j = 1; j <= fPhiAxis->GetNbins(); ++j) {
923 const Int_t bin = j * (nEta + 2) + i;
924 if (vec[bin] > threshold && rebinData.fBinData[bin] != -1) {
925 cells2D.push_back(Cell2D_t(bin, vec[bin], max_e_slice[bin]));
926 cells2D.back().SetGeom(fEtaAxis->GetBinLowEdge(i), fEtaAxis->GetBinUpEdge(i),
928 }
929 }
930 }
931}
932
933////////////////////////////////////////////////////////////////////////////////
934/// Draw cells in top view.
935
937{
938 Float_t bws = -1; //smallest bin
939 Float_t logMax = -1;
940
942
944 {
945 fM->AssertPalette();
946 UChar_t col[4];
947
948 for (vCell2D_i i = cells2D.begin(); i != cells2D.end(); ++i)
949 {
950 if (rnrCtx.SecSelection()) glLoadName(i->fId);
951 glBegin(GL_POLYGON);
952 fM->fPalette->ColorFromValue(TMath::FloorNint(i->fSumVal), col);
953 col[3] = fM->GetData()->GetSliceTransparency(i->fMaxSlice);
955 Float_t z = fM->GetHasFixedHeightIn2DMode() ? baseOffset : i->fSumVal;
956 glVertex3f(i->fX0, i->fY0, z);
957 glVertex3f(i->fX1, i->fY0, z);
958 glVertex3f(i->fX1, i->fY1, z);
959 glVertex3f(i->fX0, i->fY1, z);
960 glEnd();
961 }
962 }
963 else
964 {
965 Float_t maxv = 0;
966 bws = 1e5;
967 for (vCell2D_i i = fCells2D.begin(); i != fCells2D.end(); ++i)
968 {
969 if (i->MinSize() < bws) bws = i->MinSize();
970 if (i->fSumVal > maxv) maxv = i->fSumVal;
971 }
972 bws *= 0.5f;
973 logMax = TMath::Log10(maxv + 1);
974 fValToPixel = bws / logMax;
975
976 if (rnrCtx.SecSelection())
977 {
978 // Special draw for name stack.
979
980 for (vCell2D_i i = cells2D.begin(); i != cells2D.end(); ++i)
981 {
982 glLoadName(i->fMaxSlice);
983 glPushName(i->fId);
984
985 glBegin(GL_QUADS);
986 Float_t z = fM->GetHasFixedHeightIn2DMode() ? baseOffset : i->fSumVal;
987 glVertex3f(i->fX0, i->fY0, z);
988 glVertex3f(i->fX1, i->fY0, z);
989 glVertex3f(i->fX1, i->fY1, z);
990 glVertex3f(i->fX0, i->fY1, z);
991 glEnd();
992
993 glPopName();
994 }
995 }
996 else
997 {
998 // Optimised draw without name stack.
999
1000 if ( ! rnrCtx.Highlight())
1001 {
1002 glBegin(GL_POINTS);
1003 for (vCell2D_i i = cells2D.begin(); i != cells2D.end(); ++i)
1004 {
1006 Float_t z = fM->GetHasFixedHeightIn2DMode() ? baseOffset : i->fSumVal;
1007 glVertex3f(i->X(), i->Y() , z);
1008 }
1009 glEnd();
1010 }
1011
1012 glBegin(GL_QUADS);
1013 for (vCell2D_i i = cells2D.begin(); i != cells2D.end(); ++i)
1014 {
1016 Float_t bw = fValToPixel*TMath::Log10(i->fSumVal+1);
1017 Float_t x = i->X();
1018 Float_t y = i->Y();
1019 Float_t z = fM->GetHasFixedHeightIn2DMode() ? baseOffset : i->fSumVal;
1020 glVertex3f(x - bw, y - bw, z);
1021 glVertex3f(x + bw, y - bw, z);
1022 glVertex3f(x + bw, y + bw, z);
1023 glVertex3f(x - bw, y + bw, z);
1024 }
1025 glEnd();
1026
1028 {
1029 glPushAttrib(GL_ENABLE_BIT | GL_POLYGON_BIT);
1030 Float_t z = 0;
1031 Float_t zOff = fMaxVal*0.001 ; // avoid polygon stipling
1032 glBegin(GL_QUADS);
1033 for ( vCell2D_i i = cells2D.begin(); i != cells2D.end(); ++i) {
1034 Char_t transp = TMath::Min(100, 80 + fM->fData->GetSliceTransparency(i->fMaxSlice) / 5);
1035 TGLUtil::ColorTransparency(fM->fData->GetSliceColor(i->fMaxSlice), transp);
1036 z = fM->GetHasFixedHeightIn2DMode() ? baseOffset : i->fSumVal;
1037 z -= zOff;
1038 glVertex3f(i->fX0, i->fY0, z);
1039 glVertex3f(i->fX1, i->fY0, z);
1040 glVertex3f(i->fX1, i->fY1, z);
1041 glVertex3f(i->fX0, i->fY1, z);
1042 }
1043 glEnd();
1044
1045 glPolygonMode(GL_FRONT_AND_BACK, GL_LINE);
1046 glBegin(GL_QUADS);
1047 for ( vCell2D_i i = cells2D.begin(); i != cells2D.end(); ++i) {
1048 TGLUtil::ColorTransparency(fM->fData->GetSliceColor(i->fMaxSlice), 60);
1049 z = fM->GetHasFixedHeightIn2DMode() ? baseOffset : i->fSumVal;
1050 z += zOff;
1051 glVertex3f(i->fX0, i->fY0, z);
1052 glVertex3f(i->fX1, i->fY0, z);
1053 glVertex3f(i->fX1, i->fY1, z);
1054 glVertex3f(i->fX0, i->fY1, z);
1055 }
1056 glEnd();
1057 glPopAttrib();
1058 }
1059 }
1060 }
1061
1062 // text
1064 ! rnrCtx.Selection() && ! rnrCtx.Highlight())
1065 {
1067 TGLFont font;
1069 const char* txt;
1070 for (vCell2D_i i = cells2D.begin(); i != cells2D.end(); ++i) {
1071
1072 Float_t val = i->fSumVal;
1073 if (val > 10)
1074 txt = Form("%d", TMath::Nint(val));
1075 else if (val > 1 )
1076 txt = Form("%.1f", val);
1077 else if (val > 0.01 )
1078 txt = Form("%.2f", 0.01*TMath::Nint(val*100));
1079 else
1080 txt = Form("~1e%d", TMath::Nint(TMath::Log10(val)));
1081
1082 font.Render(txt, i->X(), i->Y(), val*1.2, TGLFont::kCenterH, TGLFont::kCenterV);
1083 }
1084 }
1085}
1086
1087////////////////////////////////////////////////////////////////////////////////
1088/// Draw highligted cells.
1089
1090void TEveCaloLegoGL::DrawHighlight(TGLRnrCtx& rnrCtx, const TGLPhysicalShape* /*pshp*/, Int_t /*lvl*/) const
1091{
1092 if (fM->fData->GetCellsSelected().empty() && fM->fData->GetCellsHighlighted().empty())
1093 {
1094 return;
1095 }
1096
1097 // modelview matrix
1098 glPushMatrix();
1099 Float_t sx, sy, sz;
1100 GetScaleForMatrix(sx, sy, sz);
1101 glScalef(sx, sy, sz);
1102 glTranslatef(-fM->GetEta(), -fM->fPhi, 0);
1103
1104 if (fCells3D)
1105 {
1106 glPushAttrib(GL_ENABLE_BIT | GL_LINE_BIT | GL_POLYGON_BIT);
1107 glDisable(GL_LIGHTING);
1108 glDisable(GL_CULL_FACE);
1109 glPolygonMode(GL_FRONT_AND_BACK, GL_LINE);
1111 }
1112
1114 if (!fM->fData->GetCellsHighlighted().empty())
1115 {
1116 glColor4ubv(rnrCtx.ColorSet().Selection(3).CArr());
1118 }
1119 if (!fM->fData->GetCellsSelected().empty())
1120 {
1121 glColor4ubv(rnrCtx.ColorSet().Selection(1).CArr());
1123 }
1125
1126 if (fCells3D)
1127 {
1128 glPopAttrib();
1129 }
1130
1131 glPopMatrix();
1132}
1133
1134
1135////////////////////////////////////////////////////////////////////////////////
1136/// Draw selected cells in highlight mode.
1137
1139{
1140 // check eta&phi range of selected cells
1141 TEveCaloData::vCellId_t cellsSelected;
1142 TEveCaloData::CellData_t cellData;
1143 for (TEveCaloData::vCellId_i i = cellsSelectedInput.begin(); i != cellsSelectedInput.end(); ++i)
1144 {
1145 fM->fData->GetCellData((*i), cellData);
1146 if (fM->CellInEtaPhiRng(cellData))
1147 cellsSelected.push_back(*i);
1148 }
1149
1150 // prepare rebin for 2D or 3D if necessary
1151 TEveCaloData::RebinData_t rebinDataSelected;
1152 if (fBinStep > 1)
1153 {
1154 fM->fData->Rebin(fEtaAxis, fPhiAxis, cellsSelected, fM->fPlotEt, rebinDataSelected);
1155 if (fM->fNormalizeRebin) {
1156 Float_t scale = 1.f / (fBinStep * fBinStep);
1157 for (std::vector<Float_t>::iterator it = rebinDataSelected.fSliceData.begin();
1158 it != rebinDataSelected.fSliceData.end(); ++it)
1159 (*it) *= scale;
1160 }
1161 }
1162
1163 if (fCells3D)
1164 {
1165 Float_t offset = 0;
1166 if (fBinStep == 1)
1167 {
1168 for (TEveCaloData::vCellId_i j = cellsSelected.begin(); j != cellsSelected.end(); ++j)
1169 {
1170 offset = 0;
1171 {
1172 Int_t orig_slice = j->fSlice;
1173 for (Int_t s = 0; s < orig_slice; ++s)
1174 {
1175 j->fSlice = s;
1176 fM->fData->GetCellData(*j, cellData);
1177 offset += cellData.Value(fM->fPlotEt);
1178 }
1179 j->fSlice = orig_slice;
1180 }
1181 fM->fData->GetCellData(*j, cellData);
1182 WrapTwoPi(cellData.fPhiMin, cellData.fPhiMax);
1183 MakeQuad(cellData.EtaMin(), cellData.PhiMin(), offset,
1184 cellData.EtaDelta(), cellData.PhiDelta(), cellData.Value(fM->fPlotEt));
1185 }
1186 }
1187 else
1188 {
1189 Float_t *vals;
1190 Float_t *valsRef;
1191 Float_t y0, y1;
1192 Int_t nSlices = fM->fData->GetNSlices();
1193 for (Int_t i = 1; i <= fEtaAxis->GetNbins(); ++i)
1194 {
1195 for (Int_t j = 1; j <= fPhiAxis->GetNbins(); ++j)
1196 {
1197 const Int_t bin = (i)+(j)*(fEtaAxis->GetNbins()+2);
1198 if (rebinDataSelected.fBinData[bin] !=-1)
1199 {
1200 offset = 0;
1201 vals = rebinDataSelected.GetSliceVals(bin);
1202 valsRef = fRebinData.GetSliceVals(bin);
1203 for (Int_t s = 0; s < nSlices; ++s)
1204 {
1205 if (vals[s] > 0)
1206 {
1207 y0 = fPhiAxis->GetBinLowEdge(j);
1208 y1 = fPhiAxis->GetBinUpEdge(j);
1209 WrapTwoPi(y0, y1);
1211 fEtaAxis->GetBinWidth(i), y1-y0, vals[s]);
1212 }
1213 offset += valsRef[s];
1214 }
1215 }
1216 }
1217 }
1218 }
1219 }
1220 else
1221 {
1222 vCell2D_t cells2DSelected;
1223 if (fBinStep == 1)
1224 {
1225 // but is confusing since top view does not draw all slices at same time
1226 TEveCaloData::vCellId_i j = cellsSelectedInput.begin();
1227 TEveCaloData::vCellId_i jEnd = cellsSelectedInput.end();
1228 std::set<Int_t> towers;
1229 while (j != jEnd)
1230 {
1231 towers.insert(j->fTower);
1232 ++j;
1233 }
1234 for (vCell2D_i i = fCells2D.begin(); i != fCells2D.end(); ++i)
1235 {
1236 TEveCaloData::CellId_t cell = fM->fCellList[i->fId];
1237 // std::set<Int_t>::iterator ti = towers.find(cell.fTower);
1238 if (towers.find(cell.fTower) != towers.end())
1239 {
1240 cells2DSelected.push_back(*i);
1241 }
1242 }
1243 }
1244 else
1245 {
1246 PrepareCell2DDataRebin(rebinDataSelected, cells2DSelected);
1247 }
1248 DrawCells2D(rnrCtx, cells2DSelected);
1249 }
1250}
1251
1252////////////////////////////////////////////////////////////////////////////////
1253/// Draw the object.
1254
1256{
1257 if (! fM->fData || ! fM->fData->GetEtaBins() || ! fM->fData->GetPhiBins())
1258 return;
1259
1260 // projection type
1262 fCells3D = (!(rnrCtx.RefCamera().IsOrthographic() && rnrCtx.RefCamera().GetCamBase().GetBaseVec(1).Z()));
1263 else if (fM->fProjection == TEveCaloLego::k2D)
1264 fCells3D = kFALSE;
1265 else if (fM->fProjection == TEveCaloLego::k3D)
1266 fCells3D = kTRUE;
1267
1268 // rebin axsis , check limits, fix TwoPi cycling
1269 Int_t new_bin_step = GetGridStep(rnrCtx);
1270
1271 // rebin data
1272 if (fM->AssertCellIdCache() || fBinStep != new_bin_step)
1273 {
1274 fBinStep = new_bin_step;
1276 fRebinData.Clear();
1277
1280
1281 if (fBinStep > 1)
1282 {
1284
1285 fMaxVal = 0;
1286 for (UInt_t i = 0; i < fRebinData.fSliceData.size(); i += fRebinData.fNSlices)
1287 {
1288 Double_t sum = 0;
1289 for (Int_t s = 0; s < fRebinData.fNSlices; s++)
1290 {
1291 sum += fRebinData.fSliceData[i+s];
1292 }
1293 if (sum > fMaxVal) fMaxVal = sum;
1294 }
1295
1296 if (fM->fNormalizeRebin)
1297 {
1298 Float_t scale = 1.f / (fBinStep * fBinStep);
1299 for (std::vector<Float_t>::iterator it = fRebinData.fSliceData.begin(); it != fRebinData.fSliceData.end();
1300 ++it) {
1301 (*it) *= scale;
1302 }
1303 fMaxVal *= scale;
1304 }
1305 }
1306 else
1307 {
1308 fMaxVal = fM->GetMaxVal();
1309 }
1310 }
1311
1312 // modelview matrix
1313 glPushMatrix();
1314 Float_t sx, sy, sz;
1315 GetScaleForMatrix(sx, sy, sz);
1316 glScalef(sx, sy, sz);
1317 glTranslatef(-fM->GetEta(), -fM->fPhi, 0);
1318
1321 if (fGridColor < 0 || fFontColor < 0)
1322 {
1325 Float_t f1, f2;
1326 if (fFontColor < 0) {
1327 f1 = 0.8; f2 = 0.2;
1328 fFontColor = TColor::GetColor(c1->GetRed() *f1 + c2->GetRed() *f2,
1329 c1->GetGreen()*f1 + c2->GetGreen()*f2,
1330 c1->GetBlue() *f1 + c2->GetBlue() *f2);
1331 }
1332 if (fGridColor < 0) {
1333 f1 = 0.3; f2 = 0.3;
1334 fGridColor = TColor::GetColor(c1->GetRed() *f1 + c2->GetRed() *f2,
1335 c1->GetGreen()*f1 + c2->GetGreen()*f2,
1336 c1->GetBlue() *f1 + c2->GetBlue() *f2);
1337 }
1338 }
1339
1340 glPushAttrib(GL_ENABLE_BIT | GL_LINE_BIT | GL_POLYGON_BIT);
1342 glEnable(GL_BLEND);
1343 glBlendFunc(GL_SRC_ALPHA, GL_ONE_MINUS_SRC_ALPHA);
1344 if (!fM->fData->Empty())
1345 {
1346 glPushName(0);
1347 if (fCells3D)
1348 {
1349 if (fDLCacheOK == kFALSE)
1350 {
1351 if (fBinStep == 1)
1353 else
1355 fDLCacheOK = kTRUE;
1356 }
1357 glEnable(GL_NORMALIZE);
1358 glEnable(GL_POLYGON_OFFSET_FILL);
1359 glPolygonOffset(0.8, 1);
1360
1361 DrawCells3D(rnrCtx);
1362 }
1363 else
1364 {
1365 glDisable(GL_LIGHTING);
1366
1367 fCells2D.clear();
1368 if (fBinStep == 1)
1370 else
1372
1373 DrawCells2D(rnrCtx, fCells2D);
1374 }
1375 glPopName();
1376 }
1377 glPopAttrib();
1378
1379 // draw histogram base
1380 if (rnrCtx.Selection() == kFALSE && rnrCtx.IsDrawPassFilled())
1381 {
1382 glPushAttrib(GL_ENABLE_BIT | GL_LINE_BIT | GL_POLYGON_BIT);
1383 glDisable(GL_LIGHTING);
1384 DrawHistBase(rnrCtx);
1385 if (fM->fDrawHPlane) {
1386 glBlendFunc(GL_ONE, GL_ONE_MINUS_SRC_ALPHA);
1387 glPolygonMode(GL_FRONT_AND_BACK, GL_FILL);
1388 glDisable(GL_CULL_FACE);
1390 Float_t zhp = fM->fHPlaneVal * fMaxVal;
1391 glBegin(GL_POLYGON);
1392 glVertex3f(fM->fEtaMin, fM->GetPhiMin(), zhp);
1393 glVertex3f(fM->fEtaMax, fM->GetPhiMin(), zhp);
1394 glVertex3f(fM->fEtaMax, fM->GetPhiMax(), zhp);
1395 glVertex3f(fM->fEtaMin, fM->GetPhiMax(), zhp);
1396 glEnd();
1397 }
1398 glPopAttrib();
1399 }
1400
1401 glPopMatrix();
1402}
1403
1404////////////////////////////////////////////////////////////////////////////////
1405/// Processes tower selection from TGLViewer.
1406
1408{
1410 if (rec.GetN() > 2)
1411 {
1412 Int_t slice = rec.GetItem(1);
1413 Int_t cell = rec.GetItem(2);
1414
1415 if (fBinStep == 1)
1416 {
1417 Int_t tower = fM->fCellList[cell].fTower;
1418 while (cell > 0 && tower == fM->fCellList[cell].fTower)
1419 {
1420 sel.push_back(fM->fCellList[cell]);
1421 if (fCells3D) break;
1422 --cell;
1423 }
1424 }
1425 else
1426 {
1427 if (cell > 0)
1428 {
1429 Int_t nEta = fEtaAxis->GetNbins();
1430 Int_t phiBin = Int_t(cell/(nEta+2));
1431 Int_t etaBin = cell - phiBin*(nEta+2);
1434 fPhiAxis->GetBinCenter(phiBin), fPhiAxis->GetBinWidth(phiBin),
1435 sl);
1436
1437 for (TEveCaloData::vCellId_i it = sl.begin(); it != sl.end(); ++it)
1438 {
1439 if (fCells3D) {
1440 if ((*it).fSlice == slice ) sel.push_back(*it);
1441 } else {
1442 if ((*it).fSlice <= slice ) sel.push_back(*it);
1443 }
1444 }
1445 }
1446 }
1447 }
1448 fM->fData->ProcessSelection(sel, rec);
1449}
#define GL_QUADS
Definition GL_glu.h:290
#define GL_LINES
Definition GL_glu.h:284
#define GL_POINTS
Definition GL_glu.h:283
#define GL_POLYGON
Definition GL_glu.h:292
#define b(i)
Definition RSha256.hxx:100
#define h(i)
Definition RSha256.hxx:106
#define e(i)
Definition RSha256.hxx:103
int Int_t
Definition RtypesCore.h:45
unsigned char UChar_t
Definition RtypesCore.h:38
char Char_t
Definition RtypesCore.h:37
float Float_t
Definition RtypesCore.h:57
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
constexpr Bool_t kTRUE
Definition RtypesCore.h:93
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:382
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t sel
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h offset
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t UChar_t len
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char y2
Option_t Option_t TPoint TPoint const char y1
#define gROOT
Definition TROOT.h:406
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2489
virtual Int_t GetNdivisions() const
Definition TAttAxis.h:36
virtual Style_t GetTitleFont() const
Definition TAttAxis.h:47
virtual Float_t GetLabelOffset() const
Definition TAttAxis.h:40
virtual void SetAxisColor(Color_t color=1, Float_t alpha=1.)
Set color of the line axis and tick marks.
Definition TAttAxis.cxx:160
virtual void SetTitleFont(Style_t font=62)
Set the title font.
Definition TAttAxis.cxx:327
virtual void SetLabelOffset(Float_t offset=0.005)
Set distance between the axis and the labels.
Definition TAttAxis.cxx:191
virtual void SetTitleSize(Float_t size=0.04)
Set size of axis title.
Definition TAttAxis.cxx:309
virtual void SetTitleColor(Color_t color=1)
Set color of axis title.
Definition TAttAxis.cxx:318
virtual Float_t GetTitleSize() const
Definition TAttAxis.h:44
virtual Float_t GetLabelSize() const
Definition TAttAxis.h:41
virtual Float_t GetTickLength() const
Definition TAttAxis.h:45
virtual void SetTickLength(Float_t length=0.03)
Set tick mark length.
Definition TAttAxis.cxx:284
virtual void SetNdivisions(Int_t n=510, Bool_t optim=kTRUE)
Set the number of divisions for this axis.
Definition TAttAxis.cxx:233
virtual void SetLabelColor(Color_t color=1, Float_t alpha=1.)
Set color of labels.
Definition TAttAxis.cxx:170
Class to manage histogram axis.
Definition TAxis.h:31
const char * GetTitle() const override
Returns title of object.
Definition TAxis.h:135
virtual Double_t GetBinCenter(Int_t bin) const
Return center of bin.
Definition TAxis.cxx:478
Double_t GetXmax() const
Definition TAxis.h:140
virtual Int_t FindBin(Double_t x)
Find bin number corresponding to abscissa x.
Definition TAxis.cxx:293
virtual Double_t GetBinLowEdge(Int_t bin) const
Return low edge of bin.
Definition TAxis.cxx:518
virtual void Set(Int_t nbins, Double_t xmin, Double_t xmax)
Initialize axis with fix bins.
Definition TAxis.cxx:794
virtual void SetLimits(Double_t xmin, Double_t xmax)
Definition TAxis.h:164
Double_t GetXmin() const
Definition TAxis.h:139
Int_t GetNbins() const
Definition TAxis.h:125
virtual Double_t GetBinWidth(Int_t bin) const
Return bin width.
Definition TAxis.cxx:540
virtual Double_t GetBinUpEdge(Int_t bin) const
Return up edge of bin.
Definition TAxis.cxx:528
The color creation and management class.
Definition TColor.h:21
static Int_t GetColor(const char *hexcolor)
Static method returning color number for color specified by hex color string of form: "#rrggbb",...
Definition TColor.cxx:1920
virtual void GetEtaLimits(Double_t &min, Double_t &max) const =0
std::vector< CellId_t > vCellId_t
Bool_t Empty() const
virtual void GetCellList(Float_t etaMin, Float_t etaMax, Float_t phi, Float_t phiRng, vCellId_t &out) const =0
void ProcessSelection(vCellId_t &sel_cells, TGLSelectRecord &rec)
Process newly selected cells with given select-record.
virtual void GetCellData(const CellId_t &id, CellData_t &data) const =0
Char_t GetSliceTransparency(Int_t slice) const
Get transparency for given slice.
virtual TAxis * GetEtaBins() const
Color_t GetSliceColor(Int_t slice) const
Get color for given slice.
Int_t GetNSlices() const
vCellId_t & GetCellsHighlighted()
virtual void GetPhiLimits(Double_t &min, Double_t &max) const =0
virtual TAxis * GetPhiBins() const
vCellId_t & GetCellsSelected()
virtual void Rebin(TAxis *ax, TAxis *ay, vCellId_t &in, Bool_t et, RebinData_t &out) const =0
std::vector< CellId_t >::iterator vCellId_i
OpenGL renderer class for TEveCaloLego.
std::vector< Cell2D_t > vCell2D_t
void MakeQuad(Float_t x, Float_t y, Float_t z, Float_t xw, Float_t yw, Float_t zh) const
Draw an axis-aligned box using quads.
TEveVector fBackPlaneYConst[2]
TEveVector fZAxisTitlePos
void WrapTwoPi(Float_t &min, Float_t &max) const
void PrepareCell2DData(TEveCaloData::vCellId_t &cellList, vCell2D_t &cells2D) const
Prepare cells 2D data non-rebinned for drawing.
void DLCacheDrop() override
Drop all display-list definitions.
Int_t fCurrentPixelsPerBin
void DrawSelectedCells(TGLRnrCtx &rnrCtx, TEveCaloData::vCellId_t cells) const
Draw selected cells in highlight mode.
void DrawCells3D(TGLRnrCtx &rnrCtx) const
Render the calo lego-plot with OpenGL.
void DrawHistBase(TGLRnrCtx &rnrCtx) const
Draw basic histogram components: x-y grid.
TEveVector fXAxisTitlePos
void DrawCells2D(TGLRnrCtx &rnrCtx, vCell2D_t &cells2D) const
Draw cells in top view.
void DrawAxis3D(TGLRnrCtx &rnrCtx) const
Draw z-axis and z-box at the appropriate grid corner-point including tick-marks and labels.
Bool_t SetModel(TObject *obj, const Option_t *opt=nullptr) override
Set model object.
void RebinAxis(TAxis *orig, TAxis *curr) const
Rebin eta, phi axis.
void ProcessSelection(TGLRnrCtx &rnrCtx, TGLSelectRecord &rec) override
Processes tower selection from TGLViewer.
Int_t GetGridStep(TGLRnrCtx &rnrCtx) const
Calculate view-dependent grid density.
void DLCachePurge() override
Unregister all display-lists.
void DirectDraw(TGLRnrCtx &rnrCtx) const override
Draw the object.
std::vector< Cell2D_t >::iterator vCell2D_i
vCell2D_t fCells2D
void GetScaleForMatrix(Float_t &sx, Float_t &sy, Float_t &sz) const
Get scale for matrix.
TEveVector fBackPlaneXConst[2]
~TEveCaloLegoGL() override
Destructor.
void SetBBox() override
Set bounding box.
void Make3DDisplayListRebin(TEveCaloData::RebinData_t &rebinData, SliceDLMap_t &map, Bool_t select) const
Create display-list that draws histogram bars for rebinned data.
void DrawHighlight(TGLRnrCtx &rnrCtx, const TGLPhysicalShape *ps, Int_t lvl=-1) const override
Draw highligted cells.
void SetAxis3DTitlePos(TGLRnrCtx &rnrCtx, Float_t x0, Float_t x1, Float_t y0, Float_t y1) const
Set the axis 3D title position.
TEveCaloLegoGL()
Constructor.
TEveCaloLego * fM
void DrawAxis2D(TGLRnrCtx &rnrCtx) const
Draw XY axis.
void Make3DDisplayList(TEveCaloData::vCellId_t &cellList, SliceDLMap_t &map, Bool_t select) const
Create display-list that draws histogram bars for non-rebinned data.
TGLAxisPainter fAxisPainter
std::map< Int_t, UInt_t > SliceDLMap_t
TEveVector fYAxisTitlePos
std::map< Int_t, UInt_t >::iterator SliceDLMap_i
void PrepareCell2DDataRebin(TEveCaloData::RebinData_t &rebinData, vCell2D_t &cells2D) const
Prepare cells 2D rebinned data for drawing.
SliceDLMap_t fDLMap
TEveCaloData::RebinData_t fRebinData
Visualization of calorimeter data as eta/phi histogram.
Definition TEveCalo.h:251
EProjection_e fProjection
Definition TEveCalo.h:279
Bool_t fNormalizeRebin
Definition TEveCalo.h:277
Char_t fPlaneTransparency
Definition TEveCalo.h:270
Color_t fPlaneColor
Definition TEveCalo.h:269
EBoxMode_e fBoxMode
Definition TEveCalo.h:281
Int_t fDrawNumberCellPixels
Definition TEveCalo.h:289
E2DMode_e f2DMode
Definition TEveCalo.h:280
Float_t fHPlaneVal
Definition TEveCalo.h:284
float GetFixedHeightValIn2DMode() const
Definition TEveCalo.h:334
TEveCaloData::vCellId_t fCellList
Definition TEveCalo.h:265
Color_t fFontColor
Definition TEveCalo.h:267
Bool_t fDrawHPlane
Definition TEveCalo.h:283
Bool_t fAutoRebin
Definition TEveCalo.h:275
bool GetHasFixedHeightIn2DMode() const
Definition TEveCalo.h:331
Int_t fCellPixelFontSize
Definition TEveCalo.h:290
Color_t fGridColor
Definition TEveCalo.h:268
Int_t fNZSteps
Definition TEveCalo.h:272
Int_t fPixelsPerBin
Definition TEveCalo.h:276
Color_t GetDataSliceColor(Int_t slice) const
Get slice color from data.
Definition TEveCalo.cxx:121
Float_t GetEtaMin() const
Definition TEveCalo.h:137
Bool_t AssertCellIdCache() const
Assert cell id cache is ok.
Definition TEveCalo.cxx:293
Float_t fPlotEt
Definition TEveCalo.h:55
Float_t GetDataSliceThreshold(Int_t slice) const
Get threshold for given slice.
Definition TEveCalo.cxx:86
Float_t GetMaxTowerH() const
Definition TEveCalo.h:113
Float_t GetPhiRng() const
Definition TEveCalo.h:147
TEveRGBAPalette * fPalette
Definition TEveCalo.h:62
Double_t fEtaMax
Definition TEveCalo.h:44
Float_t GetMaxVal() const
Definition TEveCalo.cxx:159
TEveCaloData * GetData() const
Definition TEveCalo.h:87
Float_t GetEta() const
Definition TEveCalo.h:136
TEveRGBAPalette * AssertPalette()
Make sure the TEveRGBAPalette pointer is not null.
Definition TEveCalo.cxx:378
Bool_t GetPlotEt() const
Definition TEveCalo.h:109
Float_t GetEtaMax() const
Definition TEveCalo.h:138
Bool_t fScaleAbs
Definition TEveCalo.h:58
Double_t fPhi
Definition TEveCalo.h:46
TEveCaloData * fData
Definition TEveCalo.h:40
Float_t GetPhiMax() const
Definition TEveCalo.h:146
Double_t fEtaMin
Definition TEveCalo.h:43
Float_t GetEtaRng() const
Definition TEveCalo.h:139
Float_t fMaxValAbs
Definition TEveCalo.h:59
Bool_t CellInEtaPhiRng(TEveCaloData::CellData_t &) const
Returns true if given cell is in the ceta phi range.
Definition TEveCalo.cxx:307
Float_t GetPhiMin() const
Definition TEveCalo.h:145
const UChar_t * ColorFromValue(Int_t val) const
void Set(const Float_t *v)
Definition TEveVector.h:82
TGLVector3 & RefDir()
void SetLabelAlign(TGLFont::ETextAlignH_e, TGLFont::ETextAlignV_e)
Set label align.
TGLVector3 & RefTitlePos()
void SetAttAxis(TAttAxis *a)
void SetLabelPixelFontSize(Int_t fs)
void SetTMNDim(Int_t x)
TGLVector3 & RefTMOff(Int_t i)
void SetFontMode(TGLFont::EMode m)
Int_t GetLabelPixelFontSize() const
void PaintAxis(TGLRnrCtx &ctx, TAxis *ax)
GL render TAxis.
void SetTitlePixelFontSize(Int_t fs)
Abstract base camera class - concrete classes for orthographic and perspective cameras derive from it...
Definition TGLCamera.h:44
TGLMatrix & RefLastNoPickProjM() const
Definition TGLCamera.h:174
const TGLMatrix & GetCamBase() const
Definition TGLCamera.h:166
TGLVector3 ViewportDeltaToWorld(const TGLVertex3 &worldRef, Double_t viewportXDelta, Double_t viewportYDelta, TGLMatrix *modviewMat=nullptr) const
Apply a 2D viewport delta (shift) to the projection of worldRef onto viewport, returning the resultan...
virtual Bool_t IsOrthographic() const
Definition TGLCamera.h:118
const TGLPlane & FrustumPlane(EFrustumPlane plane) const
Definition TGLCamera.h:219
TGLRect & RefViewport()
Definition TGLCamera.h:128
TGLColor & Markup()
Definition TGLUtil.h:854
TGLColor & Selection(Int_t i)
Definition TGLUtil.h:855
TGLColor & Background()
Definition TGLUtil.h:851
const UChar_t * CArr() const
Definition TGLUtil.h:800
Color_t GetColorIndex() const
Returns color-index representing the color.
Definition TGLUtil.cxx:1217
A wrapper class for FTFont.
void Render(const char *txt, Double_t x, Double_t y, Double_t angle, Double_t mgn) const
virtual void DLCachePurge()
Purge all entries for all LODs for this drawable from the display list cache, returning the reserved ...
TObject * fExternalObj
first replica
virtual void DLCacheDrop()
Drop all entries for all LODs for this drawable from the display list cache, WITHOUT returning the re...
void PurgeDLRange(UInt_t base, Int_t size) const
External object is a fake.
Bool_t fDLCache
display-list validity bit-field
16 component (4x4) transform matrix - column MAJOR as per GL.
Definition TGLUtil.h:598
Double_t * Arr()
Definition TGLUtil.h:665
TGLVector3 GetBaseVec(Int_t b) const
Definition TGLUtil.h:754
const Double_t * CArr() const
Definition TGLUtil.h:664
Base-class for direct OpenGL renderers.
Definition TGLObject.h:22
void SetAxisAlignedBBox(Float_t xmin, Float_t xmax, Float_t ymin, Float_t ymax, Float_t zmin, Float_t zmax)
Set axis-aligned bounding-box.
Definition TGLObject.cxx:86
Concrete physical shape - a GL drawable.
Double_t D() const
Definition TGLUtil.h:556
const Int_t * CArr() const
Definition TGLUtil.h:443
Int_t Diagonal() const
Return the diagonal of the rectangle.
Definition TGLUtil.cxx:286
The TGLRnrCtx class aggregates data for a given redering context as needed by various parts of the RO...
Definition TGLRnrCtx.h:41
Short_t SceneStyle() const
Definition TGLRnrCtx.h:184
Bool_t SecSelection() const
Definition TGLRnrCtx.h:224
void RegisterFontNoScale(Int_t size, Int_t file, Int_t mode, TGLFont &out)
Get font in the GL rendering context.
Bool_t IsDrawPassFilled() const
Returns true if current render-pass uses filled polygon style.
TGLColorSet & ColorSet()
Return reference to current color-set (top of the stack).
TGLCamera & RefCamera()
Definition TGLRnrCtx.h:157
Bool_t Highlight() const
Definition TGLRnrCtx.h:218
Bool_t Selection() const
Definition TGLRnrCtx.h:222
Standard selection record including information about containing scene and details ob out selected ob...
static void Color4ubv(const UChar_t *rgba)
Wrapper for glColor4ubv.
Definition TGLUtil.cxx:1779
static UInt_t LockColor()
Prevent further color changes.
Definition TGLUtil.cxx:1669
static void ColorTransparency(Color_t color_index, Char_t transparency=0)
Set color from color_index and ROOT-style transparency (default 0).
Definition TGLUtil.cxx:1741
static UInt_t UnlockColor()
Allow color changes.
Definition TGLUtil.cxx:1677
static void Color(const TGLColor &color)
Set color from TGLColor.
Definition TGLUtil.cxx:1697
static Float_t LineWidth()
Get the line-width, taking the global scaling into account.
Definition TGLUtil.cxx:1943
3 component (x/y/z) vector class.
Definition TGLUtil.h:248
3 component (x/y/z) vertex class.
Definition TGLUtil.h:84
Double_t X() const
Definition TGLUtil.h:119
Double_t Z() const
Definition TGLUtil.h:123
void Set(Double_t x, Double_t y, Double_t z)
Definition TGLUtil.h:210
Double_t Y() const
Definition TGLUtil.h:121
static void Optimize(Double_t A1, Double_t A2, Int_t nold, Double_t &BinLow, Double_t &BinHigh, Int_t &nbins, Double_t &BWID, Option_t *option="")
Static function to compute reasonable axis limits.
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:164
Mother of all ROOT objects.
Definition TObject.h:41
Double_t y[n]
Definition legend1.C:17
return c1
Definition legend1.C:41
Double_t x[n]
Definition legend1.C:17
TF1 * f1
Definition legend1.C:11
return c2
Definition legend2.C:14
Int_t Nint(T x)
Round to nearest integer. Rounds half integers to the nearest even integer.
Definition TMath.h:693
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:250
constexpr Double_t Sqrt2()
Definition TMath.h:86
Int_t FloorNint(Double_t x)
Returns the nearest integer of TMath::Floor(x).
Definition TMath.h:686
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:662
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Definition TMathBase.h:198
Double_t Hypot(Double_t x, Double_t y)
Returns sqrt(x*x + y*y)
Definition TMath.cxx:59
Double_t Log10(Double_t x)
Returns the common (base-10) logarithm of x.
Definition TMath.h:762
Cell data inner structure.
Float_t Value(Bool_t) const
Return energy value associated with the cell, usually Et.
Float_t PhiDelta() const
Float_t PhiMin() const
Float_t EtaDelta() const
Float_t EtaMin() const
std::vector< Float_t > fSliceData
std::vector< Int_t > fBinData
Float_t * GetSliceVals(Int_t bin)
TLine l
Definition textangle.C:4
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345