Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
collection_proxies.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_eve7
3///
4/// This is an example of visualization of containers
5/// with REveDataCollection and REveDataProxyBuilders.
6/// \macro_code
7///
8
9
12#include "ROOT/REveManager.hxx"
15#include <ROOT/REveGeoShape.hxx>
16#include <ROOT/REveJetCone.hxx>
17#include <ROOT/REvePointSet.hxx>
20#include <ROOT/REveScene.hxx>
23#include <ROOT/REveTrack.hxx>
25#include <ROOT/REveViewer.hxx>
27#include <ROOT/REveBoxSet.hxx>
29#include <ROOT/REveCalo.hxx>
30
31#include "TGeoTube.h"
32#include "TROOT.h"
33#include "TList.h"
34#include "TParticle.h"
35#include "TRandom.h"
36#include "TApplication.h"
37#include "TFile.h"
38#include "TH2F.h"
39#include <iostream>
40
41
42const Double_t kR_min = 299;
43const Double_t kR_max = 300;
44const Double_t kZ_d = 500;
45
46
47namespace fw3dlego {
48 const int xbins_n = 83;
49 const double xbins[xbins_n] = {
50 -5.191, -4.889, -4.716, -4.538, -4.363, -4.191, -4.013, -3.839, -3.664, -3.489, -3.314, -3.139, -2.964, -2.853,
51 -2.650, -2.500, -2.322, -2.172, -2.043, -1.930, -1.830, -1.740, -1.653, -1.566, -1.479, -1.392, -1.305, -1.218,
52 -1.131, -1.044, -0.957, -0.870, -0.783, -0.696, -0.609, -0.522, -0.435, -0.348, -0.261, -0.174, -0.087, 0.000,
53 0.087, 0.174, 0.261, 0.348, 0.435, 0.522, 0.609, 0.696, 0.783, 0.870, 0.957, 1.044, 1.131, 1.218,
54 1.305, 1.392, 1.479, 1.566, 1.653, 1.740, 1.830, 1.930, 2.043, 2.172, 2.322, 2.500, 2.650, 2.853,
55 2.964, 3.139, 3.314, 3.489, 3.664, 3.839, 4.013, 4.191, 4.363, 4.538, 4.716, 4.889, 5.191};
56} // namespace fw3dlego
57
58
61using namespace ROOT::Experimental;
62
63//==============================================================================
64//============== EMULATE FRAMEWORK CLASSES =====================================
65//==============================================================================
66
67
68// a demo class, can be provided from experiment framework
69class Jet : public TParticle
70{
71public:
72 float fEtaSize{0};
73 float fPhiSize{0};
74
75 float GetEtaSize() const { return fEtaSize; }
76 float GetPhiSize() const { return fPhiSize; }
77 void SetEtaSize(float iEtaSize) { fEtaSize = iEtaSize; }
78 void SetPhiSize(float iPhiSize) { fPhiSize = iPhiSize; }
79
80 Jet(Int_t pdg, Int_t status, Int_t mother1, Int_t mother2, Int_t daughter1, Int_t daughter2,
81 Double_t px, Double_t py, Double_t pz, Double_t etot) :
82 TParticle(pdg, status, mother1, mother2, daughter1, daughter2, px, py, pz, etot, 0, 0, 0, 0)
83 {}
84
86};
87
88class RecHit : public TObject
89{
90public:
91 float fX{0};
92 float fY{0};
93 float fZ{0};
94 float fPt{0};
95
96 RecHit(float pt, float x, float y, float z): fPt(pt), fX(x), fY(y), fZ(z) {}
98};
99
100class RCaloTower : public TObject
101{
102public:
103 float fEta{0};
104 float fPhi{0};
105 float fEt{0};
106
107 RCaloTower(float eta, float phi, float et): fEta(eta), fPhi(phi), fEt(et) {}
109};
110
112{
113private:
116
117public:
119
121 void ProcessSelection(REveCaloData::vCellId_t& sel_cells, UInt_t selectionId, Bool_t multi) override
122 {
123 std::set<int> item_set;
125 for (auto &cellId : sel_cells)
126 {
127 fCaloData->GetCellData(cellId, cd);
128
129 // loop over enire collection and check its eta/phi range
130 for (int t = 0; t < fCollection->GetNItems(); ++t)
131 {
133 if (tower->fEta > cd.fEtaMin && tower->fEta < cd.fEtaMax &&
134 tower->fPhi > cd.fPhiMin && tower->fPhi < cd.fPhiMax)
135 item_set.insert(t);
136 }
137 }
139 fCollection->GetItemList()->RefSelectedSet() = item_set;
140 sel->NewElementPicked(fCollection->GetItemList()->GetElementId(), multi, true, item_set);
141 }
142
144 void GetCellsFromSecondaryIndices(const std::set<int>& idcs, REveCaloData::vCellId_t& out) override
145 {
147 std::set<int> cbins;
148 float total = 0;
149 for( auto &i : idcs ) {
151 int bin = hist->FindBin(tower->fEta, tower->fPhi);
152 float frac = tower->fEt/hist->GetBinContent(bin);
153 bool ex = false;
154 for (size_t ci = 0; ci < out.size(); ++ci)
155 {
156 if (out[ci].fTower == bin && out[ci].fSlice == GetSliceIndex())
157 {
158 float oldv = out[ci].fFraction;
159 out[ci].fFraction = oldv + frac;
160 ex = true;
161 break;
162 }
163 }
164 if (!ex) {
165 out.push_back(REveCaloData::CellId_t(bin, GetSliceIndex(), frac));
166 }
167 }
168 }
169};
170
171class Event
172{
173public:
174 int eventId{0};
175 int N_tracks{0};
176 int N_jets{0};
177 std::vector<TList*> fListData;
178
180
182 {
183 auto baseHist = new TH2F("dummy", "dummy", fw3dlego::xbins_n - 1, fw3dlego::xbins, 72, -TMath::Pi(), TMath::Pi());
185 fCaloData->AddHistogram(baseHist);
186
187 auto selector = new REveCaloDataSelector();
188 fCaloData->SetSelector(selector);
189
191 }
192
193 void MakeJets(int N)
194 {
195 TRandom &r = *gRandom;
196 r.SetSeed(0);
197 TList* list = new TList();
198 list->SetName("Jets");
199 for (int i = 1; i <= N; ++i)
200 {
201 double pt = r.Uniform(0.5, 10);
202 double eta = r.Uniform(-2.55, 2.55);
203 double phi = r.Uniform(-TMath::Pi(), TMath::Pi());
204
205 double px = pt * std::cos(phi);
206 double py = pt * std::sin(phi);
207 double pz = pt * (1. / (std::tan(2*std::atan(std::exp(-eta)))));
208
209 auto jet = new Jet(0, 0, 0, 0, 0, 0, px, py, pz, std::sqrt(px*px + py*py + pz*pz + 80*80));
210 jet->SetEtaSize(r.Uniform(0.02, 0.2));
211 jet->SetPhiSize(r.Uniform(0.01, 0.3));
212 list->Add(jet);
213 }
214 fListData.push_back(list);
215 }
216
217 void MakeParticles(int N)
218 {
219 TRandom &r = *gRandom;
220 r.SetSeed(0);
221 TList* list = new TList();
222 list->SetName("Tracks");
223 for (int i = 1; i <= N; ++i)
224 {
225 double pt = r.Uniform(0.5, 10);
226 double eta = r.Uniform(-2.55, 2.55);
227 double phi = r.Uniform(0, TMath::TwoPi());
228
229 double px = pt * std::cos(phi);
230 double py = pt * std::sin(phi);
231 double pz = pt * (1. / (std::tan(2*std::atan(std::exp(-eta)))));
232
233 // printf("Event::MakeParticles %2d: pt=%.2f, eta=%.2f, phi=%.2f\n", i, pt, eta, phi);
234 auto particle = new TParticle(0, 0, 0, 0, 0, 0,
235 px, py, pz, std::sqrt(px*px + py*py + pz*pz + 80*80),
236 0, 0, 0, 0 );
237
238 int pdg = 11 * (r.Integer(2) > 0 ? 1 : -1);
239 particle->SetPdgCode(pdg);
240
241 list->Add(particle);
242 }
243 fListData.push_back(list);
244 }
245
246 void MakeRecHits(int N)
247 {
248 TRandom &r = *gRandom;
249 r.SetSeed(0);
250 TList* list = new TList();
251 list->SetName("RecHits");
252
253 for (int i = 1; i <= N; ++i)
254 {
255 float pt = r.Uniform(0.5, 10);
256 float x = r.Uniform(-200, 200);
257 float y = r.Uniform(-200, 200);
258 float z = r.Uniform(-500, 500);
259 auto rechit = new RecHit(pt, x, y, z);
260 list->Add(rechit);
261 }
262 fListData.push_back(list);
263 }
264
265 void Clear()
266 {
267 for (auto &l : fListData)
268 delete l;
269 fListData.clear();
270 }
271
272 void Create()
273 {
274 Clear();
275 MakeJets(4);
276 MakeParticles(100);
277 MakeRecHits(20);
278
279 // refill calo data from jet list
280 TList* jlist = fListData[0];
281 TList* elist = new TList();
282 elist->SetName("ECAL");
283 fListData.push_back(elist);
284 TList* hlist = new TList();
285 hlist->SetName("HCAL");
286 fListData.push_back(hlist);
287 for (int i = 0; i <= jlist->GetLast(); ++i) {
288 const Jet* j = (Jet*)jlist->At(i);
289 float offX = j->Eta();
290 float offY = j->Phi() > TMath::Pi() ? j->Phi() - TMath::TwoPi() : j->Phi();
291 for (int k=0; k<20; ++k) {
292 double x, y, v;
293 x = gRandom->Uniform(-j->GetEtaSize(), j->GetEtaSize());
294 y = gRandom->Uniform(-j->GetPhiSize(),j->GetPhiSize());
295 v = j->Pt();
296 auto etower = new RCaloTower(offX + x, offY + y, v + gRandom->Uniform(2,3));
297 elist->Add(etower);
298 auto htower = new RCaloTower(offX + x, offY + y, v + gRandom->Uniform(1,2));
299 hlist->Add(htower);
300 }
301 }
303 eventId++;
304 }
305};
306
307
308//==============================================================================
309//== PROXY BUILDERS ============================================================
310//==============================================================================
311
313{
314 bool HaveSingleProduct() const override { return false; }
315
317 void BuildItemViewType(const Jet& dj, int idx, REveElement* iItemHolder,
318 const std::string& viewType, const REveViewContext* context) override
319 {
320 auto jet = new REveJetCone();
321 jet->SetCylinder(context->GetMaxR(), context->GetMaxZ());
322 jet->AddEllipticCone(dj.Eta(), dj.Phi(), dj.GetEtaSize(), dj.GetPhiSize());
323 SetupAddElement(jet, iItemHolder, true);
324 jet->SetLineColor(jet->GetMainColor());
325
326 float size = 50.f * dj.Pt(); // values are saved in scale
327 double theta = dj.Theta();
328 // printf("%s jet theta = %f, phi = %f \n", iItemHolder->GetCName(), theta, dj.Phi());
329 double phi = dj.Phi();
330
331
332 if (viewType == "Projected" )
333 {
334 static const float_t offr = 6;
335 float r_ecal = context->GetMaxR() + offr;
336 float z_ecal = context->GetMaxZ() + offr;
337
338 float transAngle = abs(atan(r_ecal/z_ecal));
339 double r = 0;
340 bool debug = false;
341 if (theta < transAngle || 3.14-theta < transAngle)
342 {
343 z_ecal = context->GetMaxZ() + offr/transAngle;
344 r = z_ecal/fabs(cos(theta));
345 }
346 else
347 {
348 debug = true;
349 r = r_ecal/sin(theta);
350 }
351
352 REveVector p1(0, (phi<TMath::Pi() ? r*fabs(sin(theta)) : -r*fabs(sin(theta))), r*cos(theta));
353 REveVector p2(0, (phi<TMath::Pi() ? (r+size)*fabs(sin(theta)) : -(r+size)*fabs(sin(theta))), (r+size)*cos(theta));
354
355 auto marker = new REveScalableStraightLineSet("jetline");
356 marker->SetScaleCenter(p1.fX, p1.fY, p1.fZ);
357 marker->AddLine(p1, p2);
358 marker->SetLineWidth(4);
359 if (debug)
360 marker->AddMarker(0, 0.9);
361
362 SetupAddElement(marker, iItemHolder, true);
363 marker->SetName(Form("line %s %d", Collection()->GetCName(), idx));
364 }
365 }
366
367
369
370 void LocalModelChanges(int idx, REveElement* el, const REveViewContext* ctx) override
371 {
372 // printf("LocalModelChanges jet %s ( %s )\n", el->GetCName(), el->FirstChild()->GetCName());
373 REveJetCone* cone = dynamic_cast<REveJetCone*>(el->FirstChild());
374 cone->SetLineColor(cone->GetMainColor());
375 }
376};
377
378
380{
382
383 void BuildItem(const TParticle& p, int idx, REveElement* iItemHolder, const REveViewContext* context) override
384 {
385 const TParticle *x = &p;
386 auto track = new REveTrack((TParticle*)(x), 1, context->GetPropagator());
387 track->MakeTrack();
388 SetupAddElement(track, iItemHolder, true);
389 }
390};
391
392
394private:
395 class FWBoxSet : public REveBoxSet {
396 public:
397 using REveElement::GetSelectionMaster;
399 {
400 if (fSelectionMaster) {
403 return il;
404 }
405 return nullptr;
406 }
407 };
408
411 {
412 auto collection = Collection();
413 boxset->SetMainColor(collection->GetMainColor());
414 boxset->SetName(collection->GetCName());
415 boxset->SetPickable(true);
416 boxset->SetAlwaysSecSelect(true);
417 boxset->SetDetIdsAsSecondaryIndices(true);
418 boxset->SetSelectionMaster(((REveDataCollection *)collection)->GetItemList());
419 boxset->Reset(REveBoxSet::kBT_FreeBox, true, collection->GetNItems());
420 TRandom r(0);
421
422#define RND_BOX(x) (Float_t) r.Uniform(-(x), (x))
423 for (int h = 0; h < collection->GetNItems(); ++h) {
424 RecHit *hit = (RecHit *)collection->GetDataPtr(h);
425 const REveDataItem *item = Collection()->GetDataItem(h);
426
427 Float_t x = hit->fX;
428 Float_t y = hit->fY;
429 Float_t z = hit->fZ;
430 Float_t a = hit->fPt;
431 Float_t d = 0.05;
432 Float_t verts[24] = {x - a + RND_BOX(d), y - a + RND_BOX(d), z - a + RND_BOX(d), x - a + RND_BOX(d),
433 y + a + RND_BOX(d), z - a + RND_BOX(d), x + a + RND_BOX(d), y + a + RND_BOX(d),
434 z - a + RND_BOX(d), x + a + RND_BOX(d), y - a + RND_BOX(d), z - a + RND_BOX(d),
435 x - a + RND_BOX(d), y - a + RND_BOX(d), z + a + RND_BOX(d), x - a + RND_BOX(d),
436 y + a + RND_BOX(d), z + a + RND_BOX(d), x + a + RND_BOX(d), y + a + RND_BOX(d),
437 z + a + RND_BOX(d), x + a + RND_BOX(d), y - a + RND_BOX(d), z + a + RND_BOX(d)};
438 boxset->AddBox(verts);
439
440 boxset->DigitValue(item->GetVisible() ? 1 : 0);
441 if (item->GetVisible())
442 boxset->DigitColor(item->GetMainColor());
443 }
444 boxset->RefitPlex();
445 boxset->StampObjProps();
446 }
447
448public:
450 void BuildProduct(const REveDataCollection* collection, REveElement* product, const REveViewContext*)override
451 {
452 fBoxSet = new FWBoxSet();
454 product->AddElement(fBoxSet);
455 }
456
458 void FillImpliedSelected(REveElement::Set_t& impSet, const std::set<int>& sec_idcs, Product* p) override
459 {
460 // printf("RecHit fill implioed ----------------- !!!%zu\n", Collection()->GetItemList()->RefSelectedSet().size());
461 impSet.insert(fBoxSet);
462 }
463
465 void ModelChanges(const REveDataCollection::Ids_t& ids, Product* product) override
466 {
467 for (auto &i : ids)
468 {
469 auto digi = fBoxSet->GetDigit(i);
470 auto item = Collection()->GetDataItem(i);
472 if (item->GetVisible()) {
474 fBoxSet->DigitColor(item->GetMainColor());
475 } else {
477 }
478 }
480 }
481}; // RecHitProxyBuilder
482
484{
485private:
487 TH2F* fHist {nullptr};
488 int fSliceIndex {-1};
489
490 void assertSlice() {
491 if (!fHist) {
493
494 TH1::AddDirectory(kFALSE); //Keeps histogram from going into memory
495 fHist = new TH2F("caloHist", "caloHist", fw3dlego::xbins_n - 1, fw3dlego::xbins, 72, -M_PI, M_PI);
496 TH1::AddDirectory(status);
498
500 .Setup(Collection()->GetCName(),
501 0.,
502 Collection()->GetMainColor(),
503 Collection()->GetMainTransparency());
504
505 fCaloData->GetSelector()->AddSliceSelector(std::unique_ptr<REveCaloDataSliceSelector>
507 }
508 }
509
510public:
512
514 void BuildProduct(const REveDataCollection* collection, REveElement* product, const REveViewContext*)override
515 {
516 assertSlice();
517 fHist->Reset();
518 if (collection->GetRnrSelf())
519 {
521 .Setup(Collection()->GetCName(),
522 0.,
523 Collection()->GetMainColor(),
524 Collection()->GetMainTransparency());
525
526
527 for (int h = 0; h < collection->GetNItems(); ++h)
528 {
529 RCaloTower* tower = (RCaloTower*)collection->GetDataPtr(h);
530 const REveDataItem* item = Collection()->GetDataItem(h);
531
532 if (!item->GetVisible())
533 continue;
534 fHist->Fill(tower->fEta, tower->fPhi, tower->fEt);
535 }
536 }
538 }
539
541 void FillImpliedSelected(REveElement::Set_t& impSet, const std::set<int>& sec_idcs, Product*) override
542 {
544 impSet.insert(fCaloData);
545 fCaloData->FillImpliedSelectedSet(impSet, sec_idcs);
546 }
547
549 void ModelChanges(const REveDataCollection::Ids_t& ids, Product* product) override
550 {
551 BuildProduct(Collection(), nullptr, nullptr);
552 }
553
554}; // CaloTowerProxyBuilder
555
556//==============================================================================
557//== COLLECTION MANGER ================================================================
558//==============================================================================
559
561{
562private:
563 Event *fEvent{nullptr};
564
565 std::vector<REveScene *> m_scenes;
567
568 std::vector<REveDataProxyBuilderBase *> m_builders;
569
571 bool m_inEventLoading {false};
572
573public:
575 {
576 //view context
577 float r = 300;
578 float z = 300;
579 auto prop = new REveTrackPropagator();
580 prop->SetMagFieldObj(new REveMagFieldDuo(350, 3.5, -2.0));
581 prop->SetMaxR(r);
582 prop->SetMaxZ(z);
583 prop->SetMaxOrbs(6);
584 prop->IncRefCount();
585
589
590 // table specs
591 auto tableInfo = new REveTableViewInfo();
592
593 tableInfo->table("TParticle").
594 column("pt", 1, "i.Pt()").
595 column("eta", 3, "i.Eta()").
596 column("phi", 3, "i.Phi()");
597
598 tableInfo->table("Jet").
599 column("eta", 1, "i.Eta()").
600 column("phi", 1, "i.Phi()").
601 column("etasize", 2, "i.GetEtaSize()").
602 column("phisize", 2, "i.GetPhiSize()");
603
604 tableInfo->table("RecHit").
605 column("pt", 1, "i.fPt");
606
607 tableInfo->table("RCaloTower").
608 column("eta", 3, "i.fEta").
609 column("phi", 3, "i.fPhi").
610 column("Et", 3, "i.fEt");
611
613
614 for (auto &c : eveMng->GetScenes()->RefChildren()) {
615 if (c != eveMng->GetGlobalScene() && strncmp(c->GetCName(), "Geometry", 8) )
616 {
617 m_scenes.push_back((REveScene*)c);
618 }
619 if (!strncmp(c->GetCName(),"Table", 5))
620 c->AddElement(m_viewContext->GetTableViewInfo());
621
622 }
623
624 m_collections = eveMng->SpawnNewScene("Collections", "Collections");
625 }
626
628 {
629 for (auto &l : fEvent->fListData) {
630 TIter next(l);
631 if (collection->GetName() == std::string(l->GetName()))
632 {
633 collection->ClearItems();
634
635 for (int i = 0; i <= l->GetLast(); ++i)
636 {
637 std::string cname = collection->GetName();
638 auto len = cname.size();
639 char end = cname[len-1];
640 if (end == 's') {
641 cname = cname.substr(0, len-1);
642 }
643 TString pname(Form("%s %2d", cname.c_str(), i));
644 collection->AddItem(l->At(i), pname.Data(), "");
645 }
646 }
647 collection->ApplyFilter();
648 }
649 }
650
652 {
653 m_inEventLoading = true;
654
655 for (auto &el: m_collections->RefChildren())
656 {
657 auto c = dynamic_cast<REveDataCollection *>(el);
659 }
660
661 for (auto proxy : m_builders)
662 {
663 proxy->Build();
664 }
665
667 m_inEventLoading = false;
668 }
669
670 void addCollection(REveDataCollection* collection, REveDataProxyBuilderBase* glBuilder, bool showInTable = false)
671 {
672 m_collections->AddElement(collection);
673
674 // load data
675 SetDataItemsFromEvent(collection);
676 glBuilder->SetCollection(collection);
677 glBuilder->SetHaveAWindow(true);
678 for (auto scene : m_scenes)
679 {
680 if (strncmp(scene->GetCName(), "Tables", 5) == 0) continue;
681
682 REveElement *product = glBuilder->CreateProduct(scene->GetTitle(), m_viewContext);
683
684 if (!strncmp(scene->GetCTitle(), "Projected", 8))
685 {
686 g_projMng->ImportElements(product, scene);
687 }
688 else
689 {
690 scene->AddElement(product);
691 }
692 }
693 m_builders.push_back(glBuilder);
694 glBuilder->Build();
695
696 // Tables
697 auto tableBuilder = new REveTableProxyBuilder();
698 tableBuilder->SetHaveAWindow(true);
699 tableBuilder->SetCollection(collection);
700 REveElement* tablep = tableBuilder->CreateProduct("table-type", m_viewContext);
701 auto tableMng = m_viewContext->GetTableViewInfo();
702 if (showInTable)
703 {
704 tableMng->SetDisplayedCollection(collection->GetElementId());
705 }
706
707 for (auto s : m_scenes)
708 {
709 if (strncmp(s->GetCTitle(), "Table", 5) == 0)
710 {
711 s->AddElement(tablep);
712 tableBuilder->Build();
713 }
714 }
715 tableMng->AddDelegate([=]() { tableBuilder->ConfigChanged(); });
716 m_builders.push_back(tableBuilder);
717
718
719 // set tooltip expression for items
720 auto tableEntries = tableMng->RefTableEntries(collection->GetItemClass()->GetName());
721 int N = TMath::Min(int(tableEntries.size()), 3);
722 for (int t = 0; t < N; t++) {
723 auto te = tableEntries[t];
724 collection->GetItemList()->AddTooltipExpression(te.fName, te.fExpression);
725 }
726
727 collection->GetItemList()->SetItemsChangeDelegate([&] (REveDataItemList* collection, const REveDataCollection::Ids_t& ids)
728 {
729 this->ModelChanged( collection, ids );
730 });
731 collection->GetItemList()->SetFillImpliedSelectedDelegate([&] (REveDataItemList* collection, REveElement::Set_t& impSelSet, const std::set<int>& sec_idcs)
732 {
733 this->FillImpliedSelected( collection, impSelSet, sec_idcs);
734 });
735 }
736
738 {
739 auto mngTable = m_viewContext->GetTableViewInfo();
740 if (mngTable)
741 {
742 for (auto &el : m_collections->RefChildren())
743 {
744 if (el->GetName() == "Tracks")
745 mngTable->SetDisplayedCollection(el->GetElementId());
746 }
747 }
748 }
749
750
752 {
753 if (m_inEventLoading) return;
754
755 for (auto proxy : m_builders)
756 {
757 if (proxy->Collection()->GetItemList() == itemList)
758 {
759 // printf("Model changes check proxy %s: \n", proxy->Type().c_str());
760 proxy->ModelChanges(ids);
761 }
762 }
763 }
764
765 void FillImpliedSelected(REveDataItemList* itemList, REveElement::Set_t& impSelSet, const std::set<int>& sec_idcs)
766 {
767 if (m_inEventLoading) return;
768
769 for (auto proxy : m_builders)
770 {
771 if (proxy->Collection()->GetItemList() == itemList)
772 {
773 proxy->FillImpliedSelected(impSelSet, sec_idcs);
774 }
775 }
776 }
777
778};
779
780
781//==============================================================================
782//== Event Manager =============================================================
783//==============================================================================
784
786{
787private:
790
791public:
793
794 ~EventManager() override {}
795
796 virtual void NextEvent()
797 {
800 fEvent->Create();
801 fCMng->LoadEvent();
802 }
803};
804
806public:
808
810 bool DeviateSelection(REveSelection *selection, REveElement *el, bool multi, bool secondary,
811 const std::set<int> &secondary_idcs) override
812 {
813 if (el) {
814 auto *colItems = dynamic_cast<REveDataItemList *>(el);
815 if (colItems) {
816 // std::cout << "Deviate RefSelected=" << colItems->RefSelectedSet().size() << " passed set " << secondary_idcs.size() << "\n";
817 ExecuteNewElementPicked(selection, colItems, multi, true, colItems->RefSelectedSet());
818 return true;
819 }
820 }
821 return false;
822 }
823};
824//==============================================================================
825//== main() ====================================================================
826//==============================================================================
827
828void collection_proxies(bool proj=true)
829{
830 eveMng = REveManager::Create();
831 auto event = new Event();
832 event->Create();
833
834 // divert selection to map proxy builder products with collection
835 auto deviator = std::make_shared<FWSelectionDeviator>();
836 eveMng->GetSelection()->SetDeviator(deviator);
837 eveMng->GetHighlight()->SetDeviator(deviator);
838
839 // create scenes and views
840 REveScene* rhoZEventScene = nullptr;
841
842 auto b1 = new REveGeoShape("Barrel 1");
843 b1->SetShape(new TGeoTube(kR_min, kR_max, kZ_d));
844 b1->SetMainColor(kCyan);
846
847 rhoZEventScene = eveMng->SpawnNewScene("RhoZ Scene","Projected");
848 g_projMng = new REveProjectionManager(REveProjection::kPT_RhoZ);
850
851 auto rhoZView = eveMng->SpawnNewViewer("RhoZ View");
852 rhoZView->SetCameraType(REveViewer::kCameraOrthoXOY);
854 auto pgeoScene = eveMng->SpawnNewScene("Geometry projected");
855 rhoZView->AddScene(pgeoScene);
856 g_projMng->ImportElements(b1, pgeoScene);
857
858 auto tableScene = eveMng->SpawnNewScene ("Tables", "Tables");
859 auto tableView = eveMng->SpawnNewViewer("Table", "Table View");
860 tableView->AddScene(tableScene);
861
862 // create event data from list
863 auto collectionMng = new CollectionManager(event);
864
865 REveDataCollection* trackCollection = new REveDataCollection("Tracks");
866 trackCollection->SetItemClass(TParticle::Class());
867 trackCollection->SetMainColor(kGreen);
868 trackCollection->SetFilterExpr("i.Pt() > 4.1 && std::abs(i.Eta()) < 1");
869 collectionMng->addCollection(trackCollection, new TParticleProxyBuilder(), true);
870
871 REveDataCollection* jetCollection = new REveDataCollection("Jets");
872 jetCollection->SetItemClass(Jet::Class());
873 jetCollection->SetMainColor(kYellow);
874 jetCollection->SetFilterExpr("i.Pt() > 1");
875 collectionMng->addCollection(jetCollection, new JetProxyBuilder());
876
877 REveDataCollection* hitCollection = new REveDataCollection("RecHits");
878 hitCollection->SetItemClass(RecHit::Class());
879 hitCollection->SetMainColor(kOrange + 7);
880 hitCollection->SetFilterExpr("i.fPt > 5");
881 collectionMng->addCollection(hitCollection, new RecHitProxyBuilder(), true);
882
883 // add calorimeters
884 auto calo3d = new REveCalo3D(event->fCaloData);
885 calo3d->SetBarrelRadius(kR_max);
886 calo3d->SetEndCapPos(kZ_d);
887 calo3d->SetMaxTowerH(300);
888 eveMng->GetEventScene()->AddElement(calo3d);
890
891 REveDataCollection* ecalCollection = new REveDataCollection("ECAL");
892 ecalCollection->SetItemClass(RCaloTower::Class());
893 ecalCollection->SetMainColor(kRed);
894 collectionMng->addCollection(ecalCollection, new CaloTowerProxyBuilder(event->fCaloData));
895
896 REveDataCollection* hcalCollection = new REveDataCollection("HCAL");
897 hcalCollection->SetItemClass(RCaloTower::Class());
898 hcalCollection->SetMainColor(kBlue);
899 collectionMng->addCollection(hcalCollection, new CaloTowerProxyBuilder(event->fCaloData));
900
901 // event navigation
902 auto eventMng = new EventManager(event, collectionMng);
903 eventMng->SetName("EventManager");
904 eveMng->GetWorld()->AddElement(eventMng);
905
906 eveMng->GetWorld()->AddCommand("NextEvent", "sap-icon://step", eventMng, "NextEvent()");
907
908 eveMng->Show();
909}
#define d(i)
Definition RSha256.hxx:102
#define c(i)
Definition RSha256.hxx:101
#define a(i)
Definition RSha256.hxx:99
#define h(i)
Definition RSha256.hxx:106
#define e(i)
Definition RSha256.hxx:103
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
#define M_PI
Definition Rotated.cxx:105
constexpr Bool_t kFALSE
Definition RtypesCore.h:101
#define ClassDef(name, id)
Definition Rtypes.h:337
@ kRed
Definition Rtypes.h:66
@ kOrange
Definition Rtypes.h:67
@ kGreen
Definition Rtypes.h:66
@ kCyan
Definition Rtypes.h:66
@ kBlue
Definition Rtypes.h:66
@ kYellow
Definition Rtypes.h:66
#define N
static unsigned int total
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t sel
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h prop
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char cname
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t UChar_t len
R__EXTERN TRandom * gRandom
Definition TRandom.h:62
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2489
CaloTowerProxyBuilder(REveCaloDataHist *cd)
REveCaloDataHist * fCaloData
void FillImpliedSelected(REveElement::Set_t &impSet, const std::set< int > &sec_idcs, Product *) override
void BuildProduct(const REveDataCollection *collection, REveElement *product, const REveViewContext *) override
void ModelChanges(const REveDataCollection::Ids_t &ids, Product *product) override
void addCollection(REveDataCollection *collection, REveDataProxyBuilderBase *glBuilder, bool showInTable=false)
CollectionManager(Event *event)
std::vector< REveScene * > m_scenes
void FillImpliedSelected(REveDataItemList *itemList, REveElement::Set_t &impSelSet, const std::set< int > &sec_idcs)
void ModelChanged(REveDataItemList *itemList, const REveDataCollection::Ids_t &ids)
std::vector< REveDataProxyBuilderBase * > m_builders
REveViewContext * m_viewContext
void SetDataItemsFromEvent(REveDataCollection *collection)
CollectionManager * fCMng
~EventManager() override
EventManager(Event *e, CollectionManager *m)
virtual void NextEvent()
void MakeParticles(int N)
std::vector< TList * > fListData
void MakeRecHits(int N)
REveCaloDataHist * fCaloData
void MakeJets(int N)
bool DeviateSelection(REveSelection *selection, REveElement *el, bool multi, bool secondary, const std::set< int > &secondary_idcs) override
void LocalModelChanges(int idx, REveElement *el, const REveViewContext *ctx) override
bool HaveSingleProduct() const override
void BuildItemViewType(const Jet &dj, int idx, REveElement *iItemHolder, const std::string &viewType, const REveViewContext *context) override
void SetEtaSize(float iEtaSize)
static TClass * Class()
float GetPhiSize() const
Jet(Int_t pdg, Int_t status, Int_t mother1, Int_t mother2, Int_t daughter1, Int_t daughter2, Double_t px, Double_t py, Double_t pz, Double_t etot)
float fPhiSize
void SetPhiSize(float iPhiSize)
float fEtaSize
float GetEtaSize() const
RCaloTower(float eta, float phi, float et)
static TClass * Class()
Cell data inner structure.
void ProcessSelection(REveCaloData::vCellId_t &sel_cells, UInt_t selectionId, Bool_t multi) override
void GetCellsFromSecondaryIndices(const std::set< int > &idcs, REveCaloData::vCellId_t &out) override
REveCaloTowerSliceSelector(int s, REveDataCollection *c, REveCaloDataHist *h)
REveDataCollection * fCollection
void AddBox(const Float_t *verts)
void Reset(EBoxType_e boxType, Bool_t valIsCol, Int_t chunkSize)
Int_t AddHistogram(TH2F *hist)
Add new slice to calo tower.
void DataChanged() override
Update limits and notify data users.
TH2F * GetHist(Int_t slice) const
Get histogram in given slice.
void GetCellData(const REveCaloData::CellId_t &id, REveCaloData::CellData_t &data) const override
Get cell geometry and value from cell ID.
void AddSliceSelector(std::unique_ptr< REveCaloDataSliceSelector > s)
virtual void ProcessSelection(REveCaloData::vCellId_t &sel_cells, UInt_t selectionId, bool multi)=0
virtual void GetCellsFromSecondaryIndices(const std::set< int > &idcs, REveCaloData::vCellId_t &out)=0
void SetSelector(REveCaloDataSelector *iSelector)
REveCaloDataSelector * GetSelector()
SliceInfo_t & RefSliceInfo(Int_t s)
std::vector< CellId_t > vCellId_t
void FillImpliedSelectedSet(Set_t &impSelSet, const std::set< int > &sec_idcs) override
Populate set impSelSet with derived / dependant elements.
void AddItem(void *data_ptr, const std::string &n, const std::string &t)
void SetMainColor(Color_t) override
Set main color of the element.
const REveDataItem * GetDataItem(Int_t i) const
void AddTooltipExpression(const std::string &title, const std::string &expr, bool init=true)
void SetFillImpliedSelectedDelegate(FillImpliedSelectedFunc_t)
virtual void LocalModelChanges(int idx, REveElement *el, const REveViewContext *ctx)
virtual REveElement * CreateProduct(const std::string &viewType, const REveViewContext *)
void FillImpliedSelected(REveElement::Set_t &impSet, const std::set< int > &)
void ModelChanges(const REveDataCollection::Ids_t &)
void SetupAddElement(REveElement *el, REveElement *parent, bool set_color=true)
void SetMainColor(Color_t color) override
Override from REveElement, forward to Frame.
void DigitColor(Color_t ci)
Set color for the last digit added.
void SetCurrentDigit(Int_t idx)
Set current digit – the one that will receive calls to DigitValue/Color/Id/UserData() functions.
void RefitPlex()
Instruct underlying memory allocator to regroup itself into a contiguous memory chunk.
DigitBase_t * GetDigit(Int_t n) const
void DigitValue(Int_t value)
Set signal value for the last digit added.
const std::string & GetName() const
virtual void AddElement(REveElement *el)
Add el to the list of children.
virtual Bool_t GetRnrSelf() const
REveElement * FirstChild() const
Returns the first child element or 0 if the list is empty.
void SetSelectionMaster(REveElement *el)
std::set< REveElement * > Set_t
ElementId_t GetElementId() const
virtual Color_t GetMainColor() const
void SetName(const std::string &name)
Set name of an element.
REveMagFieldDuo Interface to magnetic field with two different values depending on radius.
REveScene * GetEventScene() const
REveSelection * GetHighlight() const
REveSceneList * GetScenes() const
REveSelection * GetSelection() const
REveElement * FindElementById(ElementId_t id) const
Lookup ElementId in element map and return corresponding REveElement*.
REveScene * GetGlobalScene() const
REveScene * SpawnNewScene(const char *name, const char *title="")
Create a new scene.
REveViewer * SpawnNewViewer(const char *name, const char *title="")
Create a new GL viewer.
void Show(const RWebDisplayArgs &args="")
Show eve manager in specified browser.
REveProjectionManager Manager class for steering of projections and managing projected objects.
virtual REveElement * ImportElements(REveElement *el, REveElement *ext_list=nullptr)
Recursively import elements and apply projection to the newly imported objects.
void AddCommand(const std::string &name, const std::string &icon, const REveElement *element, const std::string &action)
Definition REveScene.cxx:88
virtual bool DeviateSelection(REveSelection *s, REveElement *el, bool multi, bool secondary, const std::set< int > &secondary_idcs)=0
REveSelection Container for selected and highlighted elements.
void SetDeviator(std::shared_ptr< Deviator > d)
void ClearSelection()
Clear selection if not empty.
virtual void SetLineColor(Color_t c)
Definition REveShape.hxx:65
REveTrackPropagator Calculates path of a particle taking into account special path-marks and imposed ...
REveTrack Track with given vertex, momentum and optional referece-points (path-marks) along its path.
Definition REveTrack.hxx:40
REveTableViewInfo * GetTableViewInfo() const
void SetTrackPropagator(REveTrackPropagator *p)
void SetTableViewInfo(REveTableViewInfo *ti)
REveTrackPropagator * GetPropagator() const
void SetCameraType(ECameraType t)
virtual void AddScene(REveScene *scene)
Add 'scene' to the list of scenes.
REveElement * GetSelectionMaster() override
Returns the master element - that is:
void ModelChanges(const REveDataCollection::Ids_t &ids, Product *product) override
void buildBoxSet(REveBoxSet *boxset)
void FillImpliedSelected(REveElement::Set_t &impSet, const std::set< int > &sec_idcs, Product *p) override
void BuildProduct(const REveDataCollection *collection, REveElement *product, const REveViewContext *) override
static TClass * Class()
RecHit(float pt, float x, float y, float z)
void SetName(const char *name)
Cylindrical tube class.
Definition TGeoTube.h:17
static void AddDirectory(Bool_t add=kTRUE)
Sets the flag controlling the automatic add of histograms in memory.
Definition TH1.cxx:1294
virtual Int_t FindBin(Double_t x, Double_t y=0, Double_t z=0)
Return Global bin number corresponding to x,y,z.
Definition TH1.cxx:3672
static Bool_t AddDirectoryStatus()
Static function: cannot be inlined on Windows/NT.
Definition TH1.cxx:754
2-D histogram with a float per channel (see TH1 documentation)
Definition TH2.h:307
void Reset(Option_t *option="") override
Reset this histogram: contents, errors, etc.
Definition TH2.cxx:3974
Double_t GetBinContent(Int_t binx, Int_t biny) const override
Definition TH2.h:93
Int_t Fill(Double_t) override
Invalid Fill method.
Definition TH2.cxx:393
A doubly linked list.
Definition TList.h:38
void Add(TObject *obj) override
Definition TList.h:81
TObject * At(Int_t idx) const override
Returns the object at position idx. Returns 0 if idx is out of range.
Definition TList.cxx:355
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
Mother of all ROOT objects.
Definition TObject.h:41
virtual const char * GetName() const
Returns name of object.
Definition TObject.cxx:439
void BuildItem(const TParticle &p, int idx, REveElement *iItemHolder, const REveViewContext *context) override
Description of the dynamic properties of a particle.
Definition TParticle.h:26
static TClass * Class()
Double_t Pt() const
Definition TParticle.h:135
Double_t Phi() const
Definition TParticle.h:149
Double_t Eta() const
Definition TParticle.h:137
Double_t Theta(const TParticle &p)
Definition TParticle.h:115
This is the base class for the ROOT Random number generators.
Definition TRandom.h:27
virtual Double_t Uniform(Double_t x1=1)
Returns a uniform deviate on the interval (0, x1).
Definition TRandom.cxx:682
virtual Int_t GetLast() const
Returns index of last object in collection.
Basic string class.
Definition TString.h:139
const char * Data() const
Definition TString.h:376
ROOT::Experimental::REveProjectionManager * g_projMng
ROOT::Experimental::REveManager * eveMng
#define RND_BOX(x)
const Double_t kR_max
void collection_proxies(bool proj=true)
const Double_t kZ_d
const Double_t kR_min
TPaveText * pt
REX::REveManager * eveMng
Definition event_demo.C:41
REX::REveViewer * rhoZView
Definition event_demo.C:47
REX::REveScene * rhoZEventScene
Definition event_demo.C:45
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
Double_t ex[n]
Definition legend1.C:17
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Definition TMathBase.h:198
constexpr Double_t Pi()
Definition TMath.h:37
constexpr Double_t TwoPi()
Definition TMath.h:44
const int xbins_n
const double xbins[xbins_n]
void Setup(const char *name, Float_t threshold, Color_t col, Char_t transp=101)
TMarker m
Definition textangle.C:8
TLine l
Definition textangle.C:4