Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
REveCaloData.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 "ROOT/REveCaloData.hxx"
13#include "ROOT/REveCalo.hxx"
14#include "ROOT/REveUtil.hxx"
15#include "ROOT/REveManager.hxx"
17
18#include "TAxis.h"
19#include "THStack.h"
20#include "TH2.h"
21#include "TMath.h"
22#include "TList.h"
23
24#include <cassert>
25#include <algorithm>
26#include <set>
27
28#include <nlohmann/json.hpp>
29
30using namespace ROOT::Experimental;
31
32/** \class REveCaloData::CellGeom_t
33\ingroup REve
34Cell geometry inner structure.
35*/
36
37////////////////////////////////////////////////////////////////////////////////
38/// Print member data.
39
41{
42 printf("%f, %f %f, %f \n", fEtaMin, fEtaMax, fPhiMin, fPhiMax);
43}
44
45////////////////////////////////////////////////////////////////////////////////
46
48{
49 fEtaMin = etaMin;
50 fEtaMax = etaMax;
51
52 fPhiMin = phiMin;
53 fPhiMax = phiMax;
54
55 // Complain if phi is out of [-2*pi, 2*pi] range.
56 if (fPhiMin < - TMath::TwoPi() || fPhiMin > TMath::TwoPi() ||
57 fPhiMax < - TMath::TwoPi() || fPhiMax > TMath::TwoPi())
58 {
59 ::Error("REveCaloData::CellGeom_t::Configure", "phiMin and phiMax should be between -2*pi and 2*pi (min=%f, max=%f). RhoZ projection will be wrong.",
60 fPhiMin, fPhiMax);
61 }
62
63 fThetaMin = EtaToTheta(fEtaMax);
64 fThetaMax = EtaToTheta(fEtaMin);
65}
66
67/** \class REveCaloData::CellData_t
68\ingroup REve
69Cell data inner structure.
70*/
71
72////////////////////////////////////////////////////////////////////////////////
73/// Return energy value associated with the cell, usually Et.
74/// If isEt is false it is transformed into energy E.
75
77{
78 if (isEt)
79 return fValue;
80 else
81 return TMath::Abs(fValue/TMath::Sin(Theta()));
82}
83
84////////////////////////////////////////////////////////////////////////////////
85/// Print member data.
86
88{
89 printf("%f, %f %f, %f \n", fEtaMin, fEtaMax, fPhiMin, fPhiMax);
90}
91
92////////////////////////////////////////////////////////////////////////////////
93
95{
96 // printf("get val vec bin %d size %d\n", bin, fBinData.size());
97 if (fBinData[bin] == -1)
98 {
99 fBinData[bin] = fSliceData.size();
100
101 for (Int_t i=0; i<fNSlices; i++)
102 fSliceData.push_back(0.f);
103 }
104
105 return &fSliceData[fBinData[bin]];
106}
107
108/** \class REveCaloData
109\ingroup REve
110A central manager for calorimeter event data. It provides a list of
111cells within requested phi and eta range.
112*/
113
114
115////////////////////////////////////////////////////////////////////////////////
116
117REveCaloData::REveCaloData(const char* n, const char* t):
118 REveElement(),
119
121
122 fMaxValEt(0),
123 fMaxValE(0),
124
125 fEps(0)
126{
127 SetNameTitle(n,t);
128 // Constructor.
129}
130
131////////////////////////////////////////////////////////////////////////////////
132/// Process newly selected cells with given select-record.
133
134void REveCaloData::ProcessSelection(vCellId_t& sel_cells, UInt_t selectionId, Bool_t multiple)
135{
136 if (fSelector)
137 {
138 fSelector->ProcessSelection(sel_cells, selectionId, multiple);
139 }
140 else
141 {
142 REveSelection* selection = dynamic_cast<REveSelection*> (ROOT::Experimental::gEve->FindElementById(selectionId));
143
144 std::set<int> secondary_idcs;
145 for (vCellId_i i = sel_cells.begin(); i != sel_cells.end(); ++i)
146 {
147 int id = (i->fSlice << 24) + i->fTower;
148 secondary_idcs.insert(id);
149 }
150 selection->NewElementPicked(GetElementId(), multiple, true, secondary_idcs);
151 }
152}
153
154////////////////////////////////////////////////////////////////////////////////
155/// Populate set impSelSet with derived / dependant elements.
156///
157
158void REveCaloData::FillImpliedSelectedSet(Set_t& impSelSet, const std::set<int>&)
159{
160 // printf("REveCaloData::FillImpliedSelectedSet\n");
161 for (auto &n : fNieces)
162 {
163 impSelSet.insert(n);
164 }
165}
166////////////////////////////////////////////////////////////////////////////////
167/// Set threshold for given slice.
168
170{
171 fSliceInfos[slice].fThreshold = val;
173}
174
175////////////////////////////////////////////////////////////////////////////////
176/// Get threshold for given slice.
177
179{
180 return fSliceInfos[slice].fThreshold;
181}
182
183////////////////////////////////////////////////////////////////////////////////
184/// Set color for given slice.
185
187{
188 fSliceInfos[slice].fColor = col;
189 for (auto &c : fNieces)
190 {
191 c->AddStamp(REveElement::kCBObjProps);
192 }
194}
195
196////////////////////////////////////////////////////////////////////////////////
197/// Get color for given slice.
198
200{
201 return fSliceInfos[slice].fColor;
202}
203
204////////////////////////////////////////////////////////////////////////////////
205/// Set transparency for given slice.
206
208{
209 fSliceInfos[slice].fTransparency = t;
210 for (auto &c : fNieces)
211 {
212 c->AddStamp(REveElement::kCBObjProps);
213 }
215}
216
217////////////////////////////////////////////////////////////////////////////////
218/// Get transparency for given slice.
219
221{
222 return fSliceInfos[slice].fTransparency;
223}
224
225////////////////////////////////////////////////////////////////////////////////
226/// Invalidate cell ids cache on back ptr references.
227
229{
230 REveCaloViz* calo;
231 for (auto &c : fNieces)
232 {
233 calo = dynamic_cast<REveCaloViz*>(c);
234 calo->InvalidateCellIdCache();
235 calo->StampObjProps();
236 }
238}
239
240////////////////////////////////////////////////////////////////////////////////
241/// Tell users (REveCaloViz instances using this data) that data
242/// has changed and they should update the limits/scales etc.
243/// This is done by calling REveCaloViz::DataChanged().
244
246{
247 REveCaloViz* calo;
248 for (auto &c : fNieces)
249 {
250 calo = dynamic_cast<REveCaloViz*>(c);
251 calo->DataChanged();
252 calo->StampObjProps();
253 }
254}
255
256
257////////////////////////////////////////////////////////////////////////////////
258
260{
261 Int_t ret = REveElement::WriteCoreJson(j, rnr_offset);
262
263 auto sarr = nlohmann::json::array();
264 for (auto &s : fSliceInfos)
265 {
266 nlohmann::json slice = {};
267 slice["name"] = s.fName;
268 slice["threshold"] = s.fThreshold;
269 slice["color"] = s.fColor;
270 sarr.push_back(slice);
271 }
272 j["sliceInfos"] = sarr;
273 return ret;
274}
275////////////////////////////////////////////////////////////////////////////////
276
277void REveCaloData::FillExtraSelectionData(nlohmann::json& j, const std::set<int>& secondary_idcs) const
278{
279 vCellId_t cells;
280
281 if (fSelector)
282 {
283 fSelector->GetCellsFromSecondaryIndices(secondary_idcs, cells);
284 }
285 else
286 {
287 for (auto &id : secondary_idcs ) {
288
289 int s = (id >> 24);
290 int t = id & 0xffffff;
291 REveCaloData::CellId_t cell(t, s, 1.0f);
292 cells.push_back(cell);
293 }
294 }
295 for (auto &c : fNieces)
296 ((REveCaloViz*)c)->WriteCoreJsonSelection(j, cells);
297}
298
299////////////////////////////////////////////////////////////////////////////////
300
302{
303 using namespace TMath;
304
305 if (eta < 0)
306 return Pi() - 2*ATan(Exp(- Abs(eta)));
307 else
308 return 2*ATan(Exp(- Abs(eta)));
309}
310
311////////////////////////////////////////////////////////////////////////////////
312std::string REveCaloData::GetHighlightTooltip(const std::set<int>& secondary_idcs) const
313{
314 std::string s;
315 CellData_t cellData;
316
317 Bool_t single = secondary_idcs.size() == 1;
318 Float_t sum = 0;
319
320
321 for (auto &id : secondary_idcs ) {
322
323 int slice = (id >> 24);
324 int tower = id & 0xffffff;
325 REveCaloData::CellId_t cell(tower, slice, 1.0f);
326 GetCellData(cell, cellData);
327
328 s += TString::Format("%s %.2f (%.3f, %.3f)",
329 fSliceInfos[slice].fName.Data(), cellData.fValue,
330 cellData.Eta(), cellData.Phi());
331
332 if (single) return s;
333 s += "\n";
334 sum += cellData.fValue;
335 }
336 s += TString::Format("Sum = %.2f", sum);
337 return s;
338}
339
340/** \class REveCaloDataVec
341\ingroup REve
342Calo data for universal cell geometry.
343*/
344
345////////////////////////////////////////////////////////////////////////////////
346
348 REveCaloData(),
349
350 fTower(0),
351 fEtaMin( 1e3),
352 fEtaMax(-1e3),
353 fPhiMin( 1e3),
354 fPhiMax(-1e3)
355{
356 // Constructor.
357
358 fSliceInfos.assign(nslices, SliceInfo_t());
359
360 fSliceVec.assign(nslices, std::vector<Float_t> ());
361}
362
363////////////////////////////////////////////////////////////////////////////////
364/// Destructor.
365
367{
368}
369
370////////////////////////////////////////////////////////////////////////////////
371/// Add new slice.
372
374{
375 fSliceInfos.push_back(SliceInfo_t());
376 fSliceVec.push_back(std::vector<Float_t> ());
377 fSliceVec.back().resize(fGeomVec.size(), 0.f);
378
379 return fSliceInfos.size() - 1;
380}
381
382////////////////////////////////////////////////////////////////////////////////
383/// Add tower within eta/phi range.
384
386{
387 assert (etaMin < etaMax);
388 assert (phiMin < phiMax);
389
390 fGeomVec.push_back(CellGeom_t(etaMin, etaMax, phiMin, phiMax));
391
392 for (vvFloat_i it=fSliceVec.begin(); it!=fSliceVec.end(); ++it)
393 (*it).push_back(0);
394
395 if (etaMin < fEtaMin) fEtaMin = etaMin;
396 if (etaMax > fEtaMax) fEtaMax = etaMax;
397
398 if (phiMin < fPhiMin) fPhiMin = phiMin;
399 if (phiMax > fPhiMax) fPhiMax = phiMax;
400
401 fTower = fGeomVec.size() - 1;
402 return fTower;
403}
404
405////////////////////////////////////////////////////////////////////////////////
406/// Fill given slice in the current tower.
407
409{
410 fSliceVec[slice][fTower] = val;
411}
412
413////////////////////////////////////////////////////////////////////////////////
414/// Fill given slice in a given tower.
415
417{
418 fSliceVec[slice][tower] = val;
419}
420
421
422////////////////////////////////////////////////////////////////////////////////
423/// Get list of cell-ids for given eta/phi range.
424
426 Float_t phi, Float_t phiD,
427 REveCaloData::vCellId_t &out) const
428{
429 using namespace TMath;
430
431 Float_t etaMin = eta - etaD*0.5;
432 Float_t etaMax = eta + etaD*0.5;
433
434 Float_t phiMin = phi - phiD*0.5;
435 Float_t phiMax = phi + phiD*0.5;
436
437 Int_t nS = fSliceVec.size();
438
439 Int_t tower = 0;
440 Float_t fracx=0, fracy=0, frac;
441 Float_t minQ, maxQ;
442
443 for(vCellGeom_ci i=fGeomVec.begin(); i!=fGeomVec.end(); i++)
444 {
445 const CellGeom_t &cg = *i;
446 fracx = REveUtil::GetFraction(etaMin, etaMax, cg.fEtaMin, cg.fEtaMax);
447 if (fracx > 1e-3)
448 {
449 minQ = cg.fPhiMin;
450 maxQ = cg.fPhiMax;
451
452 if (fWrapTwoPi)
453 {
454 if (maxQ < phiMin)
455 {
456 minQ += TwoPi(); maxQ += TwoPi();
457 }
458 else if (minQ > phiMax)
459 {
460 minQ -= TwoPi(); maxQ -= TwoPi();
461 }
462 }
463
464 if (maxQ >= phiMin && minQ <= phiMax)
465 {
466 fracy = REveUtil::GetFraction(phiMin, phiMax, minQ, maxQ);
467 if (fracy > 1e-3)
468 {
469 frac = fracx*fracy;
470 for (Int_t s=0; s<nS; s++)
471 {
472 if (fSliceVec[s][tower] > fSliceInfos[s].fThreshold)
473 out.push_back(CellId_t(tower, s, frac));
474 }
475 }
476 }
477 }
478 tower++;
479 }
480}
481
482////////////////////////////////////////////////////////////////////////////////
483/// Rebin cells.
484
485void REveCaloDataVec::Rebin(TAxis* ax, TAxis* ay, vCellId_t &ids, Bool_t et, RebinData_t& rdata) const
486{
487 rdata.fNSlices = GetNSlices();
488 rdata.fBinData.assign((ax->GetNbins()+2)*(ay->GetNbins()+2), -1);
489
490 CellData_t cd;
491 for (vCellId_i it = ids.begin(); it != ids.end(); ++it)
492 {
493 GetCellData(*it, cd);
494 Int_t iMin = ax->FindBin(cd.EtaMin());
495 Int_t iMax = ax->FindBin(cd.EtaMax());
496 Int_t jMin = ay->FindBin(cd.PhiMin());
497 Int_t jMax = ay->FindBin(cd.PhiMax());
498 for (Int_t i = iMin; i <= iMax; ++i)
499 {
500 if (i < 0 || i > ax->GetNbins()) continue;
501 for (Int_t j = jMin; j <= jMax; ++j)
502 {
503 if (j < 0 || j > ay->GetNbins()) continue;
504
505 Double_t ratio = REveUtil::GetFraction(ax->GetBinLowEdge(i), ax->GetBinUpEdge(i), cd.EtaMin(), cd.EtaMax())
506 * REveUtil::GetFraction(ay->GetBinLowEdge(j), ay->GetBinUpEdge(j), cd.PhiMin(), cd.PhiMax());
507
508 if (ratio > 1e-6f)
509 {
510 Float_t* slices = rdata.GetSliceVals(i + j*(ax->GetNbins()+2));
511 slices[(*it).fSlice] += ratio * cd.Value(et);
512 }
513 }
514 }
515 }
516}
517
518////////////////////////////////////////////////////////////////////////////////
519/// Get cell geometry and value from cell ID.
520
522 REveCaloData::CellData_t& cellData) const
523{
524 cellData.CellGeom_t::operator=( fGeomVec[id.fTower] );
525 cellData.fValue = fSliceVec[id.fSlice][id.fTower];
526}
527
528////////////////////////////////////////////////////////////////////////////////
529/// Update limits and notify data users.
530
532{
533 using namespace TMath;
534
535 // update max E/Et values
536
537 fMaxValE = 0;
538 fMaxValEt = 0;
539 Float_t sum=0;
540 // printf("geom vec %d slices %d\n",fGeomVec.size(), fSliceVec.size() );
541
542 for (UInt_t tw=0; tw<fGeomVec.size(); tw++)
543 {
544 sum=0;
545 for (vvFloat_i it=fSliceVec.begin(); it!=fSliceVec.end(); ++it)
546 sum += (*it)[tw];
547
548 if (sum > fMaxValEt ) fMaxValEt=sum;
549
550 sum /= Abs(Sin(EtaToTheta(fGeomVec[tw].Eta())));
551
552 if (sum > fMaxValE) fMaxValE=sum;
553 }
554
556}
557
558
559////////////////////////////////////////////////////////////////////////////////
560/// Set XY axis from cells geometry.
561
563{
564 std::vector<Double_t> binX;
565 std::vector<Double_t> binY;
566
567 for(vCellGeom_ci i=fGeomVec.begin(); i!=fGeomVec.end(); i++)
568 {
569 const CellGeom_t &ch = *i;
570
571 binX.push_back(ch.EtaMin());
572 binX.push_back(ch.EtaMax());
573 binY.push_back(ch.PhiMin());
574 binY.push_back(ch.PhiMax());
575 }
576
577 std::sort(binX.begin(), binX.end());
578 std::sort(binY.begin(), binY.end());
579
580 Int_t cnt = 0;
581 Double_t sum = 0;
582 Double_t val;
583
584 // X axis
585 Double_t dx = binX.back() - binX.front();
586 epsX *= dx;
587 std::vector<Double_t> newX;
588 newX.push_back(binX.front()); // underflow
589 Int_t nX = binX.size()-1;
590 for(Int_t i=0; i<nX; i++)
591 {
592 val = (sum +binX[i])/(cnt+1);
593 if (binX[i+1] -val > epsX)
594 {
595 newX.push_back(val);
596 cnt = 0;
597 sum = 0;
598 }
599 else
600 {
601 sum += binX[i];
602 cnt++;
603 }
604 }
605 newX.push_back(binX.back()); // overflow
606
607 // Y axis
608 cnt = 0;
609 sum = 0;
610 std::vector<Double_t> newY;
611 Double_t dy = binY.back() - binY.front();
612 epsY *= dy;
613 newY.push_back(binY.front());// underflow
614 Int_t nY = binY.size()-1;
615 for(Int_t i=0 ; i<nY; i++)
616 {
617 val = (sum +binY[i])/(cnt+1);
618 if (binY[i+1] -val > epsY )
619 {
620 newY.push_back(val);
621 cnt = 0;
622 sum = 0;
623 }
624 else
625 {
626 sum += binY[i];
627 cnt++;
628 }
629
630 }
631 newY.push_back(binY.back()); // overflow
632 fEtaAxis = std::make_unique<TAxis>(newX.size()-1, &newX[0]);
633 fEtaAxis = std::make_unique<TAxis>(newY.size()-1, &newY[0]);
634 fEtaAxis->SetNdivisions(510);
635 fPhiAxis->SetNdivisions(510);
636}
637
638/** \class REveCaloDataHist
639\ingroup REve
640A central manager for calorimeter data of an event written in TH2F.
641X axis is used for eta and Y axis for phi.
642*/
643
644
645////////////////////////////////////////////////////////////////////////////////
646/// Constructor.
647
649 REveCaloData(),
650
651 fHStack(0)
652{
653 fHStack = new THStack();
654 fEps = 1e-5;
655}
656
657////////////////////////////////////////////////////////////////////////////////
658/// Destructor.
659
661{
662 delete fHStack;
663}
664
665////////////////////////////////////////////////////////////////////////////////
666/// Update limits and notify data users.
667
669{
670 using namespace TMath;
671
672 // update max E/Et values
673 fMaxValE = 0;
674 fMaxValEt = 0;
675
676 if (GetNSlices() < 1) return;
677
678 TH2* hist = GetHist(0);
679 for (Int_t ieta = 1; ieta <= fEtaAxis->GetNbins(); ++ieta)
680 {
681 Double_t eta = fEtaAxis->GetBinCenter(ieta); // conversion E/Et
682 for (Int_t iphi = 1; iphi <= fPhiAxis->GetNbins(); ++iphi)
683 {
684 Double_t value = 0;
685 for (Int_t i = 0; i < GetNSlices(); ++i)
686 {
687 hist = GetHist(i);
688 Int_t bin = hist->GetBin(ieta, iphi);
689 value += hist->GetBinContent(bin);
690 }
691
692 if (value > fMaxValEt ) fMaxValEt = value;
693
694 value /= Abs(Sin(EtaToTheta(eta)));
695
696 if (value > fMaxValE) fMaxValE = value;
697 }
698 }
700}
701
702////////////////////////////////////////////////////////////////////////////////
703/// Get list of cell IDs in given eta and phi range.
704
706 Float_t phi, Float_t phiD,
707 REveCaloData::vCellId_t &out) const
708{
709 using namespace TMath;
710
711 Float_t etaMin = eta - etaD*0.5 -fEps;
712 Float_t etaMax = eta + etaD*0.5 +fEps;
713
714 Float_t phiMin = phi - phiD*0.5 -fEps;
715 Float_t phiMax = phi + phiD*0.5 +fEps;
716
717 Int_t nEta = fEtaAxis->GetNbins();
718 Int_t nPhi = fPhiAxis->GetNbins();
719 Int_t nSlices = GetNSlices();
720
721 Int_t bin = 0;
722
723 Bool_t accept;
724 for (Int_t ieta = 1; ieta <= nEta; ++ieta)
725 {
726 if (fEtaAxis->GetBinLowEdge(ieta) >= etaMin && fEtaAxis->GetBinUpEdge(ieta) <= etaMax)
727 {
728 for (Int_t iphi = 1; iphi <= nPhi; ++iphi)
729 {
730 if (fWrapTwoPi )
731 {
733 (phiMin, phiMax, fPhiAxis->GetBinLowEdge(iphi), fPhiAxis->GetBinUpEdge(iphi));
734 }
735 else
736 {
737 accept = fPhiAxis->GetBinLowEdge(iphi) >= phiMin && fPhiAxis->GetBinUpEdge(iphi) <= phiMax &&
738 fPhiAxis->GetBinLowEdge(iphi) >= phiMin && fPhiAxis->GetBinUpEdge(iphi) <= phiMax;
739 }
740
741 if (accept)
742 {
743 for (Int_t s = 0; s < nSlices; ++s)
744 {
745 TH2F *hist = GetHist(s);
746 bin = hist->GetBin(ieta, iphi);
747 if (hist->GetBinContent(bin) > fSliceInfos[s].fThreshold)
748 out.push_back(REveCaloData::CellId_t(bin, s));
749 } // hist slices
750 }
751 } // phi bins
752 }
753 } // eta bins
754}
755
756////////////////////////////////////////////////////////////////////////////////
757/// Rebin
758
760{
761 rdata.fNSlices = GetNSlices();
762 rdata.fBinData.assign((ax->GetNbins()+2)*(ay->GetNbins()+2), -1);
764 Float_t *val;
765 Int_t i, j, w;
766 Int_t binx, biny;
767 Int_t bin;
768
769 for (vCellId_i it=ids.begin(); it!=ids.end(); ++it)
770 {
771 GetCellData(*it, cd);
772 GetHist(it->fSlice)->GetBinXYZ((*it).fTower, i, j, w);
773 binx = ax->FindBin(fEtaAxis->GetBinCenter(i));
774 biny = ay->FindBin(fPhiAxis->GetBinCenter(j));
775 bin = biny*(ax->GetNbins()+2)+binx;
776 val = rdata.GetSliceVals(bin);
777 Double_t ratio = REveUtil::GetFraction(ax->GetBinLowEdge(binx), ax->GetBinUpEdge(binx), cd.EtaMin(), cd.EtaMax())
778 * REveUtil::GetFraction(ay->GetBinLowEdge(biny), ay->GetBinUpEdge(biny), cd.PhiMin(), cd.PhiMax());
779
780 val[(*it).fSlice] += cd.Value(et)*ratio;
781 }
782}
783
784////////////////////////////////////////////////////////////////////////////////
785/// Get cell geometry and value from cell ID.
786
788 REveCaloData::CellData_t& cellData) const
789{
790 TH2F* hist = GetHist(id.fSlice);
791
792 Int_t x, y, z;
793 hist->GetBinXYZ(id.fTower, x, y, z);
794
795 cellData.fValue = hist->GetBinContent(id.fTower);
796 cellData.Configure(hist->GetXaxis()->GetBinLowEdge(x),
797 hist->GetXaxis()->GetBinUpEdge(x),
798 hist->GetYaxis()->GetBinLowEdge(y),
799 hist->GetYaxis()->GetBinUpEdge(y));
800}
801
802////////////////////////////////////////////////////////////////////////////////
803/// Add new slice to calo tower. Updates cached variables fMaxValE
804/// and fMaxValEt
805/// Return last index in the vector of slice infos.
806
808{
809 if (!fEtaAxis) {
810 fEtaAxis.reset(new TAxis(*hist->GetXaxis()));
811 fPhiAxis.reset(new TAxis(*hist->GetYaxis()));
812 }
813 fHStack->Add(hist);
814 fSliceInfos.push_back(SliceInfo_t());
815 fSliceInfos.back().fName = hist->GetName();
816 fSliceInfos.back().fColor = hist->GetLineColor();
817
818 DataChanged();
819
820 return fSliceInfos.size() - 1;
821}
822
823////////////////////////////////////////////////////////////////////////////////
824/// Get histogram in given slice.
825
827{
828 assert(slice >= 0 && slice < fHStack->GetHists()->GetSize());
829 return (TH2F*) fHStack->GetHists()->At(slice);
830}
831
832////////////////////////////////////////////////////////////////////////////////
833/// Get eta limits.
834
836{
837 min = fEtaAxis->GetXmin();
838 max = fEtaAxis->GetXmax();
839}
840
841////////////////////////////////////////////////////////////////////////////////
842/// Get phi limits.
843
845{
846 min = fPhiAxis->GetXmin();
847 max = fPhiAxis->GetXmax();
848}
849
850////////////////////////////////////////////////////////////////////////////////
851/// Process selection. Called from REveCaloViz object
852
854{
855 // only one slice can be user selected at once
856 fActiveSlice = sel_cells.front().fSlice;
857 for (auto &si : fSliceSelectors)
858 {
859 if (si->GetSliceIndex() == fActiveSlice) {
860 si->ProcessSelection(sel_cells, selectionId, multi);
861 break;
862 }
863 }
864}
865
866////////////////////////////////////////////////////////////////////////////////
867/// GetCellsFromSecondaryIndices used in implied selection
868
870{
871 for (auto &si : fSliceSelectors)
872 {
873 if (si->GetSliceIndex() == fActiveSlice) {
874 si->GetCellsFromSecondaryIndices(idcs, out);
875 break;
876 }
877 }
878}
#define c(i)
Definition RSha256.hxx:101
#define e(i)
Definition RSha256.hxx:103
short Color_t
Definition RtypesCore.h:92
char Char_t
Definition RtypesCore.h:37
constexpr Bool_t kTRUE
Definition RtypesCore.h:100
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
Definition TError.cxx:197
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
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.
TH2F * GetHist(Int_t slice) const
Get histogram in given slice.
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.
virtual void GetPhiLimits(Double_t &min, Double_t &max) const
Get phi limits.
virtual void GetCellData(const REveCaloData::CellId_t &id, REveCaloData::CellData_t &data) const
Get cell geometry and value from cell ID.
virtual void GetEtaLimits(Double_t &min, Double_t &max) const
Get eta limits.
virtual void DataChanged()
Update limits and notify data users.
void GetCellsFromSecondaryIndices(const std::set< int > &, REveCaloData::vCellId_t &out)
GetCellsFromSecondaryIndices used in implied selection.
void ProcessSelection(REveCaloData::vCellId_t &sel_cells, UInt_t selectionId, Bool_t multi)
Process selection. Called from REveCaloViz object.
std::vector< std::unique_ptr< REveCaloDataSliceSelector > > fSliceSelectors
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.
Int_t AddTower(Float_t etaMin, Float_t etaMax, Float_t phiMin, Float_t phiMax)
Add tower within eta/phi range.
virtual void Rebin(TAxis *ax, TAxis *ay, vCellId_t &in, Bool_t et, RebinData_t &out) const
Rebin cells.
void SetAxisFromBins(Double_t epsX=0.001, Double_t epsY=0.001)
Set XY axis from cells geometry.
virtual void DataChanged()
Update limits and notify data users.
std::vector< vFloat_t >::iterator vvFloat_i
virtual void GetCellData(const REveCaloData::CellId_t &id, REveCaloData::CellData_t &data) const
Get cell geometry and value from cell ID.
void FillSlice(Int_t slice, Float_t value)
Fill given slice in the current tower.
Color_t GetSliceColor(Int_t slice) const
Get color for given slice.
Char_t GetSliceTransparency(Int_t slice) const
Get transparency for given slice.
std::vector< CellId_t >::iterator vCellId_i
void SetSliceThreshold(Int_t slice, Float_t threshold)
Set threshold for given slice.
virtual void GetCellData(const CellId_t &id, CellData_t &data) const =0
std::vector< CellGeom_t >::const_iterator vCellGeom_ci
void SetSliceTransparency(Int_t slice, Char_t t)
Set transparency for given slice.
static Float_t EtaToTheta(Float_t eta)
std::unique_ptr< REveCaloDataSelector > fSelector
REveCaloData(const char *n="REveCaloData", const char *t="")
std::unique_ptr< TAxis > fEtaAxis
void FillExtraSelectionData(nlohmann::json &, const std::set< int > &) const override
std::string GetHighlightTooltip(const std::set< int > &secondary_idcs) const override
virtual void DataChanged()
Tell users (REveCaloViz instances using this data) that data has changed and they should update the l...
std::unique_ptr< TAxis > fPhiAxis
void ProcessSelection(vCellId_t &sel_cells, UInt_t selectionId, Bool_t multi)
Process newly selected cells with given select-record.
Float_t GetSliceThreshold(Int_t slice) const
Get threshold for given slice.
virtual void InvalidateUsersCellIdCache()
Invalidate cell ids cache on back ptr references.
std::vector< CellId_t > vCellId_t
void SetSliceColor(Int_t slice, Color_t col)
Set color for given slice.
Int_t WriteCoreJson(nlohmann::json &j, Int_t rnr_offset) override
Write core json.
void FillImpliedSelectedSet(Set_t &impSelSet, const std::set< int > &sec_idcs) override
Populate set impSelSet with derived / dependant elements.
void DataChanged()
Update setting and cache on data changed.
Definition REveCalo.cxx:257
void SetNameTitle(const std::string &name, const std::string &title)
Set name and title of an element.
virtual Int_t WriteCoreJson(nlohmann::json &cj, Int_t rnr_offset)
Write core json.
virtual void AddStamp(UChar_t bits)
Add (bitwise or) given stamps to fChangeBits.
std::set< REveElement * > Set_t
ElementId_t GetElementId() const
REveElement * FindElementById(ElementId_t id) const
Lookup ElementId in element map and return corresponding REveElement*.
REveSelection Container for selected and highlighted elements.
void NewElementPicked(ElementId_t id, bool multi, bool secondary, const std::set< int > &secondary_idcs={})
Called from GUI when user picks or un-picks an element.
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 REveUtil.cxx:312
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 REveUtil.cxx:272
virtual Color_t GetLineColor() const
Return the line color.
Definition TAttLine.h:33
Class to manage histogram axis.
Definition TAxis.h:30
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
Int_t GetNbins() const
Definition TAxis.h:121
virtual Double_t GetBinUpEdge(Int_t bin) const
Return up edge of bin.
Definition TAxis.cxx:528
TAxis * GetXaxis()
Definition TH1.h:322
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:4938
TAxis * GetYaxis()
Definition TH1.h:323
2-D histogram with a float per channel (see TH1 documentation)}
Definition TH2.h:257
Service class for 2-D histogram classes.
Definition TH2.h:30
Int_t GetBin(Int_t binx, Int_t biny, Int_t binz=0) const override
Return Global bin number corresponding to binx,y,z.
Definition TH2.cxx:1040
Double_t GetBinContent(Int_t binx, Int_t biny) const override
Definition TH2.h:89
The Histogram stack class.
Definition THStack.h:38
TList * GetHists() const
Definition THStack.h:70
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:364
TObject * At(Int_t idx) const override
Returns the object at position idx. Returns 0 if idx is out of range.
Definition TList.cxx:357
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
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:2356
const char * GetHist()
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
R__EXTERN REveManager * gEve
TMath.
Definition TMathBase.h:35
Double_t Sin(Double_t)
Returns the sine of an angle of x radians.
Definition TMath.h:586
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123
constexpr Double_t TwoPi()
Definition TMath.h:44
basic_json<> json
virtual void Dump() const
Print member data.
Float_t Value(Bool_t) const
Return energy value associated with the cell, usually Et.
void Configure(Float_t etaMin, Float_t etaMax, Float_t phiMin, Float_t phiMax)
virtual void Dump() const
Print member data.
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345