REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestTrackEvent.cxx
1
20
21#include "TRestTrackEvent.h"
22
23#include "TRestRun.h"
24#include "TRestTools.h"
25
26using namespace std;
27
28ClassImp(TRestTrackEvent);
29
30TRestTrackEvent::TRestTrackEvent() {
31 // TRestTrackEvent default constructor
33 fTrack.clear();
34 fXYHit = nullptr;
35 fXZHit = nullptr;
36 fYZHit = nullptr;
37 fXYZHit = nullptr;
38 fXYTrack = nullptr;
39 fXZTrack = nullptr;
40 fYZTrack = nullptr;
41 fXYZTrack = nullptr;
42 fPad = nullptr;
43 fLevels = -1;
44
45 fPrintHitsWarning = true;
46}
47
48TRestTrackEvent::~TRestTrackEvent() {
49 // TRestTrackEvent destructor
50}
51
53 fNtracks = 0;
54 fNtracksX = 0;
55 fNtracksY = 0;
56 fTrack.clear();
58}
59
60void TRestTrackEvent::AddTrack(TRestTrack* c) {
61 if (c->isXZ()) fNtracksX++;
62 if (c->isYZ()) fNtracksY++;
63 fNtracks++;
64
65 fTrack.emplace_back(*c);
66
67 SetLevels();
68}
69
70void TRestTrackEvent::RemoveTrack(int n) {
71 if (fTrack[n].isXZ()) fNtracksX--;
72 if (fTrack[n].isYZ()) fNtracksY--;
73 fNtracks--;
74
75 fTrack.erase(fTrack.begin() + n);
76
77 SetLevels();
78}
79
80Int_t TRestTrackEvent::GetNumberOfTracks(TString option) {
81 if (option == "")
82 return fNtracks;
83 else {
84 Int_t nT = 0;
85 for (int n = 0; n < GetNumberOfTracks(); n++) {
86 if (!this->isTopLevel(n)) continue;
87
88 if (option == "X" && GetTrack(n)->isXZ())
89 nT++;
90 else if (option == "Y" && GetTrack(n)->isYZ())
91 nT++;
92 else if (option == "XYZ" && GetTrack(n)->isXYZ())
93 nT++;
94 }
95
96 return nT;
97 }
98
99 return -1;
100}
101
102TRestTrack* TRestTrackEvent::GetTrackById(Int_t id) {
103 for (int i = 0; i < GetNumberOfTracks(); i++)
104 if (GetTrack(i)->GetTrackID() == id) return GetTrack(i);
105 return nullptr;
106}
107
108TRestTrack* TRestTrackEvent::GetMaxEnergyTrackInX() {
109 Int_t track = -1;
110 Double_t maxEnergy = 0;
111 for (int tck = 0; tck < GetNumberOfTracks(); tck++) {
112 if (!this->isTopLevel(tck)) continue;
113
114 TRestTrack* t = GetTrack(tck);
115 if (t->isXZ()) {
116 if (t->GetEnergy() > maxEnergy) {
117 maxEnergy = t->GetEnergy();
118 track = tck;
119 }
120 }
121 }
122
123 if (track == -1) return nullptr;
124
125 return GetTrack(track);
126}
127
128TRestTrack* TRestTrackEvent::GetMaxEnergyTrackInY() {
129 Int_t track = -1;
130 Double_t maxEnergy = 0;
131 for (int tck = 0; tck < GetNumberOfTracks(); tck++) {
132 if (!this->isTopLevel(tck)) continue;
133
134 TRestTrack* t = GetTrack(tck);
135 if (t->isYZ()) {
136 if (t->GetEnergy() > maxEnergy) {
137 maxEnergy = t->GetEnergy();
138 track = tck;
139 }
140 }
141 }
142
143 if (track == -1) return nullptr;
144
145 return GetTrack(track);
146}
147
148TRestTrack* TRestTrackEvent::GetMaxEnergyTrack(TString option) {
149 if (option == "X") return GetMaxEnergyTrackInX();
150 if (option == "Y") return GetMaxEnergyTrackInY();
151
152 Int_t track = -1;
153 Double_t maxEnergy = 0;
154 for (int tck = 0; tck < GetNumberOfTracks(); tck++) {
155 if (!this->isTopLevel(tck)) continue;
156
157 TRestTrack* t = GetTrack(tck);
158 if (t->isXYZ()) {
159 if (t->GetEnergy() > maxEnergy) {
160 maxEnergy = t->GetEnergy();
161 track = tck;
162 }
163 }
164 }
165
166 if (track == -1) return nullptr;
167
168 return GetTrack(track);
169}
170
171TRestTrack* TRestTrackEvent::GetSecondMaxEnergyTrack(TString option) {
172 if (GetMaxEnergyTrack(option) == nullptr) return nullptr;
173
174 Int_t id = GetMaxEnergyTrack(option)->GetTrackID();
175
176 Int_t track = -1;
177 Double_t maxEnergy = 0;
178 for (int tck = 0; tck < GetNumberOfTracks(); tck++) {
179 if (!this->isTopLevel(tck)) continue;
180
181 TRestTrack* t = GetTrack(tck);
182 if (t->GetTrackID() == id) continue;
183
184 Double_t en = t->GetEnergy();
185
186 if (option == "X" && (t->isXZ() || t->isXYZ())) {
187 if (en > maxEnergy) {
188 maxEnergy = t->GetEnergy();
189 track = tck;
190 }
191 } else if (option == "Y" && (t->isYZ() || t->isXYZ())) {
192 if (t->GetEnergy() > maxEnergy) {
193 maxEnergy = t->GetEnergy();
194 track = tck;
195 }
196 } else if (t->isXYZ()) {
197 if (t->GetEnergy() > maxEnergy) {
198 maxEnergy = t->GetEnergy();
199 track = tck;
200 }
201 }
202 }
203
204 if (track == -1) return nullptr;
205
206 return GetTrack(track);
207}
208
209Double_t TRestTrackEvent::GetMaxEnergyTrackVolume(TString option) {
210 if (this->GetMaxEnergyTrack(option)) return this->GetMaxEnergyTrack(option)->GetVolume();
211 return 0;
212}
213
214Double_t TRestTrackEvent::GetMaxEnergyTrackLength(TString option) {
215 if (this->GetMaxEnergyTrack(option)) return this->GetMaxEnergyTrack(option)->GetLength();
216 return 0;
217}
218
219Double_t TRestTrackEvent::GetEnergy(TString option) {
220 Double_t en = 0;
221 for (int tck = 0; tck < this->GetNumberOfTracks(); tck++) {
222 if (!this->isTopLevel(tck)) continue;
223
224 TRestTrack* t = GetTrack(tck);
225
226 if (option == "")
227 en += t->GetEnergy();
228
229 else if (option == "X" && (t->isXZ() || t->isXYZ()))
230 en += t->GetEnergy();
231
232 else if (option == "Y" && (t->isYZ() || t->isXYZ()))
233 en += t->GetEnergy();
234
235 else if (option == "XYZ" && t->isXYZ())
236 en += t->GetEnergy();
237 }
238
239 return en;
240}
241
242Bool_t TRestTrackEvent::isXYZ() {
243 for (int tck = 0; tck < GetNumberOfTracks(); tck++)
244 if (!fTrack[tck].isXYZ()) return false;
245 return true;
246}
247
248Int_t TRestTrackEvent::GetTotalHits() {
249 Int_t totHits = 0;
250 for (int tck = 0; tck < GetNumberOfTracks(); tck++) totHits += GetTrack(tck)->GetNumberOfHits();
251 return totHits;
252}
253
254Int_t TRestTrackEvent::GetLevel(Int_t tck) {
255 Int_t lvl = 1;
256 Int_t parentTrackId = GetTrack(tck)->GetParentID();
257
258 while (parentTrackId > 0) {
259 lvl++;
260 parentTrackId = GetTrackById(parentTrackId)->GetParentID();
261 }
262 return lvl;
263}
264
265Bool_t TRestTrackEvent::isTopLevel(Int_t tck) {
266 if (GetLevels() == GetLevel(tck)) return true;
267 return false;
268}
269
270Int_t TRestTrackEvent::GetOriginTrackID(Int_t tck) {
271 Int_t originTrackID = tck;
272 Int_t pID = GetTrack(tck)->GetParentID();
273
274 while (pID != 0) {
275 originTrackID = pID;
276 pID = GetTrackById(originTrackID)->GetParentID();
277 }
278
279 return originTrackID;
280}
281
282TRestTrack* TRestTrackEvent::GetOriginTrack(Int_t tck) {
283 Int_t originTrackID = GetTrack(tck)->GetTrackID();
284 Int_t pID = GetTrackById(originTrackID)->GetParentID();
285
286 while (pID != 0) {
287 originTrackID = pID;
288 pID = GetTrackById(originTrackID)->GetParentID();
289 }
290
291 return GetTrackById(originTrackID);
292}
293
294TRestTrack* TRestTrackEvent::GetOriginTrackById(Int_t tckId) {
295 Int_t originTrackID = tckId;
296 Int_t pID = GetTrackById(tckId)->GetParentID();
297
298 while (pID != 0) {
299 originTrackID = pID;
300 pID = GetTrackById(originTrackID)->GetParentID();
301 }
302
303 return GetTrackById(originTrackID);
304}
305
306void TRestTrackEvent::SetLevels() {
307 Int_t maxLevel = 0;
308
309 for (int tck = 0; tck < GetNumberOfTracks(); tck++) {
310 Int_t lvl = GetLevel(tck);
311 if (maxLevel < lvl) maxLevel = lvl;
312 }
313 fLevels = maxLevel;
314}
315
326 TRestTrack* tckX = GetMaxEnergyTrackInX();
327 TRestTrack* tckY = GetMaxEnergyTrackInY();
328
329 if (tckX == nullptr || tckY == nullptr) {
330 RESTWarning << "Track is empty, skipping" << RESTendl;
331 return {};
332 }
333
334 TRestVolumeHits hitsX = (TRestVolumeHits) * (tckX->GetVolumeHits());
335 TRestVolumeHits hitsY = (TRestVolumeHits) * (tckY->GetVolumeHits());
336
337 const int nHits = std::min(hitsX.GetNumberOfHits(), hitsY.GetNumberOfHits());
338 TRestVolumeHits best3DHits, hits3D;
339
340 for (int i = 0; i < nHits; i++) {
341 double enX = hitsX.GetEnergy(i);
342 double enY = hitsY.GetEnergy(i);
343 double posXZ = hitsX.GetZ(i);
344 double posYZ = hitsY.GetZ(i);
345 double avgZ = (enX * posXZ + enY * posYZ) / (enX + enY);
346 best3DHits.AddHit(hitsX.GetX(i), hitsY.GetY(i), avgZ, enX + enY, 0, XYZ, 0, 0, 0);
347 const int j = nHits - i - 1;
348 enY = hitsY.GetEnergy(j);
349 posYZ = hitsY.GetZ(j);
350 avgZ = (enX * posXZ + enY * posYZ) / (enX + enY);
351 hits3D.AddHit(hitsX.GetX(i), hitsY.GetY(j), avgZ, enX + enY, 0, XYZ, 0, 0, 0);
352 }
353
354 double length = (best3DHits.GetPosition(0) - best3DHits.GetPosition(nHits - 1)).Mag();
355
356 if ((hits3D.GetPosition(0) - hits3D.GetPosition(nHits - 1)).Mag() > length) {
357 best3DHits = hits3D;
358 }
359
360 double totEn = 0;
361 for (unsigned int i = 0; i < best3DHits.GetNumberOfHits(); i++) {
362 totEn += best3DHits.GetEnergy(i);
363 }
364
365 const TVector3 pos0 = best3DHits.GetPosition(0);
366 const TVector3 posE = best3DHits.GetPosition(best3DHits.GetNumberOfHits() - 1);
367
368 double integ = 0;
369 unsigned int pos = 0;
370 for (pos = 0; pos < best3DHits.GetNumberOfHits(); pos++) {
371 integ += best3DHits.GetEnergy(pos);
372 if (integ > totEn / 2.) break;
373 }
374
375 auto intPos = best3DHits.GetPosition(pos);
376 const double intToFirst = (pos0 - intPos).Mag();
377 const double intToLast = (posE - intPos).Mag();
378
379 RESTDebug << "Integ pos " << pos << " Pos to first " << intToFirst << " last " << intToLast << RESTendl;
380 if (intToFirst < intToLast) {
381 end = pos0;
382 orig = posE;
383 } else {
384 orig = pos0;
385 end = posE;
386 }
387
388 RESTDebug << "Origin " << orig.X() << " " << orig.Y() << " " << orig.Z() << RESTendl;
389 RESTDebug << "End " << end.X() << " " << end.Y() << " " << end.Z() << RESTendl;
390
391 return best3DHits;
392}
393
400void TRestTrackEvent::GetMaxTrackBoundaries(TVector3& orig, TVector3& end) {
401 TRestTrack* tckX = GetMaxEnergyTrackInX();
402 TRestTrack* tckY = GetMaxEnergyTrackInY();
403
404 if (tckX == nullptr || tckY == nullptr) {
405 RESTWarning << "Track is empty, skipping" << RESTendl;
406 return;
407 }
408
409 TVector3 origX, endX;
410 // Retreive origin and end of the track for the XZ projection
411 tckX->GetBoundaries(origX, endX);
412 TVector3 origY, endY;
413 // Retreive origin and end of the track for the YZ projection
414 tckY->GetBoundaries(origY, endY);
415
416 double originZ = (origX.Z() + origY.Z()) / 2.;
417 double endZ = (endX.Z() + endY.Z()) / 2.;
418
419 orig = TVector3(origX.X(), origY.Y(), originZ);
420 end = TVector3(endX.X(), endY.Y(), endZ);
421}
422
428 TRestTrack* tckX = GetMaxEnergyTrackInX();
429 TRestTrack* tckY = GetMaxEnergyTrackInY();
430
431 if (tckX == nullptr || tckY == nullptr) {
432 RESTWarning << "Track is empty, skipping" << RESTendl;
433 return -1;
434 }
435
436 std::vector<std::pair<double, double> > zEn;
437 double totEn = 0;
438
439 for (size_t i = 0; i < tckX->GetVolumeHits()->GetNumberOfHits(); i++) {
440 double en = tckX->GetVolumeHits()->GetEnergy(i);
441 double z = tckX->GetVolumeHits()->GetZ(i);
442 zEn.emplace_back(std::make_pair(z, en));
443 totEn += en;
444 }
445
446 for (size_t i = 0; i < tckY->GetVolumeHits()->GetNumberOfHits(); i++) {
447 double en = tckY->GetVolumeHits()->GetEnergy(i);
448 double z = tckY->GetVolumeHits()->GetZ(i);
449 zEn.emplace_back(std::make_pair(z, en));
450 totEn += en;
451 }
452
453 std::sort(zEn.begin(), zEn.end());
454
455 double integ = 0;
456 size_t pos = 0;
457 for (pos = 0; pos < zEn.size(); pos++) {
458 integ += zEn[pos].second;
459 if (integ >= totEn / 2.) break;
460 }
461
462 double length = zEn.front().first - zEn.back().first;
463 double diff = zEn.front().first - zEn[pos].first;
464
465 if (length == 0)
466 return 0;
467 else
468 return diff / length;
469}
470
471void TRestTrackEvent::PrintOnlyTracks() {
472 cout << "TrackEvent " << GetID() << endl;
473 cout << "-----------------------" << endl;
474 for (int i = 0; i < GetNumberOfTracks(); i++) {
475 cout << "Track " << i << " id : " << GetTrack(i)->GetTrackID()
476 << " parent : " << GetTrack(i)->GetParentID() << endl;
477 }
478 cout << "-----------------------" << endl;
479 cout << "Track levels : " << GetLevels() << endl;
480}
481
482void TRestTrackEvent::PrintEvent(Bool_t fullInfo) {
484
485 cout << "Number of tracks : "
486 << GetNumberOfTracks("XYZ") + GetNumberOfTracks("X") + GetNumberOfTracks("Y") << endl;
487 cout << "Number of tracks XZ " << GetNumberOfTracks("X") << endl;
488 cout << "Number of tracks YZ " << GetNumberOfTracks("Y") << endl;
489 cout << "Track levels : " << GetLevels() << endl;
490 cout << "+++++++++++++++++++++++++++++++++++" << endl;
491 for (int i = 0; i < GetNumberOfTracks(); i++) this->GetTrack(i)->PrintTrack(fullInfo);
492}
493
494TPad* TRestTrackEvent::DrawHits() {
495 if (fXZHits) {
496 delete fXZHits;
497 fXZHits = nullptr;
498 }
499 if (fYZHits) {
500 delete fYZHits;
501 fYZHits = nullptr;
502 }
503 if (fHitsPad) {
504 delete fHitsPad;
505 fHitsPad = nullptr;
506 }
507
508 std::vector<double> fX, fY, fZ;
509
510 for (int t = 0; t < GetNumberOfTracks(); t++) {
511 TRestTrack* tck = GetTrack(t);
512 if (GetLevel(t) != 1) continue;
513 TRestVolumeHits* hits = tck->GetVolumeHits();
514 for (unsigned int i = 0; i < hits->GetNumberOfHits(); i++) {
515 if (hits->GetType(i) % X == 0) fX.emplace_back(hits->GetX(i));
516 if (hits->GetType(i) % Y == 0) fY.emplace_back(hits->GetY(i));
517 if (hits->GetType(i) % Z == 0) fZ.emplace_back(hits->GetZ(i));
518 }
519 }
520
521 double maxX, minX, maxY, minY, maxZ, minZ;
522 int nBinsX, nBinsY, nBinsZ;
523 TRestHits::GetBoundaries(fX, maxX, minX, nBinsX);
524 TRestHits::GetBoundaries(fY, maxY, minY, nBinsY);
525 TRestHits::GetBoundaries(fZ, maxZ, minZ, nBinsZ);
526
527 fXZHits = new TH2F("TXZ", "TXZ", nBinsX, minX, maxX, nBinsZ, minZ, maxZ);
528 fYZHits = new TH2F("TYZ", "TYZ", nBinsY, minY, maxY, nBinsZ, minZ, maxZ);
529
530 for (int t = 0; t < GetNumberOfTracks(); t++) {
531 TRestTrack* tck = GetTrack(t);
532 if (GetLevel(t) != 1) continue;
533 TRestVolumeHits* hits = tck->GetVolumeHits();
534 for (unsigned int i = 0; i < hits->GetNumberOfHits(); i++) {
535 if (hits->GetType(i) == XZ) fXZHits->Fill(hits->GetX(i), hits->GetZ(i), hits->GetEnergy(i));
536 if (hits->GetType(i) == YZ) fYZHits->Fill(hits->GetY(i), hits->GetZ(i), hits->GetEnergy(i));
537 }
538 }
539
540 fHitsPad = new TPad("TrackHits", "TrackHits", 0, 0, 1, 1);
541 fHitsPad->Divide(2, 1);
542 fHitsPad->Draw();
543
544 fHitsPad->cd(1);
545 fXZHits->GetXaxis()->SetTitle("X-axis (mm)");
546 fXZHits->GetYaxis()->SetTitle("Z-axis (mm)");
547 fXZHits->Draw("COLZ");
548
549 fHitsPad->cd(2);
550 fYZHits->GetXaxis()->SetTitle("Y-axis (mm)");
551 fYZHits->GetYaxis()->SetTitle("Z-axis (mm)");
552 fYZHits->Draw("COLZ");
553
554 return fHitsPad;
555}
556
557// Draw current event in a Tpad
558TPad* TRestTrackEvent::DrawEvent(const TString& option) {
559 /* Not used for the moment
560 Bool_t drawXZ = false;
561 Bool_t drawYZ = false;
562 Bool_t drawXY = false;
563 Bool_t drawXYZ = false;
564 Bool_t drawLines = false;
565 */
566
567 Int_t maxLevel = 0;
568 Int_t minLevel = 0;
569
570 vector<string> optList = TRestTools::GetOptions((string)option);
571
572 for (unsigned int n = 0; n < optList.size(); n++) {
573 if (optList[n] == "print") this->PrintEvent();
574 if (optList[n] == "noWarning") fPrintHitsWarning = false;
575 }
576
577 optList.erase(std::remove(optList.begin(), optList.end(), "print"), optList.end());
578 optList.erase(std::remove(optList.begin(), optList.end(), "noWarning"), optList.end());
579
580 for (unsigned int n = 0; n < optList.size(); n++) {
581 /* Not used for the moment
582 if( optList[n] == "XZ" ) drawXZ = true;
583 if( optList[n] == "YZ" ) drawYZ = true;
584 if( optList[n] == "XY" ) drawXY = true;
585 if( optList[n] == "XYZ" ) drawXYZ = true;
586 if( optList[n] == "L" || optList[n] == "lines" ) drawLines = true;
587 */
588 string opt = optList[n];
589
590 if (opt.find("maxLevel=") != string::npos) maxLevel = stoi(opt.substr(9, opt.length()).c_str());
591
592 if (opt.find("minLevel=") != string::npos) minLevel = stoi(opt.substr(9, opt.length()).c_str());
593 }
594
595 if (fXYHit != nullptr) {
596 delete[] fXYHit;
597 fXYHit = nullptr;
598 }
599 if (fXZHit != nullptr) {
600 delete[] fXZHit;
601 fXZHit = nullptr;
602 }
603 if (fYZHit != nullptr) {
604 delete[] fYZHit;
605 fYZHit = nullptr;
606 }
607 if (fXYZHit != nullptr) {
608 delete[] fXYZHit;
609 fXYZHit = nullptr;
610 }
611 if (fXYTrack != nullptr) {
612 delete[] fXYTrack;
613 fXYTrack = nullptr;
614 }
615 if (fXZTrack != nullptr) {
616 delete[] fXZTrack;
617 fXZTrack = nullptr;
618 }
619 if (fYZTrack != nullptr) {
620 delete[] fYZTrack;
621 fYZTrack = nullptr;
622 }
623 if (fXYZTrack != nullptr) {
624 delete[] fXYZTrack;
625 fXYZTrack = nullptr;
626 }
627 if (fPad != nullptr) {
628 delete fPad;
629 fPad = nullptr;
630 }
631
632 int nTracks = this->GetNumberOfTracks();
633
634 if (nTracks == 0) {
635 cout << "Empty event " << endl;
636 return nullptr;
637 }
638
639 this->PrintEvent(false);
640
641 double maxX = -1e10, minX = 1e10, maxZ = -1e10, minZ = 1e10, maxY = -1e10, minY = 1e10;
642
643 Int_t nTotHits = GetTotalHits();
644
645 if (fPrintHitsWarning && nTotHits > 5000) {
646 cout << endl;
647 cout << " REST WARNING. TRestTrackEvent::DrawEvent. Number of hits is too "
648 "high."
649 << endl;
650 cout << " This drawing method is not properly optimized to draw events "
651 "with a high number of hits."
652 << endl;
653 cout << " To remove this warning you may use the DrawEvent method option : "
654 "noWarning "
655 << endl;
656 cout << endl;
657
658 fPrintHitsWarning = false;
659 }
660
661 fXYHit = new TGraph[nTotHits];
662 fXZHit = new TGraph[nTotHits];
663 fYZHit = new TGraph[nTotHits];
664 fXYZHit = new TGraph2D[nTotHits];
665 fXYTrack = new TGraph[nTracks];
666 fXZTrack = new TGraph[nTracks];
667 fYZTrack = new TGraph[nTracks];
668 fXYZTrack = new TGraph2D[nTracks];
669
670 vector<Int_t> drawLinesXY(nTracks, 0);
671 vector<Int_t> drawLinesXZ(nTracks, 0);
672 vector<Int_t> drawLinesYZ(nTracks, 0);
673 vector<Int_t> drawLinesXYZ(nTracks, 0);
674
675 int countXY = 0, countYZ = 0, countXZ = 0, countXYZ = 0;
676 int nTckXY = 0, nTckXZ = 0, nTckYZ = 0, nTckXYZ = 0;
677
678 Double_t minRadiusSize = 0.4;
679 Double_t maxRadiusSize = 2.;
680
681 Int_t maxTrackHits = 0;
682
683 Int_t tckColor = 1;
684
685 for (int tck = 0; tck < nTracks; tck++) {
686 TRestVolumeHits* hits = fTrack[tck].GetVolumeHits();
687
688 Double_t maxHitEnergy = hits->GetMaximumHitEnergy();
689 Double_t meanHitEnergy = hits->GetMeanHitEnergy();
690
691 /*
692 cout << "Max hit energy : " << maxHitEnergy << endl;
693 cout << "Mean hit energy : " << meanHitEnergy << endl;
694 cout << "Number of hits " << hits->GetNumberOfHits( ) <<endl;
695 */
696
697 Bool_t isTopLevel = this->isTopLevel(tck);
698 if (isTopLevel) tckColor++;
699 Int_t level = this->GetLevel(tck);
700
701 if (!isTopLevel && maxLevel > 0 && level > maxLevel) continue;
702
703 if (!isTopLevel && minLevel > 0 && level < minLevel) continue;
704
705 int tckXY = 0, tckYZ = 0, tckXZ = 0, tckXYZ = 0;
706 Double_t radius;
707
708 for (unsigned int nhit = 0; nhit < hits->GetNumberOfHits(); nhit++) {
709 if (hits->GetNumberOfHits() > (unsigned int)maxTrackHits) maxTrackHits = hits->GetNumberOfHits();
710
711 if (hits->isNaN(nhit)) {
712 cout << "REST Warning. TRestTrackEvent::Draw. Hit is not defined!!" << endl;
713 getchar();
714 continue;
715 }
716
717 Double_t x = hits->GetX(nhit);
718 Double_t y = hits->GetY(nhit);
719 Double_t z = hits->GetZ(nhit);
720 Double_t en = hits->GetEnergy(nhit);
721 auto type = hits->GetType(nhit);
722
723 // cout << x << " " << y << " " << z << " " << type << endl;
724
725 /* {{{ Hit size definition (radius) */
726 Double_t m = (maxRadiusSize) / (maxHitEnergy - meanHitEnergy);
727 Double_t n = (maxRadiusSize - minRadiusSize) - m * meanHitEnergy;
728
729 if (isTopLevel) {
730 radius = m * en + n;
731 if (radius < 0.1) radius = 0.1;
732 if (hits->GetNumberOfHits() == 1) radius = 2;
733 if (hits->GetNumberOfHits() > 100) radius = 0.8;
734
735 } else {
736 radius = 1.5 * minRadiusSize * level;
737 }
738 /* }}} */
739
740 if (this->isXYZ() && nhit > 1) {
741 if (tckXYZ == 0) nTckXYZ++;
742 fXYZTrack[nTckXYZ - 1].SetPoint(tckXYZ, x, y, z);
743
744 if (isTopLevel) drawLinesXYZ[nTckXYZ - 1] = 1;
745
746 if (isTopLevel) {
747 fXYZHit[countXYZ].SetPoint(0, x, y, z);
748 // If there is only one-point the TGraph2D does NOT draw the point!
749 fXYZHit[countXYZ].SetPoint(1, x + 0.001, y + 0.001, z + 0.001);
750
751 fXYZHit[countXYZ].SetMarkerColor(level + 11);
752
753 fXYZHit[countXYZ].SetMarkerSize(radius);
754 fXYZHit[countXYZ].SetMarkerStyle(20);
755 countXYZ++;
756 }
757 tckXYZ++;
758 }
759
760 if (type == XY) {
761 if (tckXY == 0) nTckXY++;
762 fXYTrack[nTckXY - 1].SetPoint(tckXY, x, y);
763 if (isTopLevel) drawLinesXY[nTckXY - 1] = 1;
764 fXYHit[countXY].SetPoint(0, x, y);
765
766 if (!isTopLevel)
767 fXYHit[countXY].SetMarkerColor(level + 11);
768 else
769 fXYHit[countXY].SetMarkerColor(tckColor);
770
771 fXYHit[countXY].SetMarkerSize(radius);
772 fXYHit[countXY].SetMarkerStyle(20);
773 tckXY++;
774 countXY++;
775 }
776
777 if (type == XZ) {
778 if (tckXZ == 0) nTckXZ++;
779 fXZTrack[nTckXZ - 1].SetPoint(tckXZ, x, z);
780 if (isTopLevel) drawLinesXZ[nTckXZ - 1] = 1;
781 fXZHit[countXZ].SetPoint(0, x, z);
782
783 if (!isTopLevel)
784 fXZHit[countXZ].SetMarkerColor(level + 11);
785 else
786 fXZHit[countXZ].SetMarkerColor(tckColor);
787
788 fXZHit[countXZ].SetMarkerSize(radius);
789 fXZHit[countXZ].SetMarkerStyle(20);
790 tckXZ++;
791 countXZ++;
792 }
793
794 if (type == YZ) {
795 if (tckYZ == 0) nTckYZ++;
796 fYZTrack[nTckYZ - 1].SetPoint(tckYZ, y, z);
797 if (isTopLevel) drawLinesYZ[nTckYZ - 1] = 1;
798 fYZHit[countYZ].SetPoint(0, y, z);
799
800 if (!isTopLevel)
801 fYZHit[countYZ].SetMarkerColor(level + 11);
802 else
803 fYZHit[countYZ].SetMarkerColor(tckColor);
804
805 fYZHit[countYZ].SetMarkerSize(radius);
806 fYZHit[countYZ].SetMarkerStyle(20);
807 tckYZ++;
808 countYZ++;
809 }
810
811 if (type % X == 0) {
812 if (x > maxX) maxX = x;
813 if (x < minX) minX = x;
814 }
815
816 if (type % Y == 0) {
817 if (y > maxY) maxY = y;
818 if (y < minY) minY = y;
819 }
820
821 if (type % Z == 0) {
822 if (z > maxZ) maxZ = z;
823 if (z < minZ) minZ = z;
824 }
825 }
826 }
827
828 fPad = new TPad(this->GetName(), " ", 0, 0, 1, 1);
829 if (this->isXYZ())
830 fPad->Divide(2, 2);
831 else
832 fPad->Divide(2, 1);
833 fPad->Draw();
834
835 char title[256];
836 sprintf(title, "Event ID %d", this->GetID());
837
838 TMultiGraph* mgXY = new TMultiGraph();
839 TMultiGraph* mgXZ = new TMultiGraph();
840 TMultiGraph* mgYZ = new TMultiGraph();
841
842 fPad->cd(1)->DrawFrame(minX - 10, minZ - 10, maxX + 10, maxZ + 10, title);
843 fPad->cd(2)->DrawFrame(minY - 10, minZ - 10, maxY + 10, maxZ + 10, title);
844
845 if (this->isXYZ()) fPad->cd(3)->DrawFrame(minX - 10, minY - 10, maxX + 10, maxY + 10, title);
846
847 for (int i = 0; i < countXZ; i++) mgXZ->Add(&fXZHit[i]);
848
849 fPad->cd(1);
850 mgXZ->GetXaxis()->SetTitle("X-axis (mm)");
851 mgXZ->GetYaxis()->SetTitle("Z-axis (mm)");
852 mgXZ->GetYaxis()->SetTitleOffset(1.75);
853 mgXZ->GetYaxis()->SetTitleSize(1.4 * mgXZ->GetYaxis()->GetTitleSize());
854 mgXZ->GetXaxis()->SetTitleSize(1.4 * mgXZ->GetXaxis()->GetTitleSize());
855 mgXZ->GetYaxis()->SetLabelSize(1.25 * mgXZ->GetYaxis()->GetLabelSize());
856 mgXZ->GetXaxis()->SetLabelSize(1.25 * mgXZ->GetXaxis()->GetLabelSize());
857 mgXZ->Draw("P");
858
859 for (int i = 0; i < countYZ; i++) mgYZ->Add(&fYZHit[i]);
860
861 fPad->cd(2);
862 mgYZ->GetXaxis()->SetTitle("Y-axis (mm)");
863 mgYZ->GetYaxis()->SetTitle("Z-axis (mm)");
864 mgYZ->GetYaxis()->SetTitleOffset(1.75);
865 mgYZ->GetYaxis()->SetTitleSize(1.4 * mgYZ->GetYaxis()->GetTitleSize());
866 mgYZ->GetXaxis()->SetTitleSize(1.4 * mgYZ->GetXaxis()->GetTitleSize());
867 mgYZ->GetYaxis()->SetLabelSize(1.25 * mgYZ->GetYaxis()->GetLabelSize());
868 mgYZ->GetXaxis()->SetLabelSize(1.25 * mgYZ->GetXaxis()->GetLabelSize());
869 mgYZ->Draw("P");
870
871 if (this->isXYZ()) {
872 for (int i = 0; i < countXY; i++) mgXY->Add(&fXYHit[i]);
873
874 fPad->cd(3);
875 mgXY->GetXaxis()->SetTitle("X-axis (mm)");
876 mgXY->GetYaxis()->SetTitle("Y-axis (mm)");
877 mgXY->GetYaxis()->SetTitleOffset(1.75);
878 mgXY->Draw("P");
879 mgXY->GetYaxis()->SetTitleSize(1.4 * mgXY->GetYaxis()->GetTitleSize());
880 mgXY->GetXaxis()->SetTitleSize(1.4 * mgXY->GetXaxis()->GetTitleSize());
881 mgXY->GetYaxis()->SetLabelSize(1.25 * mgXY->GetYaxis()->GetLabelSize());
882 mgXY->GetXaxis()->SetLabelSize(1.25 * mgXY->GetXaxis()->GetLabelSize());
883 }
884
885 for (int tck = 0; tck < nTckXZ; tck++) {
886 fPad->cd(1);
887 fXZTrack[tck].SetLineWidth(2.);
888 if (fXZTrack[tck].GetN() < 100 && drawLinesXZ[tck] == 1) fXZTrack[tck].Draw("L");
889 }
890
891 for (int tck = 0; tck < nTckYZ; tck++) {
892 fPad->cd(2);
893 fYZTrack[tck].SetLineWidth(2.);
894 if (fYZTrack[tck].GetN() < 100 && drawLinesYZ[tck] == 1) fYZTrack[tck].Draw("L");
895 }
896
897 if (this->isXYZ()) {
898 for (int tck = 0; tck < nTckXY; tck++) {
899 fPad->cd(3);
900 fXYTrack[tck].SetLineWidth(2.);
901 if (fXYTrack[tck].GetN() < 100 && drawLinesXY[tck] == 1) fXYTrack[tck].Draw("L");
902 }
903
904 fPad->cd(4);
905
906 TString option = "P";
907 for (int tck = 0; tck < nTckXYZ; tck++) {
908 fXYZTrack[tck].SetLineWidth(2.);
909 if (fXZTrack[tck].GetN() < 100 && drawLinesXYZ[tck] == 1) {
910 fXYZTrack[tck].Draw("LINE");
911 option = "same P";
912 }
913 }
914
915 for (int i = 0; i < countXYZ; i++) {
916 if (i > 0) option = "same P";
917 fXYZHit[i].Draw(option);
918 }
919 }
920
921 return fPad;
922}
923
928void TRestTrackEvent::GetOriginEnd(std::vector<TGraph*>& originGr, std::vector<TGraph*>& endGr,
929 std::vector<TLegend*>& leg) {
930 if (originGr.size() != 2 || endGr.size() != 2 || leg.size() != 2) return;
931
932 for (auto gr : originGr)
933 if (gr) delete gr;
934
935 for (auto gr : endGr)
936 if (gr) delete gr;
937
938 for (auto l : leg)
939 if (l) delete l;
940
941 TVector3 orig, end;
942 GetMaxTrackBoundaries(orig, end);
943
944 for (int i = 0; i < 2; i++) {
945 originGr[i] = new TGraph();
946 originGr[i]->SetPoint(0, orig[i], orig[2]);
947 originGr[i]->SetMarkerColor(kRed);
948 originGr[i]->SetMarkerStyle(20);
949 endGr[i] = new TGraph();
950 endGr[i]->SetPoint(0, end[i], end[2]);
951 endGr[i]->SetMarkerColor(kBlack);
952 endGr[i]->SetMarkerStyle(20);
953 leg[i] = new TLegend(0.7, 0.7, 0.9, 0.9);
954 leg[i]->AddEntry(originGr[i], "Origin", "p");
955 leg[i]->AddEntry(endGr[i], "End", "p");
956 }
957}
958
963void TRestTrackEvent::DrawOriginEnd(TPad* pad, std::vector<TGraph*>& originGr, std::vector<TGraph*>& endGr,
964 std::vector<TLegend*>& leg) {
965 if (originGr.size() != 2 || endGr.size() != 2 || leg.size() != 2) return;
966
967 for (int i = 0; i < 2; i++) {
968 pad->cd(i + 1);
969 if (originGr[i]) originGr[i]->Draw("LP");
970 if (endGr[i]) endGr[i]->Draw("LP");
971 if (leg[i]) leg[i]->Draw();
972 pad->Update();
973 }
974}
virtual void PrintEvent() const
Definition: TRestEvent.cxx:187
virtual void Initialize()=0
Definition: TRestEvent.cxx:73
Double_t GetMaximumHitEnergy() const
It returns the maximum hit energy.
Definition: TRestHits.cxx:426
Double_t GetMeanHitEnergy() const
It returns the mean hits energy.
Definition: TRestHits.cxx:452
static void GetBoundaries(std::vector< double > &val, double &max, double &min, int &nBins, double offset=10)
TODO This method is not using any TRestHits member. This probably means that it should be placed some...
Definition: TRestHits.cxx:738
Bool_t isNaN(Int_t n) const
It will return true only if all the 3-coordinates of hit number n are not a number,...
Definition: TRestHits.cxx:162
TVector3 GetPosition(int n) const
It returns the position of hit number n.
Definition: TRestHits.cxx:515
static std::vector< std::string > GetOptions(std::string optionsStr)
Returns all the options in an option string.
Definition: TRestTools.cxx:86
void GetOriginEnd(std::vector< TGraph * > &originGr, std::vector< TGraph * > &endGr, std::vector< TLegend * > &leg)
Retreive origin and end of the track and store in a TGraph and legend.
TRestVolumeHits GetMaxTrackBoundaries3D(TVector3 &orig, TVector3 &end)
This function retrieves the origin and the end track positions based after the reconstruction of a 3D...
TPad * DrawEvent(const TString &option="")
Draw the event.
void DrawOriginEnd(TPad *pad, std::vector< TGraph * > &originGr, std::vector< TGraph * > &endGr, std::vector< TLegend * > &leg)
Draw origin and end of the track in a pad passed to the function Note that GetOriginEnd has to be iss...
Double_t GetMaxTrackRelativeZ()
Function to calculate the relative Z of the most energetic track to crosscheck if the track is upward...
void GetMaxTrackBoundaries(TVector3 &orig, TVector3 &end)
This function retreive the origin and the end of the track based on the most energetic hit....
void GetBoundaries(TVector3 &orig, TVector3 &end)
This function retreive the origin and the end of a single track based on the most energetic hit....
Definition: TRestTrack.cxx:60