REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
src/TRestGeant4Metadata.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 
729 
730 #include "TRestGeant4Metadata.h"
731 
732 #include <TAxis.h>
733 #include <TGeoManager.h>
734 #include <TRandom.h>
735 #include <TRestGDMLParser.h>
736 #include <TRestRun.h>
737 
738 #include "TRestGeant4PrimaryGeneratorInfo.h"
739 
740 using namespace std;
741 using namespace TRestGeant4PrimaryGeneratorTypes;
742 
743 ClassImp(TRestGeant4Metadata);
744 
749 
764 TRestGeant4Metadata::TRestGeant4Metadata(const char* configFilename, const string& name)
765  : TRestMetadata(configFilename) {
766  Initialize();
767 
769 }
770 
775 
780  SetSectionName(this->ClassName());
781  SetLibraryVersion(LIBRARY_VERSION);
782 
783  fChance.clear();
784  fActiveVolumes.clear();
785  fBiasingVolumes.clear();
786 
788 }
789 
794  this->Initialize();
795 
796  fMagneticField = Get3DVectorParameterWithUnits("magneticField", TVector3(0, 0, 0));
797 
798  // Initialize the metadata members from a configfile
799  fGdmlFilename = GetParameter("gdmlFile");
800  if (fGdmlFilename == PARAMETER_NOT_FOUND_STR) {
801  fGdmlFilename = GetParameter("gdml_file"); // old name
802  }
803  if (fGdmlFilename == PARAMETER_NOT_FOUND_STR) {
804  cout << "\"gdmlFile\" parameter is not defined!" << endl;
805  exit(1);
806  }
807  if (TRestTools::isURL(fGdmlFilename.Data())) {
808  fGeometryPath = "";
810  if (fGdmlFilename == "") {
811  cout << "Error downloading geometry file from URL: " << fGdmlFilename << endl;
812  exit(1);
813  }
814  }
815 
816  fGeometryPath = GetParameter("geometryPath", "");
817 
818  string seedString = GetParameter("seed", "0");
819  if (ToUpper(seedString) == "RANDOM" || ToUpper(seedString) == "RAND" || ToUpper(seedString) == "AUTO" ||
820  seedString == "0") {
821  auto dd = new double();
822  fSeed = (uintptr_t)dd + (uintptr_t)this;
823  delete dd;
824  } else {
825  fSeed = (Long_t)StringToInteger(seedString);
826  }
827  gRandom->SetSeed(fSeed);
828 
829  // if "gdmlFile" is purely a file (without any path) and "geometryPath" is
830  // defined, we recombine them together
831  if ((((string)fGdmlFilename).find_first_not_of("./~") == 0 ||
832  ((string)fGdmlFilename).find('/') == string::npos) &&
833  fGeometryPath != "") {
834  if (fGeometryPath[fGeometryPath.Length() - 1] == '/') {
835  fGdmlFilename = fGeometryPath + GetParameter("gdmlFile");
836  } else {
837  fGdmlFilename = fGeometryPath + "/" + GetParameter("gdmlFile");
838  }
839  }
840 
841  Double_t defaultTime = 1. / REST_Units::s;
842  fSubEventTimeDelay = GetDblParameterWithUnits("subEventTimeDelay", defaultTime);
843 
844  auto nEventsString = GetParameter("nEvents");
845  if (nEventsString == PARAMETER_NOT_FOUND_STR) {
846  nEventsString = GetParameter("Nevents"); // old name
847  }
848  fNEvents = nEventsString == PARAMETER_NOT_FOUND_STR ? 0 : StringToInteger(nEventsString);
849 
850  const auto fNRequestedEntriesString = GetParameter("nRequestedEntries");
852  fNRequestedEntriesString == PARAMETER_NOT_FOUND_STR ? 0 : StringToInteger(fNRequestedEntriesString);
853 
854  fSaveAllEvents = ToUpper(GetParameter("saveAllEvents", "false")) == "TRUE" ||
855  ToUpper(GetParameter("saveAllEvents", "off")) == "ON";
856 
857  fPrintProgress = ToUpper(GetParameter("printProgress", "true")) == "TRUE" ||
858  ToUpper(GetParameter("printProgress", "on")) == "ON";
859 
860  fRegisterEmptyTracks = ToUpper(GetParameter("registerEmptyTracks", "false")) == "TRUE" ||
861  ToUpper(GetParameter("registerEmptyTracks", "off")) == "ON";
862 
863  ReadGenerator();
864  // Detector (old storage) section is processed after initializing geometry info in Detector Construction
865  // This allows to use regular expression to match logical or physical volumes etc.
866  ReadBiasing();
867 
868  fMaxTargetStepSize = GetDblParameterWithUnits("maxTargetStepSize", -1);
869  if (fMaxTargetStepSize > 0) {
870  cout << " " << endl;
871  RESTWarning << "IMPORTANT: *maxTargetStepSize* parameter is now obsolete!" << RESTendl;
872  RESTWarning << "The sensitive volume will not define any integration step limit" << RESTendl;
873  cout << " " << endl;
874  RESTWarning << "In order to avoid this warning REMOVE the *maxTargetStepSize* definition,"
875  << RESTendl;
876  RESTWarning << "and replace it by the following statement at the <storage> section" << RESTendl;
877  cout << " " << endl;
878  RESTWarning << "<activeVolume name=\"" << this->GetSensitiveVolume() << "\" maxStepSize=\""
879  << fMaxTargetStepSize << "mm\" />" << RESTendl;
880  cout << " " << endl;
881  RESTInfo << "Now, any active volume is allowed to define a maxStepSize" << RESTendl;
882  cout << " " << endl;
883  RESTInfo << "It is also possible to define a default step size for all active volumes," << RESTendl;
884  RESTInfo << "so that in case no step is defined, the default will be used." << RESTendl;
885  cout << " " << endl;
886  RESTInfo << "<storage sensitiveVolume=\"" << this->GetSensitiveVolume() << "\" maxStepSize=\""
887  << fMaxTargetStepSize << "mm\" />" << RESTendl;
888  cout << " " << endl;
889  GetChar();
890  }
891 
892  // should return success or fail
893 }
894 
913 
915  if (TRestGeant4PrimaryGeneratorTypes::StringToSpatialGeneratorTypes(
917  TRestGeant4PrimaryGeneratorTypes::SpatialGeneratorTypes::COSMIC) {
918  RESTError
919  << "TRestGeant4Metadata::GetEquivalentSimulatedTime can only be called for 'cosmic' generator"
920  << RESTendl;
921  exit(1);
922  }
923  const auto source = GetParticleSource();
924  if (TRestGeant4PrimaryGeneratorTypes::StringToEnergyDistributionTypes(
925  source->GetEnergyDistributionType().Data()) !=
926  TRestGeant4PrimaryGeneratorTypes::EnergyDistributionTypes::FORMULA2) {
927  RESTError << "TRestGeant4Metadata::GetEquivalentSimulatedTime can only be called for 'formula2' "
928  "energy distribution"
929  << RESTendl;
930  exit(1);
931  }
932  if (TRestGeant4PrimaryGeneratorTypes::StringToAngularDistributionTypes(
933  source->GetAngularDistributionType().Data()) !=
934  TRestGeant4PrimaryGeneratorTypes::AngularDistributionTypes::FORMULA2) {
935  RESTError << "TRestGeant4Metadata::GetEquivalentSimulatedTime can only be called for 'formula' "
936  "angular distribution"
937  << RESTendl;
938  exit(1);
939  }
940 
941  const auto energyRange = source->GetEnergyDistributionRange();
942  const auto angularRange = source->GetAngularDistributionRange();
943  auto function = (TF2*)source->GetEnergyAndAngularDistributionFunction()->Clone();
944  // counts per second per cm2 (distribution is already integrated over uniform phi)
945  const auto countsPerSecondPerCm2 =
946  function->Integral(energyRange.X(), energyRange.Y(), angularRange.X(), angularRange.Y(), 1E-9);
947  return countsPerSecondPerCm2;
948 }
949 
950 Double_t TRestGeant4Metadata::GetCosmicIntensityInCountsPerSecond() const {
951  const auto countsPerSecondPerCm2 = GetCosmicFluxInCountsPerCm2PerSecond();
953  const auto countsPerSecond = countsPerSecondPerCm2 * surface;
954  return countsPerSecond;
955 }
956 
957 Double_t TRestGeant4Metadata::GetEquivalentSimulatedTime() const {
958  const auto countsPerSecond = GetCosmicIntensityInCountsPerSecond();
959  const auto seconds = double(fNEvents) / countsPerSecond;
960  return seconds;
961 }
962 
963 void TRestGeant4Metadata::ReadBiasing() {
964  TiXmlElement* biasingDefinition = GetElement("biasing");
965  if (biasingDefinition == nullptr) {
966  fNBiasingVolumes = 0;
967  return;
968  }
969  TString biasEnabled = GetFieldValue("value", biasingDefinition);
970  TString biasType = GetFieldValue("type", biasingDefinition);
971 
972  RESTDebug << "Bias: " << biasEnabled << RESTendl;
973 
974  if (biasEnabled == "on" || biasEnabled == "ON" || biasEnabled == "On" || biasEnabled == "oN") {
975  RESTInfo << "Biasing is enabled" << RESTendl;
976 
977  TiXmlElement* biasVolumeDefinition = GetElement("biasingVolume", biasingDefinition);
978  Int_t n = 0;
979  while (biasVolumeDefinition) {
980  TRestGeant4BiasingVolume biasVolume;
981  RESTDebug << "Def: " << biasVolumeDefinition << RESTendl;
982 
983  biasVolume.SetBiasingVolumePosition(
984  Get3DVectorParameterWithUnits("position", biasVolumeDefinition));
985  biasVolume.SetBiasingFactor(StringToDouble(GetParameter("factor", biasVolumeDefinition)));
986  biasVolume.SetBiasingVolumeSize(GetDblParameterWithUnits("size", biasVolumeDefinition));
987  biasVolume.SetEnergyRange(Get2DVectorParameterWithUnits("energyRange", biasVolumeDefinition));
988  biasVolume.SetBiasingVolumeType(biasType); // For the moment all the volumes should be same type
989 
990  /* TODO check that values are right if not printBiasingVolume with
991  getchar() biasVolume.PrintBiasingVolume(); getchar();
992  */
993 
994  fBiasingVolumes.push_back(biasVolume);
995  biasVolumeDefinition = GetNextElement(biasVolumeDefinition);
996  n++;
997  }
998  fNBiasingVolumes = n;
999  }
1000 }
1001 
1027  // TODO if some fields are defined in the generator but not finally used
1028  // i.e. <generator type="volume" from="gasTarget" position="(100,0,-100)
1029  // here position is irrelevant since the event will be generated from the
1030  // volume defined in the geometry we should take care of these values so they
1031  // are not stored in metadata (to avoid future confusion) In the particular
1032  // case of position, the value is overwritten in DetectorConstruction by the
1033  // center of the volume (i.e. gasTarget) but if i.e rotation or side are
1034  // defined and not relevant we should set it to -1
1035 
1036  TiXmlElement* generatorDefinition = GetElement("generator");
1037 
1039  SpatialGeneratorTypesToString(StringToSpatialGeneratorTypes(GetParameter(
1040  "type", generatorDefinition, SpatialGeneratorTypesToString(SpatialGeneratorTypes::VOLUME))));
1042  SpatialGeneratorShapesToString(StringToSpatialGeneratorShapes(GetParameter(
1043  "shape", generatorDefinition, SpatialGeneratorShapesToString(SpatialGeneratorShapes::BOX))));
1044  fGeant4PrimaryGeneratorInfo.fSpatialGeneratorFrom = GetParameter("from", generatorDefinition);
1045  if (fGeant4PrimaryGeneratorInfo.fSpatialGeneratorFrom != PARAMETER_NOT_FOUND_STR) {
1047  SpatialGeneratorShapesToString(SpatialGeneratorShapes::GDML);
1048  }
1050  Get3DVectorParameterWithUnits("size", generatorDefinition, {0, 0, 0});
1052  Get3DVectorParameterWithUnits("position", generatorDefinition, {0, 0, 0});
1054  Get3DVectorParameterWithUnits("rotationAxis", generatorDefinition, {0, 0, 1});
1056  GetDblParameterWithUnits("rotationAngle", generatorDefinition, 0);
1058  GetParameter("densityFunc", generatorDefinition, "1");
1059 
1060  TiXmlElement* sourceDefinition = GetElement("source", generatorDefinition);
1061  while (sourceDefinition) {
1062  string use = GetParameter("use", sourceDefinition, "");
1063 
1064  TRestGeant4ParticleSource* source = TRestGeant4ParticleSource::instantiate(use);
1065  ReadParticleSource(source, sourceDefinition);
1066  AddParticleSource(source);
1067 
1068  sourceDefinition = GetNextElement(sourceDefinition);
1069  }
1070 
1071  // check if the generator is valid.
1072  if (GetNumberOfSources() == 0) {
1073  RESTError << "No sources are added inside generator!" << RESTendl;
1074  exit(1);
1075  }
1076 }
1077 
1082  TiXmlElement* sourceDefinition = definition;
1083 
1084  source->SetParticleName(GetParameter("particle", sourceDefinition));
1085  source->SetExcitationLevel(StringToDouble(GetParameter("excitedLevel", sourceDefinition)));
1086  source->SetParticleCharge(StringToInteger(GetParameter("charge", sourceDefinition, "0")));
1087 
1088  TString fullChain = GetParameter("fullChain", sourceDefinition);
1089  SetFullChain(false);
1090  if (fullChain.EqualTo("on", TString::ECaseCompare::kIgnoreCase)) {
1091  SetFullChain(true);
1092  }
1093 
1094  // Angular distribution parameters
1095  TiXmlElement* angularDefinition = GetElement("angular", sourceDefinition);
1096  if (angularDefinition == nullptr) {
1097  angularDefinition = GetElement("angularDist", sourceDefinition); // old name
1098  }
1099  source->SetAngularDistributionType(GetParameter(
1100  "type", angularDefinition, AngularDistributionTypesToString(AngularDistributionTypes::FLUX)));
1101  if (StringToAngularDistributionTypes(source->GetAngularDistributionType().Data()) ==
1102  AngularDistributionTypes::TH1D) {
1103  source->SetAngularDistributionFilename(SearchFile(GetParameter("file", angularDefinition)));
1104  auto name = GetParameter("name", angularDefinition);
1105  if (name == "NO_SUCH_PARA") {
1106  name = GetParameter("spctName", angularDefinition); // old name
1107  }
1108  source->SetAngularDistributionNameInFile(name);
1109  }
1110  source->SetAngularDistributionRange(
1111  Get2DVectorParameterWithUnits("range", angularDefinition, {0, TMath::Pi()}));
1112  if (StringToAngularDistributionTypes(source->GetAngularDistributionType().Data()) ==
1113  AngularDistributionTypes::FORMULA) {
1114  source->SetAngularDistributionFormula(GetParameter("name", angularDefinition));
1115  // We cannot use an angular range bigger than the range of the formula
1116  const auto function = source->GetAngularDistributionFunction();
1117  if (source->GetAngularDistributionRangeMin() < function->GetXaxis()->GetXmin()) {
1118  source->SetAngularDistributionRange(
1119  {function->GetXaxis()->GetXmin(), source->GetAngularDistributionRangeMax()});
1120  }
1121  if (source->GetAngularDistributionRangeMax() > function->GetXaxis()->GetXmax()) {
1122  source->SetAngularDistributionRange(
1123  {source->GetAngularDistributionRangeMin(), function->GetXaxis()->GetXmax()});
1124  }
1125  }
1126  if ((StringToAngularDistributionTypes(source->GetAngularDistributionType().Data()) ==
1127  AngularDistributionTypes::FORMULA) ||
1128  (StringToAngularDistributionTypes(source->GetAngularDistributionType().Data()) ==
1129  AngularDistributionTypes::FORMULA2)) {
1130  source->SetAngularDistributionFormulaNPoints(static_cast<size_t>(GetDblParameterWithUnits(
1131  "nPoints", angularDefinition, source->GetAngularDistributionFormulaNPoints())));
1132  }
1133  if (GetNumberOfSources() == 0 &&
1134  StringToAngularDistributionTypes(source->GetAngularDistributionType().Data()) ==
1135  AngularDistributionTypes::BACK_TO_BACK) {
1136  cout << "WARNING: First source cannot be backtoback. Setting it to isotropic" << endl;
1137  source->SetAngularDistributionType(
1138  AngularDistributionTypesToString(AngularDistributionTypes::ISOTROPIC));
1139  }
1140  source->SetDirection(StringTo3DVector(GetParameter("direction", angularDefinition, "(1,0,0)")));
1141 
1142  // Energy distribution parameters
1143  TiXmlElement* energyDefinition = GetElement("energy", sourceDefinition);
1144  if (energyDefinition == nullptr) {
1145  energyDefinition = GetElement("energyDist", sourceDefinition); // old name
1146  }
1147  source->SetEnergyDistributionType(GetParameter(
1148  "type", energyDefinition, EnergyDistributionTypesToString(EnergyDistributionTypes::MONO)));
1149  if (StringToEnergyDistributionTypes(source->GetEnergyDistributionType().Data()) ==
1150  EnergyDistributionTypes::TH1D) {
1151  source->SetEnergyDistributionFilename(SearchFile(GetParameter("file", energyDefinition)));
1152  auto name = GetParameter("name", energyDefinition);
1153  if (name == "NO_SUCH_PARA") {
1154  name = GetParameter("spctName", energyDefinition); // old name
1155  }
1156  source->SetEnergyDistributionNameInFile(name);
1157  }
1158  source->SetEnergyDistributionRange(Get2DVectorParameterWithUnits("range", energyDefinition, {0, 1E20}));
1159  if (StringToEnergyDistributionTypes(source->GetEnergyDistributionType().Data()) ==
1160  EnergyDistributionTypes::MONO) {
1161  Double_t energy;
1162  energy = GetDblParameterWithUnits("energy", energyDefinition, 0);
1163  source->SetEnergyDistributionRange({energy, energy});
1164  source->SetEnergy(energy);
1165  }
1166  if (StringToEnergyDistributionTypes(source->GetEnergyDistributionType().Data()) ==
1167  EnergyDistributionTypes::FORMULA) {
1168  source->SetEnergyDistributionFormula(GetParameter("name", energyDefinition));
1169  // We cannot use an energy range bigger than the range of the formula
1170  const auto function = source->GetEnergyDistributionFunction();
1171  if (source->GetEnergyDistributionRangeMin() < function->GetXaxis()->GetXmin()) {
1172  source->SetEnergyDistributionRange(
1173  {function->GetXaxis()->GetXmin(), source->GetEnergyDistributionRangeMax()});
1174  }
1175  if (source->GetEnergyDistributionRangeMax() > function->GetXaxis()->GetXmax()) {
1176  source->SetEnergyDistributionRange(
1177  {source->GetEnergyDistributionRangeMin(), function->GetXaxis()->GetXmax()});
1178  }
1179  }
1180  if ((StringToEnergyDistributionTypes(source->GetEnergyDistributionType().Data()) ==
1181  EnergyDistributionTypes::FORMULA) ||
1182  (StringToEnergyDistributionTypes(source->GetEnergyDistributionType().Data()) ==
1183  EnergyDistributionTypes::FORMULA2)) {
1184  source->SetEnergyDistributionFormulaNPoints(static_cast<size_t>(GetDblParameterWithUnits(
1185  "nPoints", energyDefinition, source->GetEnergyDistributionFormulaNPoints())));
1186  }
1187  if (StringToEnergyDistributionTypes(source->GetEnergyDistributionType().Data()) ==
1188  EnergyDistributionTypes::FORMULA2 &&
1189  StringToAngularDistributionTypes(source->GetAngularDistributionType().Data()) ==
1190  AngularDistributionTypes::FORMULA2) {
1191  const auto empty = TString("");
1192  auto energyDistName = GetParameter("name", energyDefinition, empty);
1193  auto angularDistName = GetParameter("name", energyDefinition, empty);
1194 
1195  if (energyDistName == empty && angularDistName == empty) {
1196  RESTError << "No name specified for energy and angular distribution" << RESTendl;
1197  exit(1);
1198  }
1199  if (energyDistName == empty) {
1200  angularDistName = energyDistName;
1201  RESTWarning << "No name specified for energy distribution. Using angular distribution name: "
1202  << angularDistName << RESTendl;
1203  }
1204  if (angularDistName == empty) {
1205  energyDistName = energyDistName;
1206  RESTWarning << "No name specified for angular distribution. Using energy distribution name: "
1207  << energyDistName << RESTendl;
1208  }
1209 
1210  if (energyDistName == empty || angularDistName == empty) {
1211  // we should never enter here, just leave this as a check
1212  RESTError << "When using 'formula2' the name of energy and angular dist must match" << RESTendl;
1213  exit(1);
1214  }
1215  const auto energyAndAngularDistName = energyDistName;
1216  source->SetEnergyAndAngularDistributionFormula(energyAndAngularDistName);
1217 
1218  const auto function = source->GetEnergyAndAngularDistributionFunction();
1219 
1220  // Set energy range
1221  if (source->GetEnergyDistributionRangeMin() < function->GetXaxis()->GetXmin()) {
1222  source->SetEnergyDistributionRange(
1223  {function->GetXaxis()->GetXmin(), source->GetEnergyDistributionRangeMax()});
1224  }
1225  if (source->GetEnergyDistributionRangeMax() > function->GetXaxis()->GetXmax()) {
1226  source->SetEnergyDistributionRange(
1227  {source->GetEnergyDistributionRangeMin(), function->GetXaxis()->GetXmax()});
1228  }
1229 
1230  // Set angular range
1231  if (source->GetAngularDistributionRangeMin() < function->GetYaxis()->GetXmin()) {
1232  source->SetAngularDistributionRange(
1233  {function->GetYaxis()->GetXmin(), source->GetAngularDistributionRangeMax()});
1234  }
1235  if (source->GetAngularDistributionRangeMax() > function->GetYaxis()->GetXmax()) {
1236  source->SetAngularDistributionRange(
1237  {source->GetAngularDistributionRangeMin(), function->GetYaxis()->GetXmax()});
1238  }
1239  }
1240  // allow custom configuration from the class
1241  source->LoadConfigFromElement(sourceDefinition, fElementGlobal, fVariables);
1242  // AddParticleSource(source);
1243 }
1244 
1246  for (auto source : fParticleSource) {
1247  delete source;
1248  }
1249  fParticleSource.clear();
1250 }
1251 
1253  fParticleSource.push_back(src);
1254 }
1255 
1271  RESTInfo << "TRestGeant4Metadata: Processing detector section" << RESTendl;
1272  TiXmlElement* detectorDefinition = GetElement("detector");
1273  if (detectorDefinition == nullptr) {
1274  RESTError << "Detector section (<detector>) is not present in metadata!" << RESTendl;
1275  exit(1);
1276  }
1277 
1278  const bool activateAllVolumes =
1279  StringToBool(GetParameter("activateAllVolumes", detectorDefinition, "true"));
1280  RESTInfo << "TRestGeant4Metadata: Setting 'fActivateAllVolumes' to " << fActivateAllVolumes << RESTendl;
1281  fActivateAllVolumes = activateAllVolumes;
1282 
1283  TiXmlElement* removeUnwantedTracksDefinition = GetElement("removeUnwantedTracks", detectorDefinition);
1284  if (removeUnwantedTracksDefinition != nullptr) {
1285  fRemoveUnwantedTracks = StringToBool(GetParameter("enabled", removeUnwantedTracksDefinition, "true"));
1287  StringToBool(GetParameter("keepZeroEnergyTracks", removeUnwantedTracksDefinition, "false"));
1288  RESTInfo << "TRestGeant4Metadata: Setting 'fRemoveUnwantedTracks' to " << fRemoveUnwantedTracks
1289  << " with 'keepZeroEnergyTracks' option set to " << fRemoveUnwantedTracksKeepZeroEnergyTracks
1290  << RESTendl;
1291  }
1292 
1293  Double_t defaultStep = GetDblParameterWithUnits("maxStepSize", detectorDefinition);
1294  if (defaultStep < 0) {
1295  defaultStep = 0;
1296  }
1297 
1298  TiXmlElement* volumeDefinition = GetElement("volume", detectorDefinition);
1299  while (volumeDefinition != nullptr) {
1300  string name = GetFieldValue("name", volumeDefinition);
1301  if (name.empty() || name == "Not defined") {
1302  RESTError << "volume inside detector section defined without name" << RESTendl;
1303  exit(1);
1304  }
1305  vector<string> physicalVolumes;
1306  if (!fGeant4GeometryInfo.IsValidGdmlName(name)) {
1307  if (fGeant4GeometryInfo.IsValidLogicalVolume(name)) {
1308  for (const auto& physical : fGeant4GeometryInfo.GetAllPhysicalVolumesFromLogical(name)) {
1309  physicalVolumes.emplace_back(
1310  fGeant4GeometryInfo.GetAlternativeNameFromGeant4PhysicalName(physical));
1311  }
1312  }
1313  // does not match a explicit (gdml) physical name or a logical name, perhaps its a regular exp
1314  if (physicalVolumes.empty()) {
1315  for (const auto& physical :
1316  fGeant4GeometryInfo.GetAllPhysicalVolumesMatchingExpression(name)) {
1317  physicalVolumes.emplace_back(
1318  fGeant4GeometryInfo.GetAlternativeNameFromGeant4PhysicalName(physical));
1319  }
1320  }
1321  if (physicalVolumes.empty()) {
1322  for (const auto& logical : fGeant4GeometryInfo.GetAllLogicalVolumesMatchingExpression(name)) {
1323  for (const auto& physical :
1324  fGeant4GeometryInfo.GetAllPhysicalVolumesFromLogical(logical)) {
1325  physicalVolumes.emplace_back(
1326  fGeant4GeometryInfo.GetAlternativeNameFromGeant4PhysicalName(physical));
1327  }
1328  }
1329  }
1330  } else {
1331  physicalVolumes.push_back(name);
1332  }
1333 
1334  if (physicalVolumes.empty()) {
1335  RESTError
1336  << "volume '" << name
1337  << "' is not defined in the geometry. Could not match to a physical, logical volume or regex"
1338  << RESTendl;
1339  exit(1);
1340  }
1341 
1342  bool isSensitive = false;
1343  const string isSensitiveValue = GetFieldValue("sensitive", volumeDefinition);
1344  if (isSensitiveValue != "Not defined") {
1345  isSensitive = StringToBool(isSensitiveValue);
1346  }
1347 
1348  bool isKeepTracks = isSensitive;
1349  const string isKeepTracksValue = GetFieldValue("keepTracks", volumeDefinition);
1350  if (isKeepTracksValue != "Not defined") {
1351  isKeepTracks = StringToBool(isKeepTracksValue);
1352  }
1353 
1354  bool isKill = false;
1355  const string isKillValue = GetFieldValue("kill", volumeDefinition);
1356  if (isKillValue != "Not defined") {
1357  isKill = StringToBool(isKillValue);
1358  }
1359 
1360  Double_t chance = StringToDouble(GetFieldValue("chance", volumeDefinition));
1361  if (chance == -1 || chance < 0) {
1362  chance = 1.0;
1363  }
1364  Double_t maxStep = GetDblParameterWithUnits("maxStepSize", volumeDefinition);
1365  if (maxStep < 0) {
1366  maxStep = defaultStep;
1367  }
1368 
1369  for (const auto& physical : physicalVolumes) {
1370  RESTInfo << "Adding " << (isSensitive ? "sensitive" : "active") << " volume from RML: '"
1371  << physical << (chance != 1 ? " with change: " + to_string(chance) : " ") << RESTendl;
1372  SetActiveVolume(physical, chance, maxStep);
1373  if (isSensitive) {
1374  InsertSensitiveVolume(physical);
1375  }
1376  if (fRemoveUnwantedTracks && isKeepTracks) {
1377  fRemoveUnwantedTracksVolumesToKeep.insert(physical);
1378  }
1379  if (isKill) {
1380  fKillVolumes.insert(physical);
1381  }
1382  }
1383  volumeDefinition = GetNextElement(volumeDefinition);
1384  }
1385 
1386  if (fSensitiveVolumes.empty()) {
1387  cout << "No sensitive volumes defined in detector section. Adding 'gas' as sensitive volume" << endl;
1388  InsertSensitiveVolume("gas");
1389  }
1390 
1391  fEnergyRangeStored = Get2DVectorParameterWithUnits("energyRange", detectorDefinition, fEnergyRangeStored);
1392  // TODO: Place this as a generic method
1393  if (fEnergyRangeStored.X() < 0) {
1395  }
1396  if (fEnergyRangeStored.Y() < 0) {
1398  }
1399  if (fEnergyRangeStored.Y() > fEnergyRangeStored.Y()) {
1401  }
1402  if (fEnergyRangeStored.Y() <= 0) {
1403  RESTError << "Energy range is not valid: (" << fEnergyRangeStored.X() << ", "
1404  << fEnergyRangeStored.Y() << ") keV" << RESTendl;
1405  exit(1);
1406  }
1407 
1408  auto gdml = new TRestGDMLParser();
1409  gdml->Load(GetGdmlFilename().Data());
1410 
1411  fGeant4GeometryInfo.PopulateFromGdml(gdml->GetOutputGDMLFile());
1412 
1413  // If the user did not add explicitly any volume to the storage section we understand
1414  // the user wants to register all the volumes
1415  if (fActivateAllVolumes) {
1416  for (const auto& name : fGeant4GeometryInfo.GetAllPhysicalVolumes()) {
1417  if (!IsActiveVolume(name)) {
1418  SetActiveVolume(name, 1, defaultStep);
1419  RESTInfo << "Automatically adding active volume: '" << name << "' with chance: " << 1
1420  << RESTendl;
1421  }
1422  }
1423  }
1424 
1425  if (GetNumberOfActiveVolumes() == 0) {
1426  RESTError << "No active volumes defined. Please check the detector section" << RESTendl;
1427  exit(1);
1428  }
1429 }
1430 
1437 
1438  RESTMetadata << "Geant4 version: " << GetGeant4Version() << RESTendl;
1439  RESTMetadata << "Random seed: " << GetSeed() << RESTendl;
1440  RESTMetadata << "GDML geometry: " << GetGdmlReference() << RESTendl;
1441  RESTMetadata << "GDML materials reference: " << GetMaterialsReference() << RESTendl;
1442  RESTMetadata << "Sub-event time delay: " << GetSubEventTimeDelay() << " us" << RESTendl;
1443  Double_t mx = GetMagneticField().X();
1444  Double_t my = GetMagneticField().Y();
1445  Double_t mz = GetMagneticField().Z();
1446  RESTMetadata << "Magnetic field: (" << mx << ", " << my << ", " << mz << ") T" << RESTendl;
1447  if (fSaveAllEvents) RESTMetadata << "Save all events was enabled!" << RESTendl;
1448  if (fRegisterEmptyTracks) {
1449  RESTMetadata << "Register empty tracks was enabled" << RESTendl;
1450  } else {
1451  RESTMetadata << "Register empty tracks was NOT enabled" << RESTendl;
1452  }
1453 
1454  RESTMetadata << " " << RESTendl;
1455  RESTMetadata << " ++++++++++ Generator +++++++++++ " << RESTendl;
1456  RESTMetadata << "Number of generated events: " << GetNumberOfEvents() << RESTendl;
1458 
1459  for (int i = 0; i < GetNumberOfSources(); i++) GetParticleSource(i)->PrintMetadata();
1460 
1461  RESTMetadata << " " << RESTendl;
1462  RESTMetadata << " ++++++++++ Detector +++++++++++ " << RESTendl;
1463  RESTMetadata << "Energy range (keV): (" << GetMinimumEnergyStored() << ", " << GetMaximumEnergyStored()
1464  << ")" << RESTendl;
1465  RESTMetadata << "Number of sensitive volumes: " << GetNumberOfSensitiveVolumes() << RESTendl;
1466  for (const auto& sensitiveVolume : fSensitiveVolumes) {
1467  RESTMetadata << "Sensitive volume: " << sensitiveVolume << RESTendl;
1468  }
1469  RESTMetadata << "Number of active volumes: " << GetNumberOfActiveVolumes() << RESTendl;
1470  for (unsigned int n = 0; n < GetNumberOfActiveVolumes(); n++) {
1471  const auto name = GetActiveVolumeName(n);
1472  RESTMetadata << "Name: " << name << ", ID: " << fGeant4GeometryInfo.GetIDFromVolume(name)
1473  << ", maxStep: " << GetMaxStepSize(name) << "mm "
1474  << ", chance: " << GetStorageChance(name) << RESTendl;
1475  }
1476 
1477  for (unsigned int n = 0; n < GetNumberOfBiasingVolumes(); n++) {
1478  GetBiasingVolume(n).PrintBiasingVolume();
1479  }
1480 
1481  if (GetRemoveUnwantedTracks()) {
1482  RESTMetadata << "removeUnwantedTracks is enabled "
1483  << (fRemoveUnwantedTracksKeepZeroEnergyTracks ? "with" : "without")
1484  << " keeping zero energy tracks" << RESTendl;
1485 
1486  /*
1487  RESTMetadata << "keep volumes for removeUnwantedTracks:" << RESTendl;
1488  for (const auto& volume : fRemoveUnwantedTracksVolumesToKeep) {
1489  RESTMetadata << "\t" << volume << RESTendl;
1490  }
1491  */
1492  }
1493 
1494  RESTMetadata << "+++++" << RESTendl;
1495 }
1496 
1499 Int_t TRestGeant4Metadata::GetActiveVolumeID(const TString& name) {
1500  Int_t id;
1501  for (id = 0; id < (Int_t)fActiveVolumes.size(); id++) {
1502  if (fActiveVolumes[id] == name) return id;
1503  }
1504  return -1;
1505 }
1506 
1521 void TRestGeant4Metadata::SetActiveVolume(const TString& name, Double_t chance, Double_t maxStep) {
1522  for (size_t i = 0; i < fActiveVolumes.size(); i++) {
1523  const auto activeVolumeName = fActiveVolumes[i];
1524  if (name == activeVolumeName) {
1525  fChance[i] = chance;
1526  fMaxStepSize[i] = maxStep;
1527  return;
1528  }
1529  }
1530  fActiveVolumes.push_back(name);
1531  fChance.push_back(chance);
1532  fMaxStepSize.push_back(maxStep);
1533  fActiveVolumesSet.insert(name.Data());
1534 }
1535 
1540 Bool_t TRestGeant4Metadata::isVolumeStored(const TString& volume) const {
1541  for (unsigned int n = 0; n < GetNumberOfActiveVolumes(); n++)
1542  if (GetActiveVolumeName(n) == volume) return true;
1543 
1544  return false;
1545 }
1546 
1550 Double_t TRestGeant4Metadata::GetStorageChance(const TString& volume) {
1551  Int_t id;
1552  for (id = 0; id < (Int_t)fActiveVolumes.size(); id++) {
1553  if (fActiveVolumes[id] == volume) return fChance[id];
1554  }
1555  RESTWarning << "TRestGeant4Metadata::GetStorageChance. Volume " << volume << " not found" << RESTendl;
1556 
1557  return 0;
1558 }
1559 
1563 Double_t TRestGeant4Metadata::GetMaxStepSize(const TString& volume) {
1564  for (Int_t i = 0; i < (Int_t)fActiveVolumes.size(); i++) {
1565  if (volume.EqualTo(fActiveVolumes[i], TString::kIgnoreCase)) {
1566  return fMaxStepSize[i];
1567  }
1568  }
1569  RESTWarning << "TRestGeant4Metadata::GetMaxStepSize. Volume " << volume << " not found" << RESTendl;
1570 
1571  return 0;
1572 }
1573 
1574 size_t TRestGeant4Metadata::GetGeant4VersionMajor() const {
1575  TString majorVersion = "";
1576  for (int i = 0; i < fGeant4Version.Length(); i++) {
1577  auto c = fGeant4Version[i];
1578  if (c == '.') {
1579  break;
1580  }
1581  majorVersion += c;
1582  }
1583  return std::stoi(majorVersion.Data());
1584 }
1585 
1586 void TRestGeant4Metadata::Merge(const TRestGeant4Metadata& metadata) {
1587  fIsMerge = true;
1588  fSeed = 0; // seed makes no sense in a merged file
1589 
1590  fNEvents += metadata.fNEvents;
1592  fSimulationTime += metadata.fSimulationTime;
1593 }
1594 
1595 TRestGeant4Metadata::TRestGeant4Metadata(const TRestGeant4Metadata& metadata) { *this = metadata; }
1596 
1597 TRestGeant4Metadata& TRestGeant4Metadata::operator=(const TRestGeant4Metadata& metadata) {
1598  fIsMerge = metadata.fIsMerge;
1602  fGeant4Version = metadata.fGeant4Version;
1603  fGdmlReference = metadata.fGdmlReference;
1606  fActiveVolumes = metadata.fActiveVolumes;
1607  fChance = metadata.fChance;
1608  fMaxStepSize = metadata.fMaxStepSize;
1609  // fParticleSource = metadata.fParticleSource; // segfaults (pointers!)
1611  fBiasingVolumes = metadata.fBiasingVolumes;
1614  fFullChain = metadata.fFullChain;
1616  fNEvents = metadata.fNEvents;
1619  fSeed = metadata.fSeed;
1620  fSaveAllEvents = metadata.fSaveAllEvents;
1624  fKillVolumes = metadata.fKillVolumes;
1626  fMagneticField = metadata.fMagneticField;
1627  return *this;
1628 }
The main class to store the Geant4 simulation conditions that will be used by restG4.
std::vector< TString > fSensitiveVolumes
The volume that serves as trigger for data storage. Only events that deposit a non-zero energy on thi...
void ReadGenerator()
Reads the generator section defined inside TRestGeant4Metadata.
Long_t GetSeed() const
// Used for faster lookup
std::vector< TRestGeant4BiasingVolume > fBiasingVolumes
A std::vector containing the biasing volume properties.
Bool_t fFullChain
Defines if a radioactive isotope decay is simulated in full chain (fFullChain=true),...
size_t GetNumberOfBiasingVolumes() const
Returns the number of biasing volumes defined.
void ReadDetector()
Reads the storage section defined inside TRestGeant4Metadata.
Bool_t fActivateAllVolumes
Sets all volume as active without having to explicitly list them.
TString fGeometryPath
The local path to the GDML geometry.
TString GetSensitiveVolume(int n=0) const
Returns a std::string with the name of the sensitive volume.
TVector3 GetMagneticField() const
Returns the world magnetic field in Tesla.
void AddParticleSource(TRestGeant4ParticleSource *src)
Adds a new particle source.
Int_t GetActiveVolumeID(const TString &name)
Returns the id of an active volume giving as parameter its name.
TRestGeant4Metadata()
Default constructor.
void InitFromConfigFile() override
Initialization of TRestGeant4Metadata members through a RML file.
TString fMaterialsReference
A GDML materials reference introduced in the header of the GDML of materials definition.
TString GetGdmlReference() const
Returns the reference provided at the GDML file header.
Int_t fNBiasingVolumes
The number of biasing volumes used in the simulation. If zero, no biasing technique is used.
void Initialize() override
Initialization of TRestGeant4Metadata members.
void PrintMetadata() override
Prints on screen the details about the Geant4 simulation conditions, stored in TRestGeant4Metadata.
void SetActiveVolume(const TString &name, Double_t chance, Double_t maxStep=0)
Adds a geometry volume to the list of active volumes.
Double_t GetMinimumEnergyStored() const
Returns the minimum event energy required for an event to be stored.
TRestGeant4GeometryInfo fGeant4GeometryInfo
Class used to store and retrieve geometry info.
unsigned int GetNumberOfActiveVolumes() const
Returns the number of active volumes, or geometry volumes that have been selected for data storage.
Double_t GetSubEventTimeDelay() const
Returns the time gap, in us, required to consider a Geant4 hit as a new independent event....
TString fGdmlReference
A GDML geometry reference introduced in the header of the GDML main setup.
Int_t fSimulationMaxTimeSeconds
Time before simulation is ended and saved.
Double_t GetMaxStepSize(const TString &volume)
Returns the maximum step at a particular active volume.
~TRestGeant4Metadata()
Default destructor.
std::vector< Double_t > fMaxStepSize
A std::vector to store the maximum step size at a particular volume.
Long64_t fSimulationTime
The time, in seconds, that the simulation took to complete.
std::set< std::string > fRemoveUnwantedTracksVolumesToKeep
A container related to fRemoveUnwantedTracks.
TVector2 fEnergyRangeStored
A 2d-std::vector storing the energy range, in keV, to decide if a particular event should be written ...
TVector3 fMagneticField
The world magnetic field.
TRestGeant4PhysicsInfo fGeant4PhysicsInfo
Class used to store and retrieve physics info such as process names or particle names.
TRestGeant4BiasingVolume GetBiasingVolume(int n)
Return the biasing volume with index n.
TString fGeant4Version
The version of Geant4 used to generate the data.
TString GetActiveVolumeName(Int_t n) const
Returns a std::string with the name of the active volume with index n.
std::vector< TString > fActiveVolumes
A std::vector to store the names of the active volumes.
TString GetGeant4Version() const
Returns a std::string with the version of Geant4 used on the event data simulation.
Long64_t fNEvents
The number of events simulated, or to be simulated.
Long_t fSeed
The seed value used for Geant4 random event generator. If it is zero its value will be assigned using...
void SetFullChain(Bool_t fullChain)
Enables/disables the full chain decay generation.
Long64_t GetNumberOfEvents() const
Returns the number of events to be simulated.
Int_t GetNumberOfSources() const
Returns the number of primary event sources defined in the generator.
void RemoveParticleSources()
Remove all the particle sources.
std::vector< Double_t > fChance
A std::vector to store the probability value to write to disk the hits in a particular event.
Double_t GetCosmicFluxInCountsPerCm2PerSecond() const
Reads the biasing section defined inside TRestGeant4Metadata.
Bool_t fPrintProgress
If this parameter is set to 'true' it will print out on screen every time 10k events are reached.
TRestGeant4ParticleSource * GetParticleSource(size_t n=0) const
Returns the name of the particle source with index n (Geant4 based names).
Bool_t fSaveAllEvents
If this parameter is set to 'true' it will save all events even if they leave no energy in the sensit...
Double_t fMaxTargetStepSize
The maximum target step size, in mm, allowed in Geant4 for the target volume (usually the gas volume)...
Double_t GetMaximumEnergyStored() const
Returns the maximum event energy required for an event to be stored.
std::vector< TRestGeant4ParticleSource * > fParticleSource
It the defines the primary source properties, particle type, momentum, energy, etc.
Bool_t fRegisterEmptyTracks
If this parameter is enabled it will register tracks with no hits inside. This is the default behavio...
Double_t GetStorageChance(Int_t n) const
Returns the probability per event to register (write to disk) hits in the storage volume with index n...
Double_t fSubEventTimeDelay
A time gap, in us, determining if an energy hit should be considered (and stored) as an independent e...
TString fGdmlFilename
The filename of the GDML geometry.
Bool_t fRemoveUnwantedTracks
If activated will remove tracks not present in volumes marked as "keep" or "sensitive".
void ReadParticleSource(TRestGeant4ParticleSource *src, TiXmlElement *sourceDefinition)
It reads the <source definition given by argument.
TRestGeant4PrimaryGeneratorInfo fGeant4PrimaryGeneratorInfo
Class used to store and retrieve Geant4 primary generator info.
Long64_t fNRequestedEntries
The number of events the user requested to be on the file.
TString GetGdmlFilename() const
Returns the main filename of the GDML geometry.
Bool_t isVolumeStored(const TString &volume) const
Returns true if the volume named volName has been registered for data storage.
void InsertSensitiveVolume(const TString &volume)
Sets the name of the sensitive volume.
Bool_t fRemoveUnwantedTracksKeepZeroEnergyTracks
Option for 'removeUnwantedTracks', if enabled tracks with hits in volumes will be kept event if they ...
TString GetMaterialsReference() const
Returns the reference provided at the materials file header.
virtual void PrintMetadata() override
Implemented it in the derived metadata class to print out specific metadata information.
TVector3 fSpatialGeneratorSize
The size of the generator. Stores up to three dimensions according to the shape box: length,...
TString GetSpatialGeneratorType() const
Returns a std::string specifying the generator type (volume, surface, point, virtualWall,...
TString fSpatialGeneratorFrom
The volume name where the events are generated, in case of volume or surface generator types.
TString fSpatialGeneratorSpatialDensityFunction
Defines density distribution of the generator shape. rho=F(x,y,z), in range 0~1.
Double_t fSpatialGeneratorRotationValue
degrees of rotation for generator virtual shape along the axis
TVector3 fSpatialGeneratorPosition
The position of the generator for virtual, and point, generator types.
Double_t GetSpatialGeneratorCosmicSurfaceTermCm2() const
Returns cosmic surface term (cm2) for simulation time computation.
TVector3 fSpatialGeneratorRotationAxis
A 3d-std::vector with the angles, measured in degrees, of a XYZ rotation applied to the virtual gener...
TString fSpatialGeneratorShape
Shape of spatial generator (wall, GDML, sphere, etc)
TString fSpatialGeneratorType
Type of spatial generator (point, surface, volume, custom)
A base class for any REST metadata class.
Definition: TRestMetadata.h:74
std::map< std::string, std::string > fVariables
Saving a list of rml variables. name-value std::pair.
virtual void PrintMetadata()
Implemented it in the derived metadata class to print out specific metadata information.
endl_t RESTendl
Termination flag object for TRestStringOutput.
TiXmlElement * GetElement(std::string eleDeclare, TiXmlElement *e=nullptr)
Get an xml element from a given parent element, according to its declaration.
Int_t LoadConfigFromFile(const std::string &configFilename, const std::string &sectionName="")
Give the file name, find out the corresponding section. Then call the main starter.
TiXmlElement * fElementGlobal
Saving the global element, to be passed to the resident class, if necessary.
void SetLibraryVersion(TString version)
Set the library version of this metadata class.
std::string GetFieldValue(std::string parName, TiXmlElement *e)
Returns the field value of an xml element which has the specified name.
Int_t LoadConfigFromElement(TiXmlElement *eSectional, TiXmlElement *eGlobal, std::map< std::string, std::string > envs={})
Main starter method.
void SetSectionName(std::string sName)
set the section name, clear the section content
std::string SearchFile(std::string filename)
Search files in current directory and directories specified in "searchPath" section.
std::string fConfigFileName
Full name of the rml file.
TiXmlElement * GetNextElement(TiXmlElement *e)
Get the next sibling xml element of this element, with same eleDeclare.
std::string GetParameter(std::string parName, TiXmlElement *e, TString defaultValue=PARAMETER_NOT_FOUND_STR)
Returns the value for the parameter named parName in the given section.
static std::string DownloadRemoteFile(const std::string &remoteFile, bool pidPrefix=false)
download the remote file automatically, returns the downloaded file name.
static bool isURL(const std::string &filename)
Returns true if filename is an http address.
Definition: TRestTools.cxx:770
Int_t GetChar(std::string hint="Press a KEY to continue ...")
Helps to pause the program, printing a message before pausing.
Double_t StringToDouble(std::string in)
Gets a double from a string.
std::string ToUpper(std::string in)
Convert string to its upper case. Alternative of TString::ToUpper.
Int_t StringToInteger(std::string in)
Gets an integer from a string.
TVector3 StringTo3DVector(std::string in)
Gets a 3D-vector from a string. Format should be : (X,Y,Z).