18 #include "TRestGeant4Event.h"
22 #include <TRestStringHelper.h>
23 #include <TRestTools.h>
26 #include "TRestGeant4Metadata.h"
32 TRestGeant4Event::TRestGeant4Event() {
39 TRestGeant4Event::~TRestGeant4Event() {
46 fPrimaryParticleNames.clear();
51 fXZHitGraph =
nullptr;
52 fYZHitGraph =
nullptr;
53 fXYHitGraph =
nullptr;
55 fXZMultiGraph =
nullptr;
56 fYZMultiGraph =
nullptr;
57 fXYMultiGraph =
nullptr;
73 fTotalDepositedEnergy = 0;
74 fSensitiveVolumeEnergy = 0;
89 void TRestGeant4Event::AddActiveVolume(
const string& volumeName) {
91 fVolumeStored.push_back(1);
92 fVolumeStoredNames.push_back(volumeName);
93 fVolumeDepositedEnergy.push_back(0);
96 void TRestGeant4Event::ClearVolumes() {
98 fVolumeStored.clear();
99 fVolumeStoredNames.clear();
100 fVolumeDepositedEnergy.clear();
103 void TRestGeant4Event::AddEnergyDepositToVolume(Int_t volID, Double_t eDep) {
104 fVolumeDepositedEnergy[volID] += eDep;
107 TVector3 TRestGeant4Event::GetMeanPositionInVolume(Int_t volID)
const {
111 for (
unsigned int t = 0; t < GetNumberOfTracks(); t++) {
112 const auto tck = GetTrack(t);
113 if (tck.GetEnergyInVolume(volID) > 0) {
114 pos += tck.GetMeanPositionInVolume(volID) * tck.GetEnergyInVolume(volID);
116 eDep += tck.GetEnergyInVolume(volID);
121 Double_t nan = TMath::QuietNaN();
122 return {nan, nan, nan};
125 pos = (1 / eDep) * pos;
136 TVector3 meanPos = this->GetMeanPositionInVolume(volID);
138 Double_t nan = TMath::QuietNaN();
139 if (meanPos == TVector3(nan, nan, nan))
return TVector3(nan, nan, nan);
141 TRestHits hitsInVolume = GetHitsInVolume(volID);
144 TVector3 deviation = TVector3(0, 0, 0);
146 for (
unsigned int n = 0; n < hitsInVolume.GetNumberOfHits(); n++) {
147 Double_t en = hitsInVolume.GetEnergy(n);
148 TVector3 diff = meanPos - hitsInVolume.
GetPosition(n);
149 diff.SetXYZ(diff.X() * diff.X(), diff.Y() * diff.Y(), diff.Z() * diff.Z());
151 deviation = deviation + en * diff;
155 deviation = (1. / edep) * deviation;
167 for (
unsigned int t = 0; t < GetNumberOfTracks(); t++)
168 if (GetTrack(t).GetEnergyInVolume(volID) > 0)
return GetTrack(t).GetFirstPositionInVolume(volID);
170 Double_t nan = TMath::QuietNaN();
171 return {nan, nan, nan};
184 for (
int t = GetNumberOfTracks() - 1; t >= 0; t--)
185 if (GetTrack(t).GetEnergyInVolume(volID) > 0)
return GetTrack(t).GetLastPositionInVolume(volID);
187 Double_t nan = TMath::QuietNaN();
188 return {nan, nan, nan};
193 if (fTrackIDToTrackIndex.count(trackID) > 0) {
194 result =
const_cast<TRestGeant4Track*
>(&fTracks[fTrackIDToTrackIndex.at(trackID)]);
195 if (result ==
nullptr || result->GetTrackID() != trackID) {
196 cerr <<
"TRestGeant4Event::GetTrackByID - ERROR: trackIDToTrackIndex map is corrupted" << endl;
209 size_t numberOfHits = 0;
210 for (
const auto& track : fTracks) {
211 numberOfHits += track.GetNumberOfHits(volID);
222 size_t numberOfHits = 0;
223 for (
const auto& track : fTracks) {
224 numberOfHits += track.GetNumberOfPhysicalHits(volID);
236 for (
unsigned int t = 0; t < GetNumberOfTracks(); t++) {
237 const auto& g4Hits = GetTrack(t).GetHits();
238 for (
unsigned int n = 0; n < g4Hits.GetNumberOfHits(); n++) {
239 if (volID != -1 && g4Hits.GetVolumeId(n) != volID)
continue;
241 Double_t x = g4Hits.GetX(n);
242 Double_t y = g4Hits.GetY(n);
243 Double_t z = g4Hits.GetZ(n);
244 Double_t en = g4Hits.GetEnergy(n);
246 hits.
AddHit({x, y, z}, en);
253 Int_t TRestGeant4Event::GetNumberOfTracksForParticle(
const TString& particleName)
const {
255 for (
const auto& track : fTracks) {
256 if (particleName.EqualTo(track.GetParticleName())) {
263 Double_t TRestGeant4Event::GetEnergyDepositedByParticle(
const TString& particleName)
const {
265 for (
const auto& track : fTracks) {
266 if (particleName.EqualTo(track.GetParticleName())) {
267 energy += track.GetTotalEnergy();
273 void TRestGeant4Event::SetBoundaries(Double_t xMin, Double_t xMax, Double_t yMin, Double_t yMax,
274 Double_t zMin, Double_t zMax) {
294 Double_t dX = fMaxX - fMinX;
295 Double_t dY = fMaxY - fMinY;
296 Double_t dZ = fMaxZ - fMinZ;
298 return TMath::Sqrt(dX * dX + dY * dY + dZ * dZ);
301 void TRestGeant4Event::SetBoundaries() {
302 Double_t maxX = -1e10, minX = 1e10, maxZ = -1e10, minZ = 1e10, maxY = -1e10, minY = 1e10;
303 Double_t minEnergy = 1e10, maxEnergy = -1e10;
306 for (
unsigned int ntck = 0; ntck < this->GetNumberOfTracks(); ntck++) {
307 Int_t nHits = GetTrack(ntck).GetNumberOfHits();
309 const auto& hits = GetTrack(ntck).GetHits();
311 for (
int nhit = 0; nhit < nHits; nhit++) {
312 Double_t x = hits.GetX(nhit);
313 Double_t y = hits.GetY(nhit);
314 Double_t z = hits.GetZ(nhit);
316 Double_t en = hits.GetEnergy(nhit);
318 if (en <= 0)
continue;
320 if (x > maxX) maxX = x;
321 if (x < minX) minX = x;
322 if (y > maxY) maxY = y;
323 if (y < minY) minY = y;
324 if (z > maxZ) maxZ = z;
325 if (z < minZ) minZ = z;
327 if (en > maxEnergy) maxEnergy = en;
328 if (en < minEnergy) minEnergy = en;
341 fMaxEnergy = maxEnergy;
342 fMinEnergy = minEnergy;
348 TMultiGraph* TRestGeant4Event::GetXYMultiGraph(Int_t gridElement, vector<TString> pcsList,
349 Double_t minPointSize, Double_t maxPointSize) {
351 delete[] fXYHitGraph;
352 fXYHitGraph =
nullptr;
354 fXYHitGraph =
new TGraph[fTotalHits];
356 fXYMultiGraph =
new TMultiGraph();
360 fLegend_XY =
nullptr;
363 fLegend_XY =
new TLegend(0.75, 0.75, 0.9, 0.85);
366 sprintf(title,
"Event ID %d (XY)", this->GetID());
367 fPad->cd(gridElement)->DrawFrame(fMinX - 10, fMinY - 10, fMaxX + 10, fMaxY + 10, title);
369 Double_t m = 1, n = 0;
370 if (fMaxEnergy - fMinEnergy > 0) {
371 m = (maxPointSize - minPointSize) / (fMaxEnergy - fMinEnergy);
372 n = maxPointSize - m * fMaxEnergy;
375 for (
unsigned int n = 0; n < fXYPcsMarker.size(); n++)
delete fXYPcsMarker[n];
376 fXYPcsMarker.clear();
379 for (
unsigned int p = 0; p < pcsList.size(); p++) legendAdded.push_back(0);
382 for (
unsigned int ntck = 0; ntck < GetNumberOfTracks(); ntck++) {
383 const auto hits = GetTrack(ntck).GetHits();
385 EColor color = GetTrack(ntck).GetParticleColor();
387 for (
unsigned int nhit = 0; nhit < hits.GetNumberOfHits(); nhit++) {
388 Double_t xPos = hits.GetX(nhit);
389 Double_t yPos = hits.GetY(nhit);
390 Double_t energy = hits.GetEnergy(nhit);
392 for (
unsigned int p = 0; p < pcsList.size(); p++) {
393 if (GetTrack(ntck).GetProcessName(hits.GetProcess(nhit)) == pcsList[p]) {
394 TGraph* pcsGraph =
new TGraph(1);
396 pcsGraph->SetPoint(0, xPos, yPos);
397 pcsGraph->SetMarkerColor(kBlack);
398 pcsGraph->SetMarkerSize(3);
399 pcsGraph->SetMarkerStyle(30 + p);
401 fXYPcsMarker.push_back(pcsGraph);
403 if (legendAdded[p] == 0) {
404 fLegend_XY->AddEntry(pcsGraph, GetTrack(ntck).GetProcessName(hits.GetProcess(nhit)),
411 Double_t radius = m * energy + n;
412 if (fMaxEnergy - fMinEnergy <= 0) radius = 2.5;
414 fXYHitGraph[cnt].SetPoint(0, xPos, yPos);
416 fXYHitGraph[cnt].SetMarkerColor(color);
417 fXYHitGraph[cnt].SetMarkerSize(radius);
418 fXYHitGraph[cnt].SetMarkerStyle(20);
420 fXYMultiGraph->Add(&fXYHitGraph[cnt]);
426 fXYMultiGraph->GetXaxis()->SetTitle(
"X-axis (mm)");
427 fXYMultiGraph->GetXaxis()->SetTitleSize(1.25 * fXYMultiGraph->GetXaxis()->GetTitleSize());
428 fXYMultiGraph->GetXaxis()->SetTitleOffset(1.25);
429 fXYMultiGraph->GetXaxis()->SetLabelSize(1.25 * fXYMultiGraph->GetXaxis()->GetLabelSize());
431 fXYMultiGraph->GetYaxis()->SetTitle(
"Y-axis (mm)");
432 fXYMultiGraph->GetYaxis()->SetTitleSize(1.25 * fXYMultiGraph->GetYaxis()->GetTitleSize());
433 fXYMultiGraph->GetYaxis()->SetTitleOffset(1.75);
434 fXYMultiGraph->GetYaxis()->SetLabelSize(1.25 * fXYMultiGraph->GetYaxis()->GetLabelSize());
435 fPad->cd(gridElement);
437 return fXYMultiGraph;
440 TMultiGraph* TRestGeant4Event::GetYZMultiGraph(Int_t gridElement, vector<TString> pcsList,
441 Double_t minPointSize, Double_t maxPointSize) {
443 delete[] fYZHitGraph;
444 fYZHitGraph =
nullptr;
446 fYZHitGraph =
new TGraph[fTotalHits];
448 fYZMultiGraph =
new TMultiGraph();
452 fLegend_YZ =
nullptr;
455 fLegend_YZ =
new TLegend(0.75, 0.75, 0.9, 0.85);
458 sprintf(title,
"Event ID %d (YZ)", this->GetID());
459 fPad->cd(gridElement)->DrawFrame(fMinY - 10, fMinZ - 10, fMaxY + 10, fMaxZ + 10, title);
461 for (
unsigned int n = 0; n < fYZPcsMarker.size(); n++)
delete fYZPcsMarker[n];
462 fYZPcsMarker.clear();
464 Double_t m = 1, n = 0;
465 if (fMaxEnergy - fMinEnergy > 0) {
466 m = (maxPointSize - minPointSize) / (fMaxEnergy - fMinEnergy);
467 n = maxPointSize - m * fMaxEnergy;
471 for (
unsigned int p = 0; p < pcsList.size(); p++) legendAdded.push_back(0);
474 for (
unsigned int ntck = 0; ntck < GetNumberOfTracks(); ntck++) {
475 const auto& hits = GetTrack(ntck).GetHits();
477 EColor color = GetTrack(ntck).GetParticleColor();
479 for (
unsigned int nhit = 0; nhit < hits.GetNumberOfHits(); nhit++) {
480 Double_t yPos = hits.GetY(nhit);
481 Double_t zPos = hits.GetZ(nhit);
482 Double_t energy = hits.GetEnergy(nhit);
484 for (
unsigned int p = 0; p < pcsList.size(); p++) {
485 if (GetTrack(ntck).GetProcessName(hits.GetProcess(nhit)) == pcsList[p]) {
486 TGraph* pcsGraph =
new TGraph(1);
488 pcsGraph->SetPoint(0, yPos, zPos);
489 pcsGraph->SetMarkerColor(kBlack);
490 pcsGraph->SetMarkerSize(3);
491 pcsGraph->SetMarkerStyle(30 + p);
493 fYZPcsMarker.push_back(pcsGraph);
495 if (legendAdded[p] == 0) {
496 fLegend_YZ->AddEntry(pcsGraph, GetTrack(ntck).GetProcessName(hits.GetProcess(nhit)),
503 Double_t radius = m * energy + n;
504 if (fMaxEnergy - fMinEnergy <= 0) radius = 2.5;
506 fYZHitGraph[cnt].SetPoint(0, yPos, zPos);
508 fYZHitGraph[cnt].SetMarkerColor(color);
509 fYZHitGraph[cnt].SetMarkerSize(radius);
510 fYZHitGraph[cnt].SetMarkerStyle(20);
512 fYZMultiGraph->Add(&fYZHitGraph[cnt]);
518 fYZMultiGraph->GetXaxis()->SetTitle(
"Y-axis (mm)");
519 fYZMultiGraph->GetXaxis()->SetTitleSize(1.25 * fYZMultiGraph->GetXaxis()->GetTitleSize());
520 fYZMultiGraph->GetXaxis()->SetTitleOffset(1.25);
521 fYZMultiGraph->GetXaxis()->SetLabelSize(1.25 * fYZMultiGraph->GetXaxis()->GetLabelSize());
523 fYZMultiGraph->GetYaxis()->SetTitle(
"Z-axis (mm)");
524 fYZMultiGraph->GetYaxis()->SetTitleSize(1.25 * fYZMultiGraph->GetYaxis()->GetTitleSize());
525 fYZMultiGraph->GetYaxis()->SetTitleOffset(1.75);
526 fYZMultiGraph->GetYaxis()->SetLabelSize(1.25 * fYZMultiGraph->GetYaxis()->GetLabelSize());
527 fPad->cd(gridElement);
529 return fYZMultiGraph;
532 TMultiGraph* TRestGeant4Event::GetXZMultiGraph(Int_t gridElement, vector<TString> pcsList,
533 Double_t minPointSize, Double_t maxPointSize) {
535 delete[] fXZHitGraph;
536 fXZHitGraph =
nullptr;
538 fXZHitGraph =
new TGraph[fTotalHits];
540 fXZMultiGraph =
new TMultiGraph();
544 fLegend_XZ =
nullptr;
547 fLegend_XZ =
new TLegend(0.75, 0.75, 0.9, 0.85);
550 sprintf(title,
"Event ID %d (XZ)", this->GetID());
551 fPad->cd(gridElement)->DrawFrame(fMinX - 10, fMinZ - 10, fMaxX + 10, fMaxZ + 10, title);
553 for (
unsigned int n = 0; n < fXZPcsMarker.size(); n++)
delete fXZPcsMarker[n];
554 fXZPcsMarker.clear();
556 Double_t m = 1, n = 0;
557 if (fMaxEnergy - fMinEnergy > 0) {
558 m = (maxPointSize - minPointSize) / (fMaxEnergy - fMinEnergy);
559 n = maxPointSize - m * fMaxEnergy;
563 for (
unsigned int p = 0; p < pcsList.size(); p++) legendAdded.push_back(0);
566 for (
unsigned int ntck = 0; ntck < GetNumberOfTracks(); ntck++) {
567 const auto& hits = GetTrack(ntck).GetHits();
569 EColor color = GetTrack(ntck).GetParticleColor();
571 for (
unsigned int nhit = 0; nhit < hits.GetNumberOfHits(); nhit++) {
572 Double_t xPos = hits.GetX(nhit);
573 Double_t zPos = hits.GetZ(nhit);
574 Double_t energy = hits.GetEnergy(nhit);
576 for (
unsigned int p = 0; p < pcsList.size(); p++) {
577 if (GetTrack(ntck).GetProcessName(hits.GetProcess(nhit)) == pcsList[p]) {
578 TGraph* pcsGraph =
new TGraph(1);
580 pcsGraph->SetPoint(0, xPos, zPos);
581 pcsGraph->SetMarkerColor(kBlack);
582 pcsGraph->SetMarkerSize(3);
583 pcsGraph->SetMarkerStyle(30 + p);
585 fXZPcsMarker.push_back(pcsGraph);
587 if (legendAdded[p] == 0) {
588 fLegend_XZ->AddEntry(pcsGraph, GetTrack(ntck).GetProcessName(hits.GetProcess(nhit)),
595 Double_t radius = m * energy + n;
596 if (fMaxEnergy - fMinEnergy <= 0) radius = 2.5;
598 fXZHitGraph[cnt].SetPoint(0, xPos, zPos);
600 fXZHitGraph[cnt].SetMarkerColor(color);
601 fXZHitGraph[cnt].SetMarkerSize(radius);
602 fXZHitGraph[cnt].SetMarkerStyle(20);
604 fXZMultiGraph->Add(&fXZHitGraph[cnt]);
610 fXZMultiGraph->GetXaxis()->SetTitle(
"X-axis (mm)");
611 fXZMultiGraph->GetXaxis()->SetTitleOffset(1.25);
612 fXZMultiGraph->GetXaxis()->SetTitleSize(1.25 * fXZMultiGraph->GetXaxis()->GetTitleSize());
613 fXZMultiGraph->GetXaxis()->SetLabelSize(1.25 * fXZMultiGraph->GetXaxis()->GetLabelSize());
615 fXZMultiGraph->GetYaxis()->SetTitle(
"Z-axis (mm)");
616 fXZMultiGraph->GetYaxis()->SetTitleOffset(1.75);
617 fXZMultiGraph->GetYaxis()->SetTitleSize(1.25 * fXZMultiGraph->GetYaxis()->GetTitleSize());
618 fXZMultiGraph->GetYaxis()->SetLabelSize(1.25 * fXZMultiGraph->GetYaxis()->GetLabelSize());
619 fPad->cd(gridElement);
621 return fXZMultiGraph;
626 TH2F* TRestGeant4Event::GetXYHistogram(Int_t gridElement, vector<TString> optList) {
632 Bool_t stats =
false;
634 for (
unsigned int n = 0; n < optList.size(); n++) {
635 if (optList[n].Contains(
"binSize=")) {
636 TString pitchStr = optList[n].Remove(0, 8);
637 pitch = stod((
string)pitchStr);
640 if (optList[n].Contains(
"stats")) stats =
true;
643 Int_t nBinsX = (fMaxX - fMinX + 20) / pitch;
644 Int_t nBinsY = (fMaxY - fMinY + 20) / pitch;
646 fXYHisto =
new TH2F(
"XY",
"", nBinsX, fMinX - 10, fMinX + pitch * nBinsX, nBinsY, fMinY - 10,
647 fMinY + pitch * nBinsY);
649 for (
unsigned int ntck = 0; ntck < GetNumberOfTracks(); ntck++) {
650 const auto& hits = GetTrack(ntck).GetHits();
652 for (
unsigned int nhit = 0; nhit < hits.GetNumberOfHits(); nhit++) {
653 Double_t xPos = hits.GetX(nhit);
654 Double_t yPos = hits.GetY(nhit);
655 Double_t energy = hits.GetEnergy(nhit);
657 fXYHisto->Fill(xPos, yPos, energy);
663 style.SetOptStat(1001111);
665 fXYHisto->SetStats(stats);
668 sprintf(title,
"Event ID %d (XY)", this->GetID());
669 fXYHisto->SetTitle(title);
671 fXYHisto->GetXaxis()->SetTitle(
"X-axis (mm)");
672 fXYHisto->GetXaxis()->SetTitleOffset(1.25);
673 fXYHisto->GetXaxis()->SetTitleSize(1.25 * fXYHisto->GetXaxis()->GetTitleSize());
674 fXYHisto->GetXaxis()->SetLabelSize(1.25 * fXYHisto->GetXaxis()->GetLabelSize());
676 fXYHisto->GetYaxis()->SetTitle(
"Y-axis (mm)");
677 fXYHisto->GetYaxis()->SetTitleOffset(1.75);
678 fXYHisto->GetYaxis()->SetTitleSize(1.25 * fXYHisto->GetYaxis()->GetTitleSize());
679 fXYHisto->GetYaxis()->SetLabelSize(1.25 * fXYHisto->GetYaxis()->GetLabelSize());
680 fPad->cd(gridElement);
685 TH2F* TRestGeant4Event::GetXZHistogram(Int_t gridElement, vector<TString> optList) {
691 Bool_t stats =
false;
693 for (
unsigned int n = 0; n < optList.size(); n++) {
694 if (optList[n].Contains(
"binSize=")) {
695 TString pitchStr = optList[n].Remove(0, 8);
696 pitch = stod((
string)pitchStr);
699 if (optList[n].Contains(
"stats")) stats =
true;
702 Int_t nBinsX = (fMaxX - fMinX + 20) / pitch;
703 Int_t nBinsZ = (fMaxZ - fMinZ + 20) / pitch;
705 fXZHisto =
new TH2F(
"XZ",
"", nBinsX, fMinX - 10, fMinX + pitch * nBinsX, nBinsZ, fMinZ - 10,
706 fMinZ + pitch * nBinsZ);
708 for (
unsigned int ntck = 0; ntck < GetNumberOfTracks(); ntck++) {
709 const auto& hits = GetTrack(ntck).GetHits();
711 for (
unsigned int nhit = 0; nhit < hits.GetNumberOfHits(); nhit++) {
712 Double_t xPos = hits.GetX(nhit);
713 Double_t zPos = hits.GetZ(nhit);
714 Double_t energy = hits.GetEnergy(nhit);
716 fXZHisto->Fill(xPos, zPos, energy);
722 style.SetOptStat(1001111);
724 fXZHisto->SetStats(stats);
727 sprintf(title,
"Event ID %d (XZ)", this->GetID());
728 fXZHisto->SetTitle(title);
730 fXZHisto->GetXaxis()->SetTitle(
"X-axis (mm)");
731 fXZHisto->GetXaxis()->SetTitleOffset(1.25);
732 fXZHisto->GetXaxis()->SetTitleSize(1.25 * fXZHisto->GetXaxis()->GetTitleSize());
733 fXZHisto->GetXaxis()->SetLabelSize(1.25 * fXZHisto->GetXaxis()->GetLabelSize());
735 fXZHisto->GetYaxis()->SetTitle(
"Z-axis (mm)");
736 fXZHisto->GetYaxis()->SetTitleOffset(1.75);
737 fXZHisto->GetYaxis()->SetTitleSize(1.25 * fXZHisto->GetYaxis()->GetTitleSize());
738 fXZHisto->GetYaxis()->SetLabelSize(1.25 * fXZHisto->GetYaxis()->GetLabelSize());
739 fPad->cd(gridElement);
744 TH2F* TRestGeant4Event::GetYZHistogram(Int_t gridElement, vector<TString> optList) {
752 for (
unsigned int n = 0; n < optList.size(); n++) {
753 if (optList[n].Contains(
"binSize=")) {
754 TString pitchStr = optList[n].Remove(0, 8);
755 pitch = stod((
string)pitchStr);
758 if (optList[n].Contains(
"stats")) stats =
true;
761 Int_t nBinsY = (fMaxY - fMinY + 20) / pitch;
762 Int_t nBinsZ = (fMaxZ - fMinZ + 20) / pitch;
764 fYZHisto =
new TH2F(
"YZ",
"", nBinsY, fMinY - 10, fMinY + pitch * nBinsY, nBinsZ, fMinZ - 10,
765 fMinZ + pitch * nBinsZ);
767 for (
unsigned int ntck = 0; ntck < GetNumberOfTracks(); ntck++) {
768 const auto& hits = GetTrack(ntck).GetHits();
770 for (
unsigned int nhit = 0; nhit < hits.GetNumberOfHits(); nhit++) {
771 Double_t yPos = hits.GetY(nhit);
772 Double_t zPos = hits.GetZ(nhit);
773 Double_t energy = hits.GetEnergy(nhit);
775 fYZHisto->Fill(yPos, zPos, energy);
781 style.SetOptStat(1001111);
783 fYZHisto->SetStats(stats);
786 sprintf(title,
"Event ID %d (YZ)", this->GetID());
787 fYZHisto->SetTitle(title);
789 fYZHisto->GetXaxis()->SetTitle(
"Y-axis (mm)");
790 fYZHisto->GetXaxis()->SetTitleOffset(1.25);
791 fYZHisto->GetXaxis()->SetTitleSize(1.25 * fYZHisto->GetXaxis()->GetTitleSize());
792 fYZHisto->GetXaxis()->SetLabelSize(1.25 * fYZHisto->GetXaxis()->GetLabelSize());
794 fYZHisto->GetYaxis()->SetTitle(
"Z-axis (mm)");
795 fYZHisto->GetYaxis()->SetTitleOffset(1.75);
796 fYZHisto->GetYaxis()->SetTitleSize(1.25 * fYZHisto->GetYaxis()->GetTitleSize());
797 fYZHisto->GetYaxis()->SetLabelSize(1.25 * fYZHisto->GetYaxis()->GetLabelSize());
798 fPad->cd(gridElement);
804 TH1D* TRestGeant4Event::GetXHistogram(Int_t gridElement, vector<TString> optList) {
811 Bool_t stats =
false;
812 for (
unsigned int n = 0; n < optList.size(); n++) {
813 if (optList[n].Contains(
"binSize=")) {
814 TString pitchStr = optList[n].Remove(0, 8);
815 pitch = stod((
string)pitchStr);
818 if (optList[n].Contains(
"stats")) stats =
true;
821 Int_t nBinsX = (fMaxX - fMinX + 20) / pitch;
823 fXHisto =
new TH1D(
"X",
"", nBinsX, fMinX - 10, fMinX + pitch * nBinsX);
825 for (
unsigned int ntck = 0; ntck < GetNumberOfTracks(); ntck++) {
826 const auto& hits = GetTrack(ntck).GetHits();
828 for (
unsigned int nhit = 0; nhit < hits.GetNumberOfHits(); nhit++) {
829 Double_t xPos = hits.GetX(nhit);
830 Double_t energy = hits.GetEnergy(nhit);
832 fXHisto->Fill(xPos, energy);
839 fXHisto->SetStats(stats);
842 sprintf(title,
"Event ID %d (X)", this->GetID());
843 fXHisto->SetTitle(title);
845 fXHisto->GetXaxis()->SetTitle(
"X-axis (mm)");
846 fXHisto->GetXaxis()->SetTitleOffset(1.25);
847 fXHisto->GetXaxis()->SetTitleSize(1.25 * fXHisto->GetXaxis()->GetTitleSize());
848 fXHisto->GetXaxis()->SetLabelSize(1.25 * fXHisto->GetXaxis()->GetLabelSize());
850 fPad->cd(gridElement);
855 TH1D* TRestGeant4Event::GetZHistogram(Int_t gridElement, vector<TString> optList) {
862 Bool_t stats =
false;
863 for (
unsigned int n = 0; n < optList.size(); n++) {
864 if (optList[n].Contains(
"binSize=")) {
865 TString pitchStr = optList[n].Remove(0, 8);
866 pitch = stod((
string)pitchStr);
869 if (optList[n].Contains(
"stats")) stats =
true;
872 Int_t nBinsZ = (fMaxZ - fMinZ + 20) / pitch;
874 fZHisto =
new TH1D(
"Z",
"", nBinsZ, fMinZ - 10, fMinZ + pitch * nBinsZ);
876 for (
unsigned int ntck = 0; ntck < GetNumberOfTracks(); ntck++) {
877 const auto& hits = GetTrack(ntck).GetHits();
879 for (
unsigned int nhit = 0; nhit < hits.GetNumberOfHits(); nhit++) {
880 Double_t zPos = hits.GetZ(nhit);
881 Double_t energy = hits.GetEnergy(nhit);
883 fZHisto->Fill(zPos, energy);
890 fZHisto->SetStats(stats);
893 sprintf(title,
"Event ID %d (Z)", this->GetID());
894 fZHisto->SetTitle(title);
896 fZHisto->GetXaxis()->SetTitle(
"Z-axis (mm)");
897 fZHisto->GetXaxis()->SetTitleOffset(1.25);
898 fZHisto->GetXaxis()->SetTitleSize(1.25 * fZHisto->GetXaxis()->GetTitleSize());
899 fZHisto->GetXaxis()->SetLabelSize(1.25 * fZHisto->GetXaxis()->GetLabelSize());
901 fPad->cd(gridElement);
906 TH1D* TRestGeant4Event::GetYHistogram(Int_t gridElement, vector<TString> optList) {
913 Bool_t stats =
false;
914 for (
unsigned int n = 0; n < optList.size(); n++) {
915 if (optList[n].Contains(
"binSize=")) {
916 TString pitchStr = optList[n].Remove(0, 8);
917 pitch = stod((
string)pitchStr);
920 if (optList[n].Contains(
"stats")) stats =
true;
923 Int_t nBinsY = (fMaxY - fMinY + 20) / pitch;
925 fYHisto =
new TH1D(
"Y",
"", nBinsY, fMinY - 10, fMinY + pitch * nBinsY);
927 for (
unsigned int ntck = 0; ntck < GetNumberOfTracks(); ntck++) {
928 const auto& hits = GetTrack(ntck).GetHits();
930 for (
unsigned int nhit = 0; nhit < hits.GetNumberOfHits(); nhit++) {
931 Double_t yPos = hits.GetY(nhit);
932 Double_t energy = hits.GetEnergy(nhit);
934 fYHisto->Fill(yPos, energy);
941 fYHisto->SetStats(stats);
944 sprintf(title,
"Event ID %d (Y)", this->GetID());
945 fYHisto->SetTitle(title);
947 fYHisto->GetXaxis()->SetTitle(
"Y-axis (mm)");
948 fYHisto->GetXaxis()->SetTitleOffset(1.25);
949 fYHisto->GetXaxis()->SetTitleSize(1.25 * fYHisto->GetXaxis()->GetTitleSize());
950 fYHisto->GetXaxis()->SetLabelSize(1.25 * fYHisto->GetXaxis()->GetLabelSize());
952 fPad->cd(gridElement);
1023 if (autoBoundaries) SetBoundaries();
1026 if (optList.size() == 0) {
1027 optList.push_back(
"graphXZ");
1028 optList.push_back(
"graphYZ");
1029 optList.push_back(
"histXZ(Cont0,colz)");
1030 optList.push_back(
"histYZ(Cont0,colz)");
1033 for (
unsigned int n = 0; n < optList.size(); n++) {
1034 if (optList[n] ==
"print") this->PrintEvent();
1037 optList.erase(remove(optList.begin(), optList.end(),
"print"), optList.end());
1039 unsigned int nPlots = optList.size();
1043 for (
unsigned int n = 0; n < nPlots; n++) {
1046 string optionStr = (string)optList[n];
1050 TString drawEventType = optList[n];
1051 size_t startPos = optionStr.find(
"(");
1052 if (startPos == string::npos) startPos = optionStr.find(
"[");
1054 if (startPos != string::npos) drawEventType = optList[n](0, startPos);
1059 string drawOption =
"";
1060 size_t endPos = optionStr.find(
")");
1061 if (endPos != string::npos) {
1062 TString drawOptionTmp = optList[n](startPos + 1, endPos - startPos - 1);
1064 drawOption = (string)drawOptionTmp;
1066 while ((pos = drawOption.find(
",", pos)) != string::npos) {
1067 drawOption.replace(pos, 1,
":");
1075 vector<TString> optList_2;
1077 startPos = optionStr.find(
"[");
1078 endPos = optionStr.find(
"]");
1080 if (endPos != string::npos) {
1081 TString tmpStr = optList[n](startPos + 1, endPos - startPos - 1);
1082 optList_2 = Vector_cast<string, TString>(
Split((
string)tmpStr,
","));
1086 if (drawEventType ==
"graphXZ") {
1087 GetXZMultiGraph(n + 1, optList_2)->Draw(
"P");
1091 if (fXZPcsMarker.size() > 0) fLegend_XZ->Draw(
"same");
1093 for (
unsigned int m = 0; m < fXZPcsMarker.size(); m++) fXZPcsMarker[m]->Draw(
"P");
1094 }
else if (drawEventType ==
"graphYZ") {
1095 GetYZMultiGraph(n + 1, optList_2)->Draw(
"P");
1097 if (fYZPcsMarker.size() > 0) fLegend_YZ->Draw(
"same");
1099 for (
unsigned int m = 0; m < fYZPcsMarker.size(); m++) fYZPcsMarker[m]->Draw(
"P");
1100 }
else if (drawEventType ==
"graphXY") {
1101 GetXYMultiGraph(n + 1, optList_2)->Draw(
"P");
1103 if (fXYPcsMarker.size() > 0) fLegend_XY->Draw(
"same");
1105 for (
unsigned int m = 0; m < fXYPcsMarker.size(); m++) fXYPcsMarker[m]->Draw(
"P");
1106 }
else if (drawEventType ==
"histXY") {
1107 GetXYHistogram(n + 1, optList_2)->Draw((TString)drawOption);
1108 }
else if (drawEventType ==
"histXZ") {
1109 GetXZHistogram(n + 1, optList_2)->Draw((TString)drawOption);
1110 }
else if (drawEventType ==
"histYZ") {
1111 GetYZHistogram(n + 1, optList_2)->Draw((TString)drawOption);
1112 }
else if (drawEventType ==
"histX") {
1113 GetXHistogram(n + 1, optList_2)->Draw((TString)drawOption);
1114 }
else if (drawEventType ==
"histY") {
1115 GetYHistogram(n + 1, optList_2)->Draw((TString)drawOption);
1116 }
else if (drawEventType ==
"histZ") {
1117 GetZHistogram(n + 1, optList_2)->Draw((TString)drawOption);
1124 cout <<
"Number of active volumes : " << GetNumberOfActiveVolumes() << endl;
1125 for (
unsigned int i = 0; i < fVolumeStoredNames.size(); i++) {
1126 if (isVolumeStored(i)) {
1127 cout <<
"Active volume " << i <<
":" << fVolumeStoredNames[i] <<
" has been stored." << endl;
1128 cout <<
"Total energy deposit in volume " << i <<
":" << fVolumeStoredNames[i] <<
" : "
1129 << fVolumeDepositedEnergy[i] <<
" keV" << endl;
1131 cout <<
"Active volume " << i <<
":" << fVolumeStoredNames[i] <<
" has not been stored" << endl;
1138 cout <<
"- Total deposited energy: " << ToEnergyString(fTotalDepositedEnergy) << endl;
1139 cout <<
"- Sensitive detectors total energy: " << ToEnergyString(fSensitiveVolumeEnergy) << endl;
1141 cout <<
"- Primary source position: " << VectorToString(fPrimaryPosition) <<
" mm" << endl;
1143 for (
unsigned int i = 0; i < GetNumberOfPrimaries(); i++) {
1144 const char* sourceNumberString =
1145 GetNumberOfPrimaries() <= 1 ?
""
1146 : TString::Format(
" %d of %zu", i + 1, GetNumberOfPrimaries()).Data();
1147 cout <<
" - Source" << sourceNumberString <<
" primary particle: " << fPrimaryParticleNames[i]
1149 cout <<
" Direction: " << VectorToString(fPrimaryDirections[i]) << endl;
1150 cout <<
" Energy: " << ToEnergyString(fPrimaryEnergies[i]) << endl;
1154 cout <<
"Total number of tracks: " << GetNumberOfTracks() << endl;
1156 int nTracks = GetNumberOfTracks();
1157 if (maxTracks > 0 && (
unsigned int)maxTracks < GetNumberOfTracks()) {
1158 nTracks = min(maxTracks,
int(GetNumberOfTracks()));
1159 cout <<
"Printing only the first " << nTracks <<
" tracks" << endl;
1162 for (
int i = 0; i < nTracks; i++) {
1163 GetTrack(i).PrintTrack(maxHits);
1167 Bool_t TRestGeant4Event::ContainsProcessInVolume(Int_t processID, Int_t volumeID)
const {
1168 for (
const auto& track : fTracks) {
1169 if (track.ContainsProcessInVolume(processID, volumeID)) {
1176 Bool_t TRestGeant4Event::ContainsProcessInVolume(
const TString& processName, Int_t volumeID)
const {
1178 if (metadata ==
nullptr) {
1182 for (
const auto& track : fTracks) {
1183 if (track.ContainsProcessInVolume(processID, volumeID)) {
1190 Bool_t TRestGeant4Event::ContainsParticle(
const TString& particleName)
const {
1191 for (
const auto& track : fTracks) {
1192 if (track.GetParticleName() == particleName) {
1199 Bool_t TRestGeant4Event::ContainsParticleInVolume(
const TString& particleName, Int_t volumeID)
const {
1200 for (
const auto& track : fTracks) {
1201 if (track.GetParticleName() != particleName) {
1204 if (track.GetHits().GetNumberOfHitsInVolume(volumeID) > 0) {
1212 if (fRun ==
nullptr) {
1213 RESTError <<
"TRestGeant4Event::GetGeant4Metadata: fRun is nullptr" << RESTendl;
1225 for (
auto& track : fTracks) {
1226 track.SetEvent(
this);
1227 track.fHits.SetTrack(&track);
1228 track.fHits.SetEvent(
this);
1232 set<string> TRestGeant4Event::GetUniqueParticles()
const {
1234 for (
const auto& track : fTracks) {
1235 result.insert(track.GetParticleName().Data());
1240 map<string, map<string, double>> TRestGeant4Event::GetEnergyInVolumePerProcessMap()
const {
1241 map<string, map<string, double>> result;
1242 for (
const auto& [volume, particleProcessMap] : fEnergyInVolumePerParticlePerProcess) {
1243 for (
const auto& [particle, processMap] : particleProcessMap) {
1244 for (
const auto& [process, energy] : processMap) {
1245 result[volume][process] += energy;
1252 map<string, map<string, double>> TRestGeant4Event::GetEnergyInVolumePerParticleMap()
const {
1253 map<string, map<string, double>> result;
1254 for (
const auto& [volume, particleProcessMap] : fEnergyInVolumePerParticlePerProcess) {
1255 for (
const auto& [particle, processMap] : particleProcessMap) {
1256 for (
const auto& [process, energy] : processMap) {
1257 result[volume][particle] += energy;
1264 map<string, double> TRestGeant4Event::GetEnergyPerProcessMap()
const {
1265 map<string, double> result;
1266 for (
const auto& [volume, particleProcessMap] : fEnergyInVolumePerParticlePerProcess) {
1267 for (
const auto& [particle, processMap] : particleProcessMap) {
1268 for (
const auto& [process, energy] : processMap) {
1269 result[process] += energy;
1276 map<string, double> TRestGeant4Event::GetEnergyPerParticleMap()
const {
1277 map<string, double> result;
1278 for (
const auto& [volume, particleProcessMap] : fEnergyInVolumePerParticlePerProcess) {
1279 for (
const auto& [particle, processMap] : particleProcessMap) {
1280 for (
const auto& [process, energy] : processMap) {
1281 result[particle] += energy;
1288 map<string, double> TRestGeant4Event::GetEnergyInVolumeMap()
const {
1289 map<string, double> result;
1290 for (
const auto& [volume, particleProcessMap] : fEnergyInVolumePerParticlePerProcess) {
1291 for (
const auto& [particle, processMap] : particleProcessMap) {
1292 for (
const auto& [process, energy] : processMap) {
1293 result[volume] += energy;
1300 map<string, map<string, map<string, double>>> TRestGeant4Event::GetEnergyInVolumePerParticlePerProcessMap()
1302 return fEnergyInVolumePerParticlePerProcess;
1305 void TRestGeant4Event::AddEnergyInVolumeForParticleForProcess(Double_t energy,
const string& volumeName,
1306 const string& particleName,
1307 const string& processName) {
1311 fEnergyInVolumePerParticlePerProcess[volumeName][particleName][processName] += energy;
1312 fTotalDepositedEnergy += energy;
virtual void InitializeReferences(TRestRun *run)
Initialize dynamical references when loading the event from a root file.
virtual void PrintEvent() const
virtual void Initialize()=0
An event class to store geant4 generated event information.
TVector3 GetFirstPositionInVolume(Int_t volID) const
Function to get the position (TVector3) of the first track that deposits energy in specified volume....
void InitializeReferences(TRestRun *run) override
Initialize dynamical references when loading the event from a root file.
size_t GetNumberOfHits(Int_t volID=-1) const
Function that returns the total number of hits in the Geant4 event. If a specific volume is given as ...
size_t GetNumberOfPhysicalHits(Int_t volID=-1) const
Function that returns the total number of hits with energy > 0 in the Geant4 event....
Double_t GetBoundingBoxSize()
This method returns the event size as the size of the bounding box enclosing all hits.
TVector3 GetPositionDeviationInVolume(Int_t volID) const
A method that gets the standard deviation from the hits happening at a particular volumeId inside the...
void Initialize() override
void PrintActiveVolumes() const
maxTracks : number of tracks to print, 0 = all
TVector3 GetLastPositionInVolume(Int_t volID) const
Function to get the position (TVector3) of the last track that deposits energy in specified volume....
TRestHits GetHits(Int_t volID=-1) const
Function that returns all the hit depositions in the Geant4 event. If a specific volume is given as a...
TPad * DrawEvent(const TString &option="") override
Draw the event.
It saves a 3-coordinate position and an energy for each punctual deposition.
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.
TVector3 GetPosition(int n) const
It returns the position of hit number n.
Data provider and manager in REST.
std::vector< std::string > Split(std::string in, std::string separator, bool allowBlankString=false, bool removeWhiteSpaces=false, int startPos=-1)
Split the input string according to the given separator. Returning a vector of fragments.