REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestGeant4AnalysisProcess.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 
194 
195 #include "TRestGeant4AnalysisProcess.h"
196 
197 using namespace std;
198 
200 
205 
219  Initialize();
220 
221  if (LoadConfigFromFile(configFilename)) LoadDefaultConfig();
222 }
223 
228 
232 void TRestGeant4AnalysisProcess::LoadDefaultConfig() { SetTitle("Default config"); }
233 
239  fG4Metadata = nullptr;
240  SetSectionName(this->ClassName());
241  SetLibraryVersion(LIBRARY_VERSION);
242 
243  fInputG4Event = nullptr;
244  fOutputG4Event = new TRestGeant4Event();
245 }
246 
259 void TRestGeant4AnalysisProcess::LoadConfig(const string& configFilename, const string& name) {
260  if (LoadConfigFromFile(configFilename, name)) LoadDefaultConfig();
261 }
262 
269  fG4Metadata = GetMetadata<TRestGeant4Metadata>();
270 
271  std::vector<string> fObservables;
272  fObservables = TRestEventProcess::ReadObservables();
273 
274  if (fPerProcessSensitiveEnergy) {
275  fObservables.emplace_back("PerProcessPhotoelectric");
276  fObservables.emplace_back("PerProcessCompton");
277  fObservables.emplace_back("PerProcessElectronicIoni");
278  fObservables.emplace_back("PerProcessAlphaIoni");
279  fObservables.emplace_back("PerProcessIonIoni");
280  fObservables.emplace_back("PerProcessHadronicIoni");
281  fObservables.emplace_back("PerProcessProtonIoni");
282  fObservables.emplace_back("PerProcessMsc");
283  fObservables.emplace_back("PerProcessHadronElastic");
284  fObservables.emplace_back("PerProcessNeutronElastic");
285  }
286  for (unsigned int i = 0; i < fObservables.size(); i++) {
287  if (fObservables[i].find("VolumeEDep") != string::npos) {
288  TString volName = fObservables[i].substr(0, fObservables[i].length() - 10).c_str();
289 
290  Int_t volId = fG4Metadata->GetActiveVolumeID(volName);
291  if (volId >= 0) {
292  fEnergyInObservables.push_back(fObservables[i]);
293  fVolumeID.push_back(volId);
294  fVolumeName.push_back(volName.Data());
295  }
296 
297  if (volId == -1) {
298  cout << endl;
299  cout << "??????????????????????????????????????????????????" << endl;
300  cout << "REST warning : TRestGeant4AnalysisProcess." << endl;
301  cout << "------------------------------------------" << endl;
302  cout << endl;
303  cout << " Volume " << volName << " is not an active volume" << endl;
304  cout << endl;
305  cout << "List of active volumes : " << endl;
306  cout << "------------------------ " << endl;
307 
308  for (unsigned int n = 0; n < fG4Metadata->GetNumberOfActiveVolumes(); n++)
309  cout << "Volume " << n << " : " << fG4Metadata->GetActiveVolumeName(n) << endl;
310  cout << "??????????????????????????????????????????????????" << endl;
311  cout << endl;
312  }
313  }
314 
315  if (fObservables[i].find("MeanPos") != string::npos) {
316  TString volName2 = fObservables[i].substr(0, fObservables[i].length() - 8).c_str();
317  std::string dirId = fObservables[i].substr(fObservables[i].length() - 1, 1);
318 
319  Int_t volId2 = fG4Metadata->GetActiveVolumeID(volName2);
320  if (volId2 >= 0) {
321  fMeanPosObservables.push_back(fObservables[i]);
322  fVolumeID2.push_back(volId2);
323  fDirID.push_back(dirId);
324  }
325 
326  if (volId2 == -1) {
327  cout << endl;
328  cout << "??????????????????????????????????????????????????" << endl;
329  cout << "REST warning : TRestGeant4AnalysisProcess." << endl;
330  cout << "------------------------------------------" << endl;
331  cout << endl;
332  cout << " Volume " << volName2 << " is not an active volume" << endl;
333  cout << endl;
334  cout << "List of active volumes : " << endl;
335  cout << "------------------------ " << endl;
336 
337  for (unsigned int n = 0; n < fG4Metadata->GetNumberOfActiveVolumes(); n++)
338  cout << "Volume " << n << " : " << fG4Metadata->GetActiveVolumeName(n) << endl;
339  cout << "??????????????????????????????????????????????????" << endl;
340  cout << endl;
341  }
342 
343  if ((dirId != "X") && (dirId != "Y") && (dirId != "Z")) {
344  cout << endl;
345  cout << "??????????????????????????????????????????????????" << endl;
346  cout << "REST warning : TRestGeant4AnalysisProcess." << endl;
347  cout << "------------------------------------------" << endl;
348  cout << endl;
349  cout << " Direction " << dirId << " is not valid" << endl;
350  cout << " Only X, Y or Z accepted" << endl;
351  cout << endl;
352  }
353  }
354  if (fObservables[i].find("Process") != string::npos) {
355  Int_t ls = 0;
356  if (fObservables[i].find("RadiactiveDecay") != string::npos) ls = 15;
357  if (fObservables[i].find("Photoelectric") != string::npos) ls = 13;
358  if (fObservables[i].find("PhotonNuclear") != string::npos) ls = 13;
359  if (fObservables[i].find("Bremstralung") != string::npos) ls = 12;
360  if (fObservables[i].find("NInelastic") != string::npos) ls = 10;
361  if (fObservables[i].find("HadElastic") != string::npos) ls = 10;
362  if (fObservables[i].find("NCapture") != string::npos) ls = 8;
363  if (fObservables[i].find("Compton") != string::npos) ls = 7;
364  if (fObservables[i].find("Neutron") != string::npos) ls = 7;
365  if (fObservables[i].find("Alpha") != string::npos) ls = 5;
366  if (fObservables[i].find("Argon") != string::npos) ls = 5;
367  if (fObservables[i].find("Xenon") != string::npos) ls = 5;
368  if (fObservables[i].find("Neon") != string::npos) ls = 4;
369 
370  TString processName = fObservables[i].substr(fObservables[i].length() - (ls + 7), ls).c_str();
371  TString volName3 = fObservables[i].substr(0, fObservables[i].length() - (ls + 7)).c_str();
372  Int_t volId3 = fG4Metadata->GetActiveVolumeID(volName3);
373 
374  if (volId3 >= 0) {
375  fProcessObservables.push_back(fObservables[i]);
376  fVolumeID3.push_back(volId3);
377  fProcessName.emplace_back(processName.Data());
378  }
379  }
380  if (fObservables[i].find("TracksCounter") != string::npos) {
381  TString particleName = fObservables[i].substr(0, fObservables[i].length() - 13).c_str();
382  fTrackCounterObservables.push_back(fObservables[i]);
383  fParticleTrackCounter.emplace_back(particleName.Data());
384  }
385 
386  if (fObservables[i].find("TracksEDep") != string::npos) {
387  TString particleName = fObservables[i].substr(0, fObservables[i].length() - 10).c_str();
388  fTracksEDepObservables.push_back(fObservables[i]);
389  fParticleTrackEdep.emplace_back(particleName.Data());
390  }
391  }
392 }
393 
398  fInputG4Event = (TRestGeant4Event*)inputEvent;
399  *fOutputG4Event = *((TRestGeant4Event*)inputEvent);
400 
401  const auto sensitiveVolumeName = fG4Metadata->GetSensitiveVolume();
402 
403  Double_t sensitiveVolumeEnergy = fOutputG4Event->GetEnergyInVolume(sensitiveVolumeName.Data());
404  // Get time of the first hit in the sensitive volume
405  double hitTime = std::numeric_limits<double>::infinity();
406  for (const auto& track : fOutputG4Event->GetTracks()) {
407  const auto hits = track.GetHits();
408  for (size_t hitIndex = 0; hitIndex < hits.GetNumberOfHits(); hitIndex++) {
409  const auto volumeName = hits.GetVolumeName(hitIndex);
410  if (volumeName != sensitiveVolumeName) {
411  continue;
412  }
413  const double energy = hits.GetEnergy(hitIndex);
414  if (energy <= 0) {
415  continue;
416  }
417 
418  const double time = hits.GetTime(hitIndex);
419  if (time < hitTime) {
420  hitTime = time;
421  }
422  }
423  }
424  SetObservableValue("sensitiveVolumeFirstHitTime", hitTime);
425 
426  if (GetVerboseLevel() >= TRestStringOutput::REST_Verbose_Level::REST_Debug) {
427  cout << "----------------------------" << endl;
428  cout << "TRestGeant4Event : " << fOutputG4Event->GetID() << endl;
429  cout << "Sensitive volume Energy : " << sensitiveVolumeEnergy << endl;
430  cout << "Total energy : " << fOutputG4Event->GetTotalDepositedEnergy() << endl;
431  }
432 
433  SetObservableValue("sensitiveVolumeEnergy", sensitiveVolumeEnergy);
434 
435  Double_t xOrigin = fOutputG4Event->GetPrimaryEventOrigin().X();
436  SetObservableValue("xOriginPrimary", xOrigin);
437 
438  Double_t yOrigin = fOutputG4Event->GetPrimaryEventOrigin().Y();
439  SetObservableValue("yOriginPrimary", yOrigin);
440 
441  Double_t zOrigin = fOutputG4Event->GetPrimaryEventOrigin().Z();
442  SetObservableValue("zOriginPrimary", zOrigin);
443 
444  Double_t xDirection = fOutputG4Event->GetPrimaryEventDirection(0).X();
445  SetObservableValue("xDirectionPrimary", xDirection);
446 
447  Double_t yDirection = fOutputG4Event->GetPrimaryEventDirection(0).Y();
448  SetObservableValue("yDirectionPrimary", yDirection);
449 
450  Double_t zDirection = fOutputG4Event->GetPrimaryEventDirection(0).Z();
451  SetObservableValue("zDirectionPrimary", zDirection);
452 
453  TVector3 v(xDirection, yDirection, zDirection);
454  SetObservableValue("thetaPrimary", v.Theta());
455  SetObservableValue("phiPrimary", v.Phi());
456 
457  Double_t energyPrimary = fOutputG4Event->GetPrimaryEventEnergy(0);
458  SetObservableValue("energyPrimary", energyPrimary);
459 
460  Double_t energyTotal = fOutputG4Event->GetTotalDepositedEnergy();
461  SetObservableValue("totalEdep", energyTotal);
462 
463  Double_t size = fOutputG4Event->GetBoundingBoxSize();
464  SetObservableValue("boundingSize", size);
465 
466  // process names as named by Geant4
467  // processes present here will be added to the list of observables which can be used to see if the event
468  // contains the process of interest.
469  vector<string> processNames = {"phot", "compt"};
470  for (auto& processName : processNames) {
471  Int_t containsProcess = 0;
472  if (fOutputG4Event->ContainsProcess(fG4Metadata->GetGeant4PhysicsInfo().GetProcessID(processName))) {
473  containsProcess = 1;
474  }
475 
476  if (!processName.empty()) {
477  processName[0] = toupper(processName[0]);
478  SetObservableValue("containsProcess" + processName, containsProcess);
479  }
480  }
481 
482  for (unsigned int n = 0; n < fParticleTrackCounter.size(); n++) {
483  Int_t nT = fOutputG4Event->GetNumberOfTracksForParticle(fParticleTrackCounter[n]);
484  string obsName = fTrackCounterObservables[n];
485  SetObservableValue(obsName, nT);
486  }
487 
488  for (unsigned int n = 0; n < fTracksEDepObservables.size(); n++) {
489  Double_t energy = fOutputG4Event->GetEnergyDepositedByParticle(fParticleTrackEdep[n]);
490  string obsName = fTracksEDepObservables[n];
491  SetObservableValue(obsName, energy);
492  }
493 
494  for (unsigned int n = 0; n < fEnergyInObservables.size(); n++) {
495  Double_t en = fOutputG4Event->GetEnergyInVolume(fVolumeName[n]);
496  string obsName = fEnergyInObservables[n];
497  SetObservableValue(obsName, en);
498  }
499 
500  for (unsigned int n = 0; n < fMeanPosObservables.size(); n++) {
501  string obsName = fMeanPosObservables[n];
502 
503  Double_t mpos = 0;
504  if (fDirID[n] == (TString) "X")
505  mpos = fOutputG4Event->GetMeanPositionInVolume(fVolumeID2[n]).X();
506 
507  else if (fDirID[n] == (TString) "Y")
508  mpos = fOutputG4Event->GetMeanPositionInVolume(fVolumeID2[n]).Y();
509 
510  else if (fDirID[n] == (TString) "Z")
511  mpos = fOutputG4Event->GetMeanPositionInVolume(fVolumeID2[n]).Z();
512 
513  SetObservableValue(obsName, mpos);
514  }
515 
516  if (GetVerboseLevel() >= TRestStringOutput::REST_Verbose_Level::REST_Debug) {
517  cout << "G4 Tracks : " << fOutputG4Event->GetNumberOfTracks() << endl;
518  cout << "----------------------------" << endl;
519  }
520 
521  return fOutputG4Event;
522 }
523 
A base class for any REST event.
Definition: TRestEvent.h:38
A pure analysis process to extract information from a TRestGeant4Event.
void Initialize() override
Function to initialize input/output event members and define the section name.
void InitProcess() override
Process initialization. Observable names are interpreted and auxiliar observable members,...
void LoadConfig(const std::string &configFilename, const std::string &name="")
Function to load the configuration from an external configuration file.
~TRestGeant4AnalysisProcess()
Default destructor.
void LoadDefaultConfig()
Function to load the default config in absence of RML input.
TRestEvent * ProcessEvent(TRestEvent *inputEvent) override
The main processing event function.
TRestGeant4AnalysisProcess()
Default constructor.
void EndProcess() override
Function to include required actions after all events have been processed.
An event class to store geant4 generated event information.
@ REST_Debug
+show the defined debug messages