Logo ROOT  
Reference Guide
Loading...
Searching...
No Matches
event_demo.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_eve_7
3/// This example display geometry, tracks and hits in web browser
4///
5/// \macro_code
6///
7
8#include <vector>
9#include <string>
10#include <iostream>
11
12#include "TClass.h"
13#include "TRandom.h"
14#include "TGeoTube.h"
15#include "TGeoSphere.h"
16#include "TParticle.h"
17#include "TApplication.h"
18#include "TMatrixDSym.h"
19#include "TVector.h"
20#include "TMatrixDEigen.h"
21
22#include <ROOT/REveGeoShape.hxx>
23#include <ROOT/REveScene.hxx>
24#include <ROOT/REveViewer.hxx>
25#include <ROOT/REveElement.hxx>
26#include <ROOT/REveCompound.hxx>
27#include <ROOT/REveManager.hxx>
28#include <ROOT/REveUtil.hxx>
29#include <ROOT/REveGeoShape.hxx>
32#include <ROOT/REvePointSet.hxx>
33#include <ROOT/REveJetCone.hxx>
34#include <ROOT/REveTrans.hxx>
35#include <ROOT/REveTrack.hxx>
37
38namespace REX = ROOT::Experimental;
39
40// globals
48
49const Double_t kR_min = 240;
50const Double_t kR_max = 250;
51const Double_t kZ_d = 300;
52
53REX::REvePointSet *getPointSet(int npoints = 2, float s = 2, int color = 28)
54{
55 TRandom &r = *gRandom;
56
57 auto ps = new REX::REvePointSet("fu", "", npoints);
58
59 for (Int_t i = 0; i < npoints; ++i)
60 ps->SetNextPoint(r.Uniform(-s, s), r.Uniform(-s, s), r.Uniform(-s, s));
61
62 ps->SetMarkerColor(color);
63 ps->SetMarkerSize(8);
64 ps->SetMarkerStyle(4);
65 return ps;
66}
67
69{
70 REX::REveElement *event = eveMng->GetEventScene();
71
72 auto pntHolder = new REX::REveCompound("Hits");
73
74 auto ps1 = getPointSet(20, 100);
75 ps1->SetName("Points_1");
76 ps1->SetTitle("Points_1 title"); // used as tooltip
77
78 pntHolder->AddElement(ps1);
79
80 auto ps2 = getPointSet(10, 200, 4);
81 ps2->SetName("Points_2");
82 ps2->SetTitle("Points_2 title"); // used as tooltip
83 ps2->SetAlwaysSecSelect(true);
84 pntHolder->AddElement(ps2);
85
86 event->AddElement(pntHolder);
87}
88
90{
91 TRandom &r = *gRandom;
92
93 REX::REveElement *event = eveMng->GetEventScene();
94 auto prop = new REX::REveTrackPropagator();
95 prop->SetMagFieldObj(new REX::REveMagFieldDuo(350, 3.5, -2.0));
96 prop->SetMaxR(300);
97 prop->SetMaxZ(600);
98 prop->SetMaxOrbs(6);
99
100 auto trackHolder = new REX::REveCompound("Tracks");
101
102 double v = 0.2;
103 double m = 5;
104
105 int N_Tracks = 10 + r.Integer(20);
106 for (int i = 0; i < N_Tracks; i++) {
107 TParticle p;
108
109 int pdg = 11 * (r.Integer(2) > 0 ? 1 : -1);
110 p.SetPdgCode(pdg);
111
112 p.SetProductionVertex(r.Uniform(-v, v), r.Uniform(-v, v), r.Uniform(-v, v), 1);
113 p.SetMomentum(r.Uniform(-m, m), r.Uniform(-m, m), r.Uniform(-m, m) * r.Uniform(1, 3), 1);
114 auto track = new REX::REveTrack(&p, 1, prop);
115 track->MakeTrack();
116 if (i % 4 == 3)
117 track->SetLineStyle(kDashed); // enabled dashed style for some tracks
118 track->SetMainColor(kBlue);
119 track->SetName(Form("RandomTrack_%d", i));
120 track->SetTitle(Form("RandomTrack_%d title", i)); // used as tooltip
121 trackHolder->AddElement(track);
122 }
123
124 event->AddElement(trackHolder);
125}
126
128{
129 TRandom &r = *gRandom;
130
131 REX::REveElement *event = eveMng->GetEventScene();
132 auto jetHolder = new REX::REveCompound("Jets");
133
134 int N_Jets = 5 + r.Integer(5);
135 for (int i = 0; i < N_Jets; i++) {
136 auto jet = new REX::REveJetCone(Form("Jet_%d", i));
137 jet->SetTitle(Form("Jet_%d\n pT = %.2f", i, r.Uniform(1, 40))); // used as tooltip
138 jet->SetCylinder(2 * kR_max, 2 * kZ_d);
139 jet->AddEllipticCone(r.Uniform(-3.5, 3.5), r.Uniform(0, TMath::TwoPi()), r.Uniform(0.02, 0.2),
140 r.Uniform(0.02, 0.3));
141 jet->SetFillColor(kPink - 8);
142 jet->SetLineColor(kBlack);
143
144 jetHolder->AddElement(jet);
145 }
146 event->AddElement(jetHolder);
147}
148
150{
151 addPoints();
152 addTracks();
153 addJets();
154}
155
157{
158 auto b1 = new REX::REveGeoShape("Barrel 1");
159 b1->SetShape(new TGeoTube(kR_min, kR_max, kZ_d));
160 b1->SetMainColor(kCyan);
161 b1->SetNSegments(80);
162 b1->SetMainTransparency(70);
163 eveMng->GetGlobalScene()->AddElement(b1);
164
165 // Debug of surface fill in RPhi (index buffer screwed).
166 // b1->SetNSegments(3);
167 b1->SetNSegments(40);
168
169 // an example of axis guides
170 eveMng->GetDefaultViewer()->SetAxesType(REX::REveViewer::EAxesType::kAxesOrigin);
171}
172
174{
175 // project RhoPhi
176 rPhiGeomScene = eveMng->SpawnNewScene("RPhi Geometry", "RPhi");
177 rPhiEventScene = eveMng->SpawnNewScene("RPhi Event Data", "RPhi");
178
180
181 rphiView = eveMng->SpawnNewViewer("RPhi View", "");
183 rphiView->AddScene(rPhiGeomScene);
184 rphiView->AddScene(rPhiEventScene);
185
186 // ----------------------------------------------------------------
187
188 rhoZGeomScene = eveMng->SpawnNewScene("RhoZ Geometry", "RhoZ");
189 rhoZEventScene = eveMng->SpawnNewScene("RhoZ Event Data", "RhoZ");
190
192
193 rhoZView = eveMng->SpawnNewViewer("RhoZ View", "");
195 rhoZView->AddScene(rhoZGeomScene);
196 rhoZView->AddScene(rhoZEventScene);
197}
198
199void projectScenes(bool geomp, bool eventp)
200{
201 if (geomp) {
202 for (auto &ie : eveMng->GetGlobalScene()->RefChildren()) {
203 mngRhoPhi->SetCurrentDepth(0);
204 mngRhoPhi->ImportElements(ie, rPhiGeomScene);
205 mngRhoZ->SetCurrentDepth(0);
206 mngRhoZ->ImportElements(ie, rhoZGeomScene);
207 }
208 }
209 if (eventp) {
210 int depth = 50;
211 for (auto &ie : eveMng->GetEventScene()->RefChildren()) {
212 mngRhoPhi->SetCurrentDepth(depth);
213 mngRhoPhi->ImportElements(ie, rPhiEventScene);
214 mngRhoZ->SetCurrentDepth(depth);
215 mngRhoZ->ImportElements(ie, rhoZEventScene);
216 depth -= 10;
217 }
218 }
219
220 // auto t0 = eveMng->GetEventScene()->FindChild("Tracks")->FirstChild();
221 // printf("t0=%p, %s %s\n", t0, t0->GetElementName(), t0->IsA()->GetName());
222 // dynamic_cast<REX::REveTrack*>(t0)->Print("all");
223
224 // auto t1 = rPhiEventScene->FindChild("Tracks [P]")->FirstChild();
225 // printf("t1=%p, %s %s\n", t1, t1->GetElementName(), t1->IsA()->GetName());
226 // dynamic_cast<REX::REveTrack*>(t1)->Print("all");
227}
228
229//==============================================================================
230
231class EventManager : public REX::REveElement {
232private:
233 bool fAutoplay{false};
234 int fPlayDelay{10};
235 int fCount{0};
236 std::chrono::duration<double> fDeltaTime{1};
237
238 std::thread *fTimerThread{nullptr};
239 std::mutex fMutex;
240 std::condition_variable fCV;
241
242public:
244 {
245 std::chrono::milliseconds ms(100);
246 fDeltaTime = ms;
247 }
248
249 ~EventManager() override {}
250
252 {
253 auto scene = eveMng->GetEventScene();
254 scene->DestroyElements();
256 projectScenes(false, true);
257 // if (++fCount % 10 == 0) printf("At event %d\n", fCount);
258 }
259
261 {
262 while (true) {
263 bool autoplay;
264 {
265 std::unique_lock<std::mutex> lock{fMutex};
266 if (!fAutoplay) {
267 // printf("exit thread pre wait\n");
268 return;
269 }
270 if (fCV.wait_for(lock, fDeltaTime) != std::cv_status::timeout) {
271 printf("autoplay not timed out \n");
272 if (!fAutoplay) {
273 printf("exit thread post wait\n");
274 return;
275 } else {
276 continue;
277 }
278 }
279 autoplay = fAutoplay;
280 }
281 if (autoplay) {
283 NextEvent();
284 } else {
285 return;
286 }
287 }
288 }
289
290 void Autoplay()
291 {
292 static std::mutex autoplay_mutex;
293 std::unique_lock<std::mutex> aplock{autoplay_mutex};
294 {
295 std::unique_lock<std::mutex> lock{fMutex};
297 if (fAutoplay) {
298 if (fTimerThread) {
299 fTimerThread->join();
300 delete fTimerThread;
301 fTimerThread = nullptr;
302 }
303 NextEvent();
304 fTimerThread = new std::thread{[this] { autoplay_scheduler(); }};
305 } else {
306 fCV.notify_all();
307 }
308 }
309 }
310
311 virtual void QuitRoot()
312 {
313 printf("Quit ROOT\n");
315 }
316};
317
319{
320 // disable browser cache - all scripts and html files will be loaded every time, useful for development
321 // gEnv->SetValue("WebGui.HttpMaxAge", 0);
322
323 gRandom->SetSeed(0); // make random seed
324
326
327 auto eventMng = new EventManager();
328 eventMng->SetName("EventManager");
329 eveMng->GetWorld()->AddElement(eventMng);
330
331 eveMng->GetWorld()->AddCommand("QuitRoot", "sap-icon://log", eventMng, "QuitRoot()");
332
333 eveMng->GetWorld()->AddCommand("NextEvent", "sap-icon://step", eventMng, "NextEvent()");
334
335 eveMng->GetWorld()->AddCommand("Autoplay", "sap-icon://refresh", eventMng, "Autoplay()");
336
339
340 if (true) {
342 projectScenes(true, true);
343 }
344
345 eveMng->Show();
346}
ROOT::R::TRInterface & r
Definition Object.C:4
int Int_t
Signed integer 4 bytes (int).
Definition RtypesCore.h:59
double Double_t
Double 8 bytes.
Definition RtypesCore.h:73
@ kPink
Definition Rtypes.h:68
@ kBlack
Definition Rtypes.h:66
@ kCyan
Definition Rtypes.h:67
@ kBlue
Definition Rtypes.h:67
@ kDashed
Definition TAttLine.h:54
externTRandom * gRandom
Definition TRandom.h:62
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2496
std::condition_variable fCV
Definition event_demo.C:240
std::chrono::duration< double > fDeltaTime
Definition event_demo.C:236
virtual void QuitRoot()
Definition event_demo.C:311
~EventManager() override
Definition event_demo.C:249
std::mutex fMutex
Definition event_demo.C:239
void Autoplay()
Definition event_demo.C:290
void NextEvent()
Definition event_demo.C:251
virtual void NextEvent()
std::thread * fTimerThread
Definition event_demo.C:238
void autoplay_scheduler()
Definition event_demo.C:260
REveMagFieldDuo Interface to magnetic field with two different values depending on radius.
static REveManager * Create()
If global REveManager* REX::gEve is not set initialize it.
REveProjectionManager Manager class for steering of projections and managing projected objects.
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
REveViewer Reve representation of TGLViewer.
Description of the dynamic properties of a particle.
Definition TParticle.h:26
void SetMomentum(Double_t px, Double_t py, Double_t pz, Double_t e)
Definition TParticle.h:166
void SetPdgCode(Int_t pdg)
Change the PDG code for this particle.
void SetProductionVertex(Double_t vx, Double_t vy, Double_t vz, Double_t t)
Definition TParticle.h:168
This is the base class for the ROOT Random number generators.
Definition TRandom.h:27
ROOT::Experimental::REveManager * eveMng
const Double_t kR_max
const Double_t kZ_d
const Double_t kR_min
REX::REveScene * rPhiGeomScene
Definition event_demo.C:44
void addTracks()
Definition event_demo.C:89
void addJets()
Definition event_demo.C:127
REX::REveScene * rPhiEventScene
Definition event_demo.C:44
void makeEventScene()
Definition event_demo.C:149
REX::REveProjectionManager * mngRhoZ
Definition event_demo.C:43
REX::REveManager * eveMng
Definition event_demo.C:41
REX::REveViewer * rhoZView
Definition event_demo.C:47
REX::REveScene * rhoZEventScene
Definition event_demo.C:45
void createProjectionStuff()
Definition event_demo.C:173
void addPoints()
Definition event_demo.C:68
REX::REveViewer * rphiView
Definition event_demo.C:46
void makeGeometryScene()
Definition event_demo.C:156
REX::REveProjectionManager * mngRhoPhi
Definition event_demo.C:42
void event_demo()
Definition event_demo.C:318
REX::REveScene * rhoZGeomScene
Definition event_demo.C:45
void projectScenes(bool geomp, bool eventp)
Definition event_demo.C:199
REX::REvePointSet * getPointSet(int npoints=2, float s=2, int color=28)
Definition event_demo.C:53
Namespace for ROOT features in testing.
Definition TROOT.h:100
constexpr Double_t TwoPi()
Definition TMath.h:47
TMarker m
Definition textangle.C:8