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 // 27 28 // Created by Z. Francis 29 30 #include "G4DNASancheExcitationModel.hh" 31 #include "G4SystemOfUnits.hh" 32 #include "G4DNAMolecularMaterial.hh" 33 34 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 35 36 using namespace std; 37 38 //#define SANCHE_VERBOSE 39 40 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 41 42 G4DNASancheExcitationModel::G4DNASancheExcitationModel(const G4ParticleDefinition*, 43 const G4String& nam) : 44 G4VEmModel(nam) 45 { 46 fpWaterDensity = nullptr; 47 48 SetLowEnergyLimit(2.*eV); 49 SetHighEnergyLimit(100*eV); 50 nLevels = 9; 51 52 verboseLevel = 0; 53 // Verbosity scale: 54 // 0 = nothing 55 // 1 = warning for energy non-conservation 56 // 2 = details of energy budget 57 // 3 = calculation of cross sections, file openings, sampling of atoms 58 // 4 = entering in methods 59 60 #ifdef SANCHE_VERBOSE 61 if (verboseLevel > 0) 62 { 63 G4cout << "Sanche Excitation model is constructed " 64 << G4endl 65 << "Energy range: " 66 << LowEnergyLimit() / eV << " eV - " 67 << HighEnergyLimit() / eV << " eV" 68 << G4endl; 69 } 70 #endif 71 72 fParticleChangeForGamma = nullptr; 73 fpWaterDensity = nullptr; 74 75 // Selection of stationary mode 76 77 statCode = false; 78 } 79 80 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 81 82 G4DNASancheExcitationModel::~G4DNASancheExcitationModel() 83 = default; 84 85 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 86 87 void 88 G4DNASancheExcitationModel:: 89 Initialise(const G4ParticleDefinition* /*particle*/, 90 const G4DataVector& /*cuts*/) 91 { 92 93 #ifdef SANCHE_VERBOSE 94 if (verboseLevel > 3) 95 { 96 G4cout << "Calling G4DNASancheExcitationModel::Initialise()" 97 << G4endl; 98 } 99 #endif 100 101 // Energy limits 102 103 if (LowEnergyLimit() < 2.*eV) 104 { 105 G4Exception("*** WARNING : the G4DNASancheExcitationModel class is not " 106 "validated below 2 eV !", 107 "", JustWarning, ""); 108 } 109 110 if (HighEnergyLimit() > 100.*eV) 111 { 112 G4cout << "G4DNASancheExcitationModel: high energy limit decreased from " << 113 HighEnergyLimit()/eV << " eV to " << 100. << " eV" << G4endl; 114 SetHighEnergyLimit(100.*eV); 115 } 116 117 // 118 #ifdef SANCHE_VERBOSE 119 if (verboseLevel > 0) 120 { 121 G4cout << "Sanche Excitation model is initialized " << G4endl 122 << "Energy range: " 123 << LowEnergyLimit() / eV << " eV - " 124 << HighEnergyLimit() / eV << " eV" 125 << G4endl; 126 } 127 #endif 128 129 // Initialize water density pointer 130 fpWaterDensity = G4DNAMolecularMaterial::Instance()-> 131 GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER")); 132 133 if (isInitialised) {return;} 134 135 fParticleChangeForGamma = GetParticleChangeForGamma(); 136 isInitialised = true; 137 138 const char *path = G4FindDataDir("G4LEDATA"); 139 std::ostringstream eFullFileName; 140 eFullFileName << path << "/dna/sigma_excitationvib_e_sanche.dat"; 141 std::ifstream input(eFullFileName.str().c_str()); 142 143 if (!input) 144 { 145 G4Exception("G4DNASancheExcitationModel::Initialise","em0003", 146 FatalException,"Missing data file:/dna/sigma_excitationvib_e_sanche.dat"); 147 } 148 149 // March 25th, 2014 - Vaclav Stepan, Sebastien Incerti 150 // Added clear for MT 151 tdummyVec.clear(); 152 // 153 154 G4double t; 155 G4double xs; 156 157 while(!input.eof()) 158 { 159 input>>t; 160 tdummyVec.push_back(t); 161 162 fEnergyLevelXS.emplace_back(); 163 fEnergyTotalXS.push_back(0); 164 std::vector<G4double>& levelXS = fEnergyLevelXS.back(); 165 levelXS.reserve(9); 166 167 // G4cout<<t; 168 169 for(size_t i = 0 ; i < 9 ;++i) 170 { 171 input>>xs; 172 levelXS.push_back(xs); 173 fEnergyTotalXS.back() += xs; 174 // G4cout <<" " << levelXS[i]; 175 } 176 177 // G4cout << G4endl; 178 } 179 } 180 181 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 182 183 G4double G4DNASancheExcitationModel::CrossSectionPerVolume(const G4Material* material, 184 const G4ParticleDefinition* 185 #ifdef SANCHE_VERBOSE 186 particleDefinition 187 #endif 188 , 189 G4double ekin, 190 G4double, 191 G4double) 192 { 193 #ifdef SANCHE_VERBOSE 194 if (verboseLevel > 3) 195 { 196 G4cout << "Calling CrossSectionPerVolume() of G4DNASancheExcitationModel" 197 << G4endl; 198 } 199 #endif 200 201 // Calculate total cross section for model 202 203 G4double sigma = 0.; 204 205 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()]; 206 207 if (ekin >= LowEnergyLimit() && ekin <= HighEnergyLimit()) 208 sigma = TotalCrossSection(ekin); 209 210 #ifdef SANCHE_VERBOSE 211 if (verboseLevel > 2) 212 { 213 G4cout << "__________________________________" << G4endl; 214 G4cout << "=== G4DNASancheExcitationModel - XS INFO START" << G4endl; 215 G4cout << "=== Kinetic energy(eV)=" << ekin/eV << " particle : " << particleDefinition->GetParticleName() << G4endl; 216 G4cout << "=== Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl; 217 G4cout << "=== Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl; 218 G4cout << "=== G4DNASancheExcitationModel - XS INFO END" << G4endl; 219 } 220 #endif 221 222 return sigma*2.*waterDensity; 223 // see papers for factor 2 description (liquid phase) 224 } 225 226 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 227 228 void G4DNASancheExcitationModel::SampleSecondaries(std::vector< 229 G4DynamicParticle*>*, 230 const G4MaterialCutsCouple*, 231 const G4DynamicParticle* aDynamicElectron, 232 G4double, 233 G4double) 234 { 235 #ifdef SANCHE_VERBOSE 236 if (verboseLevel > 3) 237 { 238 G4cout << "Calling SampleSecondaries() of G4DNASancheExcitationModel" 239 << G4endl; 240 } 241 #endif 242 243 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy(); 244 G4int level = RandomSelect(electronEnergy0); 245 G4double excitationEnergy = VibrationEnergy(level); // levels go from 0 to 8 246 G4double newEnergy = electronEnergy0 - excitationEnergy; 247 248 /* 249 if (electronEnergy0 < highEnergyLimit) 250 { 251 if (newEnergy >= lowEnergyLimit) 252 { 253 fParticleChangeForGamma->ProposeMomentumDirection(aDynamicElectron->GetMomentumDirection()); 254 fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy); 255 fParticleChangeForGamma->ProposeLocalEnergyDeposit(excitationEnergy); 256 } 257 258 else 259 { 260 fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill); 261 fParticleChangeForGamma->ProposeLocalEnergyDeposit(electronEnergy0); 262 } 263 } 264 */ 265 266 if (electronEnergy0 <= HighEnergyLimit() && newEnergy>0.) 267 { 268 269 if (!statCode) 270 { 271 fParticleChangeForGamma->ProposeMomentumDirection(aDynamicElectron->GetMomentumDirection()); 272 fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy); 273 fParticleChangeForGamma->ProposeLocalEnergyDeposit(excitationEnergy); 274 } 275 276 else 277 { 278 fParticleChangeForGamma->ProposeMomentumDirection(aDynamicElectron->GetMomentumDirection()); 279 fParticleChangeForGamma->SetProposedKineticEnergy(electronEnergy0); 280 fParticleChangeForGamma->ProposeLocalEnergyDeposit(excitationEnergy); 281 } 282 283 } 284 285 // 286 } 287 288 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 289 290 G4double G4DNASancheExcitationModel::PartialCrossSection(G4double t, 291 G4int level) 292 { 293 // Protection against out of boundary access 294 if (t/eV==tdummyVec.back()) t=t*(1.-1e-12); 295 // 296 297 auto t2 = std::upper_bound(tdummyVec.begin(), 298 tdummyVec.end(), t / eV); 299 auto t1 = t2 - 1; 300 301 size_t i1 = t1 - tdummyVec.begin(); 302 size_t i2 = t2 - tdummyVec.begin(); 303 304 G4double sigma = LinInterpolate((*t1), (*t2), 305 t / eV, 306 fEnergyLevelXS[i1][level], 307 fEnergyLevelXS[i2][level]); 308 309 static const G4double conv_factor = 1e-16 * cm * cm; 310 311 sigma *= conv_factor; 312 if (sigma == 0.) sigma = 1e-30; 313 return (sigma); 314 } 315 316 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 317 318 G4double G4DNASancheExcitationModel::TotalCrossSection(G4double t) 319 { 320 // Protection against out of boundary access 321 if (t/eV==tdummyVec.back()) t=t*(1.-1e-12); 322 // 323 324 auto t2 = std::upper_bound(tdummyVec.begin(), 325 tdummyVec.end(), t / eV); 326 auto t1 = t2 - 1; 327 328 size_t i1 = t1 - tdummyVec.begin(); 329 size_t i2 = t2 - tdummyVec.begin(); 330 331 G4double sigma = LinInterpolate((*t1), (*t2), 332 t / eV, 333 fEnergyTotalXS[i1], 334 fEnergyTotalXS[i2]); 335 336 static const G4double conv_factor = 1e-16 * cm * cm; 337 338 sigma *= conv_factor; 339 if (sigma == 0.) sigma = 1e-30; 340 return (sigma); 341 } 342 343 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 344 345 G4double G4DNASancheExcitationModel::VibrationEnergy(G4int level) 346 { 347 static G4double energies[9] = { 0.01, 0.024, 0.061, 0.092, 0.204, 0.417, 0.460, 348 0.500, 0.835 }; 349 return (energies[level] * eV); 350 } 351 352 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 353 354 G4int G4DNASancheExcitationModel::RandomSelect(G4double k) 355 { 356 357 // Level Selection Counting can be done here ! 358 359 G4int i = nLevels; 360 G4double value = 0.; 361 std::deque<G4double> values; 362 363 while (i > 0) 364 { 365 i--; 366 G4double partial = PartialCrossSection(k, i); 367 values.push_front(partial); 368 value += partial; 369 } 370 371 value *= G4UniformRand(); 372 373 i = nLevels; 374 375 while (i > 0) 376 { 377 i--; 378 if (values[i] > value) 379 { 380 //outcount<<i<<" "<<VibrationEnergy(i)<<G4endl; 381 return i; 382 } 383 value -= values[i]; 384 } 385 386 //outcount<<0<<" "<<VibrationEnergy(0)<<G4endl; 387 388 return 0; 389 } 390 391 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 392 393 G4double G4DNASancheExcitationModel::Sum(G4double k) 394 { 395 G4double totalCrossSection = 0.; 396 397 for (G4int i = 0; i < nLevels; i++) 398 { 399 totalCrossSection += PartialCrossSection(k, i); 400 } 401 402 return totalCrossSection; 403 } 404 405 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 406 407 G4double G4DNASancheExcitationModel::LinInterpolate(G4double e1, 408 G4double e2, 409 G4double e, 410 G4double xs1, 411 G4double xs2) 412 { 413 G4double a = (xs2 - xs1) / (e2 - e1); 414 G4double b = xs2 - a * e2; 415 G4double value = a * e + b; 416 // G4cout<<"interP >> "<<e1<<" "<<e2<<" "<<e<<" " 417 // <<xs1<<" "<<xs2<<" "<<a<<" "<<b<<" "<<value<<G4endl; 418 419 return value; 420 } 421 422