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