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
734
735#include "TRestGeant4Metadata.h"
736
737#include <TAxis.h>
738#include <TGeoManager.h>
739#include <TRandom.h>
740#include <TRestGDMLParser.h>
741#include <TRestRun.h>
742
743#include "TRestGeant4ParticleSourceCosmics.h"
744#include "TRestGeant4PrimaryGeneratorInfo.h"
745
746using namespace std;
747using namespace TRestGeant4PrimaryGeneratorTypes;
748
749ClassImp(TRestGeant4Metadata);
750
755
770TRestGeant4Metadata::TRestGeant4Metadata(const char* configFilename, const string& name)
771 : TRestMetadata(configFilename) {
772 Initialize();
773
775}
776
781
786 SetSectionName(this->ClassName());
787 SetLibraryVersion(LIBRARY_VERSION);
788
789 fChance.clear();
790 fActiveVolumes.clear();
791 fBiasingVolumes.clear();
792
794}
795
800 this->Initialize();
801
802 fMagneticField = Get3DVectorParameterWithUnits("magneticField", TVector3(0, 0, 0));
803
804 // Initialize the metadata members from a configfile
805 fGdmlFilename = GetParameter("gdmlFile");
806 if (fGdmlFilename == PARAMETER_NOT_FOUND_STR) {
807 fGdmlFilename = GetParameter("gdml_file"); // old name
808 }
809 if (fGdmlFilename == PARAMETER_NOT_FOUND_STR) {
810 cout << "\"gdmlFile\" parameter is not defined!" << endl;
811 exit(1);
812 }
813 if (TRestTools::isURL(fGdmlFilename.Data())) {
814 fGeometryPath = "";
816 if (fGdmlFilename == "") {
817 cout << "Error downloading geometry file from URL: " << fGdmlFilename << endl;
818 exit(1);
819 }
820 }
821
822 fGeometryPath = GetParameter("geometryPath", "");
823
824 fStoreHadronicTargetInfo = StringToBool(GetParameter("storeHadronicTargetInfo", "false"));
825
826 string seedString = GetParameter("seed", "0");
827 if (ToUpper(seedString) == "RANDOM" || ToUpper(seedString) == "RAND" || ToUpper(seedString) == "AUTO" ||
828 seedString == "0") {
829 auto dd = new double();
830 fSeed = (uintptr_t)dd + (uintptr_t)this;
831 delete dd;
832 } else {
833 fSeed = (Long_t)StringToInteger(seedString);
834 }
835 gRandom->SetSeed(fSeed);
836
837 // if "gdmlFile" is purely a file (without any path) and "geometryPath" is
838 // defined, we recombine them together
839 if ((((string)fGdmlFilename).find_first_not_of("./~") == 0 ||
840 ((string)fGdmlFilename).find('/') == string::npos) &&
841 fGeometryPath != "") {
842 if (fGeometryPath[fGeometryPath.Length() - 1] == '/') {
844 } else {
845 fGdmlFilename = fGeometryPath + "/" + GetParameter("gdmlFile");
846 }
847 }
848
849 Double_t defaultTime = 1. / REST_Units::s;
850 fSubEventTimeDelay = GetDblParameterWithUnits("subEventTimeDelay", defaultTime);
851
852 auto nEventsString = GetParameter("nEvents");
853 if (nEventsString == PARAMETER_NOT_FOUND_STR) {
854 nEventsString = GetParameter("Nevents"); // old name
855 }
856 fNEvents = nEventsString == PARAMETER_NOT_FOUND_STR ? 0 : StringToInteger(nEventsString);
857
858 const auto fNRequestedEntriesString = GetParameter("nRequestedEntries");
860 fNRequestedEntriesString == PARAMETER_NOT_FOUND_STR ? 0 : StringToInteger(fNRequestedEntriesString);
861
862 fSaveAllEvents = ToUpper(GetParameter("saveAllEvents", "false")) == "TRUE" ||
863 ToUpper(GetParameter("saveAllEvents", "off")) == "ON";
864
865 fPrintProgress = ToUpper(GetParameter("printProgress", "true")) == "TRUE" ||
866 ToUpper(GetParameter("printProgress", "on")) == "ON";
867
868 fRegisterEmptyTracks = ToUpper(GetParameter("registerEmptyTracks", "false")) == "TRUE" ||
869 ToUpper(GetParameter("registerEmptyTracks", "off")) == "ON";
870
872 // Detector (old storage) section is processed after initializing geometry info in Detector Construction
873 // This allows to use regular expression to match logical or physical volumes etc.
874 ReadBiasing();
875
876 fMaxTargetStepSize = GetDblParameterWithUnits("maxTargetStepSize", -1);
877 if (fMaxTargetStepSize > 0) {
878 cout << " " << endl;
879 RESTWarning << "IMPORTANT: *maxTargetStepSize* parameter is now obsolete!" << RESTendl;
880 RESTWarning << "The sensitive volume will not define any integration step limit" << RESTendl;
881 cout << " " << endl;
882 RESTWarning << "In order to avoid this warning REMOVE the *maxTargetStepSize* definition,"
883 << RESTendl;
884 RESTWarning << "and replace it by the following statement at the <storage> section" << RESTendl;
885 cout << " " << endl;
886 RESTWarning << "<activeVolume name=\"" << this->GetSensitiveVolume() << "\" maxStepSize=\""
887 << fMaxTargetStepSize << "mm\" />" << RESTendl;
888 cout << " " << endl;
889 RESTInfo << "Now, any active volume is allowed to define a maxStepSize" << RESTendl;
890 cout << " " << endl;
891 RESTInfo << "It is also possible to define a default step size for all active volumes," << RESTendl;
892 RESTInfo << "so that in case no step is defined, the default will be used." << RESTendl;
893 cout << " " << endl;
894 RESTInfo << "<storage sensitiveVolume=\"" << this->GetSensitiveVolume() << "\" maxStepSize=\""
895 << fMaxTargetStepSize << "mm\" />" << RESTendl;
896 cout << " " << endl;
897 GetChar();
898 }
899
900 // should return success or fail
901}
902
921
923 if (TRestGeant4PrimaryGeneratorTypes::StringToSpatialGeneratorTypes(
925 TRestGeant4PrimaryGeneratorTypes::SpatialGeneratorTypes::COSMIC) {
926 RESTError
927 << "TRestGeant4Metadata::GetEquivalentSimulatedTime can only be called for 'cosmic' generator"
928 << RESTendl;
929 exit(1);
930 }
931 const auto source = GetParticleSource();
932 if (TRestGeant4PrimaryGeneratorTypes::StringToEnergyDistributionTypes(
933 source->GetEnergyDistributionType().Data()) !=
934 TRestGeant4PrimaryGeneratorTypes::EnergyDistributionTypes::FORMULA2) {
935 RESTError << "TRestGeant4Metadata::GetEquivalentSimulatedTime can only be called for 'formula2' "
936 "energy distribution"
937 << RESTendl;
938 exit(1);
939 }
940 if (TRestGeant4PrimaryGeneratorTypes::StringToAngularDistributionTypes(
941 source->GetAngularDistributionType().Data()) !=
942 TRestGeant4PrimaryGeneratorTypes::AngularDistributionTypes::FORMULA2) {
943 RESTError << "TRestGeant4Metadata::GetEquivalentSimulatedTime can only be called for 'formula' "
944 "angular distribution"
945 << RESTendl;
946 exit(1);
947 }
948
949 const auto energyRange = source->GetEnergyDistributionRange();
950 const auto angularRange = source->GetAngularDistributionRange();
951 auto function = (TF2*)source->GetEnergyAndAngularDistributionFunction()->Clone();
952 // counts per second per cm2 (distribution is already integrated over uniform phi)
953 const auto countsPerSecondPerCm2 =
954 function->Integral(energyRange.X(), energyRange.Y(), angularRange.X(), angularRange.Y(), 1E-9);
955 return countsPerSecondPerCm2;
956}
957
958Double_t TRestGeant4Metadata::GetCosmicIntensityInCountsPerSecond() const {
959 const auto countsPerSecondPerCm2 = GetCosmicFluxInCountsPerCm2PerSecond();
961 const auto countsPerSecond = countsPerSecondPerCm2 * surface;
962 return countsPerSecond;
963}
964
965Double_t TRestGeant4Metadata::GetEquivalentSimulatedTime() const {
966 const auto countsPerSecond = GetCosmicIntensityInCountsPerSecond();
967 const auto seconds = double(fNEvents) / countsPerSecond;
968 return seconds;
969}
970
971void TRestGeant4Metadata::ReadBiasing() {
972 TiXmlElement* biasingDefinition = GetElement("biasing");
973 if (biasingDefinition == nullptr) {
975 return;
976 }
977 TString biasEnabled = GetFieldValue("value", biasingDefinition);
978 TString biasType = GetFieldValue("type", biasingDefinition);
979
980 RESTDebug << "Bias: " << biasEnabled << RESTendl;
981
982 if (biasEnabled == "on" || biasEnabled == "ON" || biasEnabled == "On" || biasEnabled == "oN") {
983 RESTInfo << "Biasing is enabled" << RESTendl;
984
985 TiXmlElement* biasVolumeDefinition = GetElement("biasingVolume", biasingDefinition);
986 Int_t n = 0;
987 while (biasVolumeDefinition) {
988 TRestGeant4BiasingVolume biasVolume;
989 RESTDebug << "Def: " << biasVolumeDefinition << RESTendl;
990
991 biasVolume.SetBiasingVolumePosition(
992 Get3DVectorParameterWithUnits("position", biasVolumeDefinition));
993 biasVolume.SetBiasingFactor(StringToDouble(GetParameter("factor", biasVolumeDefinition)));
994 biasVolume.SetBiasingVolumeSize(GetDblParameterWithUnits("size", biasVolumeDefinition));
995 biasVolume.SetEnergyRange(Get2DVectorParameterWithUnits("energyRange", biasVolumeDefinition));
996 biasVolume.SetBiasingVolumeType(biasType); // For the moment all the volumes should be same type
997
998 /* TODO check that values are right if not printBiasingVolume with
999 getchar() biasVolume.PrintBiasingVolume(); getchar();
1000 */
1001
1002 fBiasingVolumes.push_back(biasVolume);
1003 biasVolumeDefinition = GetNextElement(biasVolumeDefinition);
1004 n++;
1005 }
1006 fNBiasingVolumes = n;
1007 }
1008}
1009
1035 // TODO if some fields are defined in the generator but not finally used
1036 // i.e. <generator type="volume" from="gasTarget" position="(100,0,-100)
1037 // here position is irrelevant since the event will be generated from the
1038 // volume defined in the geometry we should take care of these values so they
1039 // are not stored in metadata (to avoid future confusion) In the particular
1040 // case of position, the value is overwritten in DetectorConstruction by the
1041 // center of the volume (i.e. gasTarget) but if i.e rotation or side are
1042 // defined and not relevant we should set it to -1
1043
1044 TiXmlElement* generatorDefinition = GetElement("generator");
1045
1047 SpatialGeneratorTypesToString(StringToSpatialGeneratorTypes(GetParameter(
1048 "type", generatorDefinition, SpatialGeneratorTypesToString(SpatialGeneratorTypes::VOLUME))));
1050 SpatialGeneratorShapesToString(StringToSpatialGeneratorShapes(GetParameter(
1051 "shape", generatorDefinition, SpatialGeneratorShapesToString(SpatialGeneratorShapes::BOX))));
1052 fGeant4PrimaryGeneratorInfo.fSpatialGeneratorFrom = GetParameter("from", generatorDefinition);
1053 if (fGeant4PrimaryGeneratorInfo.fSpatialGeneratorFrom != PARAMETER_NOT_FOUND_STR) {
1055 SpatialGeneratorShapesToString(SpatialGeneratorShapes::GDML);
1056 }
1058 Get3DVectorParameterWithUnits("size", generatorDefinition, {0, 0, 0});
1060 Get3DVectorParameterWithUnits("position", generatorDefinition, {0, 0, 0});
1062 Get3DVectorParameterWithUnits("rotationAxis", generatorDefinition, {0, 0, 1});
1064 GetDblParameterWithUnits("rotationAngle", generatorDefinition, 0);
1066 GetParameter("densityFunc", generatorDefinition, "1");
1067
1068 TiXmlElement* sourceDefinition = GetElement("source", generatorDefinition);
1069 while (sourceDefinition) {
1070 string use = GetParameter("use", sourceDefinition, "");
1071
1072 TRestGeant4ParticleSource* source = TRestGeant4ParticleSource::instantiate(use);
1073 ReadParticleSource(source, sourceDefinition);
1074 AddParticleSource(source);
1075
1076 if (string(source->GetName()) == "TRestGeant4ParticleSourceCosmics") {
1077 TRestGeant4ParticleSourceCosmics::SetSeed(GetSeed());
1078 }
1079
1080 sourceDefinition = GetNextElement(sourceDefinition);
1081 }
1082
1083 // check if the generator is valid.
1084 if (GetNumberOfSources() == 0) {
1085 RESTError << "No sources are added inside generator!" << RESTendl;
1086 exit(1);
1087 }
1088}
1089
1094 TiXmlElement* sourceDefinition = definition;
1095
1096 source->SetParticleName(GetParameter("particle", sourceDefinition));
1097 source->SetExcitationLevel(StringToDouble(GetParameter("excitedLevel", sourceDefinition)));
1098 source->SetParticleCharge(StringToInteger(GetParameter("charge", sourceDefinition, "0")));
1099
1100 TString fullChain = GetParameter("fullChain", sourceDefinition);
1101 SetFullChain(false);
1102 if (fullChain.EqualTo("on", TString::ECaseCompare::kIgnoreCase)) {
1103 SetFullChain(true);
1104 }
1105
1106 // Angular distribution parameters
1107 TiXmlElement* angularDefinition = GetElement("angular", sourceDefinition);
1108 if (angularDefinition == nullptr) {
1109 angularDefinition = GetElement("angularDist", sourceDefinition); // old name
1110 }
1111 source->SetAngularDistributionType(GetParameter(
1112 "type", angularDefinition, AngularDistributionTypesToString(AngularDistributionTypes::FLUX)));
1113 if ((StringToAngularDistributionTypes(source->GetAngularDistributionType().Data()) ==
1114 AngularDistributionTypes::TH1D) ||
1115 (StringToAngularDistributionTypes(source->GetAngularDistributionType().Data()) ==
1116 AngularDistributionTypes::TH2D)) {
1117 source->SetAngularDistributionFilename(SearchFile(GetParameter("file", angularDefinition)));
1118 auto name = GetParameter("name", angularDefinition);
1119 if (name == "NO_SUCH_PARA") {
1120 name = GetParameter("spctName", angularDefinition); // old name
1121 }
1122 source->SetAngularDistributionNameInFile(name);
1123 }
1124 source->SetAngularDistributionRange(
1125 Get2DVectorParameterWithUnits("range", angularDefinition, {0, TMath::Pi()}));
1126 if (StringToAngularDistributionTypes(source->GetAngularDistributionType().Data()) ==
1127 AngularDistributionTypes::FORMULA) {
1128 source->SetAngularDistributionFormula(GetParameter("name", angularDefinition));
1129 // We cannot use an angular range bigger than the range of the formula
1130 const auto function = source->GetAngularDistributionFunction();
1131 if (source->GetAngularDistributionRangeMin() < function->GetXaxis()->GetXmin()) {
1132 source->SetAngularDistributionRange(
1133 {function->GetXaxis()->GetXmin(), source->GetAngularDistributionRangeMax()});
1134 }
1135 if (source->GetAngularDistributionRangeMax() > function->GetXaxis()->GetXmax()) {
1136 source->SetAngularDistributionRange(
1137 {source->GetAngularDistributionRangeMin(), function->GetXaxis()->GetXmax()});
1138 }
1139 }
1140 if ((StringToAngularDistributionTypes(source->GetAngularDistributionType().Data()) ==
1141 AngularDistributionTypes::FORMULA) ||
1142 (StringToAngularDistributionTypes(source->GetAngularDistributionType().Data()) ==
1143 AngularDistributionTypes::FORMULA2)) {
1144 source->SetAngularDistributionFormulaNPoints(static_cast<size_t>(GetDblParameterWithUnits(
1145 "nPoints", angularDefinition, source->GetAngularDistributionFormulaNPoints())));
1146 }
1147 if (GetNumberOfSources() == 0 &&
1148 StringToAngularDistributionTypes(source->GetAngularDistributionType().Data()) ==
1149 AngularDistributionTypes::BACK_TO_BACK) {
1150 cout << "WARNING: First source cannot be backtoback. Setting it to isotropic" << endl;
1151 source->SetAngularDistributionType(
1152 AngularDistributionTypesToString(AngularDistributionTypes::ISOTROPIC));
1153 }
1154 source->SetDirection(StringTo3DVector(GetParameter("direction", angularDefinition, "(1,0,0)")));
1155
1156 if ((StringToAngularDistributionTypes(source->GetAngularDistributionType().Data()) ==
1157 AngularDistributionTypes::ISOTROPIC)) {
1158 source->SetAngularDistributionRange({0, TMath::Pi()});
1159 const double coneHalfAngle = GetDblParameterWithUnits("coneHalfAngle", angularDefinition, 0);
1160 source->SetAngularDistributionIsotropicConeHalfAngle(coneHalfAngle);
1161 }
1162 // Energy distribution parameters
1163 TiXmlElement* energyDefinition = GetElement("energy", sourceDefinition);
1164 if (energyDefinition == nullptr) {
1165 energyDefinition = GetElement("energyDist", sourceDefinition); // old name
1166 }
1167 source->SetEnergyDistributionType(GetParameter(
1168 "type", energyDefinition, EnergyDistributionTypesToString(EnergyDistributionTypes::MONO)));
1169 if ((StringToEnergyDistributionTypes(source->GetEnergyDistributionType().Data()) ==
1170 EnergyDistributionTypes::TH1D) ||
1171 (StringToEnergyDistributionTypes(source->GetEnergyDistributionType().Data()) ==
1172 EnergyDistributionTypes::TH2D)) {
1173 source->SetEnergyDistributionFilename(SearchFile(GetParameter("file", energyDefinition)));
1174 auto name = GetParameter("name", energyDefinition);
1175 if (name == "NO_SUCH_PARA") {
1176 name = GetParameter("spctName", energyDefinition); // old name
1177 }
1178 source->SetEnergyDistributionNameInFile(name);
1179 }
1180 source->SetEnergyDistributionRange(Get2DVectorParameterWithUnits("range", energyDefinition, {0, 1E20}));
1181 if (StringToEnergyDistributionTypes(source->GetEnergyDistributionType().Data()) ==
1182 EnergyDistributionTypes::MONO) {
1183 Double_t energy;
1184 energy = GetDblParameterWithUnits("energy", energyDefinition, 0);
1185 source->SetEnergyDistributionRange({energy, energy});
1186 source->SetEnergy(energy);
1187 }
1188 if (StringToEnergyDistributionTypes(source->GetEnergyDistributionType().Data()) ==
1189 EnergyDistributionTypes::FORMULA) {
1190 source->SetEnergyDistributionFormula(GetParameter("name", energyDefinition));
1191 // We cannot use an energy range bigger than the range of the formula
1192 const auto function = source->GetEnergyDistributionFunction();
1193 if (source->GetEnergyDistributionRangeMin() < function->GetXaxis()->GetXmin()) {
1194 source->SetEnergyDistributionRange(
1195 {function->GetXaxis()->GetXmin(), source->GetEnergyDistributionRangeMax()});
1196 }
1197 if (source->GetEnergyDistributionRangeMax() > function->GetXaxis()->GetXmax()) {
1198 source->SetEnergyDistributionRange(
1199 {source->GetEnergyDistributionRangeMin(), function->GetXaxis()->GetXmax()});
1200 }
1201 }
1202 if ((StringToEnergyDistributionTypes(source->GetEnergyDistributionType().Data()) ==
1203 EnergyDistributionTypes::FORMULA) ||
1204 (StringToEnergyDistributionTypes(source->GetEnergyDistributionType().Data()) ==
1205 EnergyDistributionTypes::FORMULA2)) {
1206 source->SetEnergyDistributionFormulaNPoints(static_cast<size_t>(GetDblParameterWithUnits(
1207 "nPoints", energyDefinition, source->GetEnergyDistributionFormulaNPoints())));
1208 }
1209 if (StringToEnergyDistributionTypes(source->GetEnergyDistributionType().Data()) ==
1210 EnergyDistributionTypes::FORMULA2 &&
1211 StringToAngularDistributionTypes(source->GetAngularDistributionType().Data()) ==
1212 AngularDistributionTypes::FORMULA2) {
1213 const auto empty = TString("");
1214 auto energyDistName = GetParameter("name", energyDefinition, empty);
1215 auto angularDistName = GetParameter("name", energyDefinition, empty);
1216
1217 if (energyDistName == empty && angularDistName == empty) {
1218 RESTError << "No name specified for energy and angular distribution" << RESTendl;
1219 exit(1);
1220 }
1221 if (energyDistName == empty) {
1222 angularDistName = energyDistName;
1223 RESTWarning << "No name specified for energy distribution. Using angular distribution name: "
1224 << angularDistName << RESTendl;
1225 }
1226 if (angularDistName == empty) {
1227 energyDistName = energyDistName;
1228 RESTWarning << "No name specified for angular distribution. Using energy distribution name: "
1229 << energyDistName << RESTendl;
1230 }
1231
1232 if (energyDistName == empty || angularDistName == empty) {
1233 // we should never enter here, just leave this as a check
1234 RESTError << "When using 'formula2' the name of energy and angular dist must match" << RESTendl;
1235 exit(1);
1236 }
1237 const auto energyAndAngularDistName = energyDistName;
1238 source->SetEnergyAndAngularDistributionFormula(energyAndAngularDistName);
1239
1240 const auto function = source->GetEnergyAndAngularDistributionFunction();
1241
1242 // Set energy range
1243 if (source->GetEnergyDistributionRangeMin() < function->GetXaxis()->GetXmin()) {
1244 source->SetEnergyDistributionRange(
1245 {function->GetXaxis()->GetXmin(), source->GetEnergyDistributionRangeMax()});
1246 }
1247 if (source->GetEnergyDistributionRangeMax() > function->GetXaxis()->GetXmax()) {
1248 source->SetEnergyDistributionRange(
1249 {source->GetEnergyDistributionRangeMin(), function->GetXaxis()->GetXmax()});
1250 }
1251
1252 // Set angular range
1253 if (source->GetAngularDistributionRangeMin() < function->GetYaxis()->GetXmin()) {
1254 source->SetAngularDistributionRange(
1255 {function->GetYaxis()->GetXmin(), source->GetAngularDistributionRangeMax()});
1256 }
1257 if (source->GetAngularDistributionRangeMax() > function->GetYaxis()->GetXmax()) {
1258 source->SetAngularDistributionRange(
1259 {source->GetAngularDistributionRangeMin(), function->GetYaxis()->GetXmax()});
1260 }
1261 }
1262 // allow custom configuration from the class
1263 source->LoadConfigFromElement(sourceDefinition, fElementGlobal, fVariables);
1264 // AddParticleSource(source);
1265}
1266
1268 for (auto source : fParticleSource) {
1269 delete source;
1270 }
1271 fParticleSource.clear();
1272}
1273
1275 fParticleSource.push_back(src);
1276}
1277
1293 RESTInfo << "TRestGeant4Metadata: Processing detector section" << RESTendl;
1294 TiXmlElement* detectorDefinition = GetElement("detector");
1295 if (detectorDefinition == nullptr) {
1296 RESTError << "Detector section (<detector>) is not present in metadata!" << RESTendl;
1297 exit(1);
1298 }
1299
1300 const bool activateAllVolumes =
1301 StringToBool(GetParameter("activateAllVolumes", detectorDefinition, "true"));
1302 RESTInfo << "TRestGeant4Metadata: Setting 'fActivateAllVolumes' to " << fActivateAllVolumes << RESTendl;
1303 fActivateAllVolumes = activateAllVolumes;
1304
1305 TiXmlElement* removeUnwantedTracksDefinition = GetElement("removeUnwantedTracks", detectorDefinition);
1306 if (removeUnwantedTracksDefinition != nullptr) {
1307 fRemoveUnwantedTracks = StringToBool(GetParameter("enabled", removeUnwantedTracksDefinition, "true"));
1309 StringToBool(GetParameter("keepZeroEnergyTracks", removeUnwantedTracksDefinition, "false"));
1310 RESTInfo << "TRestGeant4Metadata: Setting 'fRemoveUnwantedTracks' to " << fRemoveUnwantedTracks
1311 << " with 'keepZeroEnergyTracks' option set to " << fRemoveUnwantedTracksKeepZeroEnergyTracks
1312 << RESTendl;
1313 }
1314
1315 Double_t defaultStep = GetDblParameterWithUnits("maxStepSize", detectorDefinition);
1316 if (defaultStep < 0) {
1317 defaultStep = 0;
1318 }
1319
1320 TiXmlElement* volumeDefinition = GetElement("volume", detectorDefinition);
1321 while (volumeDefinition != nullptr) {
1322 string name = GetFieldValue("name", volumeDefinition);
1323 if (name.empty() || name == "Not defined") {
1324 RESTError << "volume inside detector section defined without name" << RESTendl;
1325 exit(1);
1326 }
1327 vector<string> physicalVolumes;
1328 if (!fGeant4GeometryInfo.IsValidGdmlName(name)) {
1329 if (fGeant4GeometryInfo.IsValidLogicalVolume(name)) {
1330 for (const auto& physical : fGeant4GeometryInfo.GetAllPhysicalVolumesFromLogical(name)) {
1331 physicalVolumes.emplace_back(
1332 fGeant4GeometryInfo.GetAlternativeNameFromGeant4PhysicalName(physical));
1333 }
1334 }
1335 // does not match a explicit (gdml) physical name or a logical name, perhaps its a regular exp
1336 if (physicalVolumes.empty()) {
1337 for (const auto& physical :
1338 fGeant4GeometryInfo.GetAllPhysicalVolumesMatchingExpression(name)) {
1339 physicalVolumes.emplace_back(
1340 fGeant4GeometryInfo.GetAlternativeNameFromGeant4PhysicalName(physical));
1341 }
1342 }
1343 if (physicalVolumes.empty()) {
1344 for (const auto& logical : fGeant4GeometryInfo.GetAllLogicalVolumesMatchingExpression(name)) {
1345 for (const auto& physical :
1346 fGeant4GeometryInfo.GetAllPhysicalVolumesFromLogical(logical)) {
1347 physicalVolumes.emplace_back(
1348 fGeant4GeometryInfo.GetAlternativeNameFromGeant4PhysicalName(physical));
1349 }
1350 }
1351 }
1352 } else {
1353 physicalVolumes.push_back(name);
1354 }
1355
1356 if (physicalVolumes.empty()) {
1357 RESTError
1358 << "volume '" << name
1359 << "' is not defined in the geometry. Could not match to a physical, logical volume or regex"
1360 << RESTendl;
1361 exit(1);
1362 }
1363
1364 bool isSensitive = false;
1365 const string isSensitiveValue = GetFieldValue("sensitive", volumeDefinition);
1366 if (isSensitiveValue != "Not defined") {
1367 isSensitive = StringToBool(isSensitiveValue);
1368 }
1369
1370 bool isKeepTracks = isSensitive;
1371 const string isKeepTracksValue = GetFieldValue("keepTracks", volumeDefinition);
1372 if (isKeepTracksValue != "Not defined") {
1373 isKeepTracks = StringToBool(isKeepTracksValue);
1374 }
1375
1376 bool isKill = false;
1377 const string isKillValue = GetFieldValue("kill", volumeDefinition);
1378 if (isKillValue != "Not defined") {
1379 isKill = StringToBool(isKillValue);
1380 }
1381
1382 Double_t chance = StringToDouble(GetFieldValue("chance", volumeDefinition));
1383 if (chance == -1 || chance < 0) {
1384 chance = 1.0;
1385 }
1386 Double_t maxStep = GetDblParameterWithUnits("maxStepSize", volumeDefinition);
1387 if (maxStep < 0) {
1388 maxStep = defaultStep;
1389 }
1390
1391 for (const auto& physical : physicalVolumes) {
1392 RESTInfo << "Adding " << (isSensitive ? "sensitive" : "active") << " volume from RML: '"
1393 << physical << (chance != 1 ? " with change: " + to_string(chance) : " ") << RESTendl;
1394 SetActiveVolume(physical, chance, maxStep);
1395 if (isSensitive) {
1396 InsertSensitiveVolume(physical);
1397 }
1398 if (fRemoveUnwantedTracks && isKeepTracks) {
1399 fRemoveUnwantedTracksVolumesToKeep.insert(physical);
1400 }
1401 if (isKill) {
1402 fKillVolumes.insert(physical);
1403 }
1404 }
1405 volumeDefinition = GetNextElement(volumeDefinition);
1406 }
1407
1408 if (fSensitiveVolumes.empty()) {
1409 cout << "No sensitive volumes defined in detector section. Adding 'gas' as sensitive volume" << endl;
1410 InsertSensitiveVolume("gas");
1411 }
1412
1413 fEnergyRangeStored = Get2DVectorParameterWithUnits("energyRange", detectorDefinition, fEnergyRangeStored);
1414 // TODO: Place this as a generic method
1415 if (fEnergyRangeStored.X() < 0) {
1417 }
1418 if (fEnergyRangeStored.Y() < 0) {
1420 }
1421 if (fEnergyRangeStored.Y() > fEnergyRangeStored.Y()) {
1423 }
1424 if (fEnergyRangeStored.Y() <= 0) {
1425 RESTError << "Energy range is not valid: (" << fEnergyRangeStored.X() << ", "
1426 << fEnergyRangeStored.Y() << ") keV" << RESTendl;
1427 exit(1);
1428 }
1429
1430 auto gdml = new TRestGDMLParser();
1431 gdml->Load(GetGdmlFilename().Data());
1432
1433 fGeant4GeometryInfo.PopulateFromGdml(gdml->GetOutputGDMLFile());
1434
1435 // If the user did not add explicitly any volume to the storage section we understand
1436 // the user wants to register all the volumes
1437 if (fActivateAllVolumes) {
1438 for (const auto& name : fGeant4GeometryInfo.GetAllPhysicalVolumes()) {
1439 if (!IsActiveVolume(name)) {
1440 SetActiveVolume(name, 1, defaultStep);
1441 RESTInfo << "Automatically adding active volume: '" << name << "' with chance: " << 1
1442 << RESTendl;
1443 }
1444 }
1445 }
1446
1447 if (GetNumberOfActiveVolumes() == 0) {
1448 RESTError << "No active volumes defined. Please check the detector section" << RESTendl;
1449 exit(1);
1450 }
1451}
1452
1459
1460 RESTMetadata << "Geant4 version: " << GetGeant4Version() << RESTendl;
1461 RESTMetadata << "Random seed: " << GetSeed() << RESTendl;
1462 RESTMetadata << "GDML geometry: " << GetGdmlReference() << RESTendl;
1463 RESTMetadata << "GDML materials reference: " << GetMaterialsReference() << RESTendl;
1464 RESTMetadata << "Sub-event time delay: " << GetSubEventTimeDelay() << " us" << RESTendl;
1465 Double_t mx = GetMagneticField().X();
1466 Double_t my = GetMagneticField().Y();
1467 Double_t mz = GetMagneticField().Z();
1468 RESTMetadata << "Magnetic field: (" << mx << ", " << my << ", " << mz << ") T" << RESTendl;
1469 if (fSaveAllEvents) RESTMetadata << "Save all events was enabled!" << RESTendl;
1471 RESTMetadata << "Register empty tracks was enabled" << RESTendl;
1472 } else {
1473 RESTMetadata << "Register empty tracks was NOT enabled" << RESTendl;
1474 }
1475
1476 RESTMetadata << " " << RESTendl;
1477 RESTMetadata << " ++++++++++ Generator +++++++++++ " << RESTendl;
1478 RESTMetadata << "Number of generated events: " << GetNumberOfEvents() << RESTendl;
1480
1481 for (int i = 0; i < GetNumberOfSources(); i++) GetParticleSource(i)->PrintMetadata();
1482
1483 RESTMetadata << " " << RESTendl;
1484 RESTMetadata << " ++++++++++ Detector +++++++++++ " << RESTendl;
1485 RESTMetadata << "Energy range (keV): (" << GetMinimumEnergyStored() << ", " << GetMaximumEnergyStored()
1486 << ")" << RESTendl;
1487 RESTMetadata << "Number of sensitive volumes: " << GetNumberOfSensitiveVolumes() << RESTendl;
1488 for (const auto& sensitiveVolume : fSensitiveVolumes) {
1489 RESTMetadata << "Sensitive volume: " << sensitiveVolume << RESTendl;
1490 }
1491 RESTMetadata << "Number of active volumes: " << GetNumberOfActiveVolumes() << RESTendl;
1492 for (unsigned int n = 0; n < GetNumberOfActiveVolumes(); n++) {
1493 const auto name = GetActiveVolumeName(n);
1494 RESTMetadata << "Name: " << name << ", ID: " << fGeant4GeometryInfo.GetIDFromVolume(name)
1495 << ", maxStep: " << GetMaxStepSize(name) << "mm "
1496 << ", chance: " << GetStorageChance(name) << RESTendl;
1497 }
1498
1499 for (unsigned int n = 0; n < GetNumberOfBiasingVolumes(); n++) {
1500 GetBiasingVolume(n).PrintBiasingVolume();
1501 }
1502
1503 if (GetRemoveUnwantedTracks()) {
1504 RESTMetadata << "removeUnwantedTracks is enabled "
1505 << (fRemoveUnwantedTracksKeepZeroEnergyTracks ? "with" : "without")
1506 << " keeping zero energy tracks" << RESTendl;
1507
1508 /*
1509 RESTMetadata << "keep volumes for removeUnwantedTracks:" << RESTendl;
1510 for (const auto& volume : fRemoveUnwantedTracksVolumesToKeep) {
1511 RESTMetadata << "\t" << volume << RESTendl;
1512 }
1513 */
1514 }
1515
1516 RESTMetadata << "+++++" << RESTendl;
1517}
1518
1521Int_t TRestGeant4Metadata::GetActiveVolumeID(const TString& name) {
1522 Int_t id;
1523 for (id = 0; id < (Int_t)fActiveVolumes.size(); id++) {
1524 if (fActiveVolumes[id] == name) return id;
1525 }
1526 return -1;
1527}
1528
1543void TRestGeant4Metadata::SetActiveVolume(const TString& name, Double_t chance, Double_t maxStep) {
1544 for (size_t i = 0; i < fActiveVolumes.size(); i++) {
1545 const auto activeVolumeName = fActiveVolumes[i];
1546 if (name == activeVolumeName) {
1547 fChance[i] = chance;
1548 fMaxStepSize[i] = maxStep;
1549 return;
1550 }
1551 }
1552 fActiveVolumes.push_back(name);
1553 fChance.push_back(chance);
1554 fMaxStepSize.push_back(maxStep);
1555 fActiveVolumesSet.insert(name.Data());
1556}
1557
1562Bool_t TRestGeant4Metadata::isVolumeStored(const TString& volume) const {
1563 for (unsigned int n = 0; n < GetNumberOfActiveVolumes(); n++)
1564 if (GetActiveVolumeName(n) == volume) return true;
1565
1566 return false;
1567}
1568
1572Double_t TRestGeant4Metadata::GetStorageChance(const TString& volume) {
1573 Int_t id;
1574 for (id = 0; id < (Int_t)fActiveVolumes.size(); id++) {
1575 if (fActiveVolumes[id] == volume) return fChance[id];
1576 }
1577 RESTWarning << "TRestGeant4Metadata::GetStorageChance. Volume " << volume << " not found" << RESTendl;
1578
1579 return 0;
1580}
1581
1585Double_t TRestGeant4Metadata::GetMaxStepSize(const TString& volume) {
1586 for (Int_t i = 0; i < (Int_t)fActiveVolumes.size(); i++) {
1587 if (volume.EqualTo(fActiveVolumes[i], TString::kIgnoreCase)) {
1588 return fMaxStepSize[i];
1589 }
1590 }
1591 RESTWarning << "TRestGeant4Metadata::GetMaxStepSize. Volume " << volume << " not found" << RESTendl;
1592
1593 return 0;
1594}
1595
1596size_t TRestGeant4Metadata::GetGeant4VersionMajor() const {
1597 TString majorVersion = "";
1598 for (int i = 0; i < fGeant4Version.Length(); i++) {
1599 auto c = fGeant4Version[i];
1600 if (c == '.') {
1601 break;
1602 }
1603 majorVersion += c;
1604 }
1605 return std::stoi(majorVersion.Data());
1606}
1607
1608void TRestGeant4Metadata::Merge(const TRestGeant4Metadata& metadata) {
1609 fIsMerge = true;
1610 fSeed = 0; // seed makes no sense in a merged file
1611
1612 fNEvents += metadata.fNEvents;
1614 fSimulationTime += metadata.fSimulationTime;
1615}
1616
1617TRestGeant4Metadata::TRestGeant4Metadata(const TRestGeant4Metadata& metadata) { *this = metadata; }
1618
1619TRestGeant4Metadata& TRestGeant4Metadata::operator=(const TRestGeant4Metadata& metadata) {
1620 fIsMerge = metadata.fIsMerge;
1624 fGeant4Version = metadata.fGeant4Version;
1625 fGdmlReference = metadata.fGdmlReference;
1628 fActiveVolumes = metadata.fActiveVolumes;
1629 fChance = metadata.fChance;
1630 fMaxStepSize = metadata.fMaxStepSize;
1631 // fParticleSource = metadata.fParticleSource; // segfaults (pointers!)
1636 fFullChain = metadata.fFullChain;
1638 fNEvents = metadata.fNEvents;
1641 fSeed = metadata.fSeed;
1642 fSaveAllEvents = metadata.fSaveAllEvents;
1646 fKillVolumes = metadata.fKillVolumes;
1648 fMagneticField = metadata.fMagneticField;
1649 return *this;
1650}
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.
Double_t GetMaxStepSize(const TString &volume)
Returns the maximum step at a particular active volume.
Bool_t fStoreHadronicTargetInfo
Option to store target isotope information on hadronic processes (disabled by default)
~TRestGeant4Metadata()
Default destructor.
std::vector< Double_t > fMaxStepSize
A std::vector to store the maximum step size at a particular volume.
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...
Double_t fSimulationMaxTimeSeconds
Time before simulation is ended and saved.
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.
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.
Double_t fSimulationTime
The time, in seconds, that the simulation took to complete.
TRestGeant4ParticleSource * GetParticleSource(size_t n=0) const
Returns the name of the particle source with index n (Geant4 based names).
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:70
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).