Logo ROOT   6.12/07
Reference Guide
TEveCaloData.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 "TEveCaloData.h"
13 #include "TEveCalo.h"
14 
15 #include "TGLSelectRecord.h"
16 
17 #include "TAxis.h"
18 #include "THStack.h"
19 #include "TH2.h"
20 #include "TMath.h"
21 #include "TList.h"
22 
23 #include <cassert>
24 #include <algorithm>
25 #include <set>
26 
27 /** \class TEveCaloData::CellGeom_t
28 \ingroup TEve
29 Cell geometry inner structure.
30 */
31 
32 ////////////////////////////////////////////////////////////////////////////////
33 /// Print member data.
34 
36 {
37  printf("%f, %f %f, %f \n", fEtaMin, fEtaMax, fPhiMin, fPhiMax);
38 }
39 
40 ////////////////////////////////////////////////////////////////////////////////
41 
43 {
44  fEtaMin = etaMin;
45  fEtaMax = etaMax;
46 
47  fPhiMin = phiMin;
48  fPhiMax = phiMax;
49 
50  // Complain if phi is out of [-2*pi, 2*pi] range.
51  if (fPhiMin < - TMath::TwoPi() || fPhiMin > TMath::TwoPi() ||
53  {
54  ::Error("TEveCaloData::CellGeom_t::Configure", "phiMin and phiMax should be between -2*pi and 2*pi (min=%f, max=%f). RhoZ projection will be wrong.",
55  fPhiMin, fPhiMax);
56  }
57 
60 }
61 
62 /** \class TEveCaloData::CellData_t
63 \ingroup TEve
64 Cell data inner structure.
65 */
66 
67 ////////////////////////////////////////////////////////////////////////////////
68 /// Return energy value associated with the cell, usually Et.
69 /// If isEt is false it is transformed into energy E.
70 
72 {
73  if (isEt)
74  return fValue;
75  else
76  return TMath::Abs(fValue/TMath::Sin(Theta()));
77 }
78 
79 ////////////////////////////////////////////////////////////////////////////////
80 /// Print member data.
81 
83 {
84  printf("%f, %f %f, %f \n", fEtaMin, fEtaMax, fPhiMin, fPhiMax);
85 }
86 
87 ////////////////////////////////////////////////////////////////////////////////
88 
90 {
91  // printf("get val vec bin %d size %d\n", bin, fBinData.size());
92  if (fBinData[bin] == -1)
93  {
94  fBinData[bin] = fSliceData.size();
95 
96  for (Int_t i=0; i<fNSlices; i++)
97  fSliceData.push_back(0.f);
98  }
99 
100  return &fSliceData[fBinData[bin]];
101 }
102 
103 /** \class TEveCaloData
104 \ingroup TEve
105 A central manager for calorimeter event data. It provides a list of
106 cells within requested phi and eta range.
107 */
108 
110 
111 ////////////////////////////////////////////////////////////////////////////////
112 
113 TEveCaloData::TEveCaloData(const char* n, const char* t):
114  TEveElement(),
115  TNamed(n, t),
116 
117  fEtaAxis(0),
118  fPhiAxis(0),
119 
120  fWrapTwoPi(kTRUE),
121 
122  fMaxValEt(0),
123  fMaxValE(0),
124 
125  fEps(0)
126 {
127  // Constructor.
128 }
129 
130 ////////////////////////////////////////////////////////////////////////////////
131 /// Virtual method TEveElement::UnSelect.
132 /// Clear selected towers when deselected.
133 
135 {
136  fCellsSelected.clear();
137 }
138 
139 ////////////////////////////////////////////////////////////////////////////////
140 /// Virtual method TEveElement::UnHighlighted.
141 
143 {
144  fCellsHighlighted.clear();
145 }
146 
147 ////////////////////////////////////////////////////////////////////////////////
148 
150 {
151  if (fCellsHighlighted.empty()) return "";
152 
153  CellData_t cellData;
154 
155  Bool_t single = fCellsHighlighted.size() == 1;
156  Float_t sum = 0;
157  TString s;
158  for (vCellId_i i = fCellsHighlighted.begin(); i!=fCellsHighlighted.end(); ++i)
159  {
160  GetCellData(*i, cellData);
161 
162  s += TString::Format("%s %.2f (%.3f, %.3f)",
163  fSliceInfos[i->fSlice].fName.Data(), cellData.fValue,
164  cellData.Eta(), cellData.Phi());
165 
166  if (single) return s;
167  s += "\n";
168  sum += cellData.fValue;
169  }
170  s += TString::Format("Sum = %.2f", sum);
171  return s;
172 }
173 
174 ////////////////////////////////////////////////////////////////////////////////
175 /// Populate set impSelSet with derived / dependant elements.
176 ///
177 
179 {
180  for (List_ci i=fChildren.begin(); i!=fChildren.end(); ++i)
181  {
182  impSelSet.insert(*i);
183  }
184 }
185 
186 ////////////////////////////////////////////////////////////////////////////////
187 /// Print selected cells info.
188 
190 {
191  printf("%d Selected selected cells:\n", (Int_t)fCellsSelected.size());
192  CellData_t cellData;
193 
194  for (vCellId_i i = fCellsSelected.begin(); i != fCellsSelected.end(); ++i)
195  {
196  GetCellData(*i, cellData);
197  printf("Tower [%d] Slice [%d] Value [%.2f] ", i->fTower, i->fSlice, cellData.fValue);
198  printf("Eta:(%f, %f) Phi(%f, %f)\n", cellData.fEtaMin, cellData.fEtaMax, cellData.fPhiMin, cellData.fPhiMax);
199  }
200 }
201 
202 ////////////////////////////////////////////////////////////////////////////////
203 /// Process newly selected cells with given select-record.
204 /// Secondary-select status is set.
205 /// CellSelectionChanged() is called if needed.
206 
208 {
209  typedef std::set<CellId_t> sCellId_t;
210  typedef std::set<CellId_t>::iterator sCellId_i;
211 
212  struct helper
213  {
214  static void fill_cell_set(sCellId_t& cset, vCellId_t& cvec)
215  {
216  for (vCellId_i i = cvec.begin(); i != cvec.end(); ++i)
217  cset.insert(*i);
218  }
219  static void fill_cell_vec(vCellId_t& cvec, sCellId_t& cset)
220  {
221  for (sCellId_i i = cset.begin(); i != cset.end(); ++i)
222  cvec.push_back(*i);
223  }
224  };
225 
227 
228  if (cells.empty())
229  {
230  if (!sel_cells.empty())
231  {
232  cells.swap(sel_cells);
234  }
235  }
236  else
237  {
238  if (!sel_cells.empty())
239  {
240  if (rec.GetMultiple())
241  {
242  sCellId_t cs;
243  helper::fill_cell_set(cs, cells);
244  for (vCellId_i i = sel_cells.begin(); i != sel_cells.end(); ++i)
245  {
246  std::set<CellId_t>::iterator csi = cs.find(*i);
247  if (csi == cs.end())
248  cs.insert(*i);
249  else
250  cs.erase(csi);
251  }
252  cells.clear();
253  if (cs.empty())
254  {
256  }
257  else
258  {
259  helper::fill_cell_vec(cells, cs);
261  }
262  }
263  else
264  {
265  Bool_t differ = kFALSE;
266  if (cells.size() == sel_cells.size())
267  {
268  sCellId_t cs;
269  helper::fill_cell_set(cs, cells);
270  for (vCellId_i i = sel_cells.begin(); i != sel_cells.end(); ++i)
271  {
272  if (cs.find(*i) == cs.end())
273  {
274  differ = kTRUE;
275  break;
276  }
277  }
278  }
279  else
280  {
281  differ = kTRUE;
282  }
283  if (differ)
284  {
285  cells.swap(sel_cells);
287  }
288  }
289  }
290  else
291  {
292  if (!rec.GetMultiple())
293  {
294  cells.clear();
296  }
297  }
298  }
299 
301  {
303  }
304 }
305 
306 ////////////////////////////////////////////////////////////////////////////////
307 /// Set threshold for given slice.
308 
310 {
311  fSliceInfos[slice].fThreshold = val;
313 }
314 
315 ////////////////////////////////////////////////////////////////////////////////
316 /// Get threshold for given slice.
317 
319 {
320  return fSliceInfos[slice].fThreshold;
321 }
322 
323 ////////////////////////////////////////////////////////////////////////////////
324 /// Set color for given slice.
325 
327 {
328  fSliceInfos[slice].fColor = col;
329  for (List_ci i=fChildren.begin(); i!=fChildren.end(); ++i)
330  {
331  (*i)->AddStamp(TEveElement::kCBObjProps);
332  }
333 }
334 
335 ////////////////////////////////////////////////////////////////////////////////
336 /// Get color for given slice.
337 
339 {
340  return fSliceInfos[slice].fColor;
341 }
342 
343 ////////////////////////////////////////////////////////////////////////////////
344 /// Set transparency for given slice.
345 
347 {
348  fSliceInfos[slice].fTransparency = t;
349  for (List_ci i=fChildren.begin(); i!=fChildren.end(); ++i)
350  {
351  (*i)->AddStamp(TEveElement::kCBObjProps);
352  }
353 }
354 
355 ////////////////////////////////////////////////////////////////////////////////
356 /// Get transparency for given slice.
357 
359 {
360  return fSliceInfos[slice].fTransparency;
361 }
362 
363 ////////////////////////////////////////////////////////////////////////////////
364 /// Invalidate cell ids cache on back ptr references.
365 
367 {
368  TEveCaloViz* calo;
369  for (List_ci i=fChildren.begin(); i!=fChildren.end(); ++i)
370  {
371  calo = dynamic_cast<TEveCaloViz*>(*i);
372  calo->InvalidateCellIdCache();
373  calo->StampObjProps();
374  }
375 }
376 
377 ////////////////////////////////////////////////////////////////////////////////
378 /// Tell users (TEveCaloViz instances using this data) that data
379 /// has changed and they should update the limits/scales etc.
380 /// This is done by calling TEveCaloViz::DataChanged().
381 
383 {
384  TEveCaloViz* calo;
385  for (List_ci i=fChildren.begin(); i!=fChildren.end(); ++i)
386  {
387  calo = dynamic_cast<TEveCaloViz*>(*i);
388  calo->DataChanged();
389  calo->StampObjProps();
390  }
391 }
392 
393 ////////////////////////////////////////////////////////////////////////////////
394 /// Tell users (TEveCaloViz instances using this data) that cell selection
395 /// has changed and they should update selection cache if necessary.
396 /// This is done by calling TEveCaloViz::CellSelectionChanged().
397 
399 {
400  TEveCaloViz* calo;
401  for (List_ci i=fChildren.begin(); i!=fChildren.end(); ++i)
402  {
403  calo = dynamic_cast<TEveCaloViz*>(*i);
404  calo->CellSelectionChanged();
405  calo->StampColorSelection();
406  }
407 }
408 
409 ////////////////////////////////////////////////////////////////////////////////
410 
412 {
413  using namespace TMath;
414 
415  if (eta < 0)
416  return Pi() - 2*ATan(Exp(- Abs(eta)));
417  else
418  return 2*ATan(Exp(- Abs(eta)));
419 }
420 
421 
422 /** \class TEveCaloDataVec
423 \ingroup TEve
424 Calo data for universal cell geometry.
425 */
426 
428 
429 ////////////////////////////////////////////////////////////////////////////////
430 
432  TEveCaloData(),
433 
434  fTower(0),
435  fEtaMin( 1e3),
436  fEtaMax(-1e3),
437  fPhiMin( 1e3),
438  fPhiMax(-1e3)
439 {
440  // Constructor.
441 
442  fSliceInfos.assign(nslices, SliceInfo_t());
443 
444  fSliceVec.assign(nslices, std::vector<Float_t> ());
445 }
446 
447 ////////////////////////////////////////////////////////////////////////////////
448 /// Destructor.
449 
451 {
452  if (fEtaAxis) delete fEtaAxis;
453  if (fPhiAxis) delete fPhiAxis;
454 }
455 
456 ////////////////////////////////////////////////////////////////////////////////
457 /// Add new slice.
458 
460 {
461  fSliceInfos.push_back(SliceInfo_t());
462  fSliceVec.push_back(std::vector<Float_t> ());
463  fSliceVec.back().resize(fGeomVec.size(), 0.f);
464 
465  return fSliceInfos.size() - 1;
466 }
467 
468 ////////////////////////////////////////////////////////////////////////////////
469 /// Add tower within eta/phi range.
470 
472 {
473  assert (etaMin < etaMax);
474  assert (phiMin < phiMax);
475 
476  fGeomVec.push_back(CellGeom_t(etaMin, etaMax, phiMin, phiMax));
477 
478  for (vvFloat_i it=fSliceVec.begin(); it!=fSliceVec.end(); ++it)
479  (*it).push_back(0);
480 
481  if (etaMin < fEtaMin) fEtaMin = etaMin;
482  if (etaMax > fEtaMax) fEtaMax = etaMax;
483 
484  if (phiMin < fPhiMin) fPhiMin = phiMin;
485  if (phiMax > fPhiMax) fPhiMax = phiMax;
486 
487  fTower = fGeomVec.size() - 1;
488  return fTower;
489 }
490 
491 ////////////////////////////////////////////////////////////////////////////////
492 /// Fill given slice in the current tower.
493 
495 {
496  fSliceVec[slice][fTower] = val;
497 }
498 
499 ////////////////////////////////////////////////////////////////////////////////
500 /// Fill given slice in a given tower.
501 
503 {
504  fSliceVec[slice][tower] = val;
505 }
506 
507 
508 ////////////////////////////////////////////////////////////////////////////////
509 /// Get list of cell-ids for given eta/phi range.
510 
512  Float_t phi, Float_t phiD,
513  TEveCaloData::vCellId_t &out) const
514 {
515  using namespace TMath;
516 
517  Float_t etaMin = eta - etaD*0.5;
518  Float_t etaMax = eta + etaD*0.5;
519 
520  Float_t phiMin = phi - phiD*0.5;
521  Float_t phiMax = phi + phiD*0.5;
522 
523  Int_t nS = fSliceVec.size();
524 
525  Int_t tower = 0;
526  Float_t fracx=0, fracy=0, frac;
527  Float_t minQ, maxQ;
528 
529  for(vCellGeom_ci i=fGeomVec.begin(); i!=fGeomVec.end(); i++)
530  {
531  const CellGeom_t &cg = *i;
532  fracx = TEveUtil::GetFraction(etaMin, etaMax, cg.fEtaMin, cg.fEtaMax);
533  if (fracx > 1e-3)
534  {
535  minQ = cg.fPhiMin;
536  maxQ = cg.fPhiMax;
537 
538  if (fWrapTwoPi)
539  {
540  if (maxQ < phiMin)
541  {
542  minQ += TwoPi(); maxQ += TwoPi();
543  }
544  else if (minQ > phiMax)
545  {
546  minQ -= TwoPi(); maxQ -= TwoPi();
547  }
548  }
549 
550  if (maxQ >= phiMin && minQ <= phiMax)
551  {
552  fracy = TEveUtil::GetFraction(phiMin, phiMax, minQ, maxQ);
553  if (fracy > 1e-3)
554  {
555  frac = fracx*fracy;
556  for (Int_t s=0; s<nS; s++)
557  {
558  if (fSliceVec[s][tower] > fSliceInfos[s].fThreshold)
559  out.push_back(CellId_t(tower, s, frac));
560  }
561  }
562  }
563  }
564  tower++;
565  }
566 }
567 
568 ////////////////////////////////////////////////////////////////////////////////
569 /// Rebin cells.
570 
571 void TEveCaloDataVec::Rebin(TAxis* ax, TAxis* ay, vCellId_t &ids, Bool_t et, RebinData_t& rdata) const
572 {
573  rdata.fNSlices = GetNSlices();
574  rdata.fBinData.assign((ax->GetNbins()+2)*(ay->GetNbins()+2), -1);
575 
576  CellData_t cd;
577  for (vCellId_i it = ids.begin(); it != ids.end(); ++it)
578  {
579  GetCellData(*it, cd);
580  Int_t iMin = ax->FindBin(cd.EtaMin());
581  Int_t iMax = ax->FindBin(cd.EtaMax());
582  Int_t jMin = ay->FindBin(cd.PhiMin());
583  Int_t jMax = ay->FindBin(cd.PhiMax());
584  for (Int_t i = iMin; i <= iMax; ++i)
585  {
586  if (i < 0 || i > ax->GetNbins()) continue;
587  for (Int_t j = jMin; j <= jMax; ++j)
588  {
589  if (j < 0 || j > ay->GetNbins()) continue;
590 
591  Double_t ratio = TEveUtil::GetFraction(ax->GetBinLowEdge(i), ax->GetBinUpEdge(i), cd.EtaMin(), cd.EtaMax())
592  * TEveUtil::GetFraction(ay->GetBinLowEdge(j), ay->GetBinUpEdge(j), cd.PhiMin(), cd.PhiMax());
593 
594  if (ratio > 1e-6f)
595  {
596  Float_t* slices = rdata.GetSliceVals(i + j*(ax->GetNbins()+2));
597  slices[(*it).fSlice] += ratio * cd.Value(et);
598  }
599  }
600  }
601  }
602 }
603 
604 ////////////////////////////////////////////////////////////////////////////////
605 /// Get cell geometry and value from cell ID.
606 
608  TEveCaloData::CellData_t& cellData) const
609 {
610  cellData.CellGeom_t::operator=( fGeomVec[id.fTower] );
611  cellData.fValue = fSliceVec[id.fSlice][id.fTower];
612 }
613 
614 ////////////////////////////////////////////////////////////////////////////////
615 /// Update limits and notify data users.
616 
618 {
619  using namespace TMath;
620 
621  // update max E/Et values
622 
623  fMaxValE = 0;
624  fMaxValEt = 0;
625  Float_t sum=0;
626  // printf("geom vec %d slices %d\n",fGeomVec.size(), fSliceVec.size() );
627 
628  for (UInt_t tw=0; tw<fGeomVec.size(); tw++)
629  {
630  sum=0;
631  for (vvFloat_i it=fSliceVec.begin(); it!=fSliceVec.end(); ++it)
632  sum += (*it)[tw];
633 
634  if (sum > fMaxValEt ) fMaxValEt=sum;
635 
636  sum /= Abs(Sin(EtaToTheta(fGeomVec[tw].Eta())));
637 
638  if (sum > fMaxValE) fMaxValE=sum;
639  }
640 
642 }
643 
644 
645 ////////////////////////////////////////////////////////////////////////////////
646 /// Set XY axis from cells geometry.
647 
649 {
650  std::vector<Double_t> binX;
651  std::vector<Double_t> binY;
652 
653  for(vCellGeom_ci i=fGeomVec.begin(); i!=fGeomVec.end(); i++)
654  {
655  const CellGeom_t &ch = *i;
656 
657  binX.push_back(ch.EtaMin());
658  binX.push_back(ch.EtaMax());
659  binY.push_back(ch.PhiMin());
660  binY.push_back(ch.PhiMax());
661  }
662 
663  std::sort(binX.begin(), binX.end());
664  std::sort(binY.begin(), binY.end());
665 
666  Int_t cnt = 0;
667  Double_t sum = 0;
668  Double_t val;
669 
670  // X axis
671  Double_t dx = binX.back() - binX.front();
672  epsX *= dx;
673  std::vector<Double_t> newX;
674  newX.push_back(binX.front()); // underflow
675  Int_t nX = binX.size()-1;
676  for(Int_t i=0; i<nX; i++)
677  {
678  val = (sum +binX[i])/(cnt+1);
679  if (binX[i+1] -val > epsX)
680  {
681  newX.push_back(val);
682  cnt = 0;
683  sum = 0;
684  }
685  else
686  {
687  sum += binX[i];
688  cnt++;
689  }
690  }
691  newX.push_back(binX.back()); // overflow
692 
693  // Y axis
694  cnt = 0;
695  sum = 0;
696  std::vector<Double_t> newY;
697  Double_t dy = binY.back() - binY.front();
698  epsY *= dy;
699  newY.push_back(binY.front());// underflow
700  Int_t nY = binY.size()-1;
701  for(Int_t i=0 ; i<nY; i++)
702  {
703  val = (sum +binY[i])/(cnt+1);
704  if (binY[i+1] -val > epsY )
705  {
706  newY.push_back(val);
707  cnt = 0;
708  sum = 0;
709  }
710  else
711  {
712  sum += binY[i];
713  cnt++;
714  }
715 
716  }
717  newY.push_back(binY.back()); // overflow
718 
719  if (fEtaAxis) delete fEtaAxis;
720  if (fPhiAxis) delete fPhiAxis;
721 
722  fEtaAxis = new TAxis(newX.size()-1, &newX[0]);
723  fPhiAxis = new TAxis(newY.size()-1, &newY[0]);
724  fEtaAxis->SetNdivisions(510);
725  fPhiAxis->SetNdivisions(510);
726 }
727 
728 /** \class TEveCaloDataHist
729 \ingroup TEve
730 A central manager for calorimeter data of an event written in TH2F.
731 X axis is used for eta and Y axis for phi.
732 */
733 
735 
736 ////////////////////////////////////////////////////////////////////////////////
737 /// Constructor.
738 
740  TEveCaloData(),
741 
742  fHStack(0)
743 {
744  fHStack = new THStack();
745  fEps = 1e-5;
746 }
747 
748 ////////////////////////////////////////////////////////////////////////////////
749 /// Destructor.
750 
752 {
753  delete fHStack;
754 }
755 
756 ////////////////////////////////////////////////////////////////////////////////
757 /// Update limits and notify data users.
758 
760 {
761  using namespace TMath;
762 
763  // update max E/Et values
764  fMaxValE = 0;
765  fMaxValEt = 0;
766 
767  if (GetNSlices() < 1) return;
768 
769  TH2* hist = GetHist(0);
770  fEtaAxis = hist->GetXaxis();
771  fPhiAxis = hist->GetYaxis();
772  for (Int_t ieta = 1; ieta <= fEtaAxis->GetNbins(); ++ieta)
773  {
774  Double_t eta = fEtaAxis->GetBinCenter(ieta); // conversion E/Et
775  for (Int_t iphi = 1; iphi <= fPhiAxis->GetNbins(); ++iphi)
776  {
777  Double_t value = 0;
778  for (Int_t i = 0; i < GetNSlices(); ++i)
779  {
780  hist = GetHist(i);
781  Int_t bin = hist->GetBin(ieta, iphi);
782  value += hist->GetBinContent(bin);
783  }
784 
785  if (value > fMaxValEt ) fMaxValEt = value;
786 
787  value /= Abs(Sin(EtaToTheta(eta)));
788 
789  if (value > fMaxValE) fMaxValE = value;
790  }
791  }
793 }
794 
795 ////////////////////////////////////////////////////////////////////////////////
796 /// Get list of cell IDs in given eta and phi range.
797 
799  Float_t phi, Float_t phiD,
800  TEveCaloData::vCellId_t &out) const
801 {
802  using namespace TMath;
803 
804  Float_t etaMin = eta - etaD*0.5 -fEps;
805  Float_t etaMax = eta + etaD*0.5 +fEps;
806 
807  Float_t phiMin = phi - phiD*0.5 -fEps;
808  Float_t phiMax = phi + phiD*0.5 +fEps;
809 
810  Int_t nEta = fEtaAxis->GetNbins();
811  Int_t nPhi = fPhiAxis->GetNbins();
812  Int_t nSlices = GetNSlices();
813 
814  Int_t bin = 0;
815 
816  Bool_t accept;
817  for (Int_t ieta = 1; ieta <= nEta; ++ieta)
818  {
819  if (fEtaAxis->GetBinLowEdge(ieta) >= etaMin && fEtaAxis->GetBinUpEdge(ieta) <= etaMax)
820  {
821  for (Int_t iphi = 1; iphi <= nPhi; ++iphi)
822  {
823  if (fWrapTwoPi )
824  {
826  (phiMin, phiMax, fPhiAxis->GetBinLowEdge(iphi), fPhiAxis->GetBinUpEdge(iphi));
827  }
828  else
829  {
830  accept = fPhiAxis->GetBinLowEdge(iphi) >= phiMin && fPhiAxis->GetBinUpEdge(iphi) <= phiMax &&
831  fPhiAxis->GetBinLowEdge(iphi) >= phiMin && fPhiAxis->GetBinUpEdge(iphi) <= phiMax;
832  }
833 
834  if (accept)
835  {
836  for (Int_t s = 0; s < nSlices; ++s)
837  {
838  TH2F *hist = GetHist(s);
839  bin = hist->GetBin(ieta, iphi);
840  if (hist->GetBinContent(bin) > fSliceInfos[s].fThreshold)
841  out.push_back(TEveCaloData::CellId_t(bin, s));
842  } // hist slices
843  }
844  } // phi bins
845  }
846  } // eta bins
847 }
848 
849 ////////////////////////////////////////////////////////////////////////////////
850 /// Rebin
851 
853 {
854  rdata.fNSlices = GetNSlices();
855  rdata.fBinData.assign((ax->GetNbins()+2)*(ay->GetNbins()+2), -1);
857  Float_t *val;
858  Int_t i, j, w;
859  Int_t binx, biny;
860  Int_t bin;
861 
862  for (vCellId_i it=ids.begin(); it!=ids.end(); ++it)
863  {
864  GetCellData(*it, cd);
865  GetHist(it->fSlice)->GetBinXYZ((*it).fTower, i, j, w);
866  binx = ax->FindBin(fEtaAxis->GetBinCenter(i));
867  biny = ay->FindBin(fPhiAxis->GetBinCenter(j));
868  bin = biny*(ax->GetNbins()+2)+binx;
869  val = rdata.GetSliceVals(bin);
870  Double_t ratio = TEveUtil::GetFraction(ax->GetBinLowEdge(binx), ax->GetBinUpEdge(binx), cd.EtaMin(), cd.EtaMax())
871  * TEveUtil::GetFraction(ay->GetBinLowEdge(biny), ay->GetBinUpEdge(biny), cd.PhiMin(), cd.PhiMax());
872 
873  val[(*it).fSlice] += cd.Value(et)*ratio;
874  }
875 }
876 
877 ////////////////////////////////////////////////////////////////////////////////
878 /// Get cell geometry and value from cell ID.
879 
881  TEveCaloData::CellData_t& cellData) const
882 {
883  TH2F* hist = GetHist(id.fSlice);
884 
885  Int_t x, y, z;
886  hist->GetBinXYZ(id.fTower, x, y, z);
887 
888  cellData.fValue = hist->GetBinContent(id.fTower);
889  cellData.Configure(hist->GetXaxis()->GetBinLowEdge(x),
890  hist->GetXaxis()->GetBinUpEdge(x),
891  hist->GetYaxis()->GetBinLowEdge(y),
892  hist->GetYaxis()->GetBinUpEdge(y));
893 }
894 
895 ////////////////////////////////////////////////////////////////////////////////
896 /// Add new slice to calo tower. Updates cached variables fMaxValE
897 /// and fMaxValEt
898 /// Return last index in the vector of slice infos.
899 
901 {
902  fHStack->Add(hist);
903  fSliceInfos.push_back(SliceInfo_t());
904  fSliceInfos.back().fName = hist->GetName();
905  fSliceInfos.back().fColor = hist->GetLineColor();
906 
907  DataChanged();
908 
909  return fSliceInfos.size() - 1;
910 }
911 
912 ////////////////////////////////////////////////////////////////////////////////
913 /// Get histogram in given slice.
914 
916 {
917  assert(slice >= 0 && slice < fHStack->GetHists()->GetSize());
918  return (TH2F*) fHStack->GetHists()->At(slice);
919 }
920 
921 ////////////////////////////////////////////////////////////////////////////////
922 /// Get eta limits.
923 
925 {
926  min = fEtaAxis->GetXmin();
927  max = fEtaAxis->GetXmax();
928 }
929 
930 ////////////////////////////////////////////////////////////////////////////////
931 /// Get phi limits.
932 
934 {
935  min = fPhiAxis->GetXmin();
936  max = fPhiAxis->GetXmax();
937 }
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
void SetSecSelResult(ESecSelResult r)
static long int sum(long int i)
Definition: Factory.cxx:2173
TAxis * fEtaAxis
Definition: TEveCaloData.h:160
A central manager for calorimeter event data.
Definition: TEveCaloData.h:26
std::vector< CellId_t >::iterator vCellId_i
Definition: TEveCaloData.h:147
The Histogram stack class.
Definition: THStack.h:31
TEveCaloDataVec(const TEveCaloDataVec &)
virtual void GetCellList(Float_t etaMin, Float_t etaMax, Float_t phi, Float_t phiRng, vCellId_t &out) const
Get list of cell IDs in given eta and phi range.
vCellId_t fCellsHighlighted
Definition: TEveCaloData.h:171
float Float_t
Definition: RtypesCore.h:53
Bool_t GetMultiple() const
std::set< TEveElement * > Set_t
Definition: TEveElement.h:73
void Configure(Float_t etaMin, Float_t etaMax, Float_t phiMin, Float_t phiMax)
virtual Double_t GetBinLowEdge(Int_t bin) const
Return low edge of bin.
Definition: TAxis.cxx:504
constexpr Double_t TwoPi()
Definition: TMath.h:44
List_t fChildren
Definition: TEveElement.h:79
std::vector< Int_t > fBinData
Definition: TEveCaloData.h:133
virtual void DataChanged()
Tell users (TEveCaloViz instances using this data) that data has changed and they should update the l...
virtual void SetNdivisions(Int_t n=510, Bool_t optim=kTRUE)
Set the number of divisions for this axis.
Definition: TAttAxis.cxx:229
virtual TString GetHighlightTooltip()
virtual void GetBinXYZ(Int_t binglobal, Int_t &binx, Int_t &biny, Int_t &binz) const
Return binx, biny, binz corresponding to the global bin number globalbin see TH1::GetBin function abo...
Definition: TH1.cxx:4678
std::vector< vFloat_t >::iterator vvFloat_i
Definition: TEveCaloData.h:251
Calo data for universal cell geometry.
Definition: TEveCaloData.h:239
Basic string class.
Definition: TString.h:125
TEveCaloData(const TEveCaloData &)
int Int_t
Definition: RtypesCore.h:41
void FillSlice(Int_t slice, Float_t value)
Fill given slice in the current tower.
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
Char_t GetSliceTransparency(Int_t slice) const
Get transparency for given slice.
Float_t PhiMin() const
Definition: TEveCaloData.h:95
void SetSliceColor(Int_t slice, Color_t col)
Set color for given slice.
Short_t Abs(Short_t d)
Definition: TMathBase.h:108
Bool_t fWrapTwoPi
Definition: TEveCaloData.h:163
virtual Double_t GetBinUpEdge(Int_t bin) const
Return up edge of bin.
Definition: TAxis.cxx:514
virtual void UnHighlighted()
Virtual method TEveElement::UnHighlighted.
virtual void FillImpliedSelectedSet(Set_t &impSelSet)
Populate set impSelSet with derived / dependant elements.
virtual void GetCellData(const TEveCaloData::CellId_t &id, TEveCaloData::CellData_t &data) const
Get cell geometry and value from cell ID.
Float_t fEps
Definition: TEveCaloData.h:168
Double_t GetXmin() const
Definition: TAxis.h:133
static Float_t GetFraction(Float_t minM, Float_t maxM, Float_t minQ, Float_t maxQ)
Get fraction of interval [minQ, maxQ] in [minM, maxM].
Definition: TEveUtil.cxx:386
Double_t x[n]
Definition: legend1.C:17
Float_t EtaMax() const
Definition: TEveCaloData.h:91
Float_t Theta() const
Definition: TEveCaloData.h:102
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString...
Definition: TString.cxx:2365
static Float_t EtaToTheta(Float_t eta)
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:29
TAxis * fPhiAxis
Definition: TEveCaloData.h:161
virtual void GetEtaLimits(Double_t &min, Double_t &max) const
Get eta limits.
virtual ~TEveCaloDataHist()
Destructor.
virtual Double_t GetBinCenter(Int_t bin) const
Return center of bin.
Definition: TAxis.cxx:464
virtual void InvalidateUsersCellIdCache()
Invalidate cell ids cache on back ptr references.
std::vector< CellId_t > vCellId_t
Definition: TEveCaloData.h:146
virtual void Rebin(TAxis *ax, TAxis *ay, vCellId_t &in, Bool_t et, RebinData_t &out) const
Rebin cells.
TH2F * GetHist(Int_t slice) const
Get histogram in given slice.
Float_t Phi() const
Definition: TEveCaloData.h:97
short Color_t
Definition: RtypesCore.h:79
Float_t PhiMax() const
Definition: TEveCaloData.h:96
virtual void UnSelected()
Virtual method TEveElement::UnSelect.
Float_t fMaxValE
Definition: TEveCaloData.h:166
Int_t AddHistogram(TH2F *hist)
Add new slice to calo tower.
virtual void Rebin(TAxis *ax, TAxis *ay, vCellId_t &in, Bool_t et, RebinData_t &out) const
Rebin.
List_t::const_iterator List_ci
Definition: TEveElement.h:71
vCellId_t fCellsSelected
Definition: TEveCaloData.h:170
Bool_t GetHighlight() const
Service class for 2-Dim histogram classes.
Definition: TH2.h:30
Class to manage histogram axis.
Definition: TAxis.h:30
void InvalidateCellIdCache()
Definition: TEveCalo.h:92
TList * GetHists() const
Definition: THStack.h:60
virtual void DataChanged()
Update limits and notify data users.
double Pi()
Mathematical constants.
Definition: Math.h:90
void ProcessSelection(vCellId_t &sel_cells, TGLSelectRecord &rec)
Process newly selected cells with given select-record.
TEveCaloDataHist()
Constructor.
PyObject * fValue
2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:249
Standard selection record including information about containing scene and details ob out selected ob...
unsigned int UInt_t
Definition: RtypesCore.h:42
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:880
void StampObjProps()
Definition: TEveElement.h:397
Cell geometry inner structure.
Definition: TEveCaloData.h:72
virtual void GetCellData(const CellId_t &id, CellData_t &data) const =0
virtual void CellSelectionChanged()
Tell users (TEveCaloViz instances using this data) that cell selection has changed and they should up...
TAxis * GetYaxis()
Definition: TH1.h:316
virtual void CellSelectionChanged()
Definition: TEveCalo.h:82
virtual TObject * At(Int_t idx) const
Returns the object at position idx. Returns 0 if idx is out of range.
Definition: TList.cxx:354
void SetAxisFromBins(Double_t epsX=0.001, Double_t epsY=0.001)
Set XY axis from cells geometry.
void StampColorSelection()
Definition: TEveElement.h:395
Float_t fMaxValEt
Definition: TEveCaloData.h:165
ESecSelResult GetSecSelResult() const
const Bool_t kFALSE
Definition: RtypesCore.h:88
virtual Color_t GetLineColor() const
Return the line color.
Definition: TAttLine.h:33
virtual Int_t FindBin(Double_t x)
Find bin number corresponding to abscissa x.
Definition: TAxis.cxx:279
Int_t AddTower(Float_t etaMin, Float_t etaMax, Float_t phiMin, Float_t phiMax)
Add tower within eta/phi range.
vSliceInfo_t fSliceInfos
Definition: TEveCaloData.h:158
Double_t Exp(Double_t x)
Definition: TMath.h:621
virtual void GetCellList(Float_t etaMin, Float_t etaMax, Float_t phi, Float_t phiRng, vCellId_t &out) const
Get list of cell-ids for given eta/phi range.
#define ClassImp(name)
Definition: Rtypes.h:359
double Double_t
Definition: RtypesCore.h:55
virtual void DataChanged()
Update limits and notify data users.
virtual void Dump() const
Print member data.
vvFloat_t fSliceVec
Definition: TEveCaloData.h:253
Color_t GetSliceColor(Int_t slice) const
Get color for given slice.
void SetSliceTransparency(Int_t slice, Char_t t)
Set transparency for given slice.
virtual void Add(TH1 *h, Option_t *option="")
add a new histogram to the list Only 1-d and 2-d histograms currently supported.
Definition: THStack.cxx:362
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
static constexpr double s
you should not use this method at all Int_t Int_t Double_t Double_t Double_t e
Definition: TRolke.cxx:630
Base class for calorimeter data visualization.
Definition: TEveCalo.h:26
void DataChanged()
Update setting and cache on data changed.
Definition: TEveCalo.cxx:256
you should not use this method at all Int_t Int_t z
Definition: TRolke.cxx:630
std::vector< CellGeom_t >::const_iterator vCellGeom_ci
Definition: TEveCaloData.h:151
char Char_t
Definition: RtypesCore.h:29
Float_t GetSliceThreshold(Int_t slice) const
Get threshold for given slice.
Float_t Eta() const
Definition: TEveCaloData.h:92
virtual void GetCellData(const TEveCaloData::CellId_t &id, TEveCaloData::CellData_t &data) const
Get cell geometry and value from cell ID.
virtual ~TEveCaloDataVec()
Destructor.
Double_t Sin(Double_t)
Definition: TMath.h:547
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition: TH2.h:84
Cell data inner structure.
Definition: TEveCaloData.h:114
vCellGeom_t fGeomVec
Definition: TEveCaloData.h:254
void PrintCellsSelected()
Print selected cells info.
Int_t GetNbins() const
Definition: TAxis.h:121
THStack * fHStack
Definition: TEveCaloData.h:304
Float_t Value(Bool_t) const
Return energy value associated with the cell, usually Et.
virtual Int_t GetBin(Int_t binx, Int_t biny, Int_t binz=0) const
Return Global bin number corresponding to binx,y,z.
Definition: TH2.cxx:966
virtual void Dump() const
Print member data.
void SetSliceThreshold(Int_t slice, Float_t threshold)
Set threshold for given slice.
const Bool_t kTRUE
Definition: RtypesCore.h:87
Base class for TEveUtil visualization elements, providing hierarchy management, rendering control and...
Definition: TEveElement.h:33
Double_t GetXmax() const
Definition: TAxis.h:134
const Int_t n
Definition: legend1.C:16
Float_t EtaMin() const
Definition: TEveCaloData.h:90
virtual void GetPhiLimits(Double_t &min, Double_t &max) const
Get phi limits.
Int_t GetNSlices() const
Definition: TEveCaloData.h:202
const char * cnt
Definition: TXMLSetup.cxx:74
TAxis * GetXaxis()
Get the behaviour adopted by the object about the statoverflows. See EStatOverflows for more informat...
Definition: TH1.h:315
Float_t * GetSliceVals(Int_t bin)
A central manager for calorimeter data of an event written in TH2F.
Definition: TEveCaloData.h:297
Int_t AddSlice()
Add new slice.
Double_t ATan(Double_t)
Definition: TMath.h:577