Geant4 Cross Reference |
1 // 2 // ******************************************************************** 3 // * License and Disclaimer * 4 // * * 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. * 10 // * * 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitation of liability. * 17 // * * 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************************************** 25 // 26 // This example is provided by the Geant4-DNA collaboration 27 // Any report or published results obtained using the Geant4-DNA software 28 // shall cite the following Geant4-DNA collaboration publications: 29 // Med. Phys. 45 (2018) e722-e739 30 // Phys. Med. 31 (2015) 861-874 31 // Med. Phys. 37 (2010) 4692-4708 32 // Int. J. Model. Simul. Sci. Comput. 1 (2010) 157–178 33 // 34 // The Geant4-DNA web site is available at http://geant4-dna.org 35 // 36 /// \file medical/dna/wvalue/src/Run.cc 37 /// \brief Implementation of the Run class 38 39 #include "Run.hh" 40 41 #include "DetectorConstruction.hh" 42 #include "HistoManager.hh" 43 #include "PrimaryGeneratorAction.hh" 44 45 #include "G4Material.hh" 46 #include "G4SystemOfUnits.hh" 47 #include "G4UnitsTable.hh" 48 49 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 50 51 Run::Run(const DetectorConstruction* detector) 52 : G4Run(), 53 fDetector(detector), 54 fParticle(0), 55 fEkin(0.), 56 fNbInelastic(0), 57 fNbInelastic2(0), 58 fEdeposit(0.), 59 fEdeposit2(0.), 60 fTrackLen(0.), 61 fTrackLen2(0.), 62 fProjRange(0.), 63 fProjRange2(0.), 64 fNbOfSteps(0), 65 fNbOfSteps2(0), 66 fStepSize(0.), 67 fStepSize2(0.) 68 {} 69 70 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 71 72 Run::~Run() {} 73 74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 75 76 void Run::SetPrimary(G4ParticleDefinition* particle, G4double energy) 77 { 78 fParticle = particle; 79 fEkin = energy; 80 } 81 82 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 83 84 void Run::AddInelastic(G4int nb) 85 { 86 fNbInelastic += nb; 87 fNbInelastic2 += nb * nb; 88 } 89 90 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 91 92 void Run::AddEdep(G4double e) 93 { 94 fEdeposit += e; 95 fEdeposit2 += e * e; 96 } 97 98 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 99 100 void Run::AddTrackLength(G4double t) 101 { 102 fTrackLen += t; 103 fTrackLen2 += t * t; 104 } 105 106 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 107 108 void Run::AddProjRange(G4double x) 109 { 110 fProjRange += x; 111 fProjRange2 += x * x; 112 } 113 114 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 115 116 void Run::AddStepSize(G4int nb, G4double st) 117 { 118 fNbOfSteps += nb; 119 fNbOfSteps2 += nb * nb; 120 fStepSize += st; 121 fStepSize2 += st * st; 122 } 123 124 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 125 126 void Run::Merge(const G4Run* run) 127 { 128 const Run* localRun = static_cast<const Run*>(run); 129 130 // Pass information about primary particle 131 132 fParticle = localRun->fParticle; 133 fEkin = localRun->fEkin; 134 135 // Accumulate sums 136 137 fNbInelastic += localRun->fNbInelastic; 138 fNbInelastic2 += localRun->fNbInelastic2; 139 fEdeposit += localRun->fEdeposit; 140 fEdeposit2 += localRun->fEdeposit2; 141 fTrackLen += localRun->fTrackLen; 142 fTrackLen2 += localRun->fTrackLen2; 143 fProjRange += localRun->fProjRange; 144 fProjRange2 += localRun->fProjRange2; 145 fNbOfSteps += localRun->fNbOfSteps; 146 fNbOfSteps2 += localRun->fNbOfSteps2; 147 fStepSize += localRun->fStepSize; 148 fStepSize2 += localRun->fStepSize2; 149 150 G4Run::Merge(run); 151 } 152 153 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 154 155 void Run::EndOfRun() 156 { 157 std::ios::fmtflags mode = G4cout.flags(); 158 G4cout.setf(std::ios::fixed, std::ios::floatfield); 159 G4int prec = G4cout.precision(2); 160 161 // Run conditions 162 163 G4Material* material = fDetector->GetAbsorMaterial(); 164 G4double density = material->GetDensity(); 165 G4String partName = fParticle->GetParticleName(); 166 167 G4cout << "\n ======================== run summary =====================\n"; 168 G4cout << "\n The run is " << numberOfEvent << " " << partName << " of " 169 << G4BestUnit(fEkin, "Energy") << " through a sphere of radius " 170 << G4BestUnit(fDetector->GetAbsorRadius(), "Length") << "of " << material->GetName() 171 << " (density: " << G4BestUnit(density, "Volumic Mass") << ")" << G4endl; 172 173 if (numberOfEvent == 0) { 174 G4cout.setf(mode, std::ios::floatfield); 175 G4cout.precision(prec); 176 return; 177 } 178 179 fNbInelastic /= numberOfEvent; 180 fNbInelastic2 /= numberOfEvent; 181 182 G4double rms = fNbInelastic2 - fNbInelastic * fNbInelastic; 183 if (rms > 0.) 184 rms = std::sqrt(rms); 185 else 186 rms = 0.; 187 188 G4cout.precision(3); 189 G4cout << "\n Nb of ionisations = " << fNbInelastic << " +- " << rms << G4endl; 190 191 G4cout.precision(3); 192 G4cout << "\n w = " << G4BestUnit((fEkin) / fNbInelastic, "Energy") << " +- " 193 << G4BestUnit((fEkin)*rms / (fNbInelastic * fNbInelastic), "Energy") << G4endl; 194 195 // Output file 196 197 if (fNbInelastic > 0.) { 198 FILE* myFile; 199 myFile = fopen("wvalue.txt", "a"); 200 fprintf(myFile, "%e %e %e %e %e \n", fEkin / eV, fNbInelastic, rms, fEkin / eV / fNbInelastic, 201 (fEkin / eV) * rms / (fNbInelastic * fNbInelastic)); 202 fclose(myFile); 203 } 204 // 205 206 fEdeposit /= numberOfEvent; 207 fEdeposit2 /= numberOfEvent; 208 rms = fEdeposit2 - fEdeposit * fEdeposit; 209 if (rms > 0.) 210 rms = std::sqrt(rms); 211 else 212 rms = 0.; 213 214 G4cout.precision(3); 215 G4cout << "\n Total Energy deposited = " << G4BestUnit(fEdeposit, "Energy") << " +- " 216 << G4BestUnit(rms, "Energy") << G4endl; 217 218 // Compute track length of primary track 219 220 fTrackLen /= numberOfEvent; 221 fTrackLen2 /= numberOfEvent; 222 rms = fTrackLen2 - fTrackLen * fTrackLen; 223 if (rms > 0.) 224 rms = std::sqrt(rms); 225 else 226 rms = 0.; 227 228 G4cout.precision(3); 229 G4cout << "\n Track length of primary track = " << G4BestUnit(fTrackLen, "Length") << " +- " 230 << G4BestUnit(rms, "Length"); 231 232 // Compute projected range of primary track 233 234 fProjRange /= numberOfEvent; 235 fProjRange2 /= numberOfEvent; 236 rms = fProjRange2 - fProjRange * fProjRange; 237 if (rms > 0.) 238 rms = std::sqrt(rms); 239 else 240 rms = 0.; 241 242 G4cout << "\n Projected range = " << G4BestUnit(fProjRange, "Length") << " +- " 243 << G4BestUnit(rms, "Length") << G4endl; 244 245 // Nb of steps and step size of primary track 246 247 G4double dNofEvents = double(numberOfEvent); 248 G4double fNbSteps = fNbOfSteps / dNofEvents, fNbSteps2 = fNbOfSteps2 / dNofEvents; 249 rms = fNbSteps2 - fNbSteps * fNbSteps; 250 if (rms > 0.) 251 rms = std::sqrt(rms); 252 else 253 rms = 0.; 254 255 G4cout.precision(2); 256 G4cout << "\n Nb of steps of primary track = " << fNbSteps << " +- " << rms << G4endl; 257 258 fStepSize /= numberOfEvent; 259 fStepSize2 /= numberOfEvent; 260 rms = fStepSize2 - fStepSize * fStepSize; 261 if (rms > 0.) 262 rms = std::sqrt(rms); 263 else 264 rms = 0.; 265 266 G4cout.precision(3); 267 G4cout << "\n Step size = " << G4BestUnit(fStepSize, "Length") << " +- " 268 << G4BestUnit(rms, "Length") << G4endl; 269 270 // Normalize histograms of longitudinal energy profile 271 272 G4AnalysisManager* analysisManager = G4AnalysisManager::Instance(); 273 G4int ih = 1; 274 G4double binWidth = analysisManager->GetH1Width(ih); 275 G4double fac = (1. / (numberOfEvent * binWidth)) * (mm / MeV); 276 analysisManager->ScaleH1(ih, fac); 277 278 // Reset default formats 279 280 G4cout.setf(mode, std::ios::floatfield); 281 G4cout.precision(prec); 282 } 283