ROOT  6.06/09
Reference Guide
TEveCalo3DGL.cxx
Go to the documentation of this file.
1 // @(#)root/eve:$Id$
2 // Author: Matevz 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 "TEveCalo3DGL.h"
13 #include "TEveCalo.h"
14 
15 #include "TMath.h"
16 #include "TAxis.h"
17 
18 #include "TGLRnrCtx.h"
19 #include "TGLSelectRecord.h"
20 #include "TGLPhysicalShape.h"
21 #include "TGLIncludes.h"
22 #include "TGLUtil.h"
23 #include "TEveRGBAPalette.h"
24 #include "TEveUtil.h"
25 
26 /** \class TEveCalo3DGL
27 \ingroup TEve
28 OpenGL renderer class for TEveCalo3D.
29 */
30 
32 
33 ////////////////////////////////////////////////////////////////////////////////
34 /// Constructor.
35 
37  TGLObject(), fM(0)
38 {
40 }
41 
42 ////////////////////////////////////////////////////////////////////////////////
43 /// Set model object.
44 
46 {
47  fM = SetModelDynCast<TEveCalo3D>(obj);
48  return kTRUE;
49 }
50 
51 ////////////////////////////////////////////////////////////////////////////////
52 /// Set bounding box.
53 
55 {
56  // !! This ok if master sub-classed from TAttBBox
57  SetAxisAlignedBBox(((TEveCalo3D*)fExternalObj)->AssertBBox());
58 }
59 
60 ////////////////////////////////////////////////////////////////////////////////
61 /// Override from TGLObject.
62 /// To account for large point-sizes we modify the projection matrix
63 /// during selection and thus we need a direct draw.
64 
66 {
67  if (rnrCtx.Highlight() || rnrCtx.Selection()) return kFALSE;
68  return TGLObject::ShouldDLCache(rnrCtx);
69 }
70 
71 ////////////////////////////////////////////////////////////////////////////////
72 /// Calculate cross-product.
73 
74 inline void TEveCalo3DGL::CrossProduct(const Float_t a[3], const Float_t b[3],
75  const Float_t c[3], Float_t out[3]) const
76 {
77  const Float_t v1[3] = { a[0] - c[0], a[1] - c[1], a[2] - c[2] };
78  const Float_t v2[3] = { b[0] - c[0], b[1] - c[1], b[2] - c[2] };
79 
80  out[0] = v1[1] * v2[2] - v1[2] * v2[1];
81  out[1] = v1[2] * v2[0] - v1[0] * v2[2];
82  out[2] = v1[0] * v2[1] - v1[1] * v2[0];
83 }
84 
85 ////////////////////////////////////////////////////////////////////////////////
86 /// Render end cap grid.
87 
89 {
90  using namespace TMath;
91 
92  Float_t rB = fM->GetBarrelRadius();
95 
96  Float_t etaMin = fM->GetEtaMin();
100  Float_t phiMin = fM->GetPhiMin();
101  Float_t phiMax = fM->GetPhiMax();
102 
103  TAxis *ax = fM->GetData()->GetEtaBins();
104  Int_t nx = ax->GetNbins();
105  TAxis *ay = fM->GetData()->GetPhiBins();
106  Int_t ny = ay->GetNbins();
107 
108 
109  Float_t r, z, theta, phiU, phiL, eta;
110 
111  // eta slices
112  for (Int_t i=0; i<=nx; ++i)
113  {
114  eta = ax->GetBinUpEdge(i);
115  if (eta >= transF && (eta > etaMin && eta < etaMax))
116  {
117  theta = TEveCaloData::EtaToTheta(eta);
118  r = Abs(zEF*Tan(theta));
119  z = Sign(zEF, ax->GetBinLowEdge(i));
120  for (Int_t j=1; j<=ny; ++j)
121  {
122  phiL = ay->GetBinLowEdge(j);
123  phiU = ay->GetBinUpEdge(j);
124  if (TEveUtil::IsU1IntervalContainedByMinMax(phiMin, phiMax, phiL, phiU))
125  {
126  glVertex3f(r*Cos(phiL), r*Sin(phiL), z);
127  glVertex3f(r*Cos(phiU), r*Sin(phiU), z);
128  }
129  }
130  } else if (eta <= transB && (eta > etaMin && eta < etaMax)) {
131  theta = TEveCaloData::EtaToTheta(eta);
132  r = Abs(zEB*Tan(theta));
133  z = Sign(zEB, ax->GetBinLowEdge(i));
134  for (Int_t j=1; j<=ny; ++j)
135  {
136  phiL = ay->GetBinLowEdge(j);
137  phiU = ay->GetBinUpEdge(j);
138  if (TEveUtil::IsU1IntervalContainedByMinMax(phiMin, phiMax, phiL, phiU))
139  {
140  glVertex3f(r*Cos(phiL), r*Sin(phiL), z);
141  glVertex3f(r*Cos(phiU), r*Sin(phiU), z);
142  }
143  }
144  }
145  }
146 
147  Float_t r1, r2;
148  // phi slices front
149  if (etaMax > transF)
150  {
151  r1 = zEF*Tan(TEveCaloData::EtaToTheta(etaMax));
152  if (etaMin < transF)
153  r2 = rB;
154  else
155  r2 = zEF*Tan(TEveCaloData::EtaToTheta(etaMin));
156 
157  for (Int_t j=1; j<=ny; ++j)
158  {
159  phiL = ay->GetBinLowEdge(j);
160  phiU = ay->GetBinUpEdge(j);
161  if (TEveUtil::IsU1IntervalContainedByMinMax(phiMin, phiMax, phiL, phiU))
162  {
163  glVertex3f( r1*Cos(phiU), r1*Sin(phiU), zEF);
164  glVertex3f( r2*Cos(phiU), r2*Sin(phiU), zEF);
165  glVertex3f( r1*Cos(phiL), r1*Sin(phiL), zEF);
166  glVertex3f( r2*Cos(phiL), r2*Sin(phiL), zEF);
167  }
168  }
169  }
170 
171  // phi slices back
172  if (etaMin < transB)
173  {
174  r1 = zEB*Tan(TEveCaloData::EtaToTheta(etaMin));
175  if (etaMax > transB)
176  r2 = rB;
177  else
178  r2 = zEB*Tan(TEveCaloData::EtaToTheta(etaMax));
179 
180  r1 = Abs(r1);
181  r2 = Abs(r2);
182  for (Int_t j=1; j<=ny; ++j)
183  {
184  phiL = ay->GetBinLowEdge(j);
185  phiU = ay->GetBinUpEdge(j);
186  if (TEveUtil::IsU1IntervalContainedByMinMax(phiMin, phiMax, phiL, phiU))
187  {
188  glVertex3f(r1*Cos(phiU), r1*Sin(phiU), zEB);
189  glVertex3f(r2*Cos(phiU), r2*Sin(phiU), zEB);
190  glVertex3f(r1*Cos(phiL), r1*Sin(phiL), zEB);
191  glVertex3f(r2*Cos(phiL), r2*Sin(phiL), zEB);
192  }
193  }
194  }
195 }
196 
197 ////////////////////////////////////////////////////////////////////////////////
198 /// Render barrel grid.
199 
201 {
202  using namespace TMath;
203 
204  Float_t etaMin = fM->GetEtaMin();
205  Float_t etaMax = fM->GetEtaMax();
206  Float_t transF = fM->GetTransitionEtaForward();
207  Float_t transB = fM->GetTransitionEtaBackward();
208  Float_t phiMin = fM->GetPhiMin();
209  Float_t phiMax = fM->GetPhiMax();
210 
211  Float_t rB = fM->GetBarrelRadius();
212  TAxis *ax = fM->GetData()->GetEtaBins();
213  Int_t nx = ax->GetNbins();
214  TAxis *ay = fM->GetData()->GetPhiBins();
215  Int_t ny = ay->GetNbins();
216 
217  Float_t z, theta, phiL, phiU, eta, x, y;
218 
219  // eta slices
220  for(Int_t i=0; i<=nx; i++)
221  {
222  eta = ax->GetBinUpEdge(i);
223  if (eta<=transF && eta>=transB && (etaMin < eta && eta < etaMax))
224  {
225  theta = TEveCaloData::EtaToTheta(eta);
226  z = rB/Tan(theta);
227  for (Int_t j=1; j<=ny; j++)
228  {
229  phiU = ay->GetBinUpEdge(j);
230  phiL = ay->GetBinLowEdge(j);
231  if (TEveUtil::IsU1IntervalContainedByMinMax(phiMin, phiMax, phiL, phiU))
232  {
233  glVertex3f(rB*Cos(phiL), rB*Sin(phiL), z);
234  glVertex3f(rB*Cos(phiU), rB*Sin(phiU), z);
235  }
236  }
237  }
238  }
239 
240  // phi slices
241  Float_t zF, zB;
242 
243  if (etaMin > transB)
244  zB = rB/Tan(TEveCaloData::EtaToTheta(etaMin));
245  else
246  zB = fM->GetBackwardEndCapPos();
247 
248 
249  if (etaMax < transF)
250  zF = rB/Tan(TEveCaloData::EtaToTheta(etaMax));
251  else
252  zF = fM->GetForwardEndCapPos();
253 
254  for (Int_t j=1; j<=ny; j++)
255  {
256  phiU = ay->GetBinUpEdge(j);
257  phiL = ay->GetBinLowEdge(j);
258  if (TEveUtil::IsU1IntervalContainedByMinMax(phiMin, phiMax, phiL, phiU))
259  {
260  x = rB * Cos(phiL);
261  y = rB * Sin(phiL);
262  glVertex3f(x, y, zB);
263  glVertex3f(x, y, zF);
264  x = rB * Cos(phiU);
265  y = rB * Sin(phiU);
266  glVertex3f(x, y, zB);
267  glVertex3f(x, y, zF);
268  }
269  }
270 }
271 
272 ////////////////////////////////////////////////////////////////////////////////
273 /// Draw frame reading eta, phi axis.
274 
276 {
277  if (rnrCtx.Highlight() || rnrCtx.Selection() || rnrCtx.IsDrawPassOutlineLine()) return;
278 
279  Bool_t transparent_p = fM->fFrameTransparency > 0;
280 
281  if (transparent_p)
282  {
283  glPushAttrib(GL_ENABLE_BIT | GL_DEPTH_BUFFER_BIT);
284 
285  glDepthMask(GL_FALSE);
286  glEnable(GL_BLEND);
287  glBlendFunc(GL_SRC_ALPHA, GL_ONE_MINUS_SRC_ALPHA);
288 
290  }
291 
292  TGLCapabilitySwitch lights_off(GL_LIGHTING, kFALSE);
293 
295  glBegin(GL_LINES);
296 
297  Float_t etaMin = fM->GetEtaMin();
298  Float_t etaMax = fM->GetEtaMax();
299 
300  Float_t transF = fM->GetTransitionEtaForward();
301  Float_t transB = fM->GetTransitionEtaBackward();
302  if (fM->GetRnrBarrelFrame() && (etaMin < transF && etaMax > transB))
303  {
305  }
306 
307  if (fM->GetRnrEndCapFrame() && (etaMax > transF || etaMin < transB))
308  {
310  }
311 
312  glEnd();
313 
314  if (transparent_p)
315  {
316  glPopAttrib();
317  }
318 }
319 
320 ////////////////////////////////////////////////////////////////////////////////
321 /// Render box with given points.
322 /// ~~~ {.cpp}
323 /// z
324 /// |
325 /// |
326 /// |________y
327 /// / 6-------7
328 /// / /| /|
329 /// x 5-------4 |
330 /// | 2-----|-3
331 /// |/ |/
332 /// 1-------0
333 ///~~~
334 
335 void TEveCalo3DGL::RenderBox(const Float_t pnts[8]) const
336 {
337  const Float_t *p = pnts;
338  Float_t cross[3];
339 
340  // bottom: 0123
341  glBegin(GL_POLYGON);
342  CrossProduct(p+3, p+9, p, cross);
343  glNormal3fv(cross);
344  glVertex3fv(p);
345  glVertex3fv(p+3);
346  glVertex3fv(p+6);
347  glVertex3fv(p+9);
348  glEnd();
349  // top: 7654
350  glBegin(GL_POLYGON);
351  CrossProduct(p+21, p+15, p+12, cross);
352  glNormal3fv(cross);
353  glVertex3fv(p+21);
354  glVertex3fv(p+18);
355  glVertex3fv(p+15);
356  glVertex3fv(p+12);
357  glEnd();
358  // back: 0451
359  glBegin(GL_POLYGON);
360  CrossProduct(p+12, p+3, p, cross);
361  glNormal3fv(cross);
362  glVertex3fv(p);
363  glVertex3fv(p+12);
364  glVertex3fv(p+15);
365  glVertex3fv(p+3);
366  glEnd();
367  //front : 3267
368  glBegin(GL_POLYGON);
369  CrossProduct(p+6, p+21, p+9, cross);
370  glNormal3fv(cross);
371  glVertex3fv(p+9);
372  glVertex3fv(p+6);
373  glVertex3fv(p+18);
374  glVertex3fv(p+21);
375  glEnd();
376  // left: 0374
377  glBegin(GL_POLYGON);
378  CrossProduct(p+21, p, p+9, cross);
379  glNormal3fv(cross);
380  glVertex3fv(p);
381  glVertex3fv(p+9);
382  glVertex3fv(p+21);
383  glVertex3fv(p+12);
384  glEnd();
385  // right: 1562
386  glBegin(GL_POLYGON);
387  CrossProduct(p+15, p+6, p+3, cross);
388  glNormal3fv(cross);
389  glVertex3fv(p+3);
390  glVertex3fv(p+15);
391  glVertex3fv(p+18);
392  glVertex3fv(p+6);
393  glEnd();
394 }
395 
396 ////////////////////////////////////////////////////////////////////////////////
397 /// Render barrel cell.
398 
399 void TEveCalo3DGL::RenderBarrelCell(const TEveCaloData::CellGeom_t &cellData, Float_t towerH, Float_t& offset ) const
400 {
401  using namespace TMath;
402 
403  Float_t r1 = fM->GetBarrelRadius() + offset;
404  Float_t r2 = r1 + towerH*Sin(cellData.ThetaMin());
405  Float_t z1In, z1Out, z2In, z2Out;
406 
407  z1In = r1/Tan(cellData.ThetaMax());
408  z1Out = r2/Tan(cellData.ThetaMax());
409  z2In = r1/Tan(cellData.ThetaMin());
410  z2Out = r2/Tan(cellData.ThetaMin());
411 
412  Float_t cos1 = Cos(cellData.PhiMin());
413  Float_t sin1 = Sin(cellData.PhiMin());
414  Float_t cos2 = Cos(cellData.PhiMax());
415  Float_t sin2 = Sin(cellData.PhiMax());
416 
417  Float_t box[24];
418  Float_t* pnts = box;
419  // 0
420  pnts[0] = r1*cos2;
421  pnts[1] = r1*sin2;
422  pnts[2] = z1In;
423  pnts += 3;
424  // 1
425  pnts[0] = r1*cos1;
426  pnts[1] = r1*sin1;
427  pnts[2] = z1In;
428  pnts += 3;
429  // 2
430  pnts[0] = r1*cos1;
431  pnts[1] = r1*sin1;
432  pnts[2] = z2In;
433  pnts += 3;
434  // 3
435  pnts[0] = r1*cos2;
436  pnts[1] = r1*sin2;
437  pnts[2] = z2In;
438  pnts += 3;
439  //---------------------------------------------------
440  // 4
441  pnts[0] = r2*cos2;
442  pnts[1] = r2*sin2;
443  pnts[2] = z1Out;
444  pnts += 3;
445  // 5
446  pnts[0] = r2*cos1;
447  pnts[1] = r2*sin1;
448  pnts[2] = z1Out;
449  pnts += 3;
450  // 6
451  pnts[0] = r2*cos1;
452  pnts[1] = r2*sin1;
453  pnts[2] = z2Out;
454  pnts += 3;
455  // 7
456  pnts[0] = r2*cos2;
457  pnts[1] = r2*sin2;
458  pnts[2] = z2Out;
459 
460  RenderBox(box);
461 
462  offset += towerH*Sin(cellData.ThetaMin());
463 
464 }// end RenderBarrelCell
465 
466 ////////////////////////////////////////////////////////////////////////////////
467 /// Render an endcap cell.
468 
469 void TEveCalo3DGL::RenderEndCapCell(const TEveCaloData::CellGeom_t &cellData, Float_t towerH, Float_t& offset ) const
470 {
471  using namespace TMath;
472  Float_t z1, r1In, r1Out, z2, r2In, r2Out;
473 
474  z1 = (cellData.EtaMin()<0) ? fM->fEndCapPosB : fM->fEndCapPosF;
475  z2 = z1 + TMath::Sign(towerH, cellData.EtaMin());
476 
477  r1In = z1*Tan(cellData.ThetaMin());
478  r2In = z2*Tan(cellData.ThetaMin());
479  r1Out = z1*Tan(cellData.ThetaMax());
480  r2Out = z2*Tan(cellData.ThetaMax());
481 
482  Float_t cos2 = Cos(cellData.PhiMin());
483  Float_t sin2 = Sin(cellData.PhiMin());
484  Float_t cos1 = Cos(cellData.PhiMax());
485  Float_t sin1 = Sin(cellData.PhiMax());
486 
487  Float_t box[24];
488  Float_t* pnts = box;
489  // 0
490  pnts[0] = r1In*cos1;
491  pnts[1] = r1In*sin1;
492  pnts[2] = z1;
493  pnts += 3;
494  // 1
495  pnts[0] = r1In*cos2;
496  pnts[1] = r1In*sin2;
497  pnts[2] = z1;
498  pnts += 3;
499  // 2
500  pnts[0] = r2In*cos2;
501  pnts[1] = r2In*sin2;
502  pnts[2] = z2;
503  pnts += 3;
504  // 3
505  pnts[0] = r2In*cos1;
506  pnts[1] = r2In*sin1;
507  pnts[2] = z2;
508  pnts += 3;
509  //---------------------------------------------------
510  // 4
511  pnts[0] = r1Out*cos1;
512  pnts[1] = r1Out*sin1;
513  pnts[2] = z1;
514  pnts += 3;
515  // 5
516  pnts[0] = r1Out*cos2;
517  pnts[1] = r1Out*sin2;
518  pnts[2] = z1;
519  pnts += 3;
520  // 6
521  pnts[0] = r2Out*cos2;
522  pnts[1] = r2Out*sin2;
523  pnts[2] = z2;
524  pnts += 3;
525  // 7
526  pnts[0] = r2Out*cos1;
527  pnts[1] = r2Out*sin1;
528  pnts[2] = z2;
529 
530  RenderBox(box);
531 
532  if (z1 > 0)
533  offset += towerH * Cos(cellData.ThetaMin());
534  else
535  offset -= towerH * Cos(cellData.ThetaMin());
536 
537 } // end RenderEndCapCell
538 
539 ////////////////////////////////////////////////////////////////////////////////
540 /// GL rendering.
541 
543 {
544  if ( fM->GetValueIsColor()) fM->AssertPalette();
545 
546  // check if eta phi range has changed
547  if (fM->fCellIdCacheOK == kFALSE)
548  fM->BuildCellIdCache();
549 
550  glPushAttrib(GL_ENABLE_BIT | GL_LINE_BIT | GL_POLYGON_BIT);
551  glEnable(GL_LIGHTING);
552  glEnable(GL_NORMALIZE);
553  glEnable(GL_BLEND);
554  glBlendFunc(GL_SRC_ALPHA, GL_ONE_MINUS_SRC_ALPHA);
555 
556  TEveCaloData::CellData_t cellData;
557  Float_t towerH = 0;
558  Int_t tower = 0;
559  Int_t prevTower = -1;
560  Float_t offset = 0;
561  Int_t cellID = 0;
562 
563  if (rnrCtx.SecSelection()) glPushName(0);
564 
565  fOffset.assign(fM->fCellList.size(), 0);
566  for (TEveCaloData::vCellId_i i = fM->fCellList.begin(); i != fM->fCellList.end(); ++i)
567  {
568  fM->fData->GetCellData((*i), cellData);
569  tower = i->fTower;
570  if (tower != prevTower)
571  {
572  offset = 0;
573  prevTower = tower;
574  }
575  fOffset[cellID] = offset;
576  fM->SetupColorHeight(cellData.Value(fM->fPlotEt), (*i).fSlice, towerH);
577 
578  if (rnrCtx.SecSelection()) glLoadName(cellID);
579 
580  if ((cellData.Eta() > 0 && cellData.Eta() < fM->GetTransitionEtaForward()) ||
581  (cellData.Eta() < 0 && cellData.Eta() > fM->GetTransitionEtaBackward()))
582  {
583  RenderBarrelCell(cellData, towerH, offset);
584  }
585  else
586  {
587  RenderEndCapCell(cellData, towerH, offset);
588  }
589  ++cellID;
590  }
591 
592  if (rnrCtx.SecSelection()) glPopName();
593 
594  RenderGrid(rnrCtx);
595 
596  glPopAttrib();
597 }
598 
599 ////////////////////////////////////////////////////////////////////////////////
600 /// Draw polygons in highlight mode.
601 
602 void TEveCalo3DGL::DrawHighlight(TGLRnrCtx & rnrCtx, const TGLPhysicalShape* /*pshp*/, Int_t /*lvl*/) const
603 {
604  if (fM->fData->GetCellsSelected().empty() && fM->fData->GetCellsHighlighted().empty())
605  {
606  return;
607  }
608 
609  glPushAttrib(GL_ENABLE_BIT | GL_LINE_BIT | GL_POLYGON_BIT);
610  glDisable(GL_LIGHTING);
611  glDisable(GL_CULL_FACE);
612  glPolygonMode(GL_FRONT_AND_BACK, GL_LINE);
613 
616 
617  if (!fM->fData->GetCellsHighlighted().empty())
618  {
619  glColor4ubv(rnrCtx.ColorSet().Selection(3).CArr());
621  }
622  if (!fM->fData->GetCellsSelected().empty())
623  {
624  Float_t dr[2];
625  glGetFloatv(GL_DEPTH_RANGE,dr);
626  glColor4ubv(rnrCtx.ColorSet().Selection(1).CArr());
627  glDepthRange(dr[0], 0.8*dr[1]);
629  glDepthRange(dr[0], dr[1]);
630  }
631 
633  glPopAttrib();
634 }
635 
636 ////////////////////////////////////////////////////////////////////////////////
637 
639 {
640  TEveCaloData::CellData_t cellData;
641  Float_t towerH = 0;
642 
643  for (TEveCaloData::vCellId_i i = cells.begin(); i != cells.end(); i++)
644  {
645  fM->fData->GetCellData(*i, cellData);
646  fM->SetupColorHeight(cellData.Value(fM->fPlotEt), (*i).fSlice, towerH);
647 
648  // find tower with offsets
649  Float_t offset = 0;
650  for (Int_t j = 0; j < (Int_t) fM->fCellList.size(); ++j)
651  {
652  if (fM->fCellList[j].fTower == i->fTower && fM->fCellList[j].fSlice == i->fSlice )
653  {
654  offset = fOffset[j];
655  break;
656  }
657  }
658 
659  if (fM->CellInEtaPhiRng(cellData))
660  {
661  if (cellData.Eta() < fM->GetTransitionEtaForward() && cellData.Eta() > fM->GetTransitionEtaBackward())
662  RenderBarrelCell(cellData, towerH, offset);
663  else
664  RenderEndCapCell(cellData, towerH, offset);
665  }
666  }
667 }
668 
669 ////////////////////////////////////////////////////////////////////////////////
670 /// Processes tower selection.
671 /// Virtual function from TGLogicalShape. Called from TGLViewer.
672 
674 {
676  if (rec.GetN() > 1)
677  {
678  sel.push_back(fM->fCellList[rec.GetItem(1)]);
679  }
680  fM->fData->ProcessSelection(sel, rec);
681 }
const int nx
Definition: kalman.C:16
TEveRGBAPalette * AssertPalette()
Make sure the TEveRGBAPalette pointer is not null.
Definition: TEveCalo.cxx:378
vCellId_t & GetCellsSelected()
Definition: TEveCaloData.h:188
const UChar_t * CArr() const
Definition: TGLUtil.h:803
The TGLRnrCtx class aggregates data for a given redering context as needed by various parts of the RO...
Definition: TGLRnrCtx.h:40
virtual TAxis * GetPhiBins() const
Definition: TEveCaloData.h:221
TEveCalo3DGL()
Constructor.
std::vector< CellId_t >::iterator vCellId_i
Definition: TEveCaloData.h:147
const Double_t * v1
Definition: TArcBall.cxx:33
T1 Sign(T1 a, T2 b)
Definition: TMathBase.h:155
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:1702
void RenderGridBarrel() const
Render barrel grid.
TVector3 cross(const TVector3 &v1, const TVector3 &v2)
Definition: CsgOps.cxx:816
float Float_t
Definition: RtypesCore.h:53
const char Option_t
Definition: RtypesCore.h:62
void DrawSelectedCells(TEveCaloData::vCellId_t cells) const
virtual Bool_t ShouldDLCache(const TGLRnrCtx &rnrCtx) const
Decide if display-list should be used for this pass rendering, as determined by rnrCtx.
Definition: TGLObject.cxx:41
Float_t GetForwardEndCapPos() const
Definition: TEveCalo.h:102
Bool_t Selection() const
Definition: TGLRnrCtx.h:222
Bool_t CellInEtaPhiRng(TEveCaloData::CellData_t &) const
Returns true if given cell is in the ceta phi range.
Definition: TEveCalo.cxx:307
Float_t PhiMax() const
Definition: TEveCaloData.h:96
virtual Double_t GetBinLowEdge(Int_t bin) const
Return low edge of bin.
Definition: TAxis.cxx:489
Bool_t GetRnrBarrelFrame() const
Definition: TEveCalo.h:188
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
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
T etaMax()
Function providing the maximum possible value of pseudorapidity for a non-zero rho, in the Scalar type with the largest dynamic range.
Definition: etaMax.h:50
TArc * a
Definition: textangle.C:12
const Bool_t kFALSE
Definition: Rtypes.h:92
TEveCaloData::vCellId_t fCellList
Definition: TEveCalo.h:164
ClassImp(TEveCalo3DGL)
void box(Int_t pat, Double_t x1, Double_t y1, Double_t x2, Double_t y2)
Definition: fillpatterns.C:1
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
Float_t fEndCapPosB
Definition: TEveCalo.h:52
Bool_t fCellIdCacheOK
Definition: TEveCalo.h:40
Concrete physical shape - a GL drawable.
void CrossProduct(const Float_t a[3], const Float_t b[3], const Float_t c[3], Float_t out[3]) const
Calculate cross-product.
Float_t PhiMin() const
Definition: TEveCaloData.h:95
Double_t x[n]
Definition: legend1.C:17
static Float_t EtaToTheta(Float_t eta)
Bool_t Highlight() const
Definition: TGLRnrCtx.h:218
virtual void DirectDraw(TGLRnrCtx &rnrCtx) const
GL rendering.
const int ny
Definition: kalman.C:17
std::vector< Float_t > fOffset
Definition: TEveCalo3DGL.h:40
Base-class for direct OpenGL renderers.
Definition: TGLObject.h:21
void RenderGridEndCap() const
Render end cap grid.
Float_t fPlotEt
Definition: TEveCalo.h:54
std::vector< CellId_t > vCellId_t
Definition: TEveCaloData.h:146
void RenderGrid(TGLRnrCtx &rnrCtx) const
Draw frame reading eta, phi axis.
Float_t EtaMin() const
Definition: TEveCaloData.h:90
static UInt_t LockColor()
Prevent further color changes.
Definition: TGLUtil.cxx:1630
static UInt_t UnlockColor()
Allow color changes.
Definition: TGLUtil.cxx:1638
virtual Bool_t SetModel(TObject *obj, const Option_t *opt=0)
Set model object.
char * out
Definition: TBase64.cxx:29
TObject * fExternalObj
first replica
Float_t GetEtaMin() const
Definition: TEveCalo.h:136
void RenderBarrelCell(const TEveCaloData::CellGeom_t &cell, Float_t towerH, Float_t &offset) const
Render barrel cell.
Float_t ThetaMin() const
Definition: TEveCaloData.h:100
Float_t ThetaMax() const
Definition: TEveCaloData.h:101
virtual void BuildCellIdCache()
Build list of drawn cell IDs. See TEveCalo3DGL::DirectDraw().
Definition: TEveCalo.cxx:459
TEveCaloData * GetData() const
Definition: TEveCalo.h:86
ROOT::R::TRInterface & r
Definition: Object.C:4
Class to manage histogram axis.
Definition: TAxis.h:36
void ProcessSelection(vCellId_t &sel_cells, TGLSelectRecord &rec)
Process newly selected cells with given select-record.
Int_t GetNbins() const
Definition: TAxis.h:125
virtual Double_t GetBinUpEdge(Int_t bin) const
Return up edge of bin.
Definition: TAxis.cxx:499
Standard selection record including information about containing scene and details ob out selected ob...
unsigned int r1[N_CITIES]
Definition: simanTSP.cxx:321
Float_t GetPhiMin() const
Definition: TEveCalo.h:144
Bool_t GetRnrEndCapFrame() const
Definition: TEveCalo.h:187
OpenGL renderer class for TEveCalo3D.
Definition: TEveCalo3DGL.h:20
void RenderBox(const Float_t pnts[8]) const
Render box with given points.
Float_t GetTransitionEtaBackward() const
Get transition eta between barrel and backward end-cap cells.
Definition: TEveCalo.cxx:227
virtual void SetBBox()
Set bounding box.
Cell geometry inner structure.
Definition: TEveCaloData.h:72
virtual void GetCellData(const CellId_t &id, CellData_t &data) const =0
Bool_t SecSelection() const
Definition: TGLRnrCtx.h:224
Bool_t IsDrawPassOutlineLine() const
Definition: TGLRnrCtx.h:207
Char_t fFrameTransparency
Definition: TEveCalo.h:171
vCellId_t & GetCellsHighlighted()
Definition: TEveCaloData.h:189
Float_t Value(Bool_t) const
Return energy value associated with the cell, usually Et.
virtual void DrawHighlight(TGLRnrCtx &rnrCtx, const TGLPhysicalShape *ps, Int_t lvl=-1) const
Draw polygons in highlight mode.
Double_t Cos(Double_t)
Definition: TMath.h:424
Float_t GetTransitionEtaForward() const
Get transition eta between barrel and forward end-cap cells.
Definition: TEveCalo.cxx:209
double Double_t
Definition: RtypesCore.h:55
UInt_t GetItem(Int_t i) const
Float_t GetPhiMax() const
Definition: TEveCalo.h:145
Double_t y[n]
Definition: legend1.C:17
static Bool_t IsU1IntervalContainedByMinMax(Float_t minM, Float_t maxM, Float_t minQ, Float_t maxQ)
Return true if interval Q is contained within interval M for U1 variables.
Definition: TEveUtil.cxx:346
Bool_t fMultiColor
Definition: TGLObject.h:28
Float_t GetEtaMax() const
Definition: TEveCalo.h:137
Mother of all ROOT objects.
Definition: TObject.h:58
TGLColor & Selection(Int_t i)
Definition: TGLUtil.h:857
TEveCalo3D * fM
Definition: TEveCalo3DGL.h:38
virtual TAxis * GetEtaBins() const
Definition: TEveCaloData.h:218
Float_t GetBarrelRadius() const
Definition: TEveCalo.h:99
Float_t GetFrameWidth() const
Definition: TEveCalo.h:181
Double_t Sin(Double_t)
Definition: TMath.h:421
Float_t fEndCapPosF
Definition: TEveCalo.h:51
void RenderEndCapCell(const TEveCaloData::CellGeom_t &cell, Float_t towerH, Float_t &offset) const
Render an endcap cell.
Cell data inner structure.
Definition: TEveCaloData.h:114
static Float_t LineWidth()
Get the line-width, taking the global scaling into account.
Definition: TGLUtil.cxx:1904
Visualization of a calorimeter event data in 3D.
Definition: TEveCalo.h:156
virtual void ProcessSelection(TGLRnrCtx &rnrCtx, TGLSelectRecord &rec)
Processes tower selection.
TGLColorSet & ColorSet()
Return reference to current color-set (top of the stack).
Definition: TGLRnrCtx.cxx:278
Float_t GetBackwardEndCapPos() const
Definition: TEveCalo.h:103
Color_t fFrameColor
Definition: TEveCalo.h:170
const Bool_t kTRUE
Definition: Rtypes.h:91
TObject * obj
TEveCaloData * fData
Definition: TEveCalo.h:39
Bool_t GetValueIsColor() const
Definition: TEveCalo.h:128
Int_t GetN() const
void SetupColorHeight(Float_t value, Int_t slice, Float_t &height) const
Set color and height for a given value and slice using slice color or TEveRGBAPalette.
Definition: TEveCalo.cxx:415
Float_t Eta() const
Definition: TEveCaloData.h:92
Double_t Tan(Double_t)
Definition: TMath.h:427
unsigned int r2[N_CITIES]
Definition: simanTSP.cxx:322
virtual Bool_t ShouldDLCache(const TGLRnrCtx &rnrCtx) const
Override from TGLObject.