File indexing completed on 2026-04-28 07:17:27
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 #include "Run.hh"
0030
0031 #include "DetectorConstruction.hh"
0032 #include "HistoManager.hh"
0033 #include "PrimaryGeneratorAction.hh"
0034
0035 #include "G4SystemOfUnits.hh"
0036 #include "G4UnitsTable.hh"
0037
0038
0039
0040 Run::Run(DetectorConstruction* det) : fDetector(det) {}
0041
0042
0043
0044 void Run::SetPrimary(G4ParticleDefinition* particle, G4double energy)
0045 {
0046 fParticle = particle;
0047 fEkin = energy;
0048 }
0049
0050
0051
0052 void Run::CountProcesses(const G4VProcess* process)
0053 {
0054 if (process == nullptr) return;
0055 G4String procName = process->GetProcessName();
0056 std::map<G4String, G4int>::iterator it = fProcCounter.find(procName);
0057 if (it == fProcCounter.end()) {
0058 fProcCounter[procName] = 1;
0059 }
0060 else {
0061 fProcCounter[procName]++;
0062 }
0063 }
0064
0065
0066
0067 void Run::ParticleCount(G4String name, G4double Ekin, G4double meanLife)
0068 {
0069 std::map<G4String, ParticleData>::iterator it = fParticleDataMap1.find(name);
0070 if (it == fParticleDataMap1.end()) {
0071 fParticleDataMap1[name] = ParticleData(1, Ekin, Ekin, Ekin, meanLife);
0072 }
0073 else {
0074 ParticleData& data = it->second;
0075 data.fCount++;
0076 data.fEmean += Ekin;
0077
0078 G4double emin = data.fEmin;
0079 if (Ekin < emin) data.fEmin = Ekin;
0080 G4double emax = data.fEmax;
0081 if (Ekin > emax) data.fEmax = Ekin;
0082 data.fTmean = meanLife;
0083 }
0084 }
0085
0086
0087
0088 void Run::SumEnergies(G4double edep, G4double eflow, G4double etot)
0089 {
0090 fEnergyDeposit += edep;
0091 fEnergyDeposit2 += edep * edep;
0092
0093 fEnergyFlow += eflow;
0094 fEnergyFlow2 += eflow * eflow;
0095
0096 fEnergyTotal += etot;
0097 fEnergyTotal2 += etot * etot;
0098 }
0099
0100
0101
0102 void Run::ParticleFlux(G4String name, G4double Ekin)
0103 {
0104 std::map<G4String, ParticleData>::iterator it = fParticleDataMap2.find(name);
0105 if (it == fParticleDataMap2.end()) {
0106 fParticleDataMap2[name] = ParticleData(1, Ekin, Ekin, Ekin, -1 * ns);
0107 }
0108 else {
0109 ParticleData& data = it->second;
0110 data.fCount++;
0111 data.fEmean += Ekin;
0112
0113 G4double emin = data.fEmin;
0114 if (Ekin < emin) data.fEmin = Ekin;
0115 G4double emax = data.fEmax;
0116 if (Ekin > emax) data.fEmax = Ekin;
0117 data.fTmean = -1 * ns;
0118 }
0119 }
0120
0121
0122
0123 void Run::Merge(const G4Run* run)
0124 {
0125 const Run* localRun = static_cast<const Run*>(run);
0126
0127
0128
0129 fParticle = localRun->fParticle;
0130 fEkin = localRun->fEkin;
0131
0132
0133
0134 fEnergyDeposit += localRun->fEnergyDeposit;
0135 fEnergyDeposit2 += localRun->fEnergyDeposit2;
0136 fEnergyFlow += localRun->fEnergyFlow;
0137 fEnergyFlow2 += localRun->fEnergyFlow2;
0138 fEnergyTotal += localRun->fEnergyTotal;
0139 fEnergyTotal2 += localRun->fEnergyTotal2;
0140
0141
0142 std::map<G4String, G4int>::const_iterator itp;
0143 for (itp = localRun->fProcCounter.begin(); itp != localRun->fProcCounter.end(); ++itp) {
0144 G4String procName = itp->first;
0145 G4int localCount = itp->second;
0146 if (fProcCounter.find(procName) == fProcCounter.end()) {
0147 fProcCounter[procName] = localCount;
0148 }
0149 else {
0150 fProcCounter[procName] += localCount;
0151 }
0152 }
0153
0154
0155 std::map<G4String, ParticleData>::const_iterator itc;
0156 for (itc = localRun->fParticleDataMap1.begin(); itc != localRun->fParticleDataMap1.end(); ++itc) {
0157 G4String name = itc->first;
0158 const ParticleData& localData = itc->second;
0159 if (fParticleDataMap1.find(name) == fParticleDataMap1.end()) {
0160 fParticleDataMap1[name] = ParticleData(localData.fCount, localData.fEmean, localData.fEmin,
0161 localData.fEmax, localData.fTmean);
0162 }
0163 else {
0164 ParticleData& data = fParticleDataMap1[name];
0165 data.fCount += localData.fCount;
0166 data.fEmean += localData.fEmean;
0167 G4double emin = localData.fEmin;
0168 if (emin < data.fEmin) data.fEmin = emin;
0169 G4double emax = localData.fEmax;
0170 if (emax > data.fEmax) data.fEmax = emax;
0171 data.fTmean = localData.fTmean;
0172 }
0173 }
0174
0175
0176 std::map<G4String, ParticleData>::const_iterator itn;
0177 for (itn = localRun->fParticleDataMap2.begin(); itn != localRun->fParticleDataMap2.end(); ++itn) {
0178 G4String name = itn->first;
0179 const ParticleData& localData = itn->second;
0180 if (fParticleDataMap2.find(name) == fParticleDataMap2.end()) {
0181 fParticleDataMap2[name] = ParticleData(localData.fCount, localData.fEmean, localData.fEmin,
0182 localData.fEmax, localData.fTmean);
0183 }
0184 else {
0185 ParticleData& data = fParticleDataMap2[name];
0186 data.fCount += localData.fCount;
0187 data.fEmean += localData.fEmean;
0188 G4double emin = localData.fEmin;
0189 if (emin < data.fEmin) data.fEmin = emin;
0190 G4double emax = localData.fEmax;
0191 if (emax > data.fEmax) data.fEmax = emax;
0192 data.fTmean = localData.fTmean;
0193 }
0194 }
0195
0196 G4Run::Merge(run);
0197 }
0198
0199
0200
0201 void Run::EndOfRun()
0202 {
0203 G4int prec = 5, wid = prec + 2;
0204 G4int dfprec = G4cout.precision(prec);
0205
0206
0207
0208 G4Material* material = fDetector->GetMaterial();
0209 G4double density = material->GetDensity();
0210
0211 G4String Particle = fParticle->GetParticleName();
0212 G4cout << "\n The run is " << numberOfEvent << " " << Particle << " of "
0213 << G4BestUnit(fEkin, "Energy") << " through "
0214 << G4BestUnit(fDetector->GetRadius(), "Length") << " of " << material->GetName()
0215 << " (density: " << G4BestUnit(density, "Volumic Mass") << ")" << G4endl;
0216
0217 if (numberOfEvent == 0) {
0218 G4cout.precision(dfprec);
0219 return;
0220 }
0221
0222
0223
0224 G4cout << "\n Process calls frequency :" << G4endl;
0225 G4int index = 0;
0226 std::map<G4String, G4int>::iterator it;
0227 for (it = fProcCounter.begin(); it != fProcCounter.end(); it++) {
0228 G4String procName = it->first;
0229 G4int count = it->second;
0230 G4String space = " ";
0231 if (++index % 3 == 0) space = "\n";
0232 G4cout << " " << std::setw(20) << procName << "=" << std::setw(7) << count << space;
0233 }
0234 G4cout << G4endl;
0235
0236
0237
0238 G4int TotNbofEvents = numberOfEvent;
0239 fEnergyDeposit /= TotNbofEvents;
0240 fEnergyDeposit2 /= TotNbofEvents;
0241 G4double rmsEdep = fEnergyDeposit2 - fEnergyDeposit * fEnergyDeposit;
0242 if (rmsEdep > 0.)
0243 rmsEdep = std::sqrt(rmsEdep);
0244 else
0245 rmsEdep = 0.;
0246
0247 G4cout << "\n Mean energy deposit per event = " << G4BestUnit(fEnergyDeposit, "Energy")
0248 << "; rms = " << G4BestUnit(rmsEdep, "Energy") << G4endl;
0249
0250
0251
0252 fEnergyFlow /= TotNbofEvents;
0253 fEnergyFlow2 /= TotNbofEvents;
0254 G4double rmsEflow = fEnergyFlow2 - fEnergyFlow * fEnergyFlow;
0255 if (rmsEflow > 0.)
0256 rmsEflow = std::sqrt(rmsEflow);
0257 else
0258 rmsEflow = 0.;
0259
0260 G4cout << " Mean energy leakage per event = " << G4BestUnit(fEnergyFlow, "Energy")
0261 << "; rms = " << G4BestUnit(rmsEflow, "Energy") << G4endl;
0262
0263
0264
0265 fEnergyTotal /= TotNbofEvents;
0266 fEnergyTotal2 /= TotNbofEvents;
0267 G4double rmsEtotal = fEnergyTotal2 - fEnergyTotal * fEnergyTotal;
0268 if (rmsEtotal > 0.)
0269 rmsEtotal = std::sqrt(rmsEtotal);
0270 else
0271 rmsEflow = 0.;
0272
0273 G4cout << "\n Mean energy total per event = " << G4BestUnit(fEnergyTotal, "Energy")
0274 << "; rms = " << G4BestUnit(rmsEtotal, "Energy") << G4endl;
0275
0276
0277
0278 if (fParticleDataMap1.size() > 0) {
0279 G4cout << "\n List of particles at creation :" << G4endl;
0280 std::map<G4String, ParticleData>::iterator itc;
0281 for (itc = fParticleDataMap1.begin(); itc != fParticleDataMap1.end(); itc++) {
0282 G4String name = itc->first;
0283 ParticleData data = itc->second;
0284 G4int count = data.fCount;
0285 G4double eMean = data.fEmean / count;
0286 G4double eMin = data.fEmin;
0287 G4double eMax = data.fEmax;
0288 G4double meanLife = data.fTmean;
0289
0290 G4cout << " " << std::setw(13) << name << ": " << std::setw(7) << count
0291 << " Emean = " << std::setw(wid) << G4BestUnit(eMean, "Energy") << "\t( "
0292 << G4BestUnit(eMin, "Energy") << " --> " << G4BestUnit(eMax, "Energy") << ")";
0293 if (meanLife >= 0.)
0294 G4cout << "\tmean life = " << G4BestUnit(meanLife, "Time") << G4endl;
0295 else
0296 G4cout << "\tstable" << G4endl;
0297 }
0298 }
0299
0300
0301
0302 G4cout << "\n List of particles emerging from the absorber :" << G4endl;
0303
0304 std::map<G4String, ParticleData>::iterator itn;
0305 for (itn = fParticleDataMap2.begin(); itn != fParticleDataMap2.end(); itn++) {
0306 G4String name = itn->first;
0307 ParticleData data = itn->second;
0308 G4int count = data.fCount;
0309 G4double eMean = data.fEmean / count;
0310 G4double eMin = data.fEmin;
0311 G4double eMax = data.fEmax;
0312 G4double Eflow = data.fEmean / TotNbofEvents;
0313
0314 G4cout << " " << std::setw(13) << name << ": " << std::setw(7) << count
0315 << " Emean = " << std::setw(wid) << G4BestUnit(eMean, "Energy") << "\t( "
0316 << G4BestUnit(eMin, "Energy") << " --> " << G4BestUnit(eMax, "Energy")
0317 << ") \tEleak/event = " << G4BestUnit(Eflow, "Energy") << G4endl;
0318 }
0319
0320
0321 G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
0322
0323 G4int ih = 2;
0324 G4double binWidth = analysisManager->GetH1Width(ih);
0325 G4double fac = (1. / (numberOfEvent * binWidth)) * (mm / MeV);
0326 analysisManager->ScaleH1(ih, fac);
0327
0328
0329 fProcCounter.clear();
0330 fParticleDataMap2.clear();
0331
0332
0333 G4cout.precision(dfprec);
0334 }
0335
0336