REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestGeant4Event.cxx
1 
18 #include "TRestGeant4Event.h"
19 
20 #include <TFrame.h>
21 #include <TRestRun.h>
22 #include <TRestStringHelper.h>
23 #include <TRestTools.h>
24 #include <TStyle.h>
25 
26 #include "TRestGeant4Metadata.h"
27 
28 using namespace std;
29 
30 ClassImp(TRestGeant4Event);
31 
32 TRestGeant4Event::TRestGeant4Event() {
33  fNVolumes = 0;
34  // TRestGeant4Event default constructor
35 
36  Initialize();
37 }
38 
39 TRestGeant4Event::~TRestGeant4Event() {
40  // TRestGeant4Event destructor
41 }
42 
45 
46  fPrimaryParticleNames.clear();
47 
48  fTracks.clear();
49 
50  // ClearVolumes();
51  fXZHitGraph = nullptr;
52  fYZHitGraph = nullptr;
53  fXYHitGraph = nullptr;
54 
55  fXZMultiGraph = nullptr;
56  fYZMultiGraph = nullptr;
57  fXYMultiGraph = nullptr;
58 
59  fXYHisto = nullptr;
60  fYZHisto = nullptr;
61  fXZHisto = nullptr;
62 
63  fXHisto = nullptr;
64  fYHisto = nullptr;
65  fZHisto = nullptr;
66 
67  fPad = nullptr;
68 
69  fLegend_XY = nullptr;
70  fLegend_XZ = nullptr;
71  fLegend_YZ = nullptr;
72 
73  fTotalDepositedEnergy = 0;
74  fSensitiveVolumeEnergy = 0;
75 
76  fMinX = 1e20;
77  fMaxX = -1e20;
78 
79  fMinY = 1e20;
80  fMaxY = -1e20;
81 
82  fMinZ = 1e20;
83  fMaxZ = -1e20;
84 
85  fMinEnergy = 1e20;
86  fMaxEnergy = -1e20;
87 }
88 
89 void TRestGeant4Event::AddActiveVolume(const string& volumeName) {
90  fNVolumes++;
91  fVolumeStored.push_back(1);
92  fVolumeStoredNames.push_back(volumeName);
93  fVolumeDepositedEnergy.push_back(0);
94 }
95 
96 void TRestGeant4Event::ClearVolumes() {
97  fNVolumes = 0;
98  fVolumeStored.clear();
99  fVolumeStoredNames.clear();
100  fVolumeDepositedEnergy.clear();
101 }
102 
103 void TRestGeant4Event::AddEnergyDepositToVolume(Int_t volID, Double_t eDep) {
104  fVolumeDepositedEnergy[volID] += eDep;
105 }
106 
107 TVector3 TRestGeant4Event::GetMeanPositionInVolume(Int_t volID) const {
108  TVector3 pos;
109  Double_t eDep = 0;
110 
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);
115 
116  eDep += tck.GetEnergyInVolume(volID);
117  }
118  }
119 
120  if (eDep == 0) {
121  Double_t nan = TMath::QuietNaN();
122  return {nan, nan, nan};
123  }
124 
125  pos = (1 / eDep) * pos;
126  return pos;
127 }
128 
136  TVector3 meanPos = this->GetMeanPositionInVolume(volID);
137 
138  Double_t nan = TMath::QuietNaN();
139  if (meanPos == TVector3(nan, nan, nan)) return TVector3(nan, nan, nan);
140 
141  TRestHits hitsInVolume = GetHitsInVolume(volID);
142 
143  Double_t edep = 0;
144  TVector3 deviation = TVector3(0, 0, 0);
145 
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());
150 
151  deviation = deviation + en * diff;
152 
153  edep += en;
154  }
155  deviation = (1. / edep) * deviation;
156  return deviation;
157 }
158 
166 TVector3 TRestGeant4Event::GetFirstPositionInVolume(Int_t volID) const {
167  for (unsigned int t = 0; t < GetNumberOfTracks(); t++)
168  if (GetTrack(t).GetEnergyInVolume(volID) > 0) return GetTrack(t).GetFirstPositionInVolume(volID);
169 
170  Double_t nan = TMath::QuietNaN();
171  return {nan, nan, nan};
172 }
173 
183 TVector3 TRestGeant4Event::GetLastPositionInVolume(Int_t volID) const {
184  for (int t = GetNumberOfTracks() - 1; t >= 0; t--)
185  if (GetTrack(t).GetEnergyInVolume(volID) > 0) return GetTrack(t).GetLastPositionInVolume(volID);
186 
187  Double_t nan = TMath::QuietNaN();
188  return {nan, nan, nan};
189 }
190 
191 TRestGeant4Track* TRestGeant4Event::GetTrackByID(Int_t trackID) const {
192  TRestGeant4Track* result = nullptr;
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;
197  exit(1);
198  }
199  }
200  return result;
201 }
202 
208 size_t TRestGeant4Event::GetNumberOfHits(Int_t volID) const {
209  size_t numberOfHits = 0;
210  for (const auto& track : fTracks) {
211  numberOfHits += track.GetNumberOfHits(volID);
212  }
213  return numberOfHits;
214 }
215 
221 size_t TRestGeant4Event::GetNumberOfPhysicalHits(Int_t volID) const {
222  size_t numberOfHits = 0;
223  for (const auto& track : fTracks) {
224  numberOfHits += track.GetNumberOfPhysicalHits(volID);
225  }
226  return numberOfHits;
227 }
228 
235  TRestHits hits;
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;
240 
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);
245 
246  hits.AddHit({x, y, z}, en);
247  }
248  }
249 
250  return hits;
251 }
252 
253 Int_t TRestGeant4Event::GetNumberOfTracksForParticle(const TString& particleName) const {
254  Int_t nTracks = 0;
255  for (const auto& track : fTracks) {
256  if (particleName.EqualTo(track.GetParticleName())) {
257  nTracks += 1;
258  }
259  }
260  return nTracks;
261 }
262 
263 Double_t TRestGeant4Event::GetEnergyDepositedByParticle(const TString& particleName) const {
264  Double_t energy = 0;
265  for (const auto& track : fTracks) {
266  if (particleName.EqualTo(track.GetParticleName())) {
267  energy += track.GetTotalEnergy();
268  }
269  }
270  return energy;
271 }
272 
273 void TRestGeant4Event::SetBoundaries(Double_t xMin, Double_t xMax, Double_t yMin, Double_t yMax,
274  Double_t zMin, Double_t zMax) {
275  SetBoundaries();
276 
277  fMinX = xMin;
278  fMaxX = xMax;
279 
280  fMinY = yMin;
281  fMaxY = yMax;
282 
283  fMinZ = zMin;
284  fMaxZ = zMax;
285 }
286 
292  SetBoundaries();
293 
294  Double_t dX = fMaxX - fMinX;
295  Double_t dY = fMaxY - fMinY;
296  Double_t dZ = fMaxZ - fMinZ;
297 
298  return TMath::Sqrt(dX * dX + dY * dY + dZ * dZ);
299 }
300 
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;
304 
305  Int_t nTHits = 0;
306  for (unsigned int ntck = 0; ntck < this->GetNumberOfTracks(); ntck++) {
307  Int_t nHits = GetTrack(ntck).GetNumberOfHits();
308  nTHits += nHits;
309  const auto& hits = GetTrack(ntck).GetHits();
310 
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);
315 
316  Double_t en = hits.GetEnergy(nhit);
317 
318  if (en <= 0) continue;
319 
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;
326 
327  if (en > maxEnergy) maxEnergy = en;
328  if (en < minEnergy) minEnergy = en;
329  }
330  }
331 
332  fMinX = minX;
333  fMaxX = maxX;
334 
335  fMinY = minY;
336  fMaxY = maxY;
337 
338  fMinZ = minZ;
339  fMaxZ = maxZ;
340 
341  fMaxEnergy = maxEnergy;
342  fMinEnergy = minEnergy;
343 
344  fTotalHits = nTHits;
345 }
346 
347 /* {{{ Get{XY,YZ,XZ}MultiGraph methods */
348 TMultiGraph* TRestGeant4Event::GetXYMultiGraph(Int_t gridElement, vector<TString> pcsList,
349  Double_t minPointSize, Double_t maxPointSize) {
350  if (fXYHitGraph) {
351  delete[] fXYHitGraph;
352  fXYHitGraph = nullptr;
353  }
354  fXYHitGraph = new TGraph[fTotalHits];
355 
356  fXYMultiGraph = new TMultiGraph();
357 
358  if (fLegend_XY) {
359  delete fLegend_XY;
360  fLegend_XY = nullptr;
361  }
362 
363  fLegend_XY = new TLegend(0.75, 0.75, 0.9, 0.85);
364 
365  char title[256];
366  sprintf(title, "Event ID %d (XY)", this->GetID());
367  fPad->cd(gridElement)->DrawFrame(fMinX - 10, fMinY - 10, fMaxX + 10, fMaxY + 10, title);
368 
369  Double_t m = 1, n = 0;
370  if (fMaxEnergy - fMinEnergy > 0) {
371  m = (maxPointSize - minPointSize) / (fMaxEnergy - fMinEnergy);
372  n = maxPointSize - m * fMaxEnergy;
373  }
374 
375  for (unsigned int n = 0; n < fXYPcsMarker.size(); n++) delete fXYPcsMarker[n];
376  fXYPcsMarker.clear();
377 
378  legendAdded.clear();
379  for (unsigned int p = 0; p < pcsList.size(); p++) legendAdded.push_back(0);
380 
381  Int_t cnt = 0;
382  for (unsigned int ntck = 0; ntck < GetNumberOfTracks(); ntck++) {
383  const auto hits = GetTrack(ntck).GetHits();
384 
385  EColor color = GetTrack(ntck).GetParticleColor();
386 
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);
391 
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);
395 
396  pcsGraph->SetPoint(0, xPos, yPos);
397  pcsGraph->SetMarkerColor(kBlack);
398  pcsGraph->SetMarkerSize(3);
399  pcsGraph->SetMarkerStyle(30 + p);
400 
401  fXYPcsMarker.push_back(pcsGraph);
402 
403  if (legendAdded[p] == 0) {
404  fLegend_XY->AddEntry(pcsGraph, GetTrack(ntck).GetProcessName(hits.GetProcess(nhit)),
405  "p");
406  legendAdded[p] = 1;
407  }
408  }
409  }
410 
411  Double_t radius = m * energy + n;
412  if (fMaxEnergy - fMinEnergy <= 0) radius = 2.5;
413 
414  fXYHitGraph[cnt].SetPoint(0, xPos, yPos);
415 
416  fXYHitGraph[cnt].SetMarkerColor(color);
417  fXYHitGraph[cnt].SetMarkerSize(radius);
418  fXYHitGraph[cnt].SetMarkerStyle(20);
419 
420  fXYMultiGraph->Add(&fXYHitGraph[cnt]);
421 
422  cnt++;
423  }
424  }
425 
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());
430 
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);
436 
437  return fXYMultiGraph;
438 }
439 
440 TMultiGraph* TRestGeant4Event::GetYZMultiGraph(Int_t gridElement, vector<TString> pcsList,
441  Double_t minPointSize, Double_t maxPointSize) {
442  if (fYZHitGraph) {
443  delete[] fYZHitGraph;
444  fYZHitGraph = nullptr;
445  }
446  fYZHitGraph = new TGraph[fTotalHits];
447 
448  fYZMultiGraph = new TMultiGraph();
449 
450  if (fLegend_YZ) {
451  delete fLegend_YZ;
452  fLegend_YZ = nullptr;
453  }
454 
455  fLegend_YZ = new TLegend(0.75, 0.75, 0.9, 0.85);
456 
457  char title[256];
458  sprintf(title, "Event ID %d (YZ)", this->GetID());
459  fPad->cd(gridElement)->DrawFrame(fMinY - 10, fMinZ - 10, fMaxY + 10, fMaxZ + 10, title);
460 
461  for (unsigned int n = 0; n < fYZPcsMarker.size(); n++) delete fYZPcsMarker[n];
462  fYZPcsMarker.clear();
463 
464  Double_t m = 1, n = 0;
465  if (fMaxEnergy - fMinEnergy > 0) {
466  m = (maxPointSize - minPointSize) / (fMaxEnergy - fMinEnergy);
467  n = maxPointSize - m * fMaxEnergy;
468  }
469 
470  legendAdded.clear();
471  for (unsigned int p = 0; p < pcsList.size(); p++) legendAdded.push_back(0);
472 
473  Int_t cnt = 0;
474  for (unsigned int ntck = 0; ntck < GetNumberOfTracks(); ntck++) {
475  const auto& hits = GetTrack(ntck).GetHits();
476 
477  EColor color = GetTrack(ntck).GetParticleColor();
478 
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);
483 
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);
487 
488  pcsGraph->SetPoint(0, yPos, zPos);
489  pcsGraph->SetMarkerColor(kBlack);
490  pcsGraph->SetMarkerSize(3);
491  pcsGraph->SetMarkerStyle(30 + p);
492 
493  fYZPcsMarker.push_back(pcsGraph);
494 
495  if (legendAdded[p] == 0) {
496  fLegend_YZ->AddEntry(pcsGraph, GetTrack(ntck).GetProcessName(hits.GetProcess(nhit)),
497  "p");
498  legendAdded[p] = 1;
499  }
500  }
501  }
502 
503  Double_t radius = m * energy + n;
504  if (fMaxEnergy - fMinEnergy <= 0) radius = 2.5;
505 
506  fYZHitGraph[cnt].SetPoint(0, yPos, zPos);
507 
508  fYZHitGraph[cnt].SetMarkerColor(color);
509  fYZHitGraph[cnt].SetMarkerSize(radius);
510  fYZHitGraph[cnt].SetMarkerStyle(20);
511 
512  fYZMultiGraph->Add(&fYZHitGraph[cnt]);
513 
514  cnt++;
515  }
516  }
517 
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());
522 
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);
528 
529  return fYZMultiGraph;
530 }
531 
532 TMultiGraph* TRestGeant4Event::GetXZMultiGraph(Int_t gridElement, vector<TString> pcsList,
533  Double_t minPointSize, Double_t maxPointSize) {
534  if (fXZHitGraph) {
535  delete[] fXZHitGraph;
536  fXZHitGraph = nullptr;
537  }
538  fXZHitGraph = new TGraph[fTotalHits];
539 
540  fXZMultiGraph = new TMultiGraph();
541 
542  if (fLegend_XZ) {
543  delete fLegend_XZ;
544  fLegend_XZ = nullptr;
545  }
546 
547  fLegend_XZ = new TLegend(0.75, 0.75, 0.9, 0.85);
548 
549  char title[256];
550  sprintf(title, "Event ID %d (XZ)", this->GetID());
551  fPad->cd(gridElement)->DrawFrame(fMinX - 10, fMinZ - 10, fMaxX + 10, fMaxZ + 10, title);
552 
553  for (unsigned int n = 0; n < fXZPcsMarker.size(); n++) delete fXZPcsMarker[n];
554  fXZPcsMarker.clear();
555 
556  Double_t m = 1, n = 0;
557  if (fMaxEnergy - fMinEnergy > 0) {
558  m = (maxPointSize - minPointSize) / (fMaxEnergy - fMinEnergy);
559  n = maxPointSize - m * fMaxEnergy;
560  }
561 
562  legendAdded.clear();
563  for (unsigned int p = 0; p < pcsList.size(); p++) legendAdded.push_back(0);
564 
565  Int_t cnt = 0;
566  for (unsigned int ntck = 0; ntck < GetNumberOfTracks(); ntck++) {
567  const auto& hits = GetTrack(ntck).GetHits();
568 
569  EColor color = GetTrack(ntck).GetParticleColor();
570 
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);
575 
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);
579 
580  pcsGraph->SetPoint(0, xPos, zPos);
581  pcsGraph->SetMarkerColor(kBlack);
582  pcsGraph->SetMarkerSize(3);
583  pcsGraph->SetMarkerStyle(30 + p);
584 
585  fXZPcsMarker.push_back(pcsGraph);
586 
587  if (legendAdded[p] == 0) {
588  fLegend_XZ->AddEntry(pcsGraph, GetTrack(ntck).GetProcessName(hits.GetProcess(nhit)),
589  "p");
590  legendAdded[p] = 1;
591  }
592  }
593  }
594 
595  Double_t radius = m * energy + n;
596  if (fMaxEnergy - fMinEnergy <= 0) radius = 2.5;
597 
598  fXZHitGraph[cnt].SetPoint(0, xPos, zPos);
599 
600  fXZHitGraph[cnt].SetMarkerColor(color);
601  fXZHitGraph[cnt].SetMarkerSize(radius);
602  fXZHitGraph[cnt].SetMarkerStyle(20);
603 
604  fXZMultiGraph->Add(&fXZHitGraph[cnt]);
605 
606  cnt++;
607  }
608  }
609 
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());
614 
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);
620 
621  return fXZMultiGraph;
622 }
623 /* }}} */
624 
625 /* {{{ Get{XY,YZ,XZ}Histogram methods */
626 TH2F* TRestGeant4Event::GetXYHistogram(Int_t gridElement, vector<TString> optList) {
627  if (fXYHisto) {
628  delete fXYHisto;
629  fXYHisto = nullptr;
630  }
631 
632  Bool_t stats = false;
633  Double_t pitch = 3;
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);
638  }
639 
640  if (optList[n].Contains("stats")) stats = true;
641  }
642 
643  Int_t nBinsX = (fMaxX - fMinX + 20) / pitch;
644  Int_t nBinsY = (fMaxY - fMinY + 20) / pitch;
645 
646  fXYHisto = new TH2F("XY", "", nBinsX, fMinX - 10, fMinX + pitch * nBinsX, nBinsY, fMinY - 10,
647  fMinY + pitch * nBinsY);
648 
649  for (unsigned int ntck = 0; ntck < GetNumberOfTracks(); ntck++) {
650  const auto& hits = GetTrack(ntck).GetHits();
651 
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);
656 
657  fXYHisto->Fill(xPos, yPos, energy);
658  }
659  }
660 
661  TStyle style;
662  style.SetPalette(1);
663  style.SetOptStat(1001111);
664 
665  fXYHisto->SetStats(stats);
666 
667  char title[256];
668  sprintf(title, "Event ID %d (XY)", this->GetID());
669  fXYHisto->SetTitle(title);
670 
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());
675 
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);
681 
682  return fXYHisto;
683 }
684 
685 TH2F* TRestGeant4Event::GetXZHistogram(Int_t gridElement, vector<TString> optList) {
686  if (fXZHisto) {
687  delete fXZHisto;
688  fXZHisto = nullptr;
689  }
690 
691  Bool_t stats = false;
692  Double_t pitch = 3;
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);
697  }
698 
699  if (optList[n].Contains("stats")) stats = true;
700  }
701 
702  Int_t nBinsX = (fMaxX - fMinX + 20) / pitch;
703  Int_t nBinsZ = (fMaxZ - fMinZ + 20) / pitch;
704 
705  fXZHisto = new TH2F("XZ", "", nBinsX, fMinX - 10, fMinX + pitch * nBinsX, nBinsZ, fMinZ - 10,
706  fMinZ + pitch * nBinsZ);
707 
708  for (unsigned int ntck = 0; ntck < GetNumberOfTracks(); ntck++) {
709  const auto& hits = GetTrack(ntck).GetHits();
710 
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);
715 
716  fXZHisto->Fill(xPos, zPos, energy);
717  }
718  }
719 
720  TStyle style;
721  style.SetPalette(1);
722  style.SetOptStat(1001111);
723 
724  fXZHisto->SetStats(stats);
725 
726  char title[256];
727  sprintf(title, "Event ID %d (XZ)", this->GetID());
728  fXZHisto->SetTitle(title);
729 
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());
734 
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);
740 
741  return fXZHisto;
742 }
743 
744 TH2F* TRestGeant4Event::GetYZHistogram(Int_t gridElement, vector<TString> optList) {
745  if (fYZHisto) {
746  delete fYZHisto;
747  fYZHisto = nullptr;
748  }
749 
750  Bool_t stats;
751  Double_t pitch = 3;
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);
756  }
757 
758  if (optList[n].Contains("stats")) stats = true;
759  }
760 
761  Int_t nBinsY = (fMaxY - fMinY + 20) / pitch;
762  Int_t nBinsZ = (fMaxZ - fMinZ + 20) / pitch;
763 
764  fYZHisto = new TH2F("YZ", "", nBinsY, fMinY - 10, fMinY + pitch * nBinsY, nBinsZ, fMinZ - 10,
765  fMinZ + pitch * nBinsZ);
766 
767  for (unsigned int ntck = 0; ntck < GetNumberOfTracks(); ntck++) {
768  const auto& hits = GetTrack(ntck).GetHits();
769 
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);
774 
775  fYZHisto->Fill(yPos, zPos, energy);
776  }
777  }
778 
779  TStyle style;
780  style.SetPalette(1);
781  style.SetOptStat(1001111); // Not Working :(
782 
783  fYZHisto->SetStats(stats);
784 
785  char title[256];
786  sprintf(title, "Event ID %d (YZ)", this->GetID());
787  fYZHisto->SetTitle(title);
788 
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());
793 
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);
799 
800  return fYZHisto;
801 }
802 /* }}} */
803 
804 TH1D* TRestGeant4Event::GetXHistogram(Int_t gridElement, vector<TString> optList) {
805  if (fXHisto) {
806  delete fXHisto;
807  fXHisto = nullptr;
808  }
809 
810  Double_t pitch = 3;
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);
816  }
817 
818  if (optList[n].Contains("stats")) stats = true;
819  }
820 
821  Int_t nBinsX = (fMaxX - fMinX + 20) / pitch;
822 
823  fXHisto = new TH1D("X", "", nBinsX, fMinX - 10, fMinX + pitch * nBinsX);
824 
825  for (unsigned int ntck = 0; ntck < GetNumberOfTracks(); ntck++) {
826  const auto& hits = GetTrack(ntck).GetHits();
827 
828  for (unsigned int nhit = 0; nhit < hits.GetNumberOfHits(); nhit++) {
829  Double_t xPos = hits.GetX(nhit);
830  Double_t energy = hits.GetEnergy(nhit);
831 
832  fXHisto->Fill(xPos, energy);
833  }
834  }
835 
836  TStyle style;
837  style.SetPalette(1);
838 
839  fXHisto->SetStats(stats);
840 
841  char title[256];
842  sprintf(title, "Event ID %d (X)", this->GetID());
843  fXHisto->SetTitle(title);
844 
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());
849 
850  fPad->cd(gridElement);
851 
852  return fXHisto;
853 }
854 
855 TH1D* TRestGeant4Event::GetZHistogram(Int_t gridElement, vector<TString> optList) {
856  if (fZHisto) {
857  delete fZHisto;
858  fZHisto = nullptr;
859  }
860 
861  Double_t pitch = 3;
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);
867  }
868 
869  if (optList[n].Contains("stats")) stats = true;
870  }
871 
872  Int_t nBinsZ = (fMaxZ - fMinZ + 20) / pitch;
873 
874  fZHisto = new TH1D("Z", "", nBinsZ, fMinZ - 10, fMinZ + pitch * nBinsZ);
875 
876  for (unsigned int ntck = 0; ntck < GetNumberOfTracks(); ntck++) {
877  const auto& hits = GetTrack(ntck).GetHits();
878 
879  for (unsigned int nhit = 0; nhit < hits.GetNumberOfHits(); nhit++) {
880  Double_t zPos = hits.GetZ(nhit);
881  Double_t energy = hits.GetEnergy(nhit);
882 
883  fZHisto->Fill(zPos, energy);
884  }
885  }
886 
887  TStyle style;
888  style.SetPalette(1);
889 
890  fZHisto->SetStats(stats);
891 
892  char title[256];
893  sprintf(title, "Event ID %d (Z)", this->GetID());
894  fZHisto->SetTitle(title);
895 
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());
900 
901  fPad->cd(gridElement);
902 
903  return fZHisto;
904 }
905 
906 TH1D* TRestGeant4Event::GetYHistogram(Int_t gridElement, vector<TString> optList) {
907  if (fYHisto) {
908  delete fYHisto;
909  fYHisto = nullptr;
910  }
911 
912  Double_t pitch = 3;
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);
918  }
919 
920  if (optList[n].Contains("stats")) stats = true;
921  }
922 
923  Int_t nBinsY = (fMaxY - fMinY + 20) / pitch;
924 
925  fYHisto = new TH1D("Y", "", nBinsY, fMinY - 10, fMinY + pitch * nBinsY);
926 
927  for (unsigned int ntck = 0; ntck < GetNumberOfTracks(); ntck++) {
928  const auto& hits = GetTrack(ntck).GetHits();
929 
930  for (unsigned int nhit = 0; nhit < hits.GetNumberOfHits(); nhit++) {
931  Double_t yPos = hits.GetY(nhit);
932  Double_t energy = hits.GetEnergy(nhit);
933 
934  fYHisto->Fill(yPos, energy);
935  }
936  }
937 
938  TStyle style;
939  style.SetPalette(1);
940 
941  fYHisto->SetStats(stats);
942 
943  char title[256];
944  sprintf(title, "Event ID %d (Y)", this->GetID());
945  fYHisto->SetTitle(title);
946 
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());
951 
952  fPad->cd(gridElement);
953 
954  return fYHisto;
955 }
956 
1020 TPad* TRestGeant4Event::DrawEvent(const TString& option, Bool_t autoBoundaries) {
1021  vector<TString> optList = Vector_cast<string, TString>(TRestTools::GetOptions((string)option));
1022 
1023  if (autoBoundaries) SetBoundaries();
1024 
1025  // If no option is given. This is the default
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)");
1031  }
1032 
1033  for (unsigned int n = 0; n < optList.size(); n++) {
1034  if (optList[n] == "print") this->PrintEvent();
1035  }
1036 
1037  optList.erase(remove(optList.begin(), optList.end(), "print"), optList.end());
1038 
1039  unsigned int nPlots = optList.size();
1040 
1041  RestartPad(nPlots);
1042 
1043  for (unsigned int n = 0; n < nPlots; n++) {
1044  fPad->cd(n + 1);
1045 
1046  string optionStr = (string)optList[n];
1047 
1048  /* {{{ Retreiving drawEventType argument. The word key before we find ( or [
1049  * character. */
1050  TString drawEventType = optList[n];
1051  size_t startPos = optionStr.find("(");
1052  if (startPos == string::npos) startPos = optionStr.find("[");
1053 
1054  if (startPos != string::npos) drawEventType = optList[n](0, startPos);
1055  /* }}} */
1056 
1057  /* {{{ Retrieving drawOption argument. Inside normal brackets ().
1058  * We do not separate the different fields we take the full string. */
1059  string drawOption = "";
1060  size_t endPos = optionStr.find(")");
1061  if (endPos != string::npos) {
1062  TString drawOptionTmp = optList[n](startPos + 1, endPos - startPos - 1);
1063 
1064  drawOption = (string)drawOptionTmp;
1065  size_t pos = 0;
1066  while ((pos = drawOption.find(",", pos)) != string::npos) {
1067  drawOption.replace(pos, 1, ":");
1068  pos = pos + 1;
1069  }
1070  }
1071  /* }}} */
1072 
1073  /* {{{ Retrieving optList. Inside squared brackets [] and separated by colon
1074  * [Field1,Field2] */
1075  vector<TString> optList_2;
1076 
1077  startPos = optionStr.find("[");
1078  endPos = optionStr.find("]");
1079 
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, ","));
1083  }
1084  /* }}} */
1085 
1086  if (drawEventType == "graphXZ") {
1087  GetXZMultiGraph(n + 1, optList_2)->Draw("P");
1088 
1089  // Markers for physical processes have been already assigned in
1090  // GetXYMultigraph method
1091  if (fXZPcsMarker.size() > 0) fLegend_XZ->Draw("same");
1092 
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");
1096 
1097  if (fYZPcsMarker.size() > 0) fLegend_YZ->Draw("same");
1098 
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");
1102 
1103  if (fXYPcsMarker.size() > 0) fLegend_XY->Draw("same");
1104 
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);
1118  }
1119  }
1120 
1121  return fPad;
1122 }
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;
1130  } else
1131  cout << "Active volume " << i << ":" << fVolumeStoredNames[i] << " has not been stored" << endl;
1132  }
1133 }
1134 
1135 void TRestGeant4Event::PrintEvent(int maxTracks, int maxHits) const {
1137 
1138  cout << "- Total deposited energy: " << ToEnergyString(fTotalDepositedEnergy) << endl;
1139  cout << "- Sensitive detectors total energy: " << ToEnergyString(fSensitiveVolumeEnergy) << endl;
1140 
1141  cout << "- Primary source position: " << VectorToString(fPrimaryPosition) << " mm" << endl;
1142 
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]
1148  << endl;
1149  cout << " Direction: " << VectorToString(fPrimaryDirections[i]) << endl;
1150  cout << " Energy: " << ToEnergyString(fPrimaryEnergies[i]) << endl;
1151  }
1152 
1153  cout << endl;
1154  cout << "Total number of tracks: " << GetNumberOfTracks() << endl;
1155 
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;
1160  }
1161 
1162  for (int i = 0; i < nTracks; i++) {
1163  GetTrack(i).PrintTrack(maxHits);
1164  }
1165 }
1166 
1167 Bool_t TRestGeant4Event::ContainsProcessInVolume(Int_t processID, Int_t volumeID) const {
1168  for (const auto& track : fTracks) {
1169  if (track.ContainsProcessInVolume(processID, volumeID)) {
1170  return true;
1171  }
1172  }
1173  return false;
1174 }
1175 
1176 Bool_t TRestGeant4Event::ContainsProcessInVolume(const TString& processName, Int_t volumeID) const {
1177  const TRestGeant4Metadata* metadata = GetGeant4Metadata();
1178  if (metadata == nullptr) {
1179  return false;
1180  }
1181  const auto& processID = metadata->GetGeant4PhysicsInfo().GetProcessID(processName);
1182  for (const auto& track : fTracks) {
1183  if (track.ContainsProcessInVolume(processID, volumeID)) {
1184  return true;
1185  }
1186  }
1187  return false;
1188 }
1189 
1190 Bool_t TRestGeant4Event::ContainsParticle(const TString& particleName) const {
1191  for (const auto& track : fTracks) {
1192  if (track.GetParticleName() == particleName) {
1193  return true;
1194  }
1195  }
1196  return false;
1197 }
1198 
1199 Bool_t TRestGeant4Event::ContainsParticleInVolume(const TString& particleName, Int_t volumeID) const {
1200  for (const auto& track : fTracks) {
1201  if (track.GetParticleName() != particleName) {
1202  continue;
1203  }
1204  if (track.GetHits().GetNumberOfHitsInVolume(volumeID) > 0) {
1205  return true;
1206  }
1207  }
1208  return false;
1209 }
1210 
1211 const TRestGeant4Metadata* TRestGeant4Event::GetGeant4Metadata() const {
1212  if (fRun == nullptr) {
1213  RESTError << "TRestGeant4Event::GetGeant4Metadata: fRun is nullptr" << RESTendl;
1214  return nullptr;
1215  }
1216  return dynamic_cast<TRestGeant4Metadata*>(fRun->GetMetadataClass("TRestGeant4Metadata"));
1217 }
1218 
1221  /*
1222  This introduces overhead to event loading, but hopefully its small enough.
1223  If this is a problem, we could rework this approach
1224  */
1225  for (auto& track : fTracks) {
1226  track.SetEvent(this);
1227  track.fHits.SetTrack(&track);
1228  track.fHits.SetEvent(this);
1229  }
1230 }
1231 
1232 set<string> TRestGeant4Event::GetUniqueParticles() const {
1233  set<string> result;
1234  for (const auto& track : fTracks) {
1235  result.insert(track.GetParticleName().Data());
1236  }
1237  return result;
1238 }
1239 
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;
1246  }
1247  }
1248  }
1249  return result;
1250 }
1251 
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;
1258  }
1259  }
1260  }
1261  return result;
1262 }
1263 
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;
1270  }
1271  }
1272  }
1273  return result;
1274 }
1275 
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;
1282  }
1283  }
1284  }
1285  return result;
1286 }
1287 
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;
1294  }
1295  }
1296  }
1297  return result;
1298 }
1299 
1300 map<string, map<string, map<string, double>>> TRestGeant4Event::GetEnergyInVolumePerParticlePerProcessMap()
1301  const {
1302  return fEnergyInVolumePerParticlePerProcess;
1303 }
1304 
1305 void TRestGeant4Event::AddEnergyInVolumeForParticleForProcess(Double_t energy, const string& volumeName,
1306  const string& particleName,
1307  const string& processName) {
1308  if (energy <= 0) {
1309  return;
1310  }
1311  fEnergyInVolumePerParticlePerProcess[volumeName][particleName][processName] += energy;
1312  fTotalDepositedEnergy += energy;
1313 }
virtual void InitializeReferences(TRestRun *run)
Initialize dynamical references when loading the event from a root file.
Definition: TRestEvent.cxx:175
virtual void PrintEvent() const
Definition: TRestEvent.cxx:187
virtual void Initialize()=0
Definition: TRestEvent.cxx:73
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.
The main class to store the Geant4 simulation conditions that will be used by restG4.
const TRestGeant4PhysicsInfo & GetGeant4PhysicsInfo() const
Returns an immutable reference to the physics info.
It saves a 3-coordinate position and an energy for each punctual deposition.
Definition: TRestHits.h:39
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.
Definition: TRestHits.cxx:345
TVector3 GetPosition(int n) const
It returns the position of hit number n.
Definition: TRestHits.cxx:515
Data provider and manager in REST.
Definition: TRestRun.h:18
static std::vector< std::string > GetOptions(std::string optionsStr)
Returns all the options in an option string.
Definition: TRestTools.cxx:86
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.