REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestGeant4Metadata.h
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
23#ifndef RestCore_TRestGeant4Metadata
24#define RestCore_TRestGeant4Metadata
25
26#include <TMath.h>
27#include <TRestMetadata.h>
28#include <TString.h>
29#include <TVector2.h>
30#include <TVector3.h>
31#include <stdio.h>
32#include <stdlib.h>
33
34#include <cstdio>
35#include <cstring>
36#include <fstream>
37#include <iostream>
38#include <vector>
39
40#include "TRestGeant4BiasingVolume.h"
41#include "TRestGeant4GeometryInfo.h"
42#include "TRestGeant4ParticleSource.h"
43#include "TRestGeant4PhysicsInfo.h"
44#include "TRestGeant4PrimaryGeneratorInfo.h"
45
48 private:
49 void Initialize() override;
50
51 void InitFromConfigFile() override;
52
53 void ReadGenerator();
54 void ReadParticleSource(TRestGeant4ParticleSource* src, TiXmlElement* sourceDefinition);
55
56 void ReadDetector();
57 void ReadBiasing();
58
59 // Metadata is the result of a merge of other metadata
60 bool fIsMerge = false;
61
64
67
70
73
75 TString fGeometryPath;
76
78 TString fGdmlFilename;
79
82
85
88 TVector2 fEnergyRangeStored = {0, 1E20};
89
91 std::vector<TString> fActiveVolumes;
92
94 std::vector<Double_t> fChance;
95
97 std::vector<Double_t> fMaxStepSize;
98
100 std::vector<TRestGeant4ParticleSource*> fParticleSource; //->
101
104
106 std::vector<TRestGeant4BiasingVolume> fBiasingVolumes;
107
110 Double_t fMaxTargetStepSize = 0;
111
114 Double_t fSubEventTimeDelay = 100;
115
118 Bool_t fFullChain = true;
119
121 std::set<std::string> fFullChainStopIsotopes;
122
125 std::vector<TString> fSensitiveVolumes;
126
128 Long64_t fNEvents = 0;
129
131 Long64_t fNRequestedEntries = 0;
132
135
137 Double_t fSimulationTime = 0;
138
141 Long_t fSeed = 0;
142
145 Bool_t fSaveAllEvents = false;
146
149
151 Bool_t fActivateAllVolumes = false;
152
154 Bool_t fRemoveUnwantedTracks = false;
155
157 Bool_t fStoreTracks = false;
158
162
165
166 std::set<std::string> fKillVolumes;
167
170 Bool_t fPrintProgress = false;
171
175 Bool_t fRegisterEmptyTracks = true;
176
178 TVector3 fMagneticField = TVector3(0, 0, 0);
179
180 public:
181 std::set<std::string> fActiveVolumesSet = {};
182
184 inline Long_t GetSeed() const { return fSeed; }
185
188
191
195 }
196
198 inline TString GetGeant4Version() const { return fGeant4Version; }
199
200 size_t GetGeant4VersionMajor() const;
201
202 bool GetStoreHadronicTargetInfo() const { return fStoreHadronicTargetInfo; }
203
205 inline TString GetGeometryPath() const { return fGeometryPath; }
206
208 inline TString GetGdmlFilename() const { return fGdmlFilename; }
209
211 inline TString GetGdmlReference() const { return fGdmlReference; }
212
214 inline TString GetMaterialsReference() const { return fMaterialsReference; }
215
217 inline Bool_t isFullChainActivated() const { return fFullChain; }
218
220 inline std::set<std::string> GetFullChainStopIsotopes() const { return fFullChainStopIsotopes; }
221
222 inline Bool_t IsIsotopeFullChainStop(const std::string& isotope) const {
223 return fFullChainStopIsotopes.count(isotope) > 0;
224 }
225
227 inline Double_t GetMaxTargetStepSize() const { return fMaxTargetStepSize; }
228
233 inline Double_t GetSubEventTimeDelay() const { return fSubEventTimeDelay; }
234
236 inline Bool_t GetSaveAllEvents() const { return fSaveAllEvents; }
237
239 inline Bool_t PrintProgress() const { return fPrintProgress; }
240
242 inline Bool_t RegisterEmptyTracks() const { return fRegisterEmptyTracks; }
243
245 inline void SetSeed(Long_t seed) { fSeed = seed; }
246
248 inline void SetSaveAllEvents(const Bool_t value) { fSaveAllEvents = value; }
249
251 inline void SetGeant4Version(const TString& versionString) { fGeant4Version = versionString; }
252
254 inline void SetFullChain(Bool_t fullChain) { fFullChain = fullChain; }
255
257 inline void SetGeometryPath(const std::string& path) { fGeometryPath = path; }
258
260 inline void SetGdmlFilename(const std::string& gdmlFile) { fGdmlFilename = gdmlFile; }
261
263 inline void SetGdmlReference(const std::string& reference) { fGdmlReference = reference; }
264
266 inline void SetMaterialsReference(const std::string& reference) { fMaterialsReference = reference; }
267
269 inline Long64_t GetNumberOfEvents() const { return fNEvents; }
270
271 inline Long64_t GetNumberOfRequestedEntries() const { return fNRequestedEntries; }
272
273 inline Double_t GetSimulationMaxTimeSeconds() const { return fSimulationMaxTimeSeconds; }
274
275 // Direct access to sources definition in primary generator
278 inline Int_t GetNumberOfSources() const { return fParticleSource.size(); }
279
281 inline TRestGeant4ParticleSource* GetParticleSource(size_t n = 0) const { return fParticleSource[n]; }
282
285
289
290 // Direct access to biasing volumes definition
293
294 inline size_t GetNumberOfBiasingVolumes() const { return fBiasingVolumes.size(); }
295
298
300 inline Int_t isBiasingActive() const { return fBiasingVolumes.size(); }
301
303 inline TString GetSensitiveVolume(int n = 0) const { return fSensitiveVolumes[n]; }
304
305 inline size_t GetNumberOfSensitiveVolumes() const { return fSensitiveVolumes.size(); }
306
307 inline const std::vector<TString>& GetSensitiveVolumes() const { return fSensitiveVolumes; }
308
310 inline void SetNumberOfEvents(Long64_t n) { fNEvents = n; }
311
312 inline void SetNumberOfRequestedEntries(Long64_t n) { fNRequestedEntries = n; }
313
314 inline void SetSimulationMaxTimeSeconds(Double_t seconds) { fSimulationMaxTimeSeconds = seconds; }
315
317 inline void InsertSensitiveVolume(const TString& volume) {
318 for (const auto& sensitiveVolume : fSensitiveVolumes) {
319 // Do not add duplicate volumes
320 if (volume == sensitiveVolume) {
321 return;
322 }
323 }
324 fSensitiveVolumes.push_back(volume);
325 }
326
329 inline Double_t GetStorageChance(Int_t n) const { return fChance[n]; }
330
333 Double_t GetStorageChance(const TString& volume);
334
335 Double_t GetMaxStepSize(const TString& volume);
336
338 inline Double_t GetMinimumEnergyStored() const { return fEnergyRangeStored.X(); }
339
341 inline Double_t GetMaximumEnergyStored() const { return fEnergyRangeStored.Y(); }
342
345 inline unsigned int GetNumberOfActiveVolumes() const { return fActiveVolumes.size(); }
346
347 inline bool IsActiveVolume(const char* volumeName) const {
348 return fActiveVolumesSet.count(volumeName) > 0;
349 }
350
351 inline bool IsKeepTracksVolume(const char* volumeName) const {
352 return fRemoveUnwantedTracksVolumesToKeep.count(volumeName) > 0;
353 }
354
355 inline bool IsKillVolume(const char* volumeName) const { return fKillVolumes.count(volumeName) > 0; }
356
357 inline std::vector<std::string> GetKillVolumes() const {
358 std::vector<std::string> result;
359 for (const auto& volume : fKillVolumes) {
360 result.emplace_back(volume);
361 }
362 return result;
363 }
364
365 inline std::vector<std::string> GetRemoveUnwantedTracksVolumesToKeep() const {
366 std::vector<std::string> result;
367 for (const auto& volume : fRemoveUnwantedTracksVolumesToKeep) {
368 result.emplace_back(volume);
369 }
370 return result;
371 }
372
373 double GetGeneratorSurfaceCm2() const;
374
376 Double_t GetCosmicIntensityInCountsPerSecond() const;
377
379 Double_t GetEquivalentSimulatedTime() const;
380
382 inline Double_t GetSimulationWallTime() const { return fSimulationTime; }
383
385 inline TString GetActiveVolumeName(Int_t n) const { return fActiveVolumes[n]; }
386
387 inline std::vector<TString> GetActiveVolumes() const { return fActiveVolumes; }
388
389 inline bool GetRemoveUnwantedTracks() const { return fRemoveUnwantedTracks; }
390
391 bool GetStoreTracks() const { return fStoreTracks; }
392
393 inline bool GetRemoveUnwantedTracksKeepZeroEnergyTracks() const {
395 }
396
398 inline TVector3 GetMagneticField() const { return fMagneticField; }
399
400 Int_t GetActiveVolumeID(const TString& name);
401
402 Bool_t isVolumeStored(const TString& volume) const;
403
404 void SetActiveVolume(const TString& name, Double_t chance, Double_t maxStep = 0);
405
406 void SetSimulationWallTime(Double_t time) { fSimulationTime = time; }
407
408 void PrintMetadata() override;
409
410 void Merge(const TRestGeant4Metadata&);
411
413 TRestGeant4Metadata(const char* configFilename, const std::string& name = "");
414
416
418 TRestGeant4Metadata& operator=(const TRestGeant4Metadata& metadata);
419
420 ClassDefOverride(TRestGeant4Metadata, 19);
421
422 // Allow modification of otherwise inaccessible / immutable members that shouldn't be modified by the user
423 friend class SteppingAction;
424 friend class DetectorConstruction;
425 friend class TRestGeant4Hits;
426};
427#endif // RestCore_TRestGeant4Metadata
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.
const TRestGeant4PhysicsInfo & GetGeant4PhysicsInfo() const
Returns an immutable reference to the physics info.
Bool_t fStoreTracks
Store event tracks (default is true)
void ReadDetector()
Reads the storage section defined inside TRestGeant4Metadata.
Bool_t fActivateAllVolumes
Sets all volume as active without having to explicitly list them.
Bool_t GetSaveAllEvents() const
It returns true if save all events is active.
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 SetSeed(Long_t seed)
Used exclusively by restG4 to set the value of the random seed used on Geant4 simulation.
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.
Bool_t PrintProgress() const
It returns true if printProgress parameter was set to true.
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.
void SetGdmlReference(const std::string &reference)
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.
const TRestGeant4PrimaryGeneratorInfo & GetGeant4PrimaryGeneratorInfo() const
Returns an immutable reference to the primary generator info.
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.
Int_t isBiasingActive() const
Returns the number of biasing volumes defined. If 0 the biasing technique is not being used.
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.
Bool_t RegisterEmptyTracks() const
It returns false if registerEmptyTracks parameter was set to false.
std::set< std::string > GetFullChainStopIsotopes() const
Returns the isotopes that will stop the full chain decay simulation.
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.
void SetNumberOfEvents(Long64_t n)
Sets the name of the sensitive volume.
Double_t GetMaxTargetStepSize() const
Returns the value of the maximum Geant4 step size in mm for the target volume.
Long64_t fNEvents
The number of events simulated, or to be simulated.
void SetGdmlFilename(const std::string &gdmlFile)
It sets the main filename to be used for the GDML geometry.
Long_t fSeed
The seed value used for Geant4 random event generator. If it is zero its value will be assigned using...
const TRestGeant4GeometryInfo & GetGeant4GeometryInfo() const
Returns an immutable reference to the geometry info.
Double_t fSimulationMaxTimeSeconds
Time before simulation is ended and saved.
void SetFullChain(Bool_t fullChain)
Enables/disables the full chain decay generation.
Bool_t isFullChainActivated() const
Returns true in case full decay chain simulation is enabled.
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.
std::set< std::string > fFullChainStopIsotopes
If defined, it will stop the full chain decay simulation when one of these isotope appears.
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 GetSimulationWallTime() const
Returns the total time of the simulation in seconds (wall time)
TString GetGeometryPath() const
Returns the local path to the GDML geometry.
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.
void SetSaveAllEvents(const Bool_t value)
Enables or disables the save all events feature.
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...
void SetMaterialsReference(const std::string &reference)
Returns the reference provided at the materials file header.
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...
void SetGeant4Version(const TString &versionString)
Sets the value of the Geant4 version.
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 SetGeometryPath(const std::string &path)
It sets the location of the geometry files.
Double_t fSimulationTime
The time, in seconds, that the simulation took to complete (wall time)
TRestGeant4ParticleSource * GetParticleSource(size_t n=0) const
Returns the name of the particle source with index n (Geant4 based names).
Double_t GetEquivalentSimulatedTime() const
Returns the equivalent simulated time in seconds (physical time simulated)
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.
A base class for any REST metadata class.
Definition: TRestMetadata.h:74