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/svalue/src/Run.cc 37 /// \brief Implementation of the Run class 38 39 #include "Run.hh" 40 41 #include "HistoManager.hh" 42 #include "MyFile.hh" 43 44 #ifdef MYFILE 45 # include "MyPrimaryGeneratorActionFromFile.hh" 46 #else 47 # include "PrimaryGeneratorAction.hh" 48 #endif 49 50 #include "G4SystemOfUnits.hh" 51 #include "G4UnitsTable.hh" 52 53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 54 55 Run::Run(const DetectorConstruction* detector) 56 : G4Run(), 57 fDetector(detector), 58 fParticle(0), 59 fEkin(0.), 60 fCytoEdeposit(0.), 61 fCytoEdeposit2(0.), 62 fNuclEdeposit(0.), 63 fNuclEdeposit2(0.), 64 fTrackLen(0.), 65 fTrackLen2(0.), 66 fProjRange(0.), 67 fProjRange2(0.), 68 fNbOfSteps(0), 69 fNbOfSteps2(0), 70 fStepSize(0.), 71 fStepSize2(0.) 72 {} 73 74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 75 76 Run::~Run() {} 77 78 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 79 80 void Run::SetPrimary(G4ParticleDefinition* particle, G4double energy) 81 { 82 fParticle = particle; 83 fEkin = energy; 84 } 85 86 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 87 88 void Run::AddCytoEdep(G4double e) 89 { 90 fCytoEdeposit += e; 91 fCytoEdeposit2 += e * e; 92 } 93 94 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 95 96 void Run::AddNuclEdep(G4double e) 97 { 98 fNuclEdeposit += e; 99 fNuclEdeposit2 += e * e; 100 } 101 102 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 103 104 void Run::AddTrackLength(G4double t) 105 { 106 fTrackLen += t; 107 fTrackLen2 += t * t; 108 } 109 110 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 111 112 void Run::AddProjRange(G4double x) 113 { 114 fProjRange += x; 115 fProjRange2 += x * x; 116 } 117 118 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 119 120 void Run::AddStepSize(G4int nb, G4double st) 121 { 122 fNbOfSteps += nb; 123 fNbOfSteps2 += nb * nb; 124 fStepSize += st; 125 fStepSize2 += st * st; 126 } 127 128 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 129 130 void Run::Merge(const G4Run* run) 131 { 132 const Run* localRun = static_cast<const Run*>(run); 133 134 // Pass information about primary particle 135 136 fParticle = localRun->fParticle; 137 fEkin = localRun->fEkin; 138 139 // Accumulate sums 140 141 fCytoEdeposit += localRun->fCytoEdeposit; 142 fCytoEdeposit2 += localRun->fCytoEdeposit2; 143 fNuclEdeposit += localRun->fNuclEdeposit; 144 fNuclEdeposit2 += localRun->fNuclEdeposit2; 145 146 fTrackLen += localRun->fTrackLen; 147 fTrackLen2 += localRun->fTrackLen2; 148 fProjRange += localRun->fProjRange; 149 fProjRange2 += localRun->fProjRange2; 150 fNbOfSteps += localRun->fNbOfSteps; 151 fNbOfSteps2 += localRun->fNbOfSteps2; 152 fStepSize += localRun->fStepSize; 153 fStepSize2 += localRun->fStepSize2; 154 155 G4Run::Merge(run); 156 } 157 158 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 159 160 void Run::EndOfRun() 161 { 162 std::ios::fmtflags mode = G4cout.flags(); 163 G4cout.setf(std::ios::fixed, std::ios::floatfield); 164 G4int prec = G4cout.precision(2); 165 166 // Run conditions 167 168 G4String partName = fParticle->GetParticleName(); 169 170 G4cout << "\n ======================== run summary =====================\n"; 171 G4cout << "\n The run is " << numberOfEvent << " " << partName << " of " 172 << G4BestUnit(fEkin, "Energy") << G4endl; 173 174 if (numberOfEvent == 0) { 175 G4cout.setf(mode, std::ios::floatfield); 176 G4cout.precision(prec); 177 return; 178 } 179 180 // Compute S-value for cytoplasm (C<-C) 181 182 fCytoEdeposit /= numberOfEvent; 183 fCytoEdeposit2 /= numberOfEvent; 184 G4double rmsCyto = fCytoEdeposit2 - fCytoEdeposit * fCytoEdeposit; 185 if (rmsCyto > 0.) 186 rmsCyto = std::sqrt(rmsCyto); 187 else 188 rmsCyto = 0.; 189 190 G4cout.precision(3); 191 G4cout << "\n Total Energy deposited in cytoplasm = " << G4BestUnit(fCytoEdeposit, "Energy") 192 << " +- " << G4BestUnit(rmsCyto, "Energy") << G4endl; 193 194 G4double sValueCyto = fCytoEdeposit / fDetector->GetCytoMass(); 195 G4double rmsSValueCyto = rmsCyto / fDetector->GetCytoMass(); 196 197 G4cout.precision(3); 198 G4cout << "\n S value for cytoplasm (C<-C) = " << sValueCyto / gray << " Gy/Bq.s " 199 << " +- " << rmsSValueCyto / gray << " Gy/Bq.s " << G4endl; 200 201 // Compute S-value for nucleus (N<-C) 202 203 fNuclEdeposit /= numberOfEvent; 204 fNuclEdeposit2 /= numberOfEvent; 205 G4double rmsNucl = fNuclEdeposit2 - fNuclEdeposit * fNuclEdeposit; 206 if (rmsNucl > 0.) 207 rmsNucl = std::sqrt(rmsNucl); 208 else 209 rmsNucl = 0.; 210 211 G4cout.precision(3); 212 G4cout << "\n Total Energy deposited in nucleus = " << G4BestUnit(fNuclEdeposit, "Energy") 213 << " +- " << G4BestUnit(rmsNucl, "Energy") << G4endl; 214 215 G4double sValueNucl = fNuclEdeposit / fDetector->GetNuclMass(); 216 G4double rmsSValueNucl = rmsNucl / fDetector->GetNuclMass(); 217 218 G4cout.precision(3); 219 G4cout << "\n S value for nucleus (N<-C) = " << sValueNucl / gray << " Gy/Bq.s " 220 << " +- " << rmsSValueNucl / gray << " Gy/Bq.s " << G4endl; 221 222 // Compute track length of primary track 223 224 fTrackLen /= numberOfEvent; 225 fTrackLen2 /= numberOfEvent; 226 G4double rms = fTrackLen2 - fTrackLen * fTrackLen; 227 if (rms > 0.) 228 rms = std::sqrt(rms); 229 else 230 rms = 0.; 231 232 G4cout.precision(3); 233 G4cout << "\n Track length of primary track = " << G4BestUnit(fTrackLen, "Length") << " +- " 234 << G4BestUnit(rms, "Length"); 235 236 // Compute projected range of primary track 237 238 fProjRange /= numberOfEvent; 239 fProjRange2 /= numberOfEvent; 240 rms = fProjRange2 - fProjRange * fProjRange; 241 if (rms > 0.) 242 rms = std::sqrt(rms); 243 else 244 rms = 0.; 245 246 G4cout << "\n Projected range = " << G4BestUnit(fProjRange, "Length") << " +- " 247 << G4BestUnit(rms, "Length") << G4endl; 248 249 // Nb of steps and step size of primary track 250 251 G4double dNofEvents = double(numberOfEvent); 252 G4double fNbSteps = fNbOfSteps / dNofEvents, fNbSteps2 = fNbOfSteps2 / dNofEvents; 253 rms = fNbSteps2 - fNbSteps * fNbSteps; 254 if (rms > 0.) 255 rms = std::sqrt(rms); 256 else 257 rms = 0.; 258 259 G4cout.precision(2); 260 G4cout << "\n Nb of steps of primary track = " << fNbSteps << " +- " << rms << G4endl; 261 262 fStepSize /= numberOfEvent; 263 fStepSize2 /= numberOfEvent; 264 rms = fStepSize2 - fStepSize * fStepSize; 265 if (rms > 0.) 266 rms = std::sqrt(rms); 267 else 268 rms = 0.; 269 270 G4cout.precision(3); 271 G4cout << "\n Step size = " << G4BestUnit(fStepSize, "Length") << " +- " 272 << G4BestUnit(rms, "Length") << G4endl; 273 274 // Normalize histograms of longitudinal energy profile 275 276 G4AnalysisManager* analysisManager = G4AnalysisManager::Instance(); 277 G4int ih = 1; 278 G4double binWidth = analysisManager->GetH1Width(ih); 279 G4double fac = (1. / (numberOfEvent * binWidth)) * (mm / MeV); 280 analysisManager->ScaleH1(ih, fac); 281 282 // Reset default formats 283 284 G4cout.setf(mode, std::ios::floatfield); 285 G4cout.precision(prec); 286 287 // Output file 288 289 FILE* myFile; 290 myFile = fopen("s.txt", "a"); 291 fprintf(myFile, "%e %e %e %e %e %e %e \n", fDetector->GetNuclRadius() / nm, 292 fDetector->GetCytoThickness() / nm, fEkin / eV, sValueCyto / gray, rmsSValueCyto / gray, 293 sValueNucl / gray, rmsSValueNucl / gray); 294 fclose(myFile); 295 } 296