Geant4 Cross Reference |
1 // 1 2 // ******************************************* 3 // * License and Disclaimer 4 // * 5 // * The Geant4 software is copyright of th 6 // * the Geant4 Collaboration. It is provided 7 // * conditions of the Geant4 Software License 8 // * LICENSE and available at http://cern.ch/ 9 // * include a list of copyright holders. 10 // * 11 // * Neither the authors of this software syst 12 // * institutes,nor the agencies providing fin 13 // * work make any representation or warran 14 // * regarding this software system or assum 15 // * use. Please see the license in the file 16 // * for the full disclaimer and the limitatio 17 // * 18 // * This code implementation is the result 19 // * technical work of the GEANT4 collaboratio 20 // * By using, copying, modifying or distri 21 // * any work based on the software) you ag 22 // * use in resulting scientific publicati 23 // * acceptance of all terms of the Geant4 Sof 24 // ******************************************* 25 // 26 /// \file Run.cc 27 /// \brief Implementation of the Run class 28 // 29 // 30 //....oooOO0OOooo........oooOO0OOooo........oo 31 //....oooOO0OOooo........oooOO0OOooo........oo 32 33 #include "Run.hh" 34 35 #include "DetectorConstruction.hh" 36 #include "HistoManager.hh" 37 #include "PrimaryGeneratorAction.hh" 38 39 #include "G4SystemOfUnits.hh" 40 #include "G4UnitsTable.hh" 41 42 //....oooOO0OOooo........oooOO0OOooo........oo 43 44 Run::Run(DetectorConstruction* det) : fDetecto 45 46 //....oooOO0OOooo........oooOO0OOooo........oo 47 48 void Run::SetPrimary(G4ParticleDefinition* par 49 { 50 fParticle = particle; 51 fEkin = energy; 52 } 53 54 //....oooOO0OOooo........oooOO0OOooo........oo 55 56 void Run::CountProcesses(const G4VProcess* pro 57 { 58 if (process == nullptr) return; 59 G4String procName = process->GetProcessName( 60 std::map<G4String, G4int>::iterator it = fPr 61 if (it == fProcCounter.end()) { 62 fProcCounter[procName] = 1; 63 } 64 else { 65 fProcCounter[procName]++; 66 } 67 } 68 69 //....oooOO0OOooo........oooOO0OOooo........oo 70 71 void Run::ParticleCount(G4String name, G4doubl 72 { 73 std::map<G4String, ParticleData>::iterator i 74 if (it == fParticleDataMap1.end()) { 75 fParticleDataMap1[name] = ParticleData(1, 76 } 77 else { 78 ParticleData& data = it->second; 79 data.fCount++; 80 data.fEmean += Ekin; 81 // update min max 82 G4double emin = data.fEmin; 83 if (Ekin < emin) data.fEmin = Ekin; 84 G4double emax = data.fEmax; 85 if (Ekin > emax) data.fEmax = Ekin; 86 data.fTmean = meanLife; 87 } 88 } 89 90 //....oooOO0OOooo........oooOO0OOooo........oo 91 92 void Run::SumEnergies(G4double edep, G4double 93 { 94 fEnergyDeposit += edep; 95 fEnergyDeposit2 += edep * edep; 96 97 fEnergyFlow += eflow; 98 fEnergyFlow2 += eflow * eflow; 99 100 fEnergyTotal += etot; 101 fEnergyTotal2 += etot * etot; 102 } 103 104 //....oooOO0OOooo........oooOO0OOooo........oo 105 106 void Run::ParticleFlux(G4String name, G4double 107 { 108 std::map<G4String, ParticleData>::iterator i 109 if (it == fParticleDataMap2.end()) { 110 fParticleDataMap2[name] = ParticleData(1, 111 } 112 else { 113 ParticleData& data = it->second; 114 data.fCount++; 115 data.fEmean += Ekin; 116 // update min max 117 G4double emin = data.fEmin; 118 if (Ekin < emin) data.fEmin = Ekin; 119 G4double emax = data.fEmax; 120 if (Ekin > emax) data.fEmax = Ekin; 121 data.fTmean = -1 * ns; 122 } 123 } 124 125 //....oooOO0OOooo........oooOO0OOooo........oo 126 127 void Run::Merge(const G4Run* run) 128 { 129 const Run* localRun = static_cast<const Run* 130 131 // primary particle info 132 // 133 fParticle = localRun->fParticle; 134 fEkin = localRun->fEkin; 135 136 // accumulate sums 137 // 138 fEnergyDeposit += localRun->fEnergyDeposit; 139 fEnergyDeposit2 += localRun->fEnergyDeposit2 140 fEnergyFlow += localRun->fEnergyFlow; 141 fEnergyFlow2 += localRun->fEnergyFlow2; 142 fEnergyTotal += localRun->fEnergyTotal; 143 fEnergyTotal2 += localRun->fEnergyTotal2; 144 145 // map: processes count 146 std::map<G4String, G4int>::const_iterator it 147 for (itp = localRun->fProcCounter.begin(); i 148 G4String procName = itp->first; 149 G4int localCount = itp->second; 150 if (fProcCounter.find(procName) == fProcCo 151 fProcCounter[procName] = localCount; 152 } 153 else { 154 fProcCounter[procName] += localCount; 155 } 156 } 157 158 // map: created particles count 159 std::map<G4String, ParticleData>::const_iter 160 for (itc = localRun->fParticleDataMap1.begin 161 G4String name = itc->first; 162 const ParticleData& localData = itc->secon 163 if (fParticleDataMap1.find(name) == fParti 164 fParticleDataMap1[name] = ParticleData(l 165 l 166 } 167 else { 168 ParticleData& data = fParticleDataMap1[n 169 data.fCount += localData.fCount; 170 data.fEmean += localData.fEmean; 171 G4double emin = localData.fEmin; 172 if (emin < data.fEmin) data.fEmin = emin 173 G4double emax = localData.fEmax; 174 if (emax > data.fEmax) data.fEmax = emax 175 data.fTmean = localData.fTmean; 176 } 177 } 178 179 // map: particles flux count 180 std::map<G4String, ParticleData>::const_iter 181 for (itn = localRun->fParticleDataMap2.begin 182 G4String name = itn->first; 183 const ParticleData& localData = itn->secon 184 if (fParticleDataMap2.find(name) == fParti 185 fParticleDataMap2[name] = ParticleData(l 186 l 187 } 188 else { 189 ParticleData& data = fParticleDataMap2[n 190 data.fCount += localData.fCount; 191 data.fEmean += localData.fEmean; 192 G4double emin = localData.fEmin; 193 if (emin < data.fEmin) data.fEmin = emin 194 G4double emax = localData.fEmax; 195 if (emax > data.fEmax) data.fEmax = emax 196 data.fTmean = localData.fTmean; 197 } 198 } 199 200 G4Run::Merge(run); 201 } 202 203 //....oooOO0OOooo........oooOO0OOooo........oo 204 205 void Run::EndOfRun() 206 { 207 G4int prec = 5, wid = prec + 2; 208 G4int dfprec = G4cout.precision(prec); 209 210 // run condition 211 // 212 G4Material* material = fDetector->GetMateria 213 G4double density = material->GetDensity(); 214 215 G4String Particle = fParticle->GetParticleNa 216 G4cout << "\n The run is " << numberOfEvent 217 << G4BestUnit(fEkin, "Energy") << " t 218 << G4BestUnit(fDetector->GetRadius(), 219 << " (density: " << G4BestUnit(densit 220 221 if (numberOfEvent == 0) { 222 G4cout.precision(dfprec); 223 return; 224 } 225 226 // frequency of processes 227 // 228 G4cout << "\n Process calls frequency :" << 229 G4int index = 0; 230 std::map<G4String, G4int>::iterator it; 231 for (it = fProcCounter.begin(); it != fProcC 232 G4String procName = it->first; 233 G4int count = it->second; 234 G4String space = " "; 235 if (++index % 3 == 0) space = "\n"; 236 G4cout << " " << std::setw(20) << procName 237 } 238 G4cout << G4endl; 239 240 // compute mean Energy deposited and rms 241 // 242 G4int TotNbofEvents = numberOfEvent; 243 fEnergyDeposit /= TotNbofEvents; 244 fEnergyDeposit2 /= TotNbofEvents; 245 G4double rmsEdep = fEnergyDeposit2 - fEnergy 246 if (rmsEdep > 0.) 247 rmsEdep = std::sqrt(rmsEdep); 248 else 249 rmsEdep = 0.; 250 251 G4cout << "\n Mean energy deposit per event 252 << "; rms = " << G4BestUnit(rmsEdep, 253 254 // compute mean Energy leakage and rms 255 // 256 fEnergyFlow /= TotNbofEvents; 257 fEnergyFlow2 /= TotNbofEvents; 258 G4double rmsEflow = fEnergyFlow2 - fEnergyFl 259 if (rmsEflow > 0.) 260 rmsEflow = std::sqrt(rmsEflow); 261 else 262 rmsEflow = 0.; 263 264 G4cout << " Mean energy leakage per event = 265 << "; rms = " << G4BestUnit(rmsEflow 266 267 // energy balance 268 // 269 fEnergyTotal /= TotNbofEvents; 270 fEnergyTotal2 /= TotNbofEvents; 271 G4double rmsEtotal = fEnergyTotal2 - fEnergy 272 if (rmsEtotal > 0.) 273 rmsEtotal = std::sqrt(rmsEtotal); 274 else 275 rmsEflow = 0.; 276 277 G4cout << "\n Mean energy total per event 278 << "; rms = " << G4BestUnit(rmsEtota 279 280 // particles at creation 281 // 282 if (fParticleDataMap1.size() > 0) { 283 G4cout << "\n List of particles at creatio 284 std::map<G4String, ParticleData>::iterator 285 for (itc = fParticleDataMap1.begin(); itc 286 G4String name = itc->first; 287 ParticleData data = itc->second; 288 G4int count = data.fCount; 289 G4double eMean = data.fEmean / count; 290 G4double eMin = data.fEmin; 291 G4double eMax = data.fEmax; 292 G4double meanLife = data.fTmean; 293 294 G4cout << " " << std::setw(13) << name 295 << " Emean = " << std::setw(wid) 296 << G4BestUnit(eMin, "Energy") << 297 if (meanLife >= 0.) 298 G4cout << "\tmean life = " << G4BestUn 299 else 300 G4cout << "\tstable" << G4endl; 301 } 302 } 303 304 // emerging particles 305 // 306 G4cout << "\n List of particles emerging fro 307 308 std::map<G4String, ParticleData>::iterator i 309 for (itn = fParticleDataMap2.begin(); itn != 310 G4String name = itn->first; 311 ParticleData data = itn->second; 312 G4int count = data.fCount; 313 G4double eMean = data.fEmean / count; 314 G4double eMin = data.fEmin; 315 G4double eMax = data.fEmax; 316 G4double Eflow = data.fEmean / TotNbofEven 317 318 G4cout << " " << std::setw(13) << name << 319 << " Emean = " << std::setw(wid) < 320 << G4BestUnit(eMin, "Energy") << " 321 << ") \tEleak/event = " << G4BestUn 322 } 323 324 // normalize histograms 325 G4AnalysisManager* analysisManager = G4Analy 326 327 G4int ih = 2; 328 G4double binWidth = analysisManager->GetH1Wi 329 G4double fac = (1. / (numberOfEvent * binWid 330 analysisManager->ScaleH1(ih, fac); 331 332 // remove all contents in fProcCounter, fCou 333 fProcCounter.clear(); 334 fParticleDataMap2.clear(); 335 336 // restore default format 337 G4cout.precision(dfprec); 338 } 339 340 //....oooOO0OOooo........oooOO0OOooo........oo 341