REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestGeant4ParticleSource.cxx
1 
18 #include "TRestGeant4ParticleSource.h"
19 
20 #include <TRestReflector.h>
21 #include <TRestStringHelper.h>
22 #include <TRestStringOutput.h>
23 
24 #include "TRestGeant4Metadata.h"
25 
26 using namespace std;
27 using namespace TRestGeant4PrimaryGeneratorTypes;
28 
30 
31 TRestGeant4ParticleSource::TRestGeant4ParticleSource() = default;
32 
33 TRestGeant4ParticleSource::~TRestGeant4ParticleSource() = default;
34 
36  RESTMetadata << " " << RESTendl;
37  if (GetParticleName() != "" && GetParticleName() != "NO_SUCH_PARA")
38  RESTMetadata << "Particle Source Name: " << GetParticleName() << RESTendl;
39  if (!fParticlesTemplate.empty() && fGenFilename != "NO_SUCH_PARA") {
40  RESTMetadata << "Generator file: " << GetGenFilename() << RESTendl;
41  RESTMetadata << "Stored templates: " << fParticlesTemplate.size() << RESTendl;
42  RESTMetadata << "Particles: ";
43  for (const auto& particle : fParticles) RESTMetadata << particle.GetParticleName() << ", ";
44  RESTMetadata << RESTendl;
45  } else {
46  if (GetParticleCharge() != 0) {
47  RESTMetadata << "Particle charge: " << GetParticleCharge() << RESTendl;
48  }
49  RESTMetadata << "Angular distribution type: " << GetAngularDistributionType() << RESTendl;
50  if (StringToAngularDistributionTypes(GetAngularDistributionType().Data()) ==
51  AngularDistributionTypes::TH1D) {
52  RESTMetadata << "Angular distribution filename: "
53  << TRestTools::GetPureFileName((string)GetAngularDistributionFilename()) << RESTendl;
54  RESTMetadata << "Angular distribution histogram name: " << GetAngularDistributionNameInFile()
55  << RESTendl;
56  }
57  RESTMetadata << "Angular distribution direction: (" << GetDirection().X() << "," << GetDirection().Y()
58  << "," << GetDirection().Z() << ")" << RESTendl;
59  if ((StringToAngularDistributionTypes(GetAngularDistributionType().Data()) ==
60  AngularDistributionTypes::FORMULA) ||
61  StringToAngularDistributionTypes(GetAngularDistributionType().Data()) ==
62  AngularDistributionTypes::FORMULA2) {
63  RESTMetadata << "Angular distribution range (deg): ("
64  << GetAngularDistributionRangeMin() * TMath::RadToDeg() << ", "
65  << GetAngularDistributionRangeMax() * TMath::RadToDeg() << ")" << RESTendl;
66  RESTMetadata << "Angular distribution random sampling grid size: "
67  << GetAngularDistributionFormulaNPoints() << RESTendl;
68  }
69  RESTMetadata << "Energy distribution type: " << GetEnergyDistributionType() << RESTendl;
70  if (StringToEnergyDistributionTypes(GetEnergyDistributionType().Data()) ==
71  EnergyDistributionTypes::TH1D) {
72  RESTMetadata << "Energy distribution filename: "
73  << TRestTools::GetPureFileName((string)GetEnergyDistributionFilename()) << RESTendl;
74  RESTMetadata << "Energy distribution histogram name: " << GetEnergyDistributionNameInFile()
75  << RESTendl;
76  }
77  if (GetEnergyDistributionRangeMin() == GetEnergyDistributionRangeMax()) {
78  RESTMetadata << "Energy distribution energy: " << GetEnergy() << " keV" << RESTendl;
79  } else {
80  RESTMetadata << "Energy distribution range (keV): (" << GetEnergyDistributionRangeMin() << ", "
81  << GetEnergyDistributionRangeMax() << ")" << RESTendl;
82  }
83  if ((StringToEnergyDistributionTypes(GetEnergyDistributionType().Data()) ==
84  EnergyDistributionTypes::FORMULA) ||
85  StringToEnergyDistributionTypes(GetEnergyDistributionType().Data()) ==
86  EnergyDistributionTypes::FORMULA2) {
87  RESTMetadata << "Energy distribution random sampling grid size: "
88  << GetEnergyDistributionFormulaNPoints() << RESTendl;
89  }
90  }
91 }
92 
93 TRestGeant4ParticleSource* TRestGeant4ParticleSource::instantiate(std::string model) {
94  if (model.empty() || model == "geant4" || model.find(".dat") != string::npos) {
95  // use default generator
96  return new TRestGeant4ParticleSource();
97  } else {
98  // use specific generator
99  // naming convention: TRestGeant4ParticleSourceXXX
100  // currently supported generators: decay0, source
101  // in future we may add wrapper of generators: pythia(for HEP), etc.
102  model[0] = *REST_StringHelper::ToUpper(std::string(&model[0], 1)).c_str();
103  TClass* c = TClass::GetClass(("TRestGeant4ParticleSource" + model).c_str());
104  if (c) // this means we have the package installed
105  {
106  return (TRestGeant4ParticleSource*)c->New();
107  } else {
108  std::cout << "REST ERROR! generator wrapper \"" << ("TRestGeant4ParticleSource" + model)
109  << "\" not found!" << std::endl;
110  }
111  }
112  return nullptr;
113 }
114 
116  TString modelUse = GetParameter("use");
117  if (((string)modelUse).find(".dat") != string::npos) {
118  string fullPathName = SearchFile((string)modelUse);
119  if (fullPathName == "") {
120  RESTError << "File not found: " << modelUse << RESTendl;
121  RESTError << "Decay0 generator file could not be found!!" << RESTendl;
122  exit(1);
123  }
124  modelUse = fullPathName;
125  }
126  SetGenFilename(modelUse);
127 
128  if (((string)fGenFilename).find(".dat") != std::string::npos) {
129  if (TRestTools::fileExists((string)fGenFilename)) {
130  ReadEventDataFile(fGenFilename);
131  }
132  }
133 }
134 
135 // base class's generator action: randomize the particle's energy/direction with distribution file
136 void TRestGeant4ParticleSource::Update() {
137  if (!fParticlesTemplate.empty()) {
138  // we use particle template to generate particles
139  Int_t rndCollection = (Int_t)(fRandomMethod() * fParticlesTemplate.size());
140  Int_t pCollectionID = rndCollection % fParticlesTemplate.size();
141  fParticles = fParticlesTemplate[pCollectionID];
142  } else {
143  TRestGeant4Particle p(*this);
144  // TODO: implement particle generation for toy simulation
145  fParticles = {p};
146  }
147 }
148 
163  if (!ReadOldDecay0File(fName)) {
164  ReadNewDecay0File(fName);
165  }
166 }
167 
175  ifstream infile;
176  infile.open(fileName);
177  if (!infile.is_open()) {
178  printf("Error when opening file %s\n", fileName.Data());
179  return false;
180  }
181 
182  int generatorEvents = 0;
183  string s;
184  // First lines are discarded.
185  for (int i = 0; i < 24; i++) {
186  getline(infile, s);
187  int pos = s.find("@nevents=");
188  if (pos != -1) generatorEvents = StringToInteger(s.substr(10, s.length() - 10));
189  }
190 
191  if (generatorEvents == 0) {
192  RESTError << "TRestG4Metadata::ReadNewDecay0File. Problem reading generator file" << RESTendl;
193  exit(1);
194  }
195 
196  TRestGeant4Particle particle;
197 
198  RESTDebug << "Reading generator file: " << fileName << RESTendl;
199  RESTDebug << "Total number of events: " << generatorEvents << RESTendl;
200  for (int n = 0; n < generatorEvents && !infile.eof(); n++) {
201  int pos = -1;
202  while (!infile.eof() && pos == -1) {
203  getline(infile, s);
204  pos = s.find("@event_start");
205  }
206 
207  // Time - nuclide is skipped
208  getline(infile, s);
209 
210  Int_t nParticles;
211  infile >> nParticles;
212  RESTDebug << "Number of particles: " << nParticles << RESTendl;
213 
214  // cout << evID <<" "<< time <<" "<< nParticles <<" "<< std::endl;
215  for (int i = 0; i < nParticles && !infile.eof(); i++) {
216  Int_t pID;
217  Double_t momx, momy, momz, mass;
218  Double_t energy = -1, momentum2;
219  TString pName;
220 
221  Double_t time;
222  infile >> pID >> time >> momx >> momy >> momz >> pName;
223 
224  RESTDebug << " ---- New particle found --- " << RESTendl;
225  RESTDebug << " Particle name: " << pName << RESTendl;
226  RESTDebug << " - pId: " << pID << RESTendl;
227  RESTDebug << " - Relative time: " << time << RESTendl;
228  RESTDebug << " - px: " << momx << " py: " << momy << " pz: " << momz << " " << RESTendl;
229 
230  if (pID == 3) {
231  momentum2 = (momx * momx) + (momy * momy) + (momz * momz);
232  mass = 0.511;
233 
234  energy = TMath::Sqrt(momentum2 + mass * mass) - mass;
235  particle.SetParticleName("e-");
236  particle.SetParticleCharge(-1);
237  particle.SetExcitationLevel(0);
238 
239  } else if (pID == 1) {
240  momentum2 = (momx * momx) + (momy * momy) + (momz * momz);
241 
242  energy = TMath::Sqrt(momentum2);
243  particle.SetParticleName("gamma");
244  particle.SetParticleCharge(0);
245  particle.SetExcitationLevel(0);
246  } else {
247  cout << "Particle id " << pID << " not recognized" << std::endl;
248  }
249 
250  TVector3 momDirection(momx, momy, momz);
251  momDirection = momDirection.Unit();
252 
253  particle.SetEnergy(1000. * energy);
254  particle.SetDirection(momDirection);
255 
256  this->AddParticle(particle);
257  }
258  this->FlushParticlesTemplate();
259  }
260 
261  return true;
262 }
263 
271  ifstream infile;
272  infile.open(fileName);
273  if (!infile.is_open()) {
274  printf("Error when opening file %s\n", fileName.Data());
275  return false;
276  }
277 
278  string s;
279  // First lines are discarded.
280  int headerFound = 0;
281  for (int i = 0; i < 30; i++) {
282  getline(infile, s);
283  if (s.find("#!bxdecay0 1.0.0") != string::npos) return false;
284  if (s.find("First event and full number of events:") != string::npos) {
285  headerFound = 1;
286  break;
287  }
288  }
289  if (!headerFound) {
290  RESTError
291  << "TRestG4Metadata::ReadOldDecay0File. Problem reading generator file: no \"First event and "
292  "full number of events:\" header.\n";
293  abort();
294  }
295  int tmpInt;
296  int fGeneratorEvents;
297  infile >> tmpInt >> fGeneratorEvents;
298 
299  // cout << "i: " << tmpInt << " fN: " << fGeneratorEvents << std::endl;
300  // TRestGeant4ParticleSource* src = TRestGeant4ParticleSource::instantiate();
301  // this->SetGenFilename(fileName);
302  // this->SetAngularDistributionType("flux");
303  // this->SetEnergyDistributionType("mono");
304 
305  TRestGeant4Particle particle;
306 
307  cout << "Reading generator file: " << fileName << std::endl;
308  cout << "Total number of events: " << fGeneratorEvents << std::endl;
309  for (int n = 0; n < fGeneratorEvents && !infile.eof(); n++) {
310  Int_t nParticles;
311  Int_t evID;
312  Double_t time;
313 
314  infile >> evID >> time >> nParticles;
315 
316  // cout << evID <<" "<< time <<" "<< nParticles <<" "<< std::endl;
317  for (int i = 0; i < nParticles && !infile.eof(); i++) {
318  Int_t pID;
319  Double_t momx, momy, momz, mass;
320  Double_t energy = -1, momentum2;
321  Double_t x = 0, y = 0, z = 0;
322 
323  infile >> pID >> momx >> momy >> momz >> time;
324  // infile >> x >> y >> z;
325 
326  // cout << momx << " " << momy << " " << momz << " " << std::endl;
327 
328  bool ise = 2 <= pID && pID <= 3, ismu = 5 <= pID && pID <= 6, isp = pID == 14, isg = pID == 1;
329  if (ise || ismu || isp || isg) {
330  momentum2 = (momx * momx) + (momy * momy) + (momz * momz);
331  if (ise) {
332  mass = 0.511;
333  particle.SetParticleName(pID == 2 ? "e+" : "e-");
334  particle.SetParticleCharge(pID == 2 ? 1 : -1);
335  } else if (ismu) {
336  mass = 105.7;
337  particle.SetParticleName(pID == 5 ? "mu+" : "mu-");
338  particle.SetParticleCharge(pID == 5 ? 1 : -1);
339  } else if (isp) {
340  mass = 938.3;
341  particle.SetParticleName("proton");
342  particle.SetParticleCharge(1);
343  } else {
344  mass = 0;
345  particle.SetParticleName("gamma");
346  particle.SetParticleCharge(0);
347  }
348 
349  energy = TMath::Sqrt(momentum2 + mass * mass) - mass;
350  particle.SetExcitationLevel(0);
351  } else {
352  cout << "Particle id " << pID << " not recognized" << std::endl;
353  }
354 
355  TVector3 momDirection(momx, momy, momz);
356  momDirection = momDirection.Unit();
357 
358  particle.SetEnergy(1000. * energy);
359  particle.SetDirection(momDirection);
360  particle.SetOrigin(TVector3(x, y, z));
361 
362  this->AddParticle(particle);
363  }
364 
365  this->FlushParticlesTemplate();
366  }
367 
368  return true;
369 }
void ReadEventDataFile(TString fName)
Reads an input file produced by Decay0.
virtual void InitFromConfigFile() override
To make settings from rml file. This method must be implemented in the derived class.
virtual void PrintMetadata() override
Implemented it in the derived metadata class to print out specific metadata information.
bool ReadOldDecay0File(TString fileName)
Reads particle information using the file format from older Decay0 versions.
bool ReadNewDecay0File(TString fileName)
Reads particle information using the file format from newer Decay0 versions.
A class used to store particle properties.
static std::string GetPureFileName(const std::string &fullPathFileName)
Removes all directories in the full path filename description given in the argument.
Definition: TRestTools.cxx:863
static bool fileExists(const std::string &filename)
Returns true if the file (or directory) with path filename exists.
Definition: TRestTools.cxx:728
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.