REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
Loading...
Searching...
No Matches
TRestDetectorHitsEvent.cxx
1/*************************************************************************
2 * This file is part of the REST software framework. *
3 * *
4 * Copyright (C) 2016 GIFNA/TREX (University of Zaragoza) *
5 * For more information see http://gifna.unizar.es/trex *
6 * *
7 * REST is free software: you can redistribute it and/or modify *
8 * it under the terms of the GNU General Public License as published by *
9 * the Free Software Foundation, either version 3 of the License, or *
10 * (at your option) any later version. *
11 * *
12 * REST is distributed in the hope that it will be useful, *
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
15 * GNU General Public License for more details. *
16 * *
17 * You should have a copy of the GNU General Public License along with *
18 * REST in $REST_PATH/LICENSE. *
19 * If not, see http://www.gnu.org/licenses/. *
20 * For the list of contributors see $REST_PATH/CREDITS. *
21 *************************************************************************/
22
40#include "TRestDetectorHitsEvent.h"
41
42#include "TCanvas.h"
43#include "TRandom.h"
44#include "TRestStringHelper.h"
45#include "TRestTools.h"
46#include "TStyle.h"
47
48using namespace std;
49using namespace TMath;
50
52
66 fHits = new TRestHits();
67
68 fPad = nullptr;
69
70 fXYHitGraph = nullptr;
71 fXZHitGraph = nullptr;
72 fYZHitGraph = nullptr;
73
74 fXYHisto = nullptr;
75 fXZHisto = nullptr;
76 fYZHisto = nullptr;
77
78 fXZHits = nullptr;
79 fYZHits = nullptr;
80 fXYZHits = nullptr;
81
82 fXHisto = nullptr;
83 fYHisto = nullptr;
84 fZHisto = nullptr;
85}
86
91
98void TRestDetectorHitsEvent::AddHit(Double_t x, Double_t y, Double_t z, Double_t en, Double_t t,
99 REST_HitType type) {
100 fHits->AddHit({x, y, z}, en, t, type);
101}
102
108void TRestDetectorHitsEvent::AddHit(const TVector3& position, Double_t energy, Double_t time,
109 REST_HitType type) {
110 fHits->AddHit(position, energy, time, type);
111}
112
118 fHits->RemoveHits();
119
120 if (fXZHits) {
121 delete fXZHits;
122 fXZHits = nullptr;
123 }
124 if (fYZHits) {
125 delete fYZHits;
126 fYZHits = nullptr;
127 }
128 if (fXYZHits) {
129 delete fXYZHits;
130 fXYZHits = nullptr;
131 }
132
133 fXZHits = new TRestHits();
134 fYZHits = new TRestHits();
135 fXYZHits = new TRestHits();
136}
137
138void TRestDetectorHitsEvent::Sort(bool(compareCondition)(const TRestHits::iterator& hit1,
139 const TRestHits::iterator& hit2)) {
140 if (compareCondition == 0) {
141 // default sort logic: z from smaller to greater
142 std::sort(fHits->begin(), fHits->end(),
143 [](const TRestHits::iterator& hit1, const TRestHits::iterator& hit2) -> bool {
144 return hit1.z() < hit2.z();
145 });
146 } else {
147 std::sort(fHits->begin(), fHits->end(), compareCondition);
148 }
149}
150
151void TRestDetectorHitsEvent::Shuffle(int NLoop) {
152 Int_t nHits = fHits->GetNumberOfHits();
153 if (nHits >= 2) {
154 for (int n = 0; n < NLoop; n++) {
155 Int_t hit1 = (Int_t)(nHits * gRandom->Uniform(0, 1));
156 Int_t hit2 = (Int_t)(nHits * gRandom->Uniform(0, 1));
157
158 fHits->SwapHits(hit1, hit2);
159 }
160 }
161}
162
172
173 for (unsigned int i = 0; i < this->GetNumberOfHits(); i++) {
174 if (GetType(i) == XZ) {
175 fXZHits->AddHit({this->GetX(i), this->GetY(i), this->GetZ(i)}, this->GetEnergy(i),
176 this->GetTime(i), XZ);
177 }
178 }
179
180 return fXZHits;
181}
182
192
193 for (unsigned int i = 0; i < this->GetNumberOfHits(); i++) {
194 if (GetType(i) == YZ) {
195 fYZHits->AddHit({this->GetX(i), this->GetY(i), this->GetZ(i)}, this->GetEnergy(i),
196 this->GetTime(i), YZ);
197 }
198 }
199
200 return fYZHits;
201}
202
212
213 for (unsigned int i = 0; i < this->GetNumberOfHits(); i++) {
214 if (GetType(i) == XYZ) {
215 {
216 fXYZHits->AddHit({this->GetX(i), this->GetY(i), this->GetZ(i)}, this->GetEnergy(i),
217 this->GetTime(i), XYZ);
218 }
219 }
220 }
221
222 return fXYZHits;
223}
224
232Bool_t TRestDetectorHitsEvent::anyHitInsideCylinder(TVector3 x0, TVector3 x1, Double_t radius) {
233 if (fHits->GetNumberOfHitsInsideCylinder(x0, x1, radius) > 0) {
234 return true;
235 }
236
237 return false;
238}
239
247Bool_t TRestDetectorHitsEvent::allHitsInsideCylinder(TVector3 x0, TVector3 x1, Double_t radius) {
248 if ((size_t)fHits->GetNumberOfHitsInsideCylinder(x0, x1, radius) == GetNumberOfHits()) return true;
249
250 return false;
251}
252
261Double_t TRestDetectorHitsEvent::GetEnergyInCylinder(TVector3 x0, TVector3 x1, Double_t radius) {
262 return fHits->GetEnergyInCylinder(x0, x1, radius);
263}
264
273Int_t TRestDetectorHitsEvent::GetNumberOfHitsInsideCylinder(TVector3 x0, TVector3 x1, Double_t radius) {
274 return fHits->GetNumberOfHitsInsideCylinder(x0, x1, radius);
275}
276
285TVector3 TRestDetectorHitsEvent::GetMeanPositionInCylinder(TVector3 x0, TVector3 x1, Double_t radius) {
286 return fHits->GetMeanPositionInCylinder(x0, x1, radius);
287}
288
298Bool_t TRestDetectorHitsEvent::anyHitInsidePrism(TVector3 x0, TVector3 x1, Double_t sizeX, Double_t sizeY,
299 Double_t theta) {
300 if (fHits->GetNumberOfHitsInsidePrism(x0, x1, sizeX, sizeY, theta) > 0) return true;
301
302 return false;
303}
304
314Bool_t TRestDetectorHitsEvent::allHitsInsidePrism(TVector3 x0, TVector3 x1, Double_t sizeX, Double_t sizeY,
315 Double_t theta) {
316 if ((size_t)fHits->GetNumberOfHitsInsidePrism(x0, x1, sizeX, sizeY, theta) == GetNumberOfHits())
317 return true;
318
319 return false;
320}
321
332Double_t TRestDetectorHitsEvent::GetEnergyInPrism(TVector3 x0, TVector3 x1, Double_t sizeX, Double_t sizeY,
333 Double_t theta) {
334 return fHits->GetEnergyInPrism(x0, x1, sizeX, sizeY, theta);
335}
336
347Int_t TRestDetectorHitsEvent::GetNumberOfHitsInsidePrism(TVector3 x0, TVector3 x1, Double_t sizeX,
348 Double_t sizeY, Double_t theta) {
349 return fHits->GetNumberOfHitsInsidePrism(x0, x1, sizeX, sizeY, theta);
350}
351
362TVector3 TRestDetectorHitsEvent::GetMeanPositionInPrism(TVector3 x0, TVector3 x1, Double_t sizeX,
363 Double_t sizeY, Double_t theta) {
364 return fHits->GetMeanPositionInPrism(x0, x1, sizeX, sizeY, theta);
365}
366
378 Double_t radius) {
379 Double_t rad2 = radius * radius;
380 Double_t hitDistance = rad2;
381
382 Double_t d2, l;
383
384 TVector3 axis = x1 - x0;
385 Double_t cylLength = axis.Mag();
386
387 Int_t nhits = 0;
388 for (unsigned int n = 0; n < GetNumberOfHits(); n++) {
389 if (fHits->isHitNInsideCylinder(n, x0, x1, radius)) {
390 l = axis.Dot(this->GetPosition(n) - x0) / cylLength;
391
392 d2 = rad2 - (this->GetPosition(n) - x0).Mag2() + l * l;
393
394 if (d2 < hitDistance) hitDistance = d2;
395
396 nhits++;
397 }
398 }
399
400 if (nhits == 0) return -1;
401
402 return TMath::Sqrt(hitDistance);
403}
404
416 Double_t radius) {
417 TVector3 axis = x1 - x0;
418 Double_t cylLength = axis.Mag();
419
420 Double_t hitDistance = cylLength;
421 Double_t d = cylLength;
422
423 Int_t nhits = 0;
424 for (unsigned int n = 0; n < GetNumberOfHits(); n++) {
425 if (fHits->isHitNInsideCylinder(n, x0, x1, radius)) {
426 d = cylLength - axis.Dot(this->GetPosition(n) - x0) / cylLength;
427
428 if (d < hitDistance) hitDistance = d;
429
430 nhits++;
431 }
432 }
433
434 if (nhits == 0) return -1;
435
436 return hitDistance;
437}
438
450 Double_t radius) {
451 TVector3 axis = x1 - x0;
452 Double_t cylLength = axis.Mag();
453
454 Double_t hitDistance = cylLength;
455 Double_t d = cylLength;
456
457 Int_t nhits = 0;
458 for (unsigned int n = 0; n < GetNumberOfHits(); n++) {
459 if (fHits->isHitNInsideCylinder(n, x0, x1, radius)) {
460 d = axis.Dot(this->GetPosition(n) - x0) / cylLength;
461
462 if (d < hitDistance) hitDistance = d;
463
464 nhits++;
465 }
466 }
467
468 if (nhits == 0) return -1;
469
470 return hitDistance;
471}
472
486 Double_t sizeX, Double_t sizeY,
487 Double_t theta) {
488 Double_t dX = sizeX / 2;
489 Double_t dY = sizeY / 2;
490
491 Double_t hitDistance = max(dX, dY);
492
493 Double_t d;
494 Int_t nhits = 0;
495 for (unsigned int n = 0; n < GetNumberOfHits(); n++) {
496 if (fHits->isHitNInsidePrism(n, x0, x1, sizeX, sizeY, theta)) {
497 dX = sizeX / 2 - TMath::Abs((this->GetPosition(n) - x0).X());
498 dY = sizeY / 2 - TMath::Abs((this->GetPosition(n) - x0).Y());
499
500 d = min(dX, dY);
501
502 if (d < hitDistance) hitDistance = d;
503
504 nhits++;
505 }
506 }
507
508 if (nhits == 0) return -1;
509
510 return hitDistance;
511}
512
526 Double_t sizeX, Double_t sizeY,
527 Double_t theta) {
528 TVector3 axis = x1 - x0;
529 Double_t prismLength = axis.Mag();
530
531 Double_t hitDistance = prismLength;
532
533 Double_t d;
534 Int_t nhits = 0;
535 for (unsigned int n = 0; n < GetNumberOfHits(); n++) {
536 if (fHits->isHitNInsidePrism(n, x0, x1, sizeX, sizeY, theta)) {
537 d = prismLength - axis.Dot(this->GetPosition(n) - x0) / prismLength;
538
539 if (d < hitDistance) hitDistance = d;
540
541 nhits++;
542 }
543 }
544
545 if (nhits == 0) return -1;
546
547 return hitDistance;
548}
549
563 const TVector3& x1, Double_t sizeX,
564 Double_t sizeY, Double_t theta) {
565 TVector3 axis = x1 - x0;
566 Double_t prismLength = axis.Mag();
567
568 Double_t hitDistance = prismLength;
569
570 Double_t d;
571 Int_t nhits = 0;
572 for (unsigned int n = 0; n < GetNumberOfHits(); n++) {
573 if (fHits->isHitNInsidePrism(n, x0, x1, sizeX, sizeY, theta)) {
574 d = axis.Dot(this->GetPosition(n) - x0) / prismLength;
575
576 if (d < hitDistance) hitDistance = d;
577
578 nhits++;
579 }
580 }
581
582 if (nhits == 0) return -1;
583
584 return hitDistance;
585}
586
611TPad* TRestDetectorHitsEvent::DrawEvent(const TString& option) {
612 vector<TString> optList = Vector_cast<string, TString>(TRestTools::GetOptions((string)option));
613
614 for (unsigned int n = 0; n < optList.size(); n++) {
615 if (optList[n] == "print") this->PrintEvent();
616 }
617
618 optList.erase(std::remove(optList.begin(), optList.end(), "print"), optList.end());
619
621 // which means that it should be extracted from the hit array
622 if (optList.size() == 0) optList.push_back("hist(Cont1,col)[3]");
623
624 if (fPad != nullptr) {
625 delete fPad;
626 fPad = nullptr;
627 }
628
629 fPad = new TPad(this->GetName(), " ", 0, 0, 1, 1);
630 fPad->Divide(3, 2 * optList.size());
631 fPad->Draw();
632
633 Int_t column = 0;
634 for (unsigned int n = 0; n < optList.size(); n++) {
635 string optionStr = (string)optList[n];
636 TString drawEventOption = optList[n];
637
638 // Generating drawOption argument
639 size_t startPos = optionStr.find("(");
640 if (startPos != string::npos) drawEventOption = optList[n](0, startPos);
641
642 // Generating histogram option argument
643 string histOption = "";
644 size_t endPos = optionStr.find(")");
645 if (endPos != string::npos) {
646 TString histOptionTmp = optList[n](startPos + 1, endPos - startPos - 1);
647
648 histOption = (string)histOptionTmp;
649 size_t pos = 0;
650 while ((pos = histOption.find(",", pos)) != string::npos) {
651 histOption.replace(pos, 1, ":");
652 pos = pos + 1;
653 }
654 }
655
656 // Generating histogram pitch argument
657 string pitchOption = "";
658
659 startPos = optionStr.find("[");
660 endPos = optionStr.find("]");
661 Double_t pitch = 0;
662 if (endPos != string::npos) {
663 TString pitchOption = optList[n](startPos + 1, endPos - startPos - 1);
664 pitch = stod((string)pitchOption);
665 }
666
667 if (drawEventOption == "graph") this->DrawGraphs(column);
668
669 if (drawEventOption == "hist") this->DrawHistograms(column, histOption, pitch);
670 }
671
672 return fPad;
673}
674
683 if (fXYHitGraph != nullptr) {
684 delete fXYHitGraph;
685 fXYHitGraph = nullptr;
686 }
687 if (fXZHitGraph != nullptr) {
688 delete fXZHitGraph;
689 fXZHitGraph = nullptr;
690 }
691 if (fYZHitGraph != nullptr) {
692 delete fYZHitGraph;
693 fYZHitGraph = nullptr;
694 }
695
696 vector<vector<Double_t>> xz(2, vector<Double_t>(this->GetNumberOfHits()));
697 vector<vector<Double_t>> yz(2, vector<Double_t>(this->GetNumberOfHits()));
698 vector<vector<Double_t>> xy(2, vector<Double_t>(this->GetNumberOfHits()));
699
700 /* {{{ Creating xz, yz, and xy arrays and initializing graphs */
701 Int_t nXZ = 0;
702 Int_t nYZ = 0;
703 Int_t nXY = 0;
704
705 for (unsigned int nhit = 0; nhit < this->GetNumberOfHits(); nhit++) {
706 Double_t x = fHits->GetX(nhit);
707 Double_t y = fHits->GetY(nhit);
708 Double_t z = fHits->GetZ(nhit);
709 int type = fHits->GetType(nhit);
710
711 if (type % XZ == 0) {
712 xz[0][nXZ] = x;
713 xz[1][nXZ] = z;
714 nXZ++;
715 }
716
717 if (type % YZ == 0) {
718 yz[0][nYZ] = y;
719 yz[1][nYZ] = z;
720 nYZ++;
721 }
722
723 if (type % XY == 0) {
724 xy[0][nXY] = x;
725 xy[1][nXY] = y;
726 nXY++;
727 }
728 }
729
730 fXZHitGraph = new TGraph(nXZ, &xz[0][0], &xz[1][0]);
731 fXZHitGraph->SetMarkerColor(kBlue);
732 fXZHitGraph->SetMarkerSize(0.3);
733 fXZHitGraph->SetMarkerStyle(20);
734
735 fYZHitGraph = new TGraph(nYZ, &yz[0][0], &yz[1][0]);
736 fYZHitGraph->SetMarkerColor(kRed);
737 fYZHitGraph->SetMarkerSize(0.3);
738 fYZHitGraph->SetMarkerStyle(20);
739
740 fXYHitGraph = new TGraph(nXY, &xy[0][0], &xy[1][0]);
741 fXYHitGraph->SetMarkerColor(kBlack);
742 fXYHitGraph->SetMarkerSize(0.3);
743 fXYHitGraph->SetMarkerStyle(20);
744 /* }}} */
745
746 char title[256];
747 sprintf(title, "Event ID %d", this->GetID());
748
749 if (nXZ > 0) {
750 fPad->cd(1 + 3 * column);
751 fXZHitGraph->SetTitle(title);
752 fXZHitGraph->Draw("AP*");
753
754 fXZHitGraph->GetXaxis()->SetTitle("X-axis (mm)");
755 fXZHitGraph->GetYaxis()->SetTitle("Z-axis (mm)");
756 }
757
758 if (nYZ > 0) {
759 fPad->cd(2 + 3 * column);
760 fYZHitGraph->SetTitle(title);
761 fYZHitGraph->Draw("AP");
762
763 fYZHitGraph->GetXaxis()->SetTitle("Y-axis (mm)");
764 fYZHitGraph->GetYaxis()->SetTitle("Z-axis (mm)");
765 }
766
767 if (nXY > 0) {
768 fPad->cd(3 + 3 * column);
769 fXYHitGraph->SetTitle(title);
770 fXYHitGraph->Draw("AP");
771
772 fXYHitGraph->GetXaxis()->SetTitle("X-axis (mm)");
773 fXYHitGraph->GetYaxis()->SetTitle("Y-axis (mm)");
774 }
775
776 column++;
777}
778
793void TRestDetectorHitsEvent::DrawHistograms(Int_t& column, const TString& histOption, double pitch) {
794 std::vector<double> fX, fY, fZ;
795 for (unsigned int i = 0; i < GetNumberOfHits(); i++) {
796 if (GetType(i) % X == 0) fX.emplace_back(GetX(i));
797 if (GetType(i) % Y == 0) fY.emplace_back(GetY(i));
798 if (GetType(i) % Z == 0) fZ.emplace_back(GetZ(i));
799 }
800
801 double maxX, minX, maxY, minY, maxZ, minZ;
802 int nBinsX, nBinsY, nBinsZ;
803 TRestHits::GetBoundaries(fX, maxX, minX, nBinsX);
804 TRestHits::GetBoundaries(fY, maxY, minY, nBinsY);
805 TRestHits::GetBoundaries(fZ, maxZ, minZ, nBinsZ);
806
807 if (pitch > 0) {
808 nBinsX = std::round((maxX - minX) / pitch);
809 nBinsY = std::round((maxY - minY) / pitch);
810 nBinsZ = std::round((maxZ - minZ) / pitch);
811 }
812
813 delete fXYHisto;
814 delete fXZHisto;
815 delete fYZHisto;
816
817 delete fXHisto;
818 delete fYHisto;
819 delete fZHisto;
820
821 fXYHisto = new TH2F("XY", "", nBinsX, minX, maxX, nBinsY, minY, maxY);
822 fXZHisto = new TH2F("XZ", "", nBinsX, minX, maxX, nBinsZ, minZ, maxZ);
823 fYZHisto = new TH2F("YZ", "", nBinsY, minY, maxY, nBinsZ, minZ, maxZ);
824
825 fXHisto = new TH1F("X", "", nBinsX, minX, maxX);
826 fYHisto = new TH1F("Y", "", nBinsY, minY, maxY);
827 fZHisto = new TH1F("Z", "", nBinsZ, minZ, maxZ);
828
829 fXYHisto->SetStats(false);
830 fXZHisto->SetStats(false);
831 fYZHisto->SetStats(false);
832
833 fXHisto->SetStats(false);
834 fYHisto->SetStats(false);
835 fZHisto->SetStats(false);
836
837 Int_t nYZ = 0, nXZ = 0, nXY = 0;
838 Int_t nX = 0, nY = 0, nZ = 0;
839
840 for (unsigned int nhit = 0; nhit < this->GetNumberOfHits(); nhit++) {
841 Double_t x = fHits->GetX(nhit);
842 Double_t y = fHits->GetY(nhit);
843 Double_t z = fHits->GetZ(nhit);
844 Double_t en = fHits->GetEnergy(nhit);
845 int type = fHits->GetType(nhit);
846
847 if (type % XZ == 0) {
848 fXZHisto->Fill(x, z, en);
849 nXZ++;
850 }
851
852 if (type % YZ == 0) {
853 fYZHisto->Fill(y, z, en);
854 nYZ++;
855 }
856
857 if (type % XY == 0) {
858 fXYHisto->Fill(x, y, en);
859 nXY++;
860 }
861
862 if (type % X == 0) {
863 fXHisto->Fill(x);
864 nX++;
865 }
866
867 if (type % Y == 0) {
868 fYHisto->Fill(y);
869 nY++;
870 }
871
872 if (type % Z == 0) {
873 fZHisto->Fill(z);
874 nZ++;
875 }
876 }
877
878 TStyle style;
879 style.SetPalette(1);
880
881 if (nXZ > 0) {
882 fPad->cd(1 + 3 * column);
883 fXZHisto->Draw(histOption);
884 fXZHisto->GetXaxis()->SetTitle("X-axis (mm)");
885 fXZHisto->GetYaxis()->SetTitle("Z-axis (mm)");
886 fXZHisto->GetYaxis()->SetTitleSize(1.4 * fXZHisto->GetYaxis()->GetTitleSize());
887 fXZHisto->GetXaxis()->SetTitleSize(1.4 * fXZHisto->GetXaxis()->GetTitleSize());
888 fXZHisto->GetYaxis()->SetLabelSize(1.25 * fXZHisto->GetYaxis()->GetLabelSize());
889 fXZHisto->GetXaxis()->SetLabelSize(1.25 * fXZHisto->GetXaxis()->GetLabelSize());
890 fXZHisto->GetYaxis()->SetTitleOffset(1);
891 }
892
893 if (nYZ > 0) {
894 fPad->cd(2 + 3 * column);
895 fYZHisto->Draw(histOption);
896 fYZHisto->GetXaxis()->SetTitle("Y-axis (mm)");
897 fYZHisto->GetYaxis()->SetTitle("Z-axis (mm)");
898 fYZHisto->GetYaxis()->SetTitleSize(1.4 * fYZHisto->GetYaxis()->GetTitleSize());
899 fYZHisto->GetXaxis()->SetTitleSize(1.4 * fYZHisto->GetXaxis()->GetTitleSize());
900 fYZHisto->GetYaxis()->SetLabelSize(1.25 * fYZHisto->GetYaxis()->GetLabelSize());
901 fYZHisto->GetXaxis()->SetLabelSize(1.25 * fYZHisto->GetXaxis()->GetLabelSize());
902 fYZHisto->GetYaxis()->SetTitleOffset(1);
903 }
904
905 if (nXY > 0) {
906 fPad->cd(3 + 3 * column);
907 fXYHisto->Draw(histOption);
908 fXYHisto->GetXaxis()->SetTitle("X-axis (mm)");
909 fXYHisto->GetYaxis()->SetTitle("Y-axis (mm)");
910 fXYHisto->GetYaxis()->SetTitleSize(1.4 * fXYHisto->GetYaxis()->GetTitleSize());
911 fXYHisto->GetXaxis()->SetTitleSize(1.4 * fXYHisto->GetXaxis()->GetTitleSize());
912 fXYHisto->GetYaxis()->SetLabelSize(1.25 * fXYHisto->GetYaxis()->GetLabelSize());
913 fXYHisto->GetXaxis()->SetLabelSize(1.25 * fXYHisto->GetXaxis()->GetLabelSize());
914 fXYHisto->GetYaxis()->SetTitleOffset(1);
915 }
916
917 column++;
918
919 if (nX > 0) {
920 fPad->cd(1 + 3 * column);
921 fXHisto->Draw(histOption);
922 fXHisto->GetXaxis()->SetTitle("X-axis (mm)");
923 fXHisto->GetYaxis()->SetTitle("Number of events");
924 fXHisto->GetYaxis()->SetTitleSize(1.4 * fXHisto->GetYaxis()->GetTitleSize());
925 fXHisto->GetXaxis()->SetTitleSize(1.4 * fXHisto->GetXaxis()->GetTitleSize());
926 fXHisto->GetYaxis()->SetLabelSize(1.25 * fXHisto->GetYaxis()->GetLabelSize());
927 fXHisto->GetXaxis()->SetLabelSize(1.25 * fXHisto->GetXaxis()->GetLabelSize());
928 fXHisto->GetYaxis()->SetTitleOffset(1);
929 }
930
931 if (nY > 0) {
932 fPad->cd(2 + 3 * column);
933 fYHisto->Draw(histOption);
934 fYHisto->GetXaxis()->SetTitle("Y-axis (mm)");
935 fYHisto->GetYaxis()->SetTitle("Number of events");
936 fYHisto->GetYaxis()->SetTitleSize(1.4 * fYHisto->GetYaxis()->GetTitleSize());
937 fYHisto->GetXaxis()->SetTitleSize(1.4 * fYHisto->GetXaxis()->GetTitleSize());
938 fYHisto->GetYaxis()->SetLabelSize(1.25 * fYHisto->GetYaxis()->GetLabelSize());
939 fYHisto->GetXaxis()->SetLabelSize(1.25 * fYHisto->GetXaxis()->GetLabelSize());
940 fYHisto->GetYaxis()->SetTitleOffset(1);
941 }
942
943 if (nZ > 0) {
944 fPad->cd(3 + 3 * column);
945 fZHisto->Draw(histOption);
946 fZHisto->GetXaxis()->SetTitle("Z-axis (mm)");
947 fZHisto->GetYaxis()->SetTitle("Number of events");
948 fZHisto->GetYaxis()->SetTitleSize(1.4 * fYHisto->GetYaxis()->GetTitleSize());
949 fZHisto->GetXaxis()->SetTitleSize(1.4 * fYHisto->GetXaxis()->GetTitleSize());
950 fZHisto->GetYaxis()->SetLabelSize(1.25 * fYHisto->GetYaxis()->GetLabelSize());
951 fZHisto->GetXaxis()->SetLabelSize(1.25 * fYHisto->GetXaxis()->GetLabelSize());
952 fZHisto->GetYaxis()->SetTitleOffset(1);
953 }
954
955 column++;
956}
957
958void TRestDetectorHitsEvent::PrintEvent(Int_t nHits) const {
960
961 cout << "Total energy : " << GetTotalEnergy() << endl;
962 cout << "Mean position : ( " << GetMeanPositionX() << " , " << GetMeanPositionY() << " , "
963 << GetMeanPositionZ() << " ) " << endl;
964 cout << "Number of hits : " << fHits->GetNumberOfHits() << endl;
965 if (nHits != -1) {
966 cout << "+++++++++++++++++++++++" << endl;
967 cout << "Printing only the first " << nHits << " hits" << endl;
968 }
969
970 fHits->PrintHits(nHits);
971}
972
984TH2F* TRestDetectorHitsEvent::GetXYHistogram(std::vector<float> ranges, Double_t pitch, Double_t border) {
985 double maxX, minX, maxY, minY;
986 int nBinsX, nBinsY;
987
988 if (ranges.size() == 6) {
989 nBinsX = ranges[0];
990 minX = ranges[1];
991 maxX = ranges[2];
992
993 nBinsY = ranges[3];
994 minY = ranges[4];
995 maxY = ranges[5];
996 } else {
997 std::vector<double> x, y, en;
998
999 for (unsigned int i = 0; i < GetNumberOfHits(); i++) {
1000 if (GetType(i) % X == 0) x.emplace_back(GetX(i));
1001 if (GetType(i) % Y == 0) y.emplace_back(GetY(i));
1002 }
1003
1004 if (x.empty() && y.empty()) {
1005 std::cout << "TRestDetectorHitsEvent::GetXYHistogram. Empty histogram!" << std::endl;
1006 return nullptr;
1007 }
1008
1009 TRestHits::GetBoundaries(x, maxX, minX, nBinsX);
1010 TRestHits::GetBoundaries(y, maxY, minY, nBinsY);
1011
1012 maxX += border;
1013 minX -= border;
1014
1015 maxY += border;
1016 minY -= border;
1017
1018 if (pitch > 0) {
1019 nBinsX = std::round((maxX - minX) / pitch);
1020 nBinsY = std::round((maxY - minY) / pitch);
1021 }
1022 }
1023
1024 delete fXYHisto;
1025 fXYHisto = new TH2F("XY", "", nBinsX, minX, maxX, nBinsY, minY, maxY);
1026 fXYHisto->SetStats(false);
1027
1028 Int_t nXY = 0;
1029
1030 for (unsigned int nhit = 0; nhit < this->GetNumberOfHits(); nhit++) {
1031 Double_t x = fHits->GetX(nhit);
1032 Double_t y = fHits->GetY(nhit);
1033 Double_t en = fHits->GetEnergy(nhit);
1034 int type = fHits->GetType(nhit);
1035
1036 if (type % XY == 0 || type % XYZ == 0) {
1037 fXYHisto->Fill(x, y, en);
1038 nXY++;
1039 }
1040 }
1041
1042 TStyle style;
1043 style.SetPalette(1);
1044
1045 if (nXY > 0) {
1046 fXYHisto->Draw("col");
1047 fXYHisto->Draw("cont3 same");
1048 fXYHisto->GetXaxis()->SetTitle("X-axis (mm)");
1049 fXYHisto->GetYaxis()->SetTitle("Y-axis (mm)");
1050 fXYHisto->GetYaxis()->SetTitleSize(1.4 * fXYHisto->GetYaxis()->GetTitleSize());
1051 fXYHisto->GetXaxis()->SetTitleSize(1.4 * fXYHisto->GetXaxis()->GetTitleSize());
1052 fXYHisto->GetYaxis()->SetLabelSize(1.25 * fXYHisto->GetYaxis()->GetLabelSize());
1053 fXYHisto->GetXaxis()->SetLabelSize(1.25 * fXYHisto->GetXaxis()->GetLabelSize());
1054 fXYHisto->GetYaxis()->SetTitleOffset(1);
1055 }
1056
1057 return fXYHisto;
1058}
1059
1071TH2F* TRestDetectorHitsEvent::GetXZHistogram(std::vector<float> ranges, Double_t pitch, Double_t border) {
1072 double maxX, minX, maxZ, minZ;
1073 int nBinsX, nBinsZ;
1074
1075 if (ranges.size() == 6) {
1076 nBinsX = ranges[0];
1077 minX = ranges[1];
1078 maxX = ranges[2];
1079
1080 nBinsZ = ranges[3];
1081 minZ = ranges[4];
1082 maxZ = ranges[5];
1083 } else {
1084 std::vector<double> x, z;
1085
1086 for (unsigned int i = 0; i < GetNumberOfHits(); i++) {
1087 if (GetType(i) % X == 0) x.emplace_back(GetX(i));
1088 if (GetType(i) % Z == 0) z.emplace_back(GetZ(i));
1089 }
1090
1091 if (x.empty() || z.empty()) {
1092 std::cout << "TRestDetectorHitsEvent::GetXZHistogram. Empty histogram!" << std::endl;
1093 return nullptr;
1094 }
1095
1096 TRestHits::GetBoundaries(x, maxX, minX, nBinsX);
1097 TRestHits::GetBoundaries(z, maxZ, minZ, nBinsZ);
1098
1099 maxX += border;
1100 minX -= border;
1101
1102 maxZ += border;
1103 minZ -= border;
1104
1105 if (pitch > 0) {
1106 nBinsX = std::round((maxX - minX) / pitch);
1107 nBinsZ = std::round((maxZ - minZ) / pitch);
1108 }
1109 }
1110
1111 delete fXZHisto;
1112 fXZHisto = new TH2F("XZ", "", nBinsX, minX, maxX, nBinsZ, minZ, maxZ);
1113 fXZHisto->SetStats(false);
1114
1115 Int_t nXZ = 0;
1116
1117 for (unsigned int nhit = 0; nhit < this->GetNumberOfHits(); nhit++) {
1118 Double_t x = fHits->GetX(nhit);
1119 Double_t z = fHits->GetZ(nhit);
1120 Double_t en = fHits->GetEnergy(nhit);
1121 int type = fHits->GetType(nhit);
1122
1123 if (type % XZ == 0 || type % XYZ == 0) {
1124 fXZHisto->Fill(x, z, en);
1125 nXZ++;
1126 }
1127 }
1128
1129 TStyle style;
1130 style.SetPalette(1);
1131
1132 if (nXZ > 0) {
1133 fXZHisto->Draw("col");
1134 fXZHisto->Draw("cont3 same");
1135 fXZHisto->GetXaxis()->SetTitle("X-axis (mm)");
1136 fXZHisto->GetYaxis()->SetTitle("Z-axis (mm)");
1137 fXZHisto->GetYaxis()->SetTitleSize(1.4 * fXZHisto->GetYaxis()->GetTitleSize());
1138 fXZHisto->GetXaxis()->SetTitleSize(1.4 * fXYHisto->GetXaxis()->GetTitleSize());
1139 fXZHisto->GetYaxis()->SetLabelSize(1.25 * fXYHisto->GetYaxis()->GetLabelSize());
1140 fXZHisto->GetXaxis()->SetLabelSize(1.25 * fXYHisto->GetXaxis()->GetLabelSize());
1141 fXZHisto->GetYaxis()->SetTitleOffset(1);
1142 }
1143
1144 return fXZHisto;
1145}
1146
1158TH2F* TRestDetectorHitsEvent::GetYZHistogram(std::vector<float> ranges, Double_t pitch, Double_t border) {
1159 double maxY, minY, maxZ, minZ;
1160 int nBinsY, nBinsZ;
1161
1162 if (ranges.size() == 6) {
1163 nBinsY = ranges[0];
1164 minY = ranges[1];
1165 maxY = ranges[2];
1166
1167 nBinsZ = ranges[3];
1168 minZ = ranges[4];
1169 maxZ = ranges[5];
1170 } else {
1171 std::vector<double> y, z, en;
1172
1173 for (unsigned int i = 0; i < GetNumberOfHits(); i++) {
1174 if (GetType(i) % Y == 0) y.emplace_back(GetY(i));
1175 if (GetType(i) % Z == 0) z.emplace_back(GetZ(i));
1176 }
1177
1178 if (y.empty() || z.empty()) {
1179 std::cout << "TRestDetectorHitsEvent::GetYZHistogram. Empty histogram!" << std::endl;
1180 return nullptr;
1181 }
1182
1183 TRestHits::GetBoundaries(y, maxY, minY, nBinsY);
1184 TRestHits::GetBoundaries(z, maxZ, minZ, nBinsZ);
1185
1186 maxY += border;
1187 minY -= border;
1188
1189 maxZ += border;
1190 minZ -= border;
1191
1192 if (pitch > 0) {
1193 nBinsY = std::round((maxY - minZ) / pitch);
1194 nBinsZ = std::round((maxZ - minZ) / pitch);
1195 }
1196 }
1197
1198 delete fYZHisto;
1199 fYZHisto = new TH2F("YZ", "", nBinsY, minY, maxY, nBinsZ, minZ, maxZ);
1200 fYZHisto->SetStats(false);
1201
1202 Int_t nYZ = 0;
1203
1204 for (unsigned int nhit = 0; nhit < this->GetNumberOfHits(); nhit++) {
1205 Double_t y = fHits->GetY(nhit);
1206 Double_t z = fHits->GetZ(nhit);
1207 Double_t en = fHits->GetEnergy(nhit);
1208 int type = fHits->GetType(nhit);
1209
1210 if (type % YZ == 0 || type % XYZ == 0) {
1211 fYZHisto->Fill(y, z, en);
1212 nYZ++;
1213 }
1214 }
1215
1216 TStyle style;
1217 style.SetPalette(1);
1218
1219 if (nYZ > 0) {
1220 fYZHisto->Draw("col");
1221 fYZHisto->Draw("cont3 same");
1222 fYZHisto->GetXaxis()->SetTitle("Y-axis (mm)");
1223 fYZHisto->GetYaxis()->SetTitle("Z-axis (mm)");
1224 fYZHisto->GetYaxis()->SetTitleSize(1.4 * fYZHisto->GetYaxis()->GetTitleSize());
1225 fYZHisto->GetXaxis()->SetTitleSize(1.4 * fYZHisto->GetXaxis()->GetTitleSize());
1226 fYZHisto->GetYaxis()->SetLabelSize(1.25 * fYZHisto->GetYaxis()->GetLabelSize());
1227 fYZHisto->GetXaxis()->SetLabelSize(1.25 * fYZHisto->GetXaxis()->GetLabelSize());
1228 fYZHisto->GetYaxis()->SetTitleOffset(1);
1229 }
1230
1231 return fYZHisto;
1232}
TH2F * GetXZHistogram(std::vector< float > ranges, Double_t pitch=3, Double_t border=5)
This method draws the hits found on XY as a TH2F and it returns the generated histogram.
TRestHits * GetXYZHits()
This method collects all hits which are compatible with a XYZ hit.
virtual void PrintEvent() const
TVector3 GetMeanPositionInCylinder(TVector3 x0, TVector3 x1, Double_t radius)
This method returns the mean position of the hits found inside the cylinder volume given by argument.
TRestHits * GetXZHits()
This method collects all hits which are compatible with a XZ-projected hit.
Double_t GetClosestHitInsideDistanceToCylinderTop(TVector3 x0, TVector3 x1, Double_t radius)
This method returns the distance to the cylinder top face from the closest hit contained inside the c...
void DrawGraphs(Int_t &column)
This method draw the hits events as a graph.
TH1F * fZHisto
An auxiliary TH1F histogram to visualize hits on Z-projection.
TPad * DrawEvent(const TString &option="")
This method draws the hits event structure into a TPad.
Double_t GetX(int n) const
Returns the X-coordinate of hit entry n in mm.
TH2F * fXYHisto
An auxiliary TH2F histogram to visualize hits on XY-projection.
Double_t GetClosestHitInsideDistanceToPrismBottom(const TVector3 &x0, const TVector3 &x1, Double_t sizeX, Double_t sizeY, Double_t theta)
This method returns the distance to the prism bottom face from the closest hit contained inside the p...
Int_t GetNumberOfHitsInsideCylinder(TVector3 x0, TVector3 x1, Double_t radius)
This method returns the total number hits found inside the cylinder volume given by argument.
TH2F * GetXYHistogram(std::vector< float > ranges, Double_t pitch=3, Double_t border=5)
This method draws the hits found on XY as a TH2F and it returns the generated histogram.
Bool_t anyHitInsideCylinder(TVector3 x0, TVector3 x1, Double_t radius)
This method returns true if at least 1 hit is found inside the cylinder volume given by argument.
TH2F * fYZHisto
An auxiliary TH2F histogram to visualize hits on YZ-projection.
Bool_t allHitsInsidePrism(TVector3 x0, TVector3 x1, Double_t sX, Double_t sY, Double_t theta)
This method returns true if all hits are found inside the prism volume given by argument.
void DrawHistograms(Int_t &column, const TString &histOption="", double pitch=0)
This method draw the hits events as an histogram.
Double_t GetEnergyInPrism(TVector3 x0, TVector3 x1, Double_t sizeX, Double_t sizeY, Double_t theta)
This method returns the total integrated energy of all hits found inside the prism volume given by ar...
Double_t GetClosestHitInsideDistanceToCylinderWall(TVector3 x0, TVector3 x1, Double_t radius)
This method returns the distance to the cylinder wall from the closest hit contained inside the cylin...
TGraph * fYZHitGraph
An auxiliary TGraph pointer to visualize hits on YZ-projection.
Double_t GetY(int n) const
Returns the Y-coordinate of hit entry n in mm.
Double_t GetZ(int n) const
Returns the Z-coordinate of hit entry n in mm.
TRestHits * fHits
The hits structure that is is saved to disk.
TRestHits * fXZHits
An auxiliary TRestHits structure to register hits on XZ projection.
TH2F * fXZHisto
An auxiliary TH2F histogram to visualize hits on XZ-projection.
TRestDetectorHitsEvent()
TRestDetectorHitsEvent default constructor.
Int_t GetNumberOfHitsInsidePrism(TVector3 x0, TVector3 x1, Double_t sizeX, Double_t sizeY, Double_t theta)
This method returns the total number of hits found inside the prism volume given by argument.
Double_t GetClosestHitInsideDistanceToCylinderBottom(TVector3 x0, TVector3 x1, Double_t radius)
This method returns the distance to the cylinder bottom face from the closest hit contained inside th...
TH1F * fYHisto
An auxiliary TH1F histogram to visualize hits on Y-projection.
TRestHits * fXYZHits
An auxiliary TRestHits structure to register hits on XYZ projection.
TVector3 GetMeanPositionInPrism(TVector3 x0, TVector3 x1, Double_t sizeX, Double_t sizeY, Double_t theta)
This method returns the mean position of all hits found inside the prism volume given by argument.
Bool_t allHitsInsideCylinder(TVector3 x0, TVector3 x1, Double_t radius)
This method returns true if all hits are contained inside the cylinder volume given by argument.
Double_t GetEnergyInCylinder(TVector3 x0, TVector3 x1, Double_t radius)
This method returns the total integrated energy of all hits found inside the cylinder volume given by...
void AddHit(Double_t x, Double_t y, Double_t z, Double_t en, Double_t t=0, REST_HitType type=XYZ)
Adds a new hit to this event.
TRestHits * GetYZHits()
This method collects all hits which are compatible with a YZ-projected hit.
TGraph * fXZHitGraph
An auxiliary TGraph pointer to visualize hits on XZ-projection.
TH2F * GetYZHistogram(std::vector< float > ranges, Double_t pitch=3, Double_t border=5)
This method draws the hits found on YZ as a TH2F and it returns the generated histogram.
TH1F * fXHisto
An auxiliary TH1F histogram to visualize hits on X-projection.
Double_t GetClosestHitInsideDistanceToPrismTop(TVector3 x0, TVector3 x1, Double_t sizeX, Double_t sizeY, Double_t theta)
This method returns the distance to the prism top face from the closest hit contained inside the pris...
Bool_t anyHitInsidePrism(TVector3 x0, TVector3 x1, Double_t sizeX, Double_t sizeY, Double_t theta)
This method returns true if at least 1 hit is found inside the prism volume given by argument.
virtual void Initialize()
Removes all hits from this event, and clears all auxiliar variables.
TGraph * fXYHitGraph
An auxiliary TGraph pointer to visualize hits on XY-projection.
Double_t GetClosestHitInsideDistanceToPrismWall(TVector3 x0, TVector3 x1, Double_t sizeX, Double_t sizeY, Double_t theta)
This method returns the distance to the prism wall from the closest hit contained inside the prism vo...
~TRestDetectorHitsEvent()
TRestDetectorHitsEvent default destructor.
TRestHits * fYZHits
An auxiliary TRestHits structure to register hits on YZ projection.
virtual void PrintEvent() const
virtual void Initialize()=0
It saves a 3-coordinate position and an energy for each punctual deposition.
Definition TRestHits.h:39
Int_t GetNumberOfHitsInsidePrism(const TVector3 &x0, const TVector3 &x1, Double_t sizeX, Double_t sizeY, Double_t theta) const
It determines the total number of hits contained inside a prisma delimited between x0 and y0 vertex,...
virtual void RemoveHits()
It removes all hits inside the class.
TVector3 GetMeanPositionInCylinder(const TVector3 &x0, const TVector3 &x1, Double_t radius) const
It determines the mean position using the hits contained inside a cylinder with a given radius and de...
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...
Double_t GetEnergyInCylinder(const TVector3 &x0, const TVector3 &x1, Double_t radius) const
It determines the total energy contained inside a cylinder with a given radius and delimited between ...
Bool_t isHitNInsidePrism(Int_t n, const TVector3 &x0, const TVector3 &x1, Double_t sizeX, Double_t sizeY, Double_t theta) const
It determines if hit n is contained inside a prisma delimited between x0 and y0 vertex,...
Double_t GetEnergyInPrism(const TVector3 &x0, const TVector3 &x1, Double_t sizeX, Double_t sizeY, Double_t theta) const
It determines the total hit energy contained inside a prisma delimited between x0 and y0 vertex,...
void AddHit(const TVector3 &position, Double_t energy, Double_t time=0, REST_HitType type=XYZ)
Adds a new hit to the list of hits using a TVector3.
Int_t GetNumberOfHitsInsideCylinder(const TVector3 &x0, const TVector3 &x1, Double_t radius) const
It determines the total number of hits contained inside a cylinder with a given radius and delimited ...
Bool_t isHitNInsideCylinder(Int_t n, const TVector3 &x0, const TVector3 &x1, Double_t radius) const
It determines if hit n is contained inside a cylinder with a given radius and delimited between x0 an...
TVector3 GetMeanPositionInPrism(const TVector3 &x0, const TVector3 &x1, Double_t sizeX, Double_t sizeY, Double_t theta) const
It determines the mean position of hits contained inside a prisma delimited between x0 and x1 vertex,...
virtual void SwapHits(Int_t i, Int_t j)
It exchanges hits n and m affecting to the ordering of the hits inside the list of hits.
virtual void PrintHits(Int_t nHits=-1) const
It prints on screen the first nHits from the list.
static std::vector< std::string > GetOptions(std::string optionsStr)
Returns all the options in an option string.