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
11#include "ROOT/REveManager.hxx"
14#include <ROOT/REveGeoShape.hxx>
15#include <ROOT/REveJetCone.hxx>
16#include <ROOT/REvePointSet.hxx>
19#include <ROOT/REveScene.hxx>
22#include <ROOT/REveTrack.hxx>
24#include <ROOT/REveViewer.hxx>
26#include <ROOT/REveBoxSet.hxx>
28#include <ROOT/REveCalo.hxx>
29
30#include "TGeoTube.h"
31#include "TROOT.h"
32#include "TList.h"
33#include "TParticle.h"
34#include "TRandom.h"
35#include "TApplication.h"
36#include "TFile.h"
37#include "TH2F.h"
38#include <iostream>
39
40const Double_t kR_min = 299;
41const Double_t kR_max = 300;
42const Double_t kZ_d = 500;
43
44namespace fw3dlego {
45const int xbins_n = 83;
46const double xbins[xbins_n] = {
47 -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,
48 -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,
49 -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,
50 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,
51 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,
52 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};
53} // namespace fw3dlego
54
57using namespace ROOT::Experimental;
58
59//==============================================================================
60//============== EMULATE FRAMEWORK CLASSES =====================================
61//==============================================================================
62
63// a demo class, can be provided from experiment framework
64class Jet : public TParticle {
65public:
66 float fEtaSize{0};
67 float fPhiSize{0};
68
69 float GetEtaSize() const { return fEtaSize; }
70 float GetPhiSize() const { return fPhiSize; }
71 void SetEtaSize(float iEtaSize) { fEtaSize = iEtaSize; }
72 void SetPhiSize(float iPhiSize) { fPhiSize = iPhiSize; }
73
74 Jet(Int_t pdg, Int_t status, Int_t mother1, Int_t mother2, Int_t daughter1, Int_t daughter2, Double_t px,
75 Double_t py, Double_t pz, Double_t etot)
76 : TParticle(pdg, status, mother1, mother2, daughter1, daughter2, px, py, pz, etot, 0, 0, 0, 0)
77 {
78 }
79
81};
82
83class RecHit : public TObject {
84public:
85 float fX{0};
86 float fY{0};
87 float fZ{0};
88 float fPt{0};
89
90 RecHit(float pt, float x, float y, float z) : fPt(pt), fX(x), fY(y), fZ(z) {}
92};
93
94class RCaloTower : public TObject {
95public:
96 float fEta{0};
97 float fPhi{0};
98 float fEt{0};
99
100 RCaloTower(float eta, float phi, float et) : fEta(eta), fPhi(phi), fEt(et) {}
102};
103
105private:
108
109public:
112 {
113 }
114
116 void ProcessSelection(REveCaloData::vCellId_t &sel_cells, UInt_t selectionId, Bool_t multi) override
117 {
118 std::set<int> item_set;
120 for (auto &cellId : sel_cells) {
121 fCaloData->GetCellData(cellId, cd);
122
123 // loop over enire collection and check its eta/phi range
124 for (int t = 0; t < fCollection->GetNItems(); ++t) {
126 if (tower->fEta > cd.fEtaMin && tower->fEta < cd.fEtaMax && tower->fPhi > cd.fPhiMin &&
127 tower->fPhi < cd.fPhiMax && fCollection->GetDataItem(t)->GetVisible()) {
128 item_set.insert(t);
129 }
130 }
131 }
133 fCollection->GetItemList()->RefSelectedSet() = item_set;
134 sel->NewElementPicked(fCollection->GetItemList()->GetElementId(), multi, true, item_set);
135 }
136
138 void GetCellsFromSecondaryIndices(const std::set<int> &idcs, REveCaloData::vCellId_t &out) override
139 {
141 std::set<int> cbins;
142 float total = 0;
143 for (auto &i : idcs) {
145 int bin = hist->FindBin(tower->fEta, tower->fPhi);
146 float frac = tower->fEt / hist->GetBinContent(bin);
147 bool ex = false;
148 for (size_t ci = 0; ci < out.size(); ++ci) {
149 if (out[ci].fTower == bin && out[ci].fSlice == GetSliceIndex()) {
150 float oldv = out[ci].fFraction;
151 out[ci].fFraction = oldv + frac;
152 ex = true;
153 break;
154 }
155 }
156 if (!ex) {
157 out.push_back(REveCaloData::CellId_t(bin, GetSliceIndex(), frac));
158 }
159 }
160 }
161};
162
163class Event {
164public:
165 int eventId{0};
166 int N_tracks{0};
167 int N_jets{0};
168 std::vector<TList *> fListData;
169
171
173 {
174 auto baseHist = new TH2F("dummy", "dummy", fw3dlego::xbins_n - 1, fw3dlego::xbins, 72, -TMath::Pi(), TMath::Pi());
176 fCaloData->AddHistogram(baseHist);
177
178 auto selector = new REveCaloDataSelector();
179 fCaloData->SetSelector(selector);
180
182 }
183
184 void MakeJets(int N)
185 {
186 TRandom &r = *gRandom;
187 r.SetSeed(0);
188 TList *list = new TList();
189 list->SetName("Jets");
190 for (int i = 1; i <= N; ++i) {
191 double pt = r.Uniform(0.5, 10);
192 double eta = r.Uniform(-2.55, 2.55);
193 double phi = r.Uniform(-TMath::Pi(), TMath::Pi());
194
195 double px = pt * std::cos(phi);
196 double py = pt * std::sin(phi);
197 double pz = pt * (1. / (std::tan(2 * std::atan(std::exp(-eta)))));
198
199 auto jet = new Jet(0, 0, 0, 0, 0, 0, px, py, pz, std::sqrt(px * px + py * py + pz * pz + 80 * 80));
200 jet->SetEtaSize(r.Uniform(0.02, 0.2));
201 jet->SetPhiSize(r.Uniform(0.01, 0.3));
202 list->Add(jet);
203 }
204 fListData.push_back(list);
205 }
206
207 void MakeParticles(int N)
208 {
209 TRandom &r = *gRandom;
210 r.SetSeed(0);
211 TList *list = new TList();
212 list->SetName("Tracks");
213 for (int i = 1; i <= N; ++i) {
214 double pt = r.Uniform(0.5, 10);
215 double eta = r.Uniform(-2.55, 2.55);
216 double phi = r.Uniform(0, TMath::TwoPi());
217
218 double px = pt * std::cos(phi);
219 double py = pt * std::sin(phi);
220 double pz = pt * (1. / (std::tan(2 * std::atan(std::exp(-eta)))));
221
222 // printf("Event::MakeParticles %2d: pt=%.2f, eta=%.2f, phi=%.2f\n", i, pt, eta, phi);
223 auto particle =
224 new TParticle(0, 0, 0, 0, 0, 0, px, py, pz, std::sqrt(px * px + py * py + pz * pz + 80 * 80), 0, 0, 0, 0);
225
226 int pdg = 11 * (r.Integer(2) > 0 ? 1 : -1);
227 particle->SetPdgCode(pdg);
228
229 list->Add(particle);
230 }
231 fListData.push_back(list);
232 }
233
234 void MakeRecHits(int N)
235 {
236 TRandom &r = *gRandom;
237 r.SetSeed(0);
238 TList *list = new TList();
239 list->SetName("RecHits");
240
241 for (int i = 1; i <= N; ++i) {
242 float pt = r.Uniform(0.5, 10);
243 float x = r.Uniform(-200, 200);
244 float y = r.Uniform(-200, 200);
245 float z = r.Uniform(-500, 500);
246 auto rechit = new RecHit(pt, x, y, z);
247 list->Add(rechit);
248 }
249 fListData.push_back(list);
250 }
251
252 void Clear()
253 {
254 for (auto &l : fListData)
255 delete l;
256 fListData.clear();
257 }
258
259 void Create()
260 {
261 Clear();
262 MakeJets(4);
263 MakeParticles(100);
264 MakeRecHits(20);
265
266 // refill calo data from jet list
267 TList *jlist = fListData[0];
268 TList *elist = new TList();
269 elist->SetName("ECAL");
270 fListData.push_back(elist);
271 TList *hlist = new TList();
272 hlist->SetName("HCAL");
273 fListData.push_back(hlist);
274 for (int i = 0; i <= jlist->GetLast(); ++i) {
275 const Jet *j = (Jet *)jlist->At(i);
276 float offX = j->Eta();
277 float offY = j->Phi() > TMath::Pi() ? j->Phi() - TMath::TwoPi() : j->Phi();
278 for (int k = 0; k < 20; ++k) {
279 double x, y, v;
280 x = gRandom->Uniform(-j->GetEtaSize(), j->GetEtaSize());
281 y = gRandom->Uniform(-j->GetPhiSize(), j->GetPhiSize());
282 v = j->Pt();
283 auto etower = new RCaloTower(offX + x, offY + y, v + gRandom->Uniform(2, 3));
284 elist->Add(etower);
285 auto htower = new RCaloTower(offX + x, offY + y, v + gRandom->Uniform(1, 2));
286 hlist->Add(htower);
287 }
288 }
290 eventId++;
291 }
292};
293
294//==============================================================================
295//== PROXY BUILDERS ============================================================
296//==============================================================================
297
299 bool HaveSingleProduct() const override { return false; }
300
302 void BuildItemViewType(const Jet &dj, int idx, REveElement *iItemHolder, const std::string &viewType,
303 const REveViewContext *context) override
304 {
305 auto jet = new REveJetCone();
306 jet->SetCylinder(context->GetMaxR(), context->GetMaxZ());
307 jet->AddEllipticCone(dj.Eta(), dj.Phi(), dj.GetEtaSize(), dj.GetPhiSize());
308 SetupAddElement(jet, iItemHolder, true);
309 jet->SetLineColor(jet->GetMainColor());
310
311 float size = 50.f * dj.Pt(); // values are saved in scale
312 double theta = dj.Theta();
313 // printf("%s jet theta = %f, phi = %f \n", iItemHolder->GetCName(), theta, dj.Phi());
314 double phi = dj.Phi();
315
316 if (viewType == "Projected") {
317 static const float_t offr = 6;
318 float r_ecal = context->GetMaxR() + offr;
319 float z_ecal = context->GetMaxZ() + offr;
320
321 float transAngle = abs(atan(r_ecal / z_ecal));
322 double r = 0;
323 bool debug = false;
324 if (theta < transAngle || 3.14 - theta < transAngle) {
325 z_ecal = context->GetMaxZ() + offr / transAngle;
326 r = z_ecal / fabs(cos(theta));
327 } else {
328 debug = true;
329 r = r_ecal / sin(theta);
330 }
331
332 REveVector p1(0, (phi < TMath::Pi() ? r * fabs(sin(theta)) : -r * fabs(sin(theta))), r * cos(theta));
333 REveVector p2(0, (phi < TMath::Pi() ? (r + size) * fabs(sin(theta)) : -(r + size) * fabs(sin(theta))),
334 (r + size) * cos(theta));
335
336 auto marker = new REveScalableStraightLineSet("jetline");
337 marker->SetScaleCenter(p1.fX, p1.fY, p1.fZ);
338 marker->AddLine(p1, p2);
339 marker->SetLineWidth(4);
340 if (debug)
341 marker->AddMarker(0, 0.9);
342
343 SetupAddElement(marker, iItemHolder, true);
344 marker->SetName(Form("line %s %d", Collection()->GetCName(), idx));
345 }
346 }
347
349
350 void LocalModelChanges(int idx, REveElement *el, const REveViewContext *ctx) override
351 {
352 // printf("LocalModelChanges jet %s ( %s )\n", el->GetCName(), el->FirstChild()->GetCName());
353 REveJetCone *cone = dynamic_cast<REveJetCone *>(el->FirstChild());
354 cone->SetLineColor(cone->GetMainColor());
355 }
356};
357
360
361 void BuildItem(const TParticle &p, int idx, REveElement *iItemHolder, const REveViewContext *context) override
362 {
363 const TParticle *x = &p;
364 auto track = new REveTrack((TParticle *)(x), 1, context->GetPropagator());
365 track->MakeTrack();
366 SetupAddElement(track, iItemHolder, true);
367 }
368};
369
371private:
372 class FWBoxSet : public REveBoxSet {
373 public:
374 using REveElement::GetSelectionMaster;
376 {
377 if (fSelectionMaster) {
380 return il;
381 }
382 return nullptr;
383 }
384 };
385
388 {
389 auto collection = Collection();
390 boxset->SetMainColor(collection->GetMainColor());
391 boxset->SetName(collection->GetCName());
392 boxset->SetPickable(true);
393 boxset->SetAlwaysSecSelect(true);
394 boxset->SetDetIdsAsSecondaryIndices(true);
395 boxset->SetSelectionMaster(((REveDataCollection *)collection)->GetItemList());
396 boxset->Reset(REveBoxSet::kBT_FreeBox, true, collection->GetNItems());
397 TRandom r(0);
398
399#define RND_BOX(x) (Float_t) r.Uniform(-(x), (x))
400 for (int h = 0; h < collection->GetNItems(); ++h) {
401 RecHit *hit = (RecHit *)collection->GetDataPtr(h);
402 const REveDataItem *item = Collection()->GetDataItem(h);
403
404 Float_t x = hit->fX;
405 Float_t y = hit->fY;
406 Float_t z = hit->fZ;
407 Float_t a = hit->fPt;
408 Float_t d = 0.05;
409 Float_t verts[24] = {x - a + RND_BOX(d), y - a + RND_BOX(d), z - a + RND_BOX(d), x - a + RND_BOX(d),
410 y + a + RND_BOX(d), z - a + RND_BOX(d), x + a + RND_BOX(d), y + a + RND_BOX(d),
411 z - a + RND_BOX(d), x + a + RND_BOX(d), y - a + RND_BOX(d), z - a + RND_BOX(d),
412 x - a + RND_BOX(d), y - a + RND_BOX(d), z + a + RND_BOX(d), x - a + RND_BOX(d),
413 y + a + RND_BOX(d), z + a + RND_BOX(d), x + a + RND_BOX(d), y + a + RND_BOX(d),
414 z + a + RND_BOX(d), x + a + RND_BOX(d), y - a + RND_BOX(d), z + a + RND_BOX(d)};
415 boxset->AddBox(verts);
416
417 boxset->DigitValue(item->GetVisible() ? 1 : 0);
418 if (item->GetVisible())
419 boxset->DigitColor(item->GetMainColor());
420 }
421 boxset->RefitPlex();
422 boxset->StampObjProps();
423 }
424
425public:
427 void BuildProduct(const REveDataCollection *collection, REveElement *product, const REveViewContext *) override
428 {
429 fBoxSet = new FWBoxSet();
431 product->AddElement(fBoxSet);
432 }
433
435 void FillImpliedSelected(REveElement::Set_t &impSet, const std::set<int> &sec_idcs, Product *p) override
436 {
437 // printf("RecHit fill implioed ----------------- !!!%zu\n",
438 // Collection()->GetItemList()->RefSelectedSet().size());
439 impSet.insert(fBoxSet);
440 }
441
443 void ModelChanges(const REveDataCollection::Ids_t &ids, Product *product) override
444 {
445 for (auto &i : ids) {
446 auto digi = fBoxSet->GetDigit(i);
447 auto item = Collection()->GetDataItem(i);
449 if (item->GetVisible()) {
451 fBoxSet->DigitColor(item->GetMainColor());
452 } else {
454 }
455 }
457 }
458}; // RecHitProxyBuilder
459
461private:
463 TH2F *fHist{nullptr};
464 int fSliceIndex{-1};
465
467 {
468 if (!fHist) {
470
471 TH1::AddDirectory(kFALSE); // Keeps histogram from going into memory
472 fHist = new TH2F("caloHist", "caloHist", fw3dlego::xbins_n - 1, fw3dlego::xbins, 72, -M_PI, M_PI);
473 TH1::AddDirectory(status);
475
477 .Setup(Collection()->GetCName(), 0., Collection()->GetMainColor(), Collection()->GetMainTransparency());
478
479 fCaloData->GetSelector()->AddSliceSelector(std::unique_ptr<REveCaloDataSliceSelector>(
481 }
482 }
483
484public:
486
488 void BuildProduct(const REveDataCollection *collection, REveElement *product, const REveViewContext *) override
489 {
490 assertSlice();
491 fHist->Reset();
492 if (collection->GetRnrSelf()) {
494 .Setup(Collection()->GetCName(), 0., Collection()->GetMainColor(), Collection()->GetMainTransparency());
495
496 for (int h = 0; h < collection->GetNItems(); ++h) {
497 RCaloTower *tower = (RCaloTower *)collection->GetDataPtr(h);
498 const REveDataItem *item = Collection()->GetDataItem(h);
499
500 if (!item->GetVisible())
501 continue;
502 fHist->Fill(tower->fEta, tower->fPhi, tower->fEt);
503 }
504 }
506 }
507
509 void FillImpliedSelected(REveElement::Set_t &impSet, const std::set<int> &sec_idcs, Product *) override
510 {
512 impSet.insert(fCaloData);
513 fCaloData->FillImpliedSelectedSet(impSet, sec_idcs);
514 }
515
517 void ModelChanges(const REveDataCollection::Ids_t &ids, Product *product) override
518 {
519 BuildProduct(Collection(), nullptr, nullptr);
520 }
521
522}; // CaloTowerProxyBuilder
523
524//==============================================================================
525//== COLLECTION MANGER ================================================================
526//==============================================================================
527
529private:
530 Event *fEvent{nullptr};
531
532 std::vector<REveScene *> m_scenes;
534
535 std::vector<REveDataProxyBuilderBase *> m_builders;
536
538 bool m_inEventLoading{false};
539
540public:
542 {
543 // view context
544 float r = 300;
545 float z = 300;
546 auto prop = new REveTrackPropagator();
547 prop->SetMagFieldObj(new REveMagFieldDuo(350, 3.5, -2.0));
548 prop->SetMaxR(r);
549 prop->SetMaxZ(z);
550 prop->SetMaxOrbs(6);
551 prop->IncRefCount();
552
556
557 // table specs
558 auto tableInfo = new REveTableViewInfo();
559
560 tableInfo->table("TParticle").column("pt", 1, "i.Pt()").column("eta", 3, "i.Eta()").column("phi", 3, "i.Phi()");
561
562 tableInfo->table("Jet")
563 .column("eta", 1, "i.Eta()")
564 .column("phi", 1, "i.Phi()")
565 .column("etasize", 2, "i.GetEtaSize()")
566 .column("phisize", 2, "i.GetPhiSize()");
567
568 tableInfo->table("RecHit").column("pt", 1, "i.fPt");
569
570 tableInfo->table("RCaloTower").column("eta", 3, "i.fEta").column("phi", 3, "i.fPhi").column("Et", 3, "i.fEt");
571
573
574 for (auto &c : eveMng->GetScenes()->RefChildren()) {
575 if (c != eveMng->GetGlobalScene() && strncmp(c->GetCName(), "Geometry", 8)) {
576 m_scenes.push_back((REveScene *)c);
577 }
578 if (!strncmp(c->GetCName(), "Table", 5))
579 c->AddElement(m_viewContext->GetTableViewInfo());
580 }
581
582 m_collections = eveMng->SpawnNewScene("Collections", "Collections");
583 }
584
586 {
587 for (auto &l : fEvent->fListData) {
588 TIter next(l);
589 if (collection->GetName() == std::string(l->GetName())) {
590 collection->ClearItems();
591
592 for (int i = 0; i <= l->GetLast(); ++i) {
593 std::string cname = collection->GetName();
594 auto len = cname.size();
595 char end = cname[len - 1];
596 if (end == 's') {
597 cname = cname.substr(0, len - 1);
598 }
599 TString pname(Form("%s %2d", cname.c_str(), i));
600 collection->AddItem(l->At(i), pname.Data(), "");
601 }
602 }
603 collection->ApplyFilter();
604 }
605 }
606
608 {
609 m_inEventLoading = true;
610
611 for (auto &el : m_collections->RefChildren()) {
612 auto c = dynamic_cast<REveDataCollection *>(el);
614 }
615
616 for (auto proxy : m_builders) {
617 proxy->Build();
618 }
619
621 m_inEventLoading = false;
622 }
623
624 void addCollection(REveDataCollection *collection, REveDataProxyBuilderBase *glBuilder, bool showInTable = false)
625 {
626 m_collections->AddElement(collection);
627
628 // load data
629 SetDataItemsFromEvent(collection);
630 glBuilder->SetCollection(collection);
631 glBuilder->SetHaveAWindow(true);
632 for (auto scene : m_scenes) {
633 if (strncmp(scene->GetCName(), "Tables", 5) == 0)
634 continue;
635
636 REveElement *product = glBuilder->CreateProduct(scene->GetTitle(), m_viewContext);
637
638 if (!strncmp(scene->GetCTitle(), "Projected", 8)) {
639 g_projMng->ImportElements(product, scene);
640 } else {
641 scene->AddElement(product);
642 }
643 }
644 m_builders.push_back(glBuilder);
645 glBuilder->Build();
646
647 // Tables
648 auto tableBuilder = new REveTableProxyBuilder();
649 tableBuilder->SetHaveAWindow(true);
650 tableBuilder->SetCollection(collection);
651 REveElement *tablep = tableBuilder->CreateProduct("table-type", m_viewContext);
652 auto tableMng = m_viewContext->GetTableViewInfo();
653 if (showInTable) {
654 tableMng->SetDisplayedCollection(collection->GetElementId());
655 }
656
657 for (auto s : m_scenes) {
658 if (strncmp(s->GetCTitle(), "Table", 5) == 0) {
659 s->AddElement(tablep);
660 tableBuilder->Build();
661 }
662 }
663 tableMng->AddDelegate([=]() { tableBuilder->ConfigChanged(); });
664 m_builders.push_back(tableBuilder);
665
666 // set tooltip expression for items
667 auto tableEntries = tableMng->RefTableEntries(collection->GetItemClass()->GetName());
668 int N = TMath::Min(int(tableEntries.size()), 3);
669 for (int t = 0; t < N; t++) {
670 auto te = tableEntries[t];
671 collection->GetItemList()->AddTooltipExpression(te.fName, te.fExpression);
672 }
673
674 collection->GetItemList()->SetItemsChangeDelegate(
675 [&](REveDataItemList *collection, const REveDataCollection::Ids_t &ids) {
676 this->ModelChanged(collection, ids);
677 });
679 [&](REveDataItemList *collection, REveElement::Set_t &impSelSet, const std::set<int> &sec_idcs) {
680 this->FillImpliedSelected(collection, impSelSet, sec_idcs);
681 });
682 }
683
685 {
686 auto mngTable = m_viewContext->GetTableViewInfo();
687 if (mngTable) {
688 for (auto &el : m_collections->RefChildren()) {
689 if (el->GetName() == "Tracks")
690 mngTable->SetDisplayedCollection(el->GetElementId());
691 }
692 }
693 }
694
696 {
698 return;
699
700 for (auto proxy : m_builders) {
701 if (proxy->Collection()->GetItemList() == itemList) {
702 // printf("Model changes check proxy %s: \n", proxy->Type().c_str());
703 proxy->ModelChanges(ids);
704 }
705 }
706 }
707
708 void FillImpliedSelected(REveDataItemList *itemList, REveElement::Set_t &impSelSet, const std::set<int> &sec_idcs)
709 {
711 return;
712
713 for (auto proxy : m_builders) {
714 if (proxy->Collection()->GetItemList() == itemList) {
715 proxy->FillImpliedSelected(impSelSet, sec_idcs);
716 }
717 }
718 }
719};
720
721//==============================================================================
722//== Event Manager =============================================================
723//==============================================================================
724
725class EventManager : public REveElement {
726private:
729
730public:
732
733 ~EventManager() override {}
734
735 virtual void NextEvent()
736 {
739 fEvent->Create();
740 fCMng->LoadEvent();
741 }
742};
743
745public:
747
749 bool DeviateSelection(REveSelection *selection, REveElement *el, bool multi, bool secondary,
750 const std::set<int> &secondary_idcs) override
751 {
752 if (el) {
753 auto *colItems = dynamic_cast<REveDataItemList *>(el);
754 if (colItems) {
755 // std::cout << "Deviate RefSelected=" << colItems->RefSelectedSet().size() << " passed set " <<
756 // secondary_idcs.size() << "\n";
757 ExecuteNewElementPicked(selection, colItems, multi, true, colItems->RefSelectedSet());
758 return true;
759 }
760 }
761 return false;
762 }
763};
764//==============================================================================
765//== main() ====================================================================
766//==============================================================================
767
768void collection_proxies(bool proj = true)
769{
770 eveMng = REveManager::Create();
771 auto event = new Event();
772 event->Create();
773
774 // divert selection to map proxy builder products with collection
775 auto deviator = std::make_shared<FWSelectionDeviator>();
776 eveMng->GetSelection()->SetDeviator(deviator);
777 eveMng->GetHighlight()->SetDeviator(deviator);
778
779 // create scenes and views
780 REveScene *rhoZEventScene = nullptr;
781
782 auto b1 = new REveGeoShape("Barrel 1");
783 b1->SetShape(new TGeoTube(kR_min, kR_max, kZ_d));
784 b1->SetMainColor(kCyan);
785 b1->SetMainTransparency(90);
787
788 rhoZEventScene = eveMng->SpawnNewScene("RhoZ Scene", "Projected");
789 g_projMng = new REveProjectionManager(REveProjection::kPT_RhoZ);
791
792 auto rhoZView = eveMng->SpawnNewViewer("RhoZ View");
793 rhoZView->SetCameraType(REveViewer::kCameraOrthoXOY);
795 auto pgeoScene = eveMng->SpawnNewScene("Geometry projected");
796 rhoZView->AddScene(pgeoScene);
797 g_projMng->ImportElements(b1, pgeoScene);
798
799 auto tableScene = eveMng->SpawnNewScene("Tables", "Tables");
800 auto tableView = eveMng->SpawnNewViewer("Table", "Table View");
801 tableView->AddScene(tableScene);
802
803 // create event data from list
804 auto collectionMng = new CollectionManager(event);
805
806 REveDataCollection *trackCollection = new REveDataCollection("Tracks");
807 trackCollection->SetItemClass(TParticle::Class());
808 trackCollection->SetMainColor(kGreen);
809 trackCollection->SetFilterExpr("i.Pt() > 4.1 && std::abs(i.Eta()) < 1");
810 collectionMng->addCollection(trackCollection, new TParticleProxyBuilder(), true);
811
812 REveDataCollection *jetCollection = new REveDataCollection("Jets");
813 jetCollection->SetItemClass(Jet::Class());
814 jetCollection->SetMainColor(kYellow);
815 jetCollection->SetFilterExpr("i.Pt() > 1");
816 collectionMng->addCollection(jetCollection, new JetProxyBuilder());
817
818 REveDataCollection *hitCollection = new REveDataCollection("RecHits");
819 hitCollection->SetItemClass(RecHit::Class());
820 hitCollection->SetMainColor(kOrange + 7);
821 hitCollection->SetFilterExpr("i.fPt > 5");
822 collectionMng->addCollection(hitCollection, new RecHitProxyBuilder(), true);
823
824 // add calorimeters
825 auto calo3d = new REveCalo3D(event->fCaloData);
826 calo3d->SetBarrelRadius(kR_max);
827 calo3d->SetEndCapPos(kZ_d);
828 calo3d->SetMaxTowerH(300);
829 eveMng->GetEventScene()->AddElement(calo3d);
831
832 REveDataCollection *ecalCollection = new REveDataCollection("ECAL");
833 ecalCollection->SetItemClass(RCaloTower::Class());
834 ecalCollection->SetMainColor(kRed);
835 collectionMng->addCollection(ecalCollection, new CaloTowerProxyBuilder(event->fCaloData));
836
837 REveDataCollection *hcalCollection = new REveDataCollection("HCAL");
838 hcalCollection->SetItemClass(RCaloTower::Class());
839 hcalCollection->SetMainColor(kBlue);
840 collectionMng->addCollection(hcalCollection, new CaloTowerProxyBuilder(event->fCaloData));
841
842 // event navigation
843 auto eventMng = new EventManager(event, collectionMng);
844 eventMng->SetName("EventManager");
845 eveMng->GetWorld()->AddElement(eventMng);
846
847 eveMng->GetWorld()->AddCommand("NextEvent", "sap-icon://step", eventMng, "NextEvent()");
848
849 eveMng->Show();
850}
#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:94
#define ClassDef(name, id)
Definition Rtypes.h:342
@ 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
Definition JetEvent.h:49
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:90
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.
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)
static void AddDirectory(Bool_t add=kTRUE)
Sets the flag controlling the automatic add of histograms in memory.
Definition TH1.cxx:1296
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:3680
static Bool_t AddDirectoryStatus()
Static function: cannot be inlined on Windows/NT.
Definition TH1.cxx:756
2-D histogram with a float per channel (see TH1 documentation)
Definition TH2.h:308
void Reset(Option_t *option="") override
Reset this histogram: contents, errors, etc.
Definition TH2.cxx:3969
Double_t GetBinContent(Int_t binx, Int_t biny) const override
Definition TH2.h:94
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:456
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