File indexing completed on 2025-02-23 09:22:40
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033 #include "Run.hh"
0034
0035 #include "DetectorConstruction.hh"
0036 #include "HistoManager.hh"
0037 #include "PrimaryGeneratorAction.hh"
0038
0039 #include "G4AutoLock.hh"
0040 #include "G4SystemOfUnits.hh"
0041 #include "G4Threading.hh"
0042 #include "G4UnitsTable.hh"
0043
0044
0045
0046 namespace
0047 {
0048
0049 G4Mutex ionIdMapMutex = G4MUTEX_INITIALIZER;
0050 }
0051
0052
0053
0054 std::map<G4String, G4int> Run::fgIonMap;
0055 G4int Run::fgIonId = kMaxHisto1;
0056
0057
0058
0059 Run::Run(DetectorConstruction* det) : fDetector(det) {}
0060
0061
0062
0063 void Run::Merge(std::map<G4String, ParticleData>& destinationMap,
0064 const std::map<G4String, ParticleData>& sourceMap) const
0065 {
0066 for (const auto& particleData : sourceMap) {
0067 G4String name = particleData.first;
0068 const ParticleData& localData = particleData.second;
0069 if (destinationMap.find(name) == destinationMap.end()) {
0070 destinationMap[name] = ParticleData(localData.fCount, localData.fEmean, localData.fEmin,
0071 localData.fEmax, localData.fTmean);
0072 }
0073 else {
0074 ParticleData& data = destinationMap[name];
0075 data.fCount += localData.fCount;
0076 data.fEmean += localData.fEmean;
0077 G4double emin = localData.fEmin;
0078 if (emin < data.fEmin) data.fEmin = emin;
0079 G4double emax = localData.fEmax;
0080 if (emax > data.fEmax) data.fEmax = emax;
0081 data.fTmean = localData.fTmean;
0082 }
0083 }
0084 }
0085
0086
0087
0088 void Run::SetPrimary(G4ParticleDefinition* particle, G4double energy)
0089 {
0090 fParticle = particle;
0091 fEkin = energy;
0092 }
0093
0094
0095
0096 void Run::CountProcesses(const G4VProcess* process)
0097 {
0098 if (process == nullptr) return;
0099 G4String procName = process->GetProcessName();
0100 std::map<G4String, G4int>::iterator it = fProcCounter.find(procName);
0101 if (it == fProcCounter.end()) {
0102 fProcCounter[procName] = 1;
0103 }
0104 else {
0105 fProcCounter[procName]++;
0106 }
0107 }
0108
0109
0110
0111 void Run::ParticleCount(G4String name, G4double Ekin, G4double meanLife)
0112 {
0113 std::map<G4String, ParticleData>::iterator it = fParticleDataMap1.find(name);
0114 if (it == fParticleDataMap1.end()) {
0115 fParticleDataMap1[name] = ParticleData(1, Ekin, Ekin, Ekin, meanLife);
0116 }
0117 else {
0118 ParticleData& data = it->second;
0119 data.fCount++;
0120 data.fEmean += Ekin;
0121
0122 G4double emin = data.fEmin;
0123 if (Ekin < emin) data.fEmin = Ekin;
0124 G4double emax = data.fEmax;
0125 if (Ekin > emax) data.fEmax = Ekin;
0126 data.fTmean = meanLife;
0127 }
0128 }
0129
0130
0131
0132 void Run::AddEdep(G4double edep)
0133 {
0134 fEnergyDeposit += edep;
0135 fEnergyDeposit2 += edep * edep;
0136 }
0137
0138
0139
0140 void Run::AddEflow(G4double eflow)
0141 {
0142 fEnergyFlow += eflow;
0143 fEnergyFlow2 += eflow * eflow;
0144 }
0145
0146
0147 void Run::ParticleFlux(G4String name, G4double Ekin)
0148 {
0149 std::map<G4String, ParticleData>::iterator it = fParticleDataMap2.find(name);
0150 if (it == fParticleDataMap2.end()) {
0151 fParticleDataMap2[name] = ParticleData(1, Ekin, Ekin, Ekin, -1 * ns);
0152 }
0153 else {
0154 ParticleData& data = it->second;
0155 data.fCount++;
0156 data.fEmean += Ekin;
0157
0158 G4double emin = data.fEmin;
0159 if (Ekin < emin) data.fEmin = Ekin;
0160 G4double emax = data.fEmax;
0161 if (Ekin > emax) data.fEmax = Ekin;
0162 data.fTmean = -1 * ns;
0163 }
0164 }
0165
0166
0167
0168 G4int Run::GetIonId(G4String ionName)
0169 {
0170 G4AutoLock lock(&ionIdMapMutex);
0171
0172
0173 std::map<G4String, G4int>::const_iterator it = fgIonMap.find(ionName);
0174 if (it == fgIonMap.end()) {
0175 fgIonMap[ionName] = fgIonId;
0176 if (fgIonId < (kMaxHisto2 - 1)) fgIonId++;
0177 }
0178 return fgIonMap[ionName];
0179 }
0180
0181
0182
0183 void Run::Merge(const G4Run* run)
0184 {
0185 const Run* localRun = static_cast<const Run*>(run);
0186
0187
0188
0189 fParticle = localRun->fParticle;
0190 fEkin = localRun->fEkin;
0191
0192
0193
0194 fEnergyDeposit += localRun->fEnergyDeposit;
0195 fEnergyDeposit2 += localRun->fEnergyDeposit2;
0196 fEnergyFlow += localRun->fEnergyFlow;
0197 fEnergyFlow2 += localRun->fEnergyFlow2;
0198
0199
0200 for (const auto& procCounter : localRun->fProcCounter) {
0201 G4String procName = procCounter.first;
0202 G4int localCount = procCounter.second;
0203 if (fProcCounter.find(procName) == fProcCounter.end()) {
0204 fProcCounter[procName] = localCount;
0205 }
0206 else {
0207 fProcCounter[procName] += localCount;
0208 }
0209 }
0210
0211
0212 Merge(fParticleDataMap1, localRun->fParticleDataMap1);
0213
0214
0215 Merge(fParticleDataMap2, localRun->fParticleDataMap2);
0216
0217 G4Run::Merge(run);
0218 }
0219
0220
0221
0222 void Run::EndOfRun()
0223 {
0224 G4int prec = 5, wid = prec + 2;
0225 G4int dfprec = G4cout.precision(prec);
0226
0227
0228
0229 G4Material* material = fDetector->GetAbsorMaterial();
0230 G4double density = material->GetDensity();
0231
0232 G4String Particle = fParticle->GetParticleName();
0233 G4cout << "\n The run is " << numberOfEvent << " " << Particle << " of "
0234 << G4BestUnit(fEkin, "Energy") << " through "
0235 << G4BestUnit(fDetector->GetAbsorThickness(), "Length") << " of " << material->GetName()
0236 << " (density: " << G4BestUnit(density, "Volumic Mass") << ")" << G4endl;
0237
0238 if (numberOfEvent == 0) {
0239 G4cout.precision(dfprec);
0240 return;
0241 }
0242
0243
0244
0245 G4cout << "\n Process calls frequency :" << G4endl;
0246 G4int index = 0;
0247 for (const auto& procCounter : fProcCounter) {
0248 G4String procName = procCounter.first;
0249 G4int count = procCounter.second;
0250 G4String space = " ";
0251 if (++index % 3 == 0) space = "\n";
0252 G4cout << " " << std::setw(20) << procName << "=" << std::setw(7) << count << space;
0253 }
0254 G4cout << G4endl;
0255
0256
0257
0258 G4cout << "\n List of generated particles (with meanLife != 0):" << G4endl;
0259
0260 for (const auto& particleData : fParticleDataMap1) {
0261 G4String name = particleData.first;
0262 ParticleData data = particleData.second;
0263 G4int count = data.fCount;
0264 G4double eMean = data.fEmean / count;
0265 G4double eMin = data.fEmin;
0266 G4double eMax = data.fEmax;
0267 G4double meanLife = data.fTmean;
0268
0269 G4cout << " " << std::setw(13) << name << ": " << std::setw(7) << count
0270 << " Emean = " << std::setw(wid) << G4BestUnit(eMean, "Energy") << "\t( "
0271 << G4BestUnit(eMin, "Energy") << " --> " << G4BestUnit(eMax, "Energy") << ")";
0272 if (meanLife >= 0.)
0273 G4cout << "\tmean life = " << G4BestUnit(meanLife, "Time") << G4endl;
0274 else
0275 G4cout << "\tstable" << G4endl;
0276 }
0277
0278
0279
0280 G4int TotNbofEvents = numberOfEvent;
0281 fEnergyDeposit /= TotNbofEvents;
0282 fEnergyDeposit2 /= TotNbofEvents;
0283 G4double rmsEdep = fEnergyDeposit2 - fEnergyDeposit * fEnergyDeposit;
0284 if (rmsEdep > 0.)
0285 rmsEdep = std::sqrt(rmsEdep);
0286 else
0287 rmsEdep = 0.;
0288
0289 G4cout << "\n Mean energy deposit per event = " << G4BestUnit(fEnergyDeposit, "Energy")
0290 << "; rms = " << G4BestUnit(rmsEdep, "Energy") << G4endl;
0291
0292
0293
0294 fEnergyFlow /= TotNbofEvents;
0295 fEnergyFlow2 /= TotNbofEvents;
0296 G4double rmsEflow = fEnergyFlow2 - fEnergyFlow * fEnergyFlow;
0297 if (rmsEflow > 0.)
0298 rmsEflow = std::sqrt(rmsEflow);
0299 else
0300 rmsEflow = 0.;
0301
0302 G4cout << " Mean energy flow per event = " << G4BestUnit(fEnergyFlow, "Energy")
0303 << "; rms = " << G4BestUnit(rmsEflow, "Energy") << G4endl;
0304
0305
0306
0307 G4cout << "\n List of particles emerging from the target :" << G4endl;
0308
0309 for (const auto& particleData : fParticleDataMap2) {
0310 G4String name = particleData.first;
0311 ParticleData data = particleData.second;
0312 G4int count = data.fCount;
0313 G4double eMean = data.fEmean / count;
0314 G4double eMin = data.fEmin;
0315 G4double eMax = data.fEmax;
0316 G4double Eflow = data.fEmean / TotNbofEvents;
0317
0318 G4cout << " " << std::setw(13) << name << ": " << std::setw(7) << count
0319 << " Emean = " << std::setw(wid) << G4BestUnit(eMean, "Energy") << "\t( "
0320 << G4BestUnit(eMin, "Energy") << " --> " << G4BestUnit(eMax, "Energy")
0321 << ") \tEflow/event = " << G4BestUnit(Eflow, "Energy") << G4endl;
0322 }
0323
0324
0325
0326 G4cout << "\n histo Id for populations :" << G4endl;
0327
0328
0329
0330 G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
0331 for (const auto& ionMapElement : fgIonMap) {
0332 G4String ionName = ionMapElement.first;
0333 G4int h1Id = ionMapElement.second;
0334
0335 G4cout << " " << std::setw(20) << ionName << " id = " << std::setw(3) << h1Id << G4endl;
0336
0337
0338 if (!analysisManager->GetH1(h1Id)) continue;
0339
0340
0341 G4String title = analysisManager->GetH1Title(h1Id);
0342 title = ionName + title;
0343 analysisManager->SetH1Title(h1Id, title);
0344 }
0345 G4cout << G4endl;
0346
0347
0348 G4int ih = 2;
0349 G4double binWidth = analysisManager->GetH1Width(ih);
0350 G4double fac = (1. / (numberOfEvent * binWidth)) * (mm / MeV);
0351 analysisManager->ScaleH1(ih, fac);
0352
0353 for (ih = 14; ih < 24; ih++) {
0354 binWidth = analysisManager->GetH1Width(ih);
0355 G4double unit = analysisManager->GetH1Unit(ih);
0356 fac = (second / (binWidth * unit));
0357 analysisManager->ScaleH1(ih, fac);
0358 }
0359
0360
0361 fProcCounter.clear();
0362 fParticleDataMap1.clear();
0363 fParticleDataMap2.clear();
0364 fgIonMap.clear();
0365
0366
0367 G4cout.precision(dfprec);
0368 }
0369
0370