18 #include "TRestGeant4ParticleSource.h"
20 #include <TRestReflector.h>
21 #include <TRestStringHelper.h>
22 #include <TRestStringOutput.h>
24 #include "TRestGeant4Metadata.h"
27 using namespace TRestGeant4PrimaryGeneratorTypes;
31 TRestGeant4ParticleSource::TRestGeant4ParticleSource() =
default;
33 TRestGeant4ParticleSource::~TRestGeant4ParticleSource() =
default;
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;
46 if (GetParticleCharge() != 0) {
47 RESTMetadata <<
"Particle charge: " << GetParticleCharge() << RESTendl;
49 RESTMetadata <<
"Angular distribution type: " << GetAngularDistributionType() << RESTendl;
50 if (StringToAngularDistributionTypes(GetAngularDistributionType().Data()) ==
51 AngularDistributionTypes::TH1D) {
52 RESTMetadata <<
"Angular distribution filename: "
54 RESTMetadata <<
"Angular distribution histogram name: " << GetAngularDistributionNameInFile()
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;
69 RESTMetadata <<
"Energy distribution type: " << GetEnergyDistributionType() << RESTendl;
70 if (StringToEnergyDistributionTypes(GetEnergyDistributionType().Data()) ==
71 EnergyDistributionTypes::TH1D) {
72 RESTMetadata <<
"Energy distribution filename: "
74 RESTMetadata <<
"Energy distribution histogram name: " << GetEnergyDistributionNameInFile()
77 if (GetEnergyDistributionRangeMin() == GetEnergyDistributionRangeMax()) {
78 RESTMetadata <<
"Energy distribution energy: " << GetEnergy() <<
" keV" << RESTendl;
80 RESTMetadata <<
"Energy distribution range (keV): (" << GetEnergyDistributionRangeMin() <<
", "
81 << GetEnergyDistributionRangeMax() <<
")" << RESTendl;
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;
94 if (model.empty() || model ==
"geant4" || model.find(
".dat") != string::npos) {
103 TClass* c = TClass::GetClass((
"TRestGeant4ParticleSource" + model).c_str());
108 std::cout <<
"REST ERROR! generator wrapper \"" << (
"TRestGeant4ParticleSource" + model)
109 <<
"\" not found!" << std::endl;
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;
124 modelUse = fullPathName;
126 SetGenFilename(modelUse);
128 if (((
string)fGenFilename).find(
".dat") != std::string::npos) {
130 ReadEventDataFile(fGenFilename);
136 void TRestGeant4ParticleSource::Update() {
137 if (!fParticlesTemplate.empty()) {
139 Int_t rndCollection = (Int_t)(fRandomMethod() * fParticlesTemplate.size());
140 Int_t pCollectionID = rndCollection % fParticlesTemplate.size();
141 fParticles = fParticlesTemplate[pCollectionID];
163 if (!ReadOldDecay0File(fName)) {
164 ReadNewDecay0File(fName);
176 infile.open(fileName);
177 if (!infile.is_open()) {
178 printf(
"Error when opening file %s\n", fileName.Data());
182 int generatorEvents = 0;
185 for (
int i = 0; i < 24; i++) {
187 int pos = s.find(
"@nevents=");
188 if (pos != -1) generatorEvents =
StringToInteger(s.substr(10, s.length() - 10));
191 if (generatorEvents == 0) {
192 RESTError <<
"TRestG4Metadata::ReadNewDecay0File. Problem reading generator file" << RESTendl;
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++) {
202 while (!infile.eof() && pos == -1) {
204 pos = s.find(
"@event_start");
211 infile >> nParticles;
212 RESTDebug <<
"Number of particles: " << nParticles << RESTendl;
215 for (
int i = 0; i < nParticles && !infile.eof(); i++) {
217 Double_t momx, momy, momz, mass;
218 Double_t energy = -1, momentum2;
222 infile >> pID >> time >> momx >> momy >> momz >> pName;
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;
231 momentum2 = (momx * momx) + (momy * momy) + (momz * momz);
234 energy = TMath::Sqrt(momentum2 + mass * mass) - mass;
235 particle.SetParticleName(
"e-");
236 particle.SetParticleCharge(-1);
237 particle.SetExcitationLevel(0);
239 }
else if (pID == 1) {
240 momentum2 = (momx * momx) + (momy * momy) + (momz * momz);
242 energy = TMath::Sqrt(momentum2);
243 particle.SetParticleName(
"gamma");
244 particle.SetParticleCharge(0);
245 particle.SetExcitationLevel(0);
247 cout <<
"Particle id " << pID <<
" not recognized" << std::endl;
250 TVector3 momDirection(momx, momy, momz);
251 momDirection = momDirection.Unit();
253 particle.SetEnergy(1000. * energy);
254 particle.SetDirection(momDirection);
256 this->AddParticle(particle);
258 this->FlushParticlesTemplate();
272 infile.open(fileName);
273 if (!infile.is_open()) {
274 printf(
"Error when opening file %s\n", fileName.Data());
281 for (
int i = 0; i < 30; i++) {
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) {
291 <<
"TRestG4Metadata::ReadOldDecay0File. Problem reading generator file: no \"First event and "
292 "full number of events:\" header.\n";
296 int fGeneratorEvents;
297 infile >> tmpInt >> fGeneratorEvents;
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++) {
314 infile >> evID >> time >> nParticles;
317 for (
int i = 0; i < nParticles && !infile.eof(); i++) {
319 Double_t momx, momy, momz, mass;
320 Double_t energy = -1, momentum2;
321 Double_t x = 0, y = 0, z = 0;
323 infile >> pID >> momx >> momy >> momz >> time;
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);
333 particle.SetParticleName(pID == 2 ?
"e+" :
"e-");
334 particle.SetParticleCharge(pID == 2 ? 1 : -1);
337 particle.SetParticleName(pID == 5 ?
"mu+" :
"mu-");
338 particle.SetParticleCharge(pID == 5 ? 1 : -1);
341 particle.SetParticleName(
"proton");
342 particle.SetParticleCharge(1);
345 particle.SetParticleName(
"gamma");
346 particle.SetParticleCharge(0);
349 energy = TMath::Sqrt(momentum2 + mass * mass) - mass;
350 particle.SetExcitationLevel(0);
352 cout <<
"Particle id " << pID <<
" not recognized" << std::endl;
355 TVector3 momDirection(momx, momy, momz);
356 momDirection = momDirection.Unit();
358 particle.SetEnergy(1000. * energy);
359 particle.SetDirection(momDirection);
360 particle.SetOrigin(TVector3(x, y, z));
362 this->AddParticle(particle);
365 this->FlushParticlesTemplate();
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.
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.