REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
Loading...
Searching...
No Matches
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
122 std::vector<TString> fSensitiveVolumes;
123
125 Long64_t fNEvents = 0;
126
128 Long64_t fNRequestedEntries = 0;
129
132
134 Long64_t fSimulationTime = 0;
135
138 Long_t fSeed = 0;
139
142 Bool_t fSaveAllEvents = false;
143
145 Bool_t fActivateAllVolumes = false;
146
148 Bool_t fRemoveUnwantedTracks = false;
149
153
156
157 std::set<std::string> fKillVolumes;
158
161 Bool_t fPrintProgress = false;
162
166 Bool_t fRegisterEmptyTracks = true;
167
169 TVector3 fMagneticField = TVector3(0, 0, 0);
170
171 public:
172 std::set<std::string> fActiveVolumesSet = {};
173
175 inline Long_t GetSeed() const { return fSeed; }
176
179
182
187
189 inline TString GetGeant4Version() const { return fGeant4Version; }
190
191 size_t GetGeant4VersionMajor() const;
192
194 inline TString GetGeometryPath() const { return fGeometryPath; }
195
197 inline TString GetGdmlFilename() const { return fGdmlFilename; }
198
200 inline TString GetGdmlReference() const { return fGdmlReference; }
201
203 inline TString GetMaterialsReference() const { return fMaterialsReference; }
204
206 inline Bool_t isFullChainActivated() const { return fFullChain; }
207
210 inline Double_t GetMaxTargetStepSize() const { return fMaxTargetStepSize; }
211
216 inline Double_t GetSubEventTimeDelay() const { return fSubEventTimeDelay; }
217
219 inline Bool_t GetSaveAllEvents() const { return fSaveAllEvents; }
220
222 inline Bool_t PrintProgress() const { return fPrintProgress; }
223
225 inline Bool_t RegisterEmptyTracks() const { return fRegisterEmptyTracks; }
226
228 inline void SetSeed(Long_t seed) { fSeed = seed; }
229
231 inline void SetSaveAllEvents(const Bool_t value) { fSaveAllEvents = value; }
232
234 inline void SetGeant4Version(const TString& versionString) { fGeant4Version = versionString; }
235
237 inline void SetFullChain(Bool_t fullChain) { fFullChain = fullChain; }
238
240 inline void SetGeometryPath(std::string path) { fGeometryPath = path; }
241
243 inline void SetGdmlFilename(std::string gdmlFile) { fGdmlFilename = gdmlFile; }
244
246 inline void SetGdmlReference(std::string reference) { fGdmlReference = reference; }
247
249 inline void SetMaterialsReference(std::string reference) { fMaterialsReference = reference; }
250
252 inline Long64_t GetNumberOfEvents() const { return fNEvents; }
253
254 inline Long64_t GetNumberOfRequestedEntries() const { return fNRequestedEntries; }
255
256 inline Int_t GetSimulationMaxTimeSeconds() const { return fSimulationMaxTimeSeconds; }
257
258 inline Long64_t GetSimulationTime() const { return fSimulationTime; }
259
260 // Direct access to sources definition in primary generator
263 inline Int_t GetNumberOfSources() const { return fParticleSource.size(); }
264
266 inline TRestGeant4ParticleSource* GetParticleSource(size_t n = 0) const { return fParticleSource[n]; }
267
270
274
275 // Direct access to biasing volumes definition
278
279 inline size_t GetNumberOfBiasingVolumes() const { return fBiasingVolumes.size(); }
280
283
285 inline Int_t isBiasingActive() const { return fBiasingVolumes.size(); }
286
288 inline TString GetSensitiveVolume(int n = 0) const { return fSensitiveVolumes[n]; }
289
290 inline size_t GetNumberOfSensitiveVolumes() const { return fSensitiveVolumes.size(); }
291
292 inline const std::vector<TString>& GetSensitiveVolumes() const { return fSensitiveVolumes; }
293
295 inline void SetNumberOfEvents(Long64_t n) { fNEvents = n; }
296
297 inline void SetNumberOfRequestedEntries(Long64_t n) { fNRequestedEntries = n; }
298
299 inline void SetSimulationMaxTimeSeconds(Int_t seconds) { fSimulationMaxTimeSeconds = seconds; }
300
302 inline void InsertSensitiveVolume(const TString& volume) {
303 for (const auto& sensitiveVolume : fSensitiveVolumes) {
304 // Do not add duplicate volumes
305 if (volume == sensitiveVolume) {
306 return;
307 }
308 }
309 fSensitiveVolumes.push_back(volume);
310 }
311
314 inline Double_t GetStorageChance(Int_t n) const { return fChance[n]; }
315
318 Double_t GetStorageChance(const TString& volume);
319
320 Double_t GetMaxStepSize(const TString& volume);
321
323 inline Double_t GetMinimumEnergyStored() const { return fEnergyRangeStored.X(); }
324
326 inline Double_t GetMaximumEnergyStored() const { return fEnergyRangeStored.Y(); }
327
330 inline unsigned int GetNumberOfActiveVolumes() const { return fActiveVolumes.size(); }
331
332 inline bool IsActiveVolume(const char* volumeName) const {
333 return fActiveVolumesSet.count(volumeName) > 0;
334 }
335
336 inline bool IsKeepTracksVolume(const char* volumeName) const {
337 return fRemoveUnwantedTracksVolumesToKeep.count(volumeName) > 0;
338 }
339
340 inline bool IsKillVolume(const char* volumeName) const { return fKillVolumes.count(volumeName) > 0; }
341
342 inline std::vector<std::string> GetKillVolumes() const {
343 std::vector<std::string> result;
344 for (const auto& volume : fKillVolumes) {
345 result.emplace_back(volume);
346 }
347 return result;
348 }
349
350 inline std::vector<std::string> GetRemoveUnwantedTracksVolumesToKeep() const {
351 std::vector<std::string> result;
352 for (const auto& volume : fRemoveUnwantedTracksVolumesToKeep) {
353 result.emplace_back(volume);
354 }
355 return result;
356 }
357
359 Double_t GetCosmicIntensityInCountsPerSecond() const;
360 Double_t GetEquivalentSimulatedTime() const;
361
363 inline TString GetActiveVolumeName(Int_t n) const { return fActiveVolumes[n]; }
364
365 inline std::vector<TString> GetActiveVolumes() const { return fActiveVolumes; }
366
367 inline bool GetRemoveUnwantedTracks() const { return fRemoveUnwantedTracks; }
368
369 inline bool GetRemoveUnwantedTracksKeepZeroEnergyTracks() const {
371 }
372
374 inline TVector3 GetMagneticField() const { return fMagneticField; }
375
376 Int_t GetActiveVolumeID(const TString& name);
377
378 Bool_t isVolumeStored(const TString& volume) const;
379
380 void SetActiveVolume(const TString& name, Double_t chance, Double_t maxStep = 0);
381
382 void SetSimulationTime(Long64_t time) { fSimulationTime = time; }
383
384 void PrintMetadata() override;
385
386 void Merge(const TRestGeant4Metadata&);
387
389 TRestGeant4Metadata(const char* configFilename, const std::string& name = "");
390
392
394 TRestGeant4Metadata& operator=(const TRestGeant4Metadata& metadata);
395
396 ClassDefOverride(TRestGeant4Metadata, 15);
397
398 // Allow modification of otherwise inaccessible / immutable members that shouldn't be modified by the user
399 friend class SteppingAction;
400 friend class DetectorConstruction;
401 friend class TRestGeant4Hits;
402};
403#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.
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.
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.
void SetGdmlReference(std::string reference)
Returns the reference provided at the GDML file header.
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.
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.
Int_t isBiasingActive() const
Returns the number of biasing volumes defined. If 0 the biasing technique is not being used.
~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.
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.
void SetGdmlFilename(std::string gdmlFile)
It sets the main filename to be used for the GDML geometry.
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.
void SetGeometryPath(std::string path)
It sets the location of the geometry files.
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.
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.
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.
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...
TString GetGeometryPath() const
Returns the local path to the GDML geometry.
void SetMaterialsReference(std::string reference)
Returns the reference provided at the materials file header.
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...
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.
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.
A base class for any REST metadata class.