Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
TEveCalo.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 "TEveCalo.h"
13 #include "TEveCaloData.h"
14 #include "TEveProjections.h"
15 #include "TEveProjectionManager.h"
16 #include "TEveRGBAPalette.h"
17 #include "TEveText.h"
18 #include "TEveTrans.h"
19 
20 #include "TClass.h"
21 #include "TMathBase.h"
22 #include "TMath.h"
23 #include "TAxis.h"
24 
25 #include "TGLUtil.h"
26 
27 #include <cassert>
28 
29 /** \class TEveCaloViz
30 \ingroup TEve
31 Base class for calorimeter data visualization.
32 See TEveCalo2D and TEveCalo3D for concrete implementations.
33 */
34 
35 ClassImp(TEveCaloViz);
36 
37 ////////////////////////////////////////////////////////////////////////////////
38 
39 TEveCaloViz::TEveCaloViz(TEveCaloData* data, const char* n, const char* t) :
40  TEveElement(),
41  TNamed(n, t),
42  TEveProjectable(),
43 
44  fData(0),
45  fCellIdCacheOK(kFALSE),
46 
47  fEtaMin(-10),
48  fEtaMax(10),
49 
50  fPhi(0.),
51  fPhiOffset(TMath::Pi()),
52 
53  fAutoRange(kTRUE),
54 
55  fBarrelRadius(-1.f),
56  fEndCapPosF(-1.f),
57  fEndCapPosB(-1.f),
58 
59  fPlotEt(kTRUE),
60 
61  fMaxTowerH(100),
62  fScaleAbs(kFALSE),
63  fMaxValAbs(100),
64 
65  fValueIsColor(kFALSE),
66  fPalette(0)
67 {
68  // Constructor.
69 
70  fPickable = kTRUE;
71  SetElementNameTitle(n, t);
72  SetData(data);
73 }
74 
75 ////////////////////////////////////////////////////////////////////////////////
76 /// Destructor.
77 
78 TEveCaloViz::~TEveCaloViz()
79 {
80  if (fPalette) fPalette->DecRefCount();
81 }
82 
83 ////////////////////////////////////////////////////////////////////////////////
84 /// Get threshold for given slice.
85 
86 Float_t TEveCaloViz::GetDataSliceThreshold(Int_t slice) const
87 {
88  return fData->RefSliceInfo(slice).fThreshold;
89 }
90 
91 ////////////////////////////////////////////////////////////////////////////////
92 /// Management of selection state and ownership of selected cell list
93 /// is done in TEveCaloData. This is a reason selection is forwarded to it.
94 
95 TEveElement* TEveCaloViz::ForwardSelection()
96 {
97  return fData;
98 }
99 
100 ////////////////////////////////////////////////////////////////////////////////
101 /// Management of selection state and ownership of selected cell list
102 /// is done in TEveCaloData. We still want GUI editor to display
103 /// concrete calo-viz object.
104 
105 TEveElement* TEveCaloViz::ForwardEdit()
106 {
107  return this;
108 }
109 
110 ////////////////////////////////////////////////////////////////////////////////
111 /// Set threshold for given slice.
112 
113 void TEveCaloViz::SetDataSliceThreshold(Int_t slice, Float_t val)
114 {
115  fData->SetSliceThreshold(slice, val);
116 }
117 
118 ////////////////////////////////////////////////////////////////////////////////
119 /// Get slice color from data.
120 
121 Color_t TEveCaloViz::GetDataSliceColor(Int_t slice) const
122 {
123  return fData->RefSliceInfo(slice).fColor;
124 }
125 
126 ////////////////////////////////////////////////////////////////////////////////
127 /// Set slice color in data.
128 
129 void TEveCaloViz::SetDataSliceColor(Int_t slice, Color_t col)
130 {
131  fData->SetSliceColor(slice, col);
132 }
133 
134 ////////////////////////////////////////////////////////////////////////////////
135 /// Set eta range.
136 
137 void TEveCaloViz::SetEta(Float_t l, Float_t u)
138 {
139  fEtaMin=l;
140  fEtaMax=u;
141 
142  InvalidateCellIdCache();
143 }
144 
145 ////////////////////////////////////////////////////////////////////////////////
146 /// Set E/Et plot.
147 
148 void TEveCaloViz::SetPlotEt(Bool_t isEt)
149 {
150  fPlotEt=isEt;
151  if (fPalette)
152  fPalette->SetLimits(0, TMath::CeilNint(GetMaxVal()));
153 
154  InvalidateCellIdCache();
155 }
156 
157 ////////////////////////////////////////////////////////////////////////////////
158 
159 Float_t TEveCaloViz::GetMaxVal() const
160 {
161  // Get maximum plotted value.
162 
163  return fData->GetMaxVal(fPlotEt);
164 
165 }
166 
167 ////////////////////////////////////////////////////////////////////////////////
168 /// Set phi range.
169 
170 void TEveCaloViz::SetPhiWithRng(Float_t phi, Float_t rng)
171 {
172  using namespace TMath;
173 
174  fPhi = phi;
175  fPhiOffset = rng;
176 
177  InvalidateCellIdCache();
178 }
179 
180 ////////////////////////////////////////////////////////////////////////////////
181 /// Get transition angle between barrel and end-cap cells, assuming fEndCapPosF = -fEndCapPosB.
182 
183 Float_t TEveCaloViz::GetTransitionTheta() const
184 {
185  return TMath::ATan(fBarrelRadius/fEndCapPosF);
186 }
187 
188 ////////////////////////////////////////////////////////////////////////////////
189 /// Get transition eta between barrel and end-cap cells, assuming fEndCapPosF = -fEndCapPosB.
190 
191 Float_t TEveCaloViz::GetTransitionEta() const
192 {
193  using namespace TMath;
194  Float_t t = GetTransitionTheta()*0.5f;
195  return -Log(Tan(t));
196 }
197 
198 ////////////////////////////////////////////////////////////////////////////////
199 /// Get transition angle between barrel and forward end-cap cells.
200 
201 Float_t TEveCaloViz::GetTransitionThetaForward() const
202 {
203  return TMath::ATan(fBarrelRadius/fEndCapPosF);
204 }
205 
206 ////////////////////////////////////////////////////////////////////////////////
207 /// Get transition eta between barrel and forward end-cap cells.
208 
209 Float_t TEveCaloViz::GetTransitionEtaForward() const
210 {
211  using namespace TMath;
212  Float_t t = GetTransitionThetaForward()*0.5f;
213  return -Log(Tan(t));
214 }
215 
216 ////////////////////////////////////////////////////////////////////////////////
217 /// Get transition angle between barrel and backward end-cap cells.
218 
219 Float_t TEveCaloViz::GetTransitionThetaBackward() const
220 {
221  return TMath::ATan(fBarrelRadius/fEndCapPosB);
222 }
223 
224 ////////////////////////////////////////////////////////////////////////////////
225 /// Get transition eta between barrel and backward end-cap cells.
226 
227 Float_t TEveCaloViz::GetTransitionEtaBackward() const
228 {
229  using namespace TMath;
230  Float_t t = GetTransitionThetaBackward()*0.5f;
231  //negative theta means negative eta
232  return Log(-Tan(t));
233 }
234 
235 
236 ////////////////////////////////////////////////////////////////////////////////
237 /// Set calorimeter event data.
238 
239 void TEveCaloViz::SetData(TEveCaloData* data)
240 {
241 
242  if (data == fData) return;
243  if (fData) fData->RemoveElement(this);
244  fData = data;
245  if (fData)
246  {
247  fData->AddElement(this);
248  DataChanged();
249  }
250 }
251 
252 ////////////////////////////////////////////////////////////////////////////////
253 /// Update setting and cache on data changed.
254 /// Called from TEvecaloData::BroadcastDataChange()
255 
256 void TEveCaloViz::DataChanged()
257 {
258  Double_t min, max, delta;
259 
260  fData->GetEtaLimits(min, max);
261  if (fAutoRange) {
262  fEtaMin = min;
263  fEtaMax = max;
264  } else {
265  if (fEtaMin < min) fEtaMin = min;
266  if (fEtaMax > max) fEtaMax = max;
267  }
268 
269  fData->GetPhiLimits(min, max);
270  delta = 0.5*(max - min);
271  if (fAutoRange || fPhi < min || fPhi > max) {
272  fPhi = 0.5*(max + min);
273  fPhiOffset = delta;
274  } else {
275  if (fPhiOffset > delta) fPhiOffset = delta;
276  }
277 
278  if (fPalette)
279  {
280  Int_t hlimit = TMath::CeilNint(GetMaxVal());
281  fPalette->SetLimits(0, hlimit);
282  fPalette->SetMin(0);
283  fPalette->SetMax(hlimit);
284  }
285 
286  InvalidateCellIdCache();
287 }
288 
289 ////////////////////////////////////////////////////////////////////////////////
290 /// Assert cell id cache is ok.
291 /// Returns true if the cache has been updated.
292 
293 Bool_t TEveCaloViz::AssertCellIdCache() const
294 {
295  TEveCaloViz* cv = const_cast<TEveCaloViz*>(this);
296  if (!fCellIdCacheOK) {
297  cv->BuildCellIdCache();
298  return kTRUE;
299  } else {
300  return kFALSE;
301  }
302 }
303 
304 ////////////////////////////////////////////////////////////////////////////////
305 /// Returns true if given cell is in the ceta phi range.
306 
307 Bool_t TEveCaloViz::CellInEtaPhiRng(TEveCaloData::CellData_t& cellData) const
308 {
309  if (cellData.EtaMin() >= fEtaMin && cellData.EtaMax() <= fEtaMax)
310  {
311  if (TEveUtil::IsU1IntervalContainedByMinMax
312  (fPhi-fPhiOffset, fPhi+fPhiOffset, cellData.PhiMin(), cellData.PhiMax()))
313  return kTRUE;
314  }
315  return kFALSE;
316 }
317 
318 ////////////////////////////////////////////////////////////////////////////////
319 /// Assign parameters from given model.
320 
321 void TEveCaloViz::AssignCaloVizParameters(TEveCaloViz* m)
322 {
323  SetData(m->fData);
324 
325  fEtaMin = m->fEtaMin;
326  fEtaMax = m->fEtaMax;
327 
328  fPhi = m->fPhi;
329  fPhiOffset = m->fPhiOffset;
330 
331  fBarrelRadius = m->fBarrelRadius;
332  fEndCapPosF = m->fEndCapPosF;
333  fEndCapPosB = m->fEndCapPosB;
334 
335  if (m->fPalette)
336  {
337  TEveRGBAPalette& mp = * m->fPalette;
338  if (fPalette) fPalette->DecRefCount();
339  fPalette = new TEveRGBAPalette(mp.GetMinVal(), mp.GetMaxVal(), mp.GetInterpolate());
340  fPalette->SetDefaultColor(mp.GetDefaultColor());
341  }
342 }
343 
344 ////////////////////////////////////////////////////////////////////////////////
345 /// Set TEveRGBAPalette object pointer.
346 
347 void TEveCaloViz::SetPalette(TEveRGBAPalette* p)
348 {
349  if ( fPalette == p) return;
350  if (fPalette) fPalette->DecRefCount();
351  fPalette = p;
352  if (fPalette) fPalette->IncRefCount();
353 }
354 
355 ////////////////////////////////////////////////////////////////////////////////
356 /// Get transformation factor from E/Et to height
357 
358 Float_t TEveCaloViz::GetValToHeight() const
359 {
360  if (fScaleAbs)
361  {
362  return fMaxTowerH/fMaxValAbs;
363  }
364  else
365  {
366  if (fData->Empty())
367  return 1;
368 
369  return fMaxTowerH/fData->GetMaxVal(fPlotEt);
370  }
371 }
372 
373 ////////////////////////////////////////////////////////////////////////////////
374 /// Make sure the TEveRGBAPalette pointer is not null.
375 /// If it is not set, a new one is instantiated and the range is set
376 /// to current min/max signal values.
377 
378 TEveRGBAPalette* TEveCaloViz::AssertPalette()
379 {
380  if (fPalette == 0) {
381  fPalette = new TEveRGBAPalette;
382  fPalette->SetDefaultColor((Color_t)4);
383 
384  Int_t hlimit = TMath::CeilNint(GetMaxVal());
385  fPalette->SetLimits(0, hlimit);
386  fPalette->SetMin(0);
387  fPalette->SetMax(hlimit);
388 
389  }
390  return fPalette;
391 }
392 
393 ////////////////////////////////////////////////////////////////////////////////
394 /// Paint this object. Only direct rendering is supported.
395 
396 void TEveCaloViz::Paint(Option_t* /*option*/)
397 {
398  if (fData)
399  {
400  PaintStandard(this);
401  }
402 }
403 
404 ////////////////////////////////////////////////////////////////////////////////
405 /// Virtual from TEveProjectable, returns TEveCalo2D class.
406 
407 TClass* TEveCaloViz::ProjectedClass(const TEveProjection*) const
408 {
409  return TEveCalo2D::Class();
410 }
411 
412 ////////////////////////////////////////////////////////////////////////////////
413 /// Set color and height for a given value and slice using slice color or TEveRGBAPalette.
414 
415 void TEveCaloViz::SetupColorHeight(Float_t value, Int_t slice, Float_t& outH) const
416 {
417  if (fValueIsColor)
418  {
419  outH = GetValToHeight()*fData->GetMaxVal(fPlotEt);
420  UChar_t c[4];
421  fPalette->ColorFromValue((Int_t)value, c);
422  c[3] = fData->GetSliceTransparency(slice);
423  TGLUtil::Color4ubv(c);
424  }
425  else
426  {
427  TGLUtil::ColorTransparency(fData->GetSliceColor(slice), fData->GetSliceTransparency(slice));
428  outH = GetValToHeight()*value;
429  }
430 }
431 
432 /** \class TEveCalo3D
433 \ingroup TEve
434 Visualization of a calorimeter event data in 3D.
435 */
436 
437 ClassImp(TEveCalo3D);
438 
439 ////////////////////////////////////////////////////////////////////////////////
440 /// Constructor.
441 
442 TEveCalo3D::TEveCalo3D(TEveCaloData* d, const char* n, const char* t):
443  TEveCaloViz(d, n, t),
444 
445  fRnrEndCapFrame (kTRUE),
446  fRnrBarrelFrame (kTRUE),
447  fFrameWidth (0.5),
448  fFrameColor (kGray+1),
449  fFrameTransparency (80)
450 {
451  fCanEditMainColor = kTRUE;
452  fCanEditMainTransparency = kTRUE;
453  fMainColorPtr = &fFrameColor;
454 }
455 
456 ////////////////////////////////////////////////////////////////////////////////
457 /// Build list of drawn cell IDs. See TEveCalo3DGL::DirectDraw().
458 
459 void TEveCalo3D::BuildCellIdCache()
460 {
461  fCellList.clear();
462 
463  fData->GetCellList(GetEta(), GetEtaRng(), GetPhi(), GetPhiRng(), fCellList);
464  fCellIdCacheOK = kTRUE;
465 }
466 
467 ////////////////////////////////////////////////////////////////////////////////
468 /// Fill bounding-box information of the base-class TAttBBox (virtual method).
469 /// If member 'TEveFrameBox* fFrame' is set, frame's corners are used as bbox.
470 
471 void TEveCalo3D::ComputeBBox()
472 {
473  BBoxInit();
474 
475  Float_t th = (fData) ? GetValToHeight() * fData->GetMaxVal(fPlotEt) : 0;
476 
477  fBBox[0] = -fBarrelRadius - th;
478  fBBox[1] = fBarrelRadius + th;
479  fBBox[2] = fBBox[0];
480  fBBox[3] = fBBox[1];
481  fBBox[4] = fEndCapPosB - th;
482  fBBox[5] = fEndCapPosF + th;
483 }
484 
485 /** \class TEveCalo2D
486 \ingroup TEve
487 Visualization of a calorimeter event data in 2D.
488 */
489 
490 ClassImp(TEveCalo2D);
491 
492 ////////////////////////////////////////////////////////////////////////////////
493 /// Constructor.
494 
495 TEveCalo2D::TEveCalo2D(const char* n, const char* t):
496  TEveCaloViz(0, n, t),
497  TEveProjected(),
498  fOldProjectionType(TEveProjection::kPT_Unknown),
499  fMaxESumBin( 0),
500  fMaxEtSumBin(0)
501 {
502 }
503 
504 ////////////////////////////////////////////////////////////////////////////////
505 /// Destructor.
506 
507 TEveCalo2D::~TEveCalo2D()
508 {
509  TEveCaloData::vCellId_t* cids;
510  UInt_t n;
511 
512  // clear selected cell ids
513  n = fCellListsSelected.size();
514  for(UInt_t i = 0; i < n; ++i) {
515  cids = fCellListsSelected[i];
516  if (cids) {
517  cids->clear(); delete cids;
518  }
519  }
520  fCellListsSelected.clear();
521 
522  // clear all cell dds
523  n = fCellLists.size();
524  for(UInt_t i = 0; i < n; ++i) {
525  cids = fCellLists[i];
526  if (cids) {
527  cids->clear(); delete cids;
528  }
529  }
530  fCellLists.clear();
531 }
532 
533 ////////////////////////////////////////////////////////////////////////////////
534 /// This is virtual method from base-class TEveProjected.
535 
536 void TEveCalo2D::UpdateProjection()
537 {
538  if (fManager->GetProjection()->GetType() != fOldProjectionType)
539  {
540  fCellIdCacheOK=kFALSE;
541  fOldProjectionType = fManager->GetProjection()->GetType();
542  }
543  ComputeBBox();
544 }
545 
546 ////////////////////////////////////////////////////////////////////////////////
547 /// Set projection manager and model object.
548 
549 void TEveCalo2D::SetProjection(TEveProjectionManager* mng, TEveProjectable* model)
550 {
551  TEveProjected::SetProjection(mng, model);
552  TEveCaloViz* viz = dynamic_cast<TEveCaloViz*>(model);
553  AssignCaloVizParameters(viz);
554 }
555 
556 ////////////////////////////////////////////////////////////////////////////////
557 /// Build lists of drawn cell IDs. See TEveCalo2DGL::DirecDraw().
558 
559 void TEveCalo2D::BuildCellIdCache()
560 {
561  // clear old cache
562  for (vBinCells_i it = fCellLists.begin(); it != fCellLists.end(); it++)
563  {
564  if (*it)
565  {
566  (*it)->clear();
567  delete *it;
568  }
569  }
570  fCellLists.clear();
571  fCellLists.push_back(0);
572 
573  TEveProjection::EPType_e pt = fManager->GetProjection()->GetType();
574  TEveCaloData::vCellId_t* clv; // ids per phi bin in r-phi projection else ids per eta bins in rho-z projection
575 
576  Bool_t isRPhi = (pt == TEveProjection::kPT_RPhi);
577 
578  const TAxis* axis = isRPhi ? fData->GetPhiBins() : fData->GetEtaBins();
579  Int_t nBins = axis->GetNbins();
580 
581  Float_t min, max;
582  if (isRPhi)
583  {
584  min = GetPhiMin() - fData->GetEps();
585  max = GetPhiMax() + fData->GetEps();
586  for (Int_t ibin = 1; ibin <= nBins; ++ibin) {
587  clv = 0;
588  if ( TEveUtil::IsU1IntervalOverlappingByMinMax
589  (min, max, axis->GetBinLowEdge(ibin), axis->GetBinUpEdge(ibin)))
590  {
591  clv = new TEveCaloData::vCellId_t();
592  fData->GetCellList(GetEta(), GetEtaRng(), axis->GetBinCenter(ibin), axis->GetBinWidth(ibin), *clv);
593  if (!clv->size()) {
594  delete clv; clv = 0;
595  }
596  }
597  fCellLists.push_back(clv);
598  }
599  }
600  else
601  {
602  min = GetEtaMin() - fData->GetEps();
603  max = GetEtaMax() + fData->GetEps();
604  for (Int_t ibin = 1; ibin <= nBins; ++ibin) {
605  clv = 0;
606  Float_t low = axis->GetBinLowEdge(ibin);
607  Float_t up = axis->GetBinUpEdge(ibin) ;
608  if (low >= min && up <= max)
609  {
610  clv = new TEveCaloData::vCellId_t();
611  fData->GetCellList(axis->GetBinCenter(ibin), axis->GetBinWidth(ibin), fPhi, GetPhiRng(), *clv);
612  if (!clv->size()) {
613  delete clv; clv = 0;
614  }
615  }
616  fCellLists.push_back(clv);
617  }
618  }
619 
620  // cache max bin sum for auto scale
621  if (!fScaleAbs)
622  {
623  fMaxESumBin = 0;
624  fMaxEtSumBin = 0;
625  Float_t sumE = 0;
626  Float_t sumEt = 0;
627  TEveCaloData::CellData_t cellData;
628  for (Int_t ibin = 1; ibin <= nBins; ++ibin) {
629  TEveCaloData::vCellId_t* cids = fCellLists[ibin];
630  if (cids)
631  {
632  sumE = 0; sumEt = 0;
633  for (TEveCaloData::vCellId_i it = cids->begin(); it != cids->end(); it++)
634  {
635  fData->GetCellData(*it, cellData);
636  sumE += cellData.Value(kFALSE);
637  sumEt += cellData.Value(kTRUE);
638  }
639  fMaxESumBin = TMath::Max(fMaxESumBin, sumE);
640  fMaxEtSumBin = TMath::Max(fMaxEtSumBin, sumEt);
641  }
642  }
643  ComputeBBox();
644  }
645 
646  fCellIdCacheOK= kTRUE;
647 }
648 
649 ////////////////////////////////////////////////////////////////////////////////
650 /// Sort selected cells in eta or phi bins for selection and highlight.
651 
652 void TEveCalo2D::CellSelectionChanged()
653 {
654  CellSelectionChangedInternal(fData->GetCellsSelected(), fCellListsSelected);
655  CellSelectionChangedInternal(fData->GetCellsHighlighted(), fCellListsHighlighted);
656 }
657 
658 ////////////////////////////////////////////////////////////////////////////////
659 /// Sort selected cells in eta or phi bins.
660 
661 void TEveCalo2D::CellSelectionChangedInternal(TEveCaloData::vCellId_t& inputCells, std::vector<TEveCaloData::vCellId_t*>& outputCellLists)
662 {
663  Bool_t isRPhi = (fManager->GetProjection()->GetType() == TEveProjection::kPT_RPhi);
664  const TAxis* axis = isRPhi ? fData->GetPhiBins() : fData->GetEtaBins();
665 
666  // clear old cache
667  for (vBinCells_i it = outputCellLists.begin(); it != outputCellLists.end(); it++)
668  {
669  if (*it)
670  {
671  (*it)->clear();
672  delete *it;
673  }
674  }
675  outputCellLists.clear();
676  UInt_t nBins = axis->GetNbins();
677  outputCellLists.resize(nBins+1);
678  for (UInt_t b = 0; b <= nBins; ++b)
679  outputCellLists[b] = 0;
680 
681  for(UInt_t bin = 1; bin <= nBins; ++bin)
682  {
683  TEveCaloData::vCellId_t* idsInBin = fCellLists[bin];
684  if (!idsInBin)
685  continue;
686 
687  for (TEveCaloData::vCellId_i i = idsInBin->begin(); i != idsInBin->end(); i++)
688  {
689  for (TEveCaloData::vCellId_i j = inputCells.begin(); j != inputCells.end(); j++)
690  {
691  if( (*i).fTower == (*j).fTower && (*i).fSlice == (*j).fSlice)
692  {
693  if (!outputCellLists[bin])
694  outputCellLists[bin] = new TEveCaloData::vCellId_t();
695 
696  outputCellLists[bin]->push_back(TEveCaloData::CellId_t((*i).fTower, (*i).fSlice, (*i).fFraction));
697  }
698  }
699  }
700  }
701 }
702 
703 ////////////////////////////////////////////////////////////////////////////////
704 /// Set absolute scale in projected calorimeter.
705 
706 void TEveCalo2D::SetScaleAbs(Bool_t sa)
707 {
708  TEveCaloViz::SetScaleAbs(sa);
709  BuildCellIdCache();
710 }
711 
712 ////////////////////////////////////////////////////////////////////////////////
713 /// Virtual function of TEveCaloViz.
714 /// Get transformation factor from E/Et to height.
715 
716 Float_t TEveCalo2D::GetValToHeight() const
717 {
718  AssertCellIdCache();
719 
720  if (fScaleAbs)
721  {
722  return fMaxTowerH/fMaxValAbs;
723  }
724  else
725  {
726  if (fData->Empty())
727  return 1;
728 
729  if (fPlotEt)
730  return fMaxTowerH/fMaxEtSumBin;
731  else
732  return fMaxTowerH/fMaxESumBin;
733  }
734 }
735 
736 ////////////////////////////////////////////////////////////////////////////////
737 /// Fill bounding-box information of the base-class TAttBBox (virtual method).
738 /// If member 'TEveFrameBox* fFrame' is set, frame's corners are used as bbox.
739 
740 void TEveCalo2D::ComputeBBox()
741 {
742  BBoxZero();
743 
744  Float_t x, y, z;
745  Float_t th = fMaxTowerH ;
746  Float_t r = fBarrelRadius + th;
747 
748  x = r, y = 0, z = 0;
749  fManager->GetProjection()->ProjectPoint(x, y, z, fDepth);
750  BBoxCheckPoint(x, y, z);
751  x = -r, y = 0, z = 0;
752  fManager->GetProjection()->ProjectPoint(x, y, z, fDepth);
753  BBoxCheckPoint(x, y, z);
754 
755  x = 0, y = 0, z = fEndCapPosF + th;
756  fManager->GetProjection()->ProjectPoint(x, y, z, fDepth);
757  BBoxCheckPoint(x, y, z);
758  x = 0, y = 0, z = fEndCapPosB - th;
759  fManager->GetProjection()->ProjectPoint(x, y, z, fDepth);
760  BBoxCheckPoint(x, y, z);
761 
762  x = 0, y = r, z = 0;
763  fManager->GetProjection()->ProjectPoint(x, y, z, fDepth);
764  BBoxCheckPoint(x, y, z);
765  x = 0, y = -r, z = 0;
766  fManager->GetProjection()->ProjectPoint(x, y, z, fDepth);
767  BBoxCheckPoint(x, y, z);
768 }
769 
770 
771 /** \class TEveCaloLego
772 \ingroup TEve
773 Visualization of calorimeter data as eta/phi histogram.
774 */
775 
776 ClassImp(TEveCaloLego);
777 
778 ////////////////////////////////////////////////////////////////////////////////
779 /// Constructor.
780 
781 TEveCaloLego::TEveCaloLego(TEveCaloData* d, const char* n, const char* t):
782  TEveCaloViz(d, n, t),
783 
784  fFontColor(-1),
785  fGridColor(-1),
786  fPlaneColor(kRed-5),
787  fPlaneTransparency(60),
788 
789  fNZSteps(6),
790  fZAxisStep(0.f),
791 
792  fAutoRebin(kTRUE),
793 
794  fPixelsPerBin(12),
795  fNormalizeRebin(kFALSE),
796 
797  fProjection(kAuto),
798  f2DMode(kValSize),
799  fBoxMode(kBack),
800 
801  fDrawHPlane(kFALSE),
802  fHPlaneVal(0),
803 
804  fHasFixedHeightIn2DMode(kFALSE),
805  fFixedHeightValIn2DMode(0.f),
806 
807  fDrawNumberCellPixels(18), // draw numbers on cell above 30 pixels
808  fCellPixelFontSize(12) // size of cell fonts in pixels
809 {
810  fMaxTowerH = 1;
811  SetElementNameTitle("TEveCaloLego", "TEveCaloLego");
812 }
813 
814 ////////////////////////////////////////////////////////////////////////////////
815 // Set data.
816 
817 void TEveCaloLego::SetData(TEveCaloData* data)
818 {
819  TEveCaloViz::SetData(data);
820 }
821 
822 ////////////////////////////////////////////////////////////////////////////////
823 /// Build list of drawn cell IDs. For more information see TEveCaloLegoGL:DirectDraw().
824 
825 void TEveCaloLego::BuildCellIdCache()
826 {
827  fCellList.clear();
828 
829  fData->GetCellList(GetEta(), GetEtaRng(), GetPhi(), GetPhiRng(), fCellList);
830  fCellIdCacheOK = kTRUE;
831 }
832 
833 ////////////////////////////////////////////////////////////////////////////////
834 /// Fill bounding-box information of the base-class TAttBBox (virtual method).
835 /// If member 'TEveFrameBox* fFrame' is set, frame's corners are used as bbox.
836 
837 void TEveCaloLego::ComputeBBox()
838 {
839 
840  // fBBox = Float_t[6] X(min,max), Y(min,max), Z(min,max)
841 
842  BBoxZero();
843 
844  Float_t ex = 1.2; // 20% offset for axis labels
845 
846  Float_t a = 0.5*ex;
847 
848  fBBox[0] = -a;
849  fBBox[1] = a;
850  fBBox[2] = -a;
851  fBBox[3] = a;
852 
853  // scaling is relative to shortest XY axis
854  Double_t em, eM, pm, pM;
855  fData->GetEtaLimits(em, eM);
856  fData->GetPhiLimits(pm, pM);
857  Double_t r = (eM-em)/(pM-pm);
858  if (r<1)
859  {
860  fBBox[2] /= r;
861  fBBox[3] /= r;
862  }
863  else
864  {
865  fBBox[0] *= r;
866  fBBox[1] *= r;
867  }
868 
869  fBBox[4] = 0;
870  if (fScaleAbs && !fData->Empty())
871  fBBox[5] = GetMaxVal()*GetValToHeight();
872  else
873  fBBox[5] = fMaxTowerH;
874 }