[1058] | 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 | // |
---|
[1192] | 26 | // $Id: G4DNADingfelderChargeDecreaseModel.cc,v 1.6 2009/08/13 11:32:47 sincerti Exp $ |
---|
[1228] | 27 | // GEANT4 tag $Name: geant4-09-03 $ |
---|
[1058] | 28 | // |
---|
| 29 | |
---|
| 30 | #include "G4DNADingfelderChargeDecreaseModel.hh" |
---|
| 31 | |
---|
| 32 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 33 | |
---|
| 34 | using namespace std; |
---|
| 35 | |
---|
| 36 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 37 | |
---|
| 38 | G4DNADingfelderChargeDecreaseModel::G4DNADingfelderChargeDecreaseModel(const G4ParticleDefinition*, |
---|
| 39 | const G4String& nam) |
---|
| 40 | :G4VEmModel(nam),isInitialised(false) |
---|
| 41 | { |
---|
| 42 | |
---|
| 43 | verboseLevel= 0; |
---|
| 44 | // Verbosity scale: |
---|
| 45 | // 0 = nothing |
---|
| 46 | // 1 = warning for energy non-conservation |
---|
| 47 | // 2 = details of energy budget |
---|
| 48 | // 3 = calculation of cross sections, file openings, sampling of atoms |
---|
| 49 | // 4 = entering in methods |
---|
| 50 | |
---|
[1192] | 51 | if( verboseLevel>0 ) |
---|
| 52 | { |
---|
| 53 | G4cout << "Dingfelder charge decrease model is constructed " << G4endl; |
---|
| 54 | } |
---|
[1058] | 55 | } |
---|
| 56 | |
---|
| 57 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 58 | |
---|
| 59 | G4DNADingfelderChargeDecreaseModel::~G4DNADingfelderChargeDecreaseModel() |
---|
| 60 | {} |
---|
| 61 | |
---|
| 62 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 63 | |
---|
| 64 | void G4DNADingfelderChargeDecreaseModel::Initialise(const G4ParticleDefinition* particle, |
---|
| 65 | const G4DataVector& /*cuts*/) |
---|
| 66 | { |
---|
| 67 | |
---|
| 68 | if (verboseLevel > 3) |
---|
| 69 | G4cout << "Calling G4DNADingfelderChargeDecreaseModel::Initialise()" << G4endl; |
---|
| 70 | |
---|
| 71 | // Energy limits |
---|
| 72 | |
---|
| 73 | G4DNAGenericIonsManager *instance; |
---|
| 74 | instance = G4DNAGenericIonsManager::Instance(); |
---|
| 75 | G4ParticleDefinition* protonDef = G4Proton::ProtonDefinition(); |
---|
| 76 | G4ParticleDefinition* alphaPlusPlusDef = instance->GetIon("alpha++"); |
---|
| 77 | G4ParticleDefinition* alphaPlusDef = instance->GetIon("alpha+"); |
---|
| 78 | |
---|
| 79 | G4String proton; |
---|
| 80 | G4String alphaPlusPlus; |
---|
| 81 | G4String alphaPlus; |
---|
| 82 | |
---|
| 83 | if (protonDef != 0) |
---|
| 84 | { |
---|
| 85 | proton = protonDef->GetParticleName(); |
---|
| 86 | lowEnergyLimit[proton] = 1. * keV; |
---|
| 87 | highEnergyLimit[proton] = 10. * MeV; |
---|
| 88 | } |
---|
| 89 | else |
---|
| 90 | { |
---|
| 91 | G4Exception("G4DNADingfelderChargeDecreaseModel::Initialise: proton is not defined"); |
---|
| 92 | } |
---|
| 93 | |
---|
| 94 | if (alphaPlusPlusDef != 0) |
---|
| 95 | { |
---|
| 96 | alphaPlusPlus = alphaPlusPlusDef->GetParticleName(); |
---|
| 97 | lowEnergyLimit[alphaPlusPlus] = 1. * keV; |
---|
| 98 | highEnergyLimit[alphaPlusPlus] = 10. * MeV; |
---|
| 99 | } |
---|
| 100 | else |
---|
| 101 | { |
---|
| 102 | G4Exception("G4DNADingfelderChargeDecreaseModel::Initialise: alphaPlusPlus is not defined"); |
---|
| 103 | } |
---|
| 104 | |
---|
| 105 | if (alphaPlusDef != 0) |
---|
| 106 | { |
---|
| 107 | alphaPlus = alphaPlusDef->GetParticleName(); |
---|
| 108 | lowEnergyLimit[alphaPlus] = 1. * keV; |
---|
| 109 | highEnergyLimit[alphaPlus] = 10. * MeV; |
---|
| 110 | } |
---|
| 111 | else |
---|
| 112 | { |
---|
| 113 | G4Exception("G4DNADingfelderChargeDecreaseModel::Initialise: alphaPlus is not defined"); |
---|
| 114 | } |
---|
| 115 | |
---|
| 116 | if (particle==protonDef) |
---|
| 117 | { |
---|
| 118 | SetLowEnergyLimit(lowEnergyLimit[proton]); |
---|
| 119 | SetHighEnergyLimit(highEnergyLimit[proton]); |
---|
| 120 | } |
---|
| 121 | |
---|
| 122 | if (particle==alphaPlusPlusDef) |
---|
| 123 | { |
---|
| 124 | SetLowEnergyLimit(lowEnergyLimit[alphaPlusPlus]); |
---|
| 125 | SetHighEnergyLimit(highEnergyLimit[alphaPlusPlus]); |
---|
| 126 | } |
---|
| 127 | |
---|
| 128 | if (particle==alphaPlusDef) |
---|
| 129 | { |
---|
| 130 | SetLowEnergyLimit(lowEnergyLimit[alphaPlus]); |
---|
| 131 | SetHighEnergyLimit(highEnergyLimit[alphaPlus]); |
---|
| 132 | } |
---|
| 133 | |
---|
| 134 | // Final state |
---|
| 135 | |
---|
| 136 | //PROTON |
---|
| 137 | f0[0][0]=1.; |
---|
| 138 | a0[0][0]=-0.180; |
---|
| 139 | a1[0][0]=-3.600; |
---|
| 140 | b0[0][0]=-18.22; |
---|
| 141 | b1[0][0]=-1.997; |
---|
| 142 | c0[0][0]=0.215; |
---|
| 143 | d0[0][0]=3.550; |
---|
| 144 | x0[0][0]=3.450; |
---|
| 145 | x1[0][0]=5.251; |
---|
| 146 | |
---|
| 147 | numberOfPartialCrossSections[0] = 1; |
---|
| 148 | |
---|
| 149 | //ALPHA++ |
---|
| 150 | f0[0][1]=1.; a0[0][1]=0.95; |
---|
| 151 | a1[0][1]=-2.75; |
---|
| 152 | b0[0][1]=-23.00; |
---|
| 153 | c0[0][1]=0.215; |
---|
| 154 | d0[0][1]=2.95; |
---|
| 155 | x0[0][1]=3.50; |
---|
| 156 | |
---|
| 157 | f0[1][1]=1.; |
---|
| 158 | a0[1][1]=0.95; |
---|
| 159 | a1[1][1]=-2.75; |
---|
| 160 | b0[1][1]=-23.73; |
---|
| 161 | c0[1][1]=0.250; |
---|
| 162 | d0[1][1]=3.55; |
---|
| 163 | x0[1][1]=3.72; |
---|
| 164 | |
---|
| 165 | x1[0][1]=-1.; |
---|
| 166 | b1[0][1]=-1.; |
---|
| 167 | |
---|
| 168 | x1[1][1]=-1.; |
---|
| 169 | b1[1][1]=-1.; |
---|
| 170 | |
---|
| 171 | numberOfPartialCrossSections[1] = 2; |
---|
| 172 | |
---|
| 173 | // ALPHA+ |
---|
| 174 | f0[0][2]=1.; |
---|
| 175 | a0[0][2]=0.65; |
---|
| 176 | a1[0][2]=-2.75; |
---|
| 177 | b0[0][2]=-21.81; |
---|
| 178 | c0[0][2]=0.232; |
---|
| 179 | d0[0][2]=2.95; |
---|
| 180 | x0[0][2]=3.53; |
---|
| 181 | |
---|
| 182 | x1[0][2]=-1.; |
---|
| 183 | b1[0][2]=-1.; |
---|
| 184 | |
---|
| 185 | numberOfPartialCrossSections[2] = 1; |
---|
| 186 | |
---|
| 187 | // |
---|
| 188 | |
---|
[1192] | 189 | if( verboseLevel>0 ) |
---|
| 190 | { |
---|
| 191 | G4cout << "Dingfelder charge decrease model is initialized " << G4endl |
---|
| 192 | << "Energy range: " |
---|
| 193 | << LowEnergyLimit() / keV << " keV - " |
---|
| 194 | << HighEnergyLimit() / MeV << " MeV for " |
---|
| 195 | << particle->GetParticleName() |
---|
| 196 | << G4endl; |
---|
| 197 | } |
---|
| 198 | |
---|
[1058] | 199 | if(!isInitialised) |
---|
| 200 | { |
---|
| 201 | isInitialised = true; |
---|
| 202 | |
---|
| 203 | if(pParticleChange) |
---|
| 204 | fParticleChangeForGamma = reinterpret_cast<G4ParticleChangeForGamma*>(pParticleChange); |
---|
| 205 | else |
---|
| 206 | fParticleChangeForGamma = new G4ParticleChangeForGamma(); |
---|
| 207 | } |
---|
| 208 | |
---|
| 209 | // InitialiseElementSelectors(particle,cuts); |
---|
| 210 | |
---|
| 211 | // Test if water material |
---|
| 212 | |
---|
| 213 | flagMaterialIsWater= false; |
---|
| 214 | densityWater = 0; |
---|
| 215 | |
---|
| 216 | const G4ProductionCutsTable* theCoupleTable = G4ProductionCutsTable::GetProductionCutsTable(); |
---|
| 217 | |
---|
| 218 | if(theCoupleTable) |
---|
| 219 | { |
---|
| 220 | G4int numOfCouples = theCoupleTable->GetTableSize(); |
---|
| 221 | |
---|
| 222 | if(numOfCouples>0) |
---|
| 223 | { |
---|
| 224 | for (G4int i=0; i<numOfCouples; i++) |
---|
| 225 | { |
---|
| 226 | const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(i); |
---|
| 227 | const G4Material* material = couple->GetMaterial(); |
---|
| 228 | |
---|
[1192] | 229 | if (material->GetName() == "G4_WATER") |
---|
[1058] | 230 | { |
---|
[1192] | 231 | G4double density = material->GetAtomicNumDensityVector()[1]; |
---|
| 232 | flagMaterialIsWater = true; |
---|
| 233 | densityWater = density; |
---|
| 234 | |
---|
| 235 | if (verboseLevel > 3) |
---|
| 236 | G4cout << "****** Water material is found with density(cm^-3)=" << density/(cm*cm*cm) << G4endl; |
---|
[1058] | 237 | } |
---|
| 238 | |
---|
| 239 | } |
---|
| 240 | |
---|
[1192] | 241 | } // if(numOfCouples>0) |
---|
| 242 | |
---|
[1058] | 243 | } // if (theCoupleTable) |
---|
| 244 | |
---|
| 245 | } |
---|
| 246 | |
---|
| 247 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 248 | |
---|
| 249 | G4double G4DNADingfelderChargeDecreaseModel::CrossSectionPerVolume(const G4Material*, |
---|
| 250 | const G4ParticleDefinition* particleDefinition, |
---|
| 251 | G4double k, |
---|
| 252 | G4double, |
---|
| 253 | G4double) |
---|
| 254 | { |
---|
| 255 | if (verboseLevel > 3) |
---|
| 256 | G4cout << "Calling CrossSectionPerVolume() of G4DNADingfelderChargeDecreaseModel" << G4endl; |
---|
| 257 | |
---|
| 258 | // Calculate total cross section for model |
---|
| 259 | |
---|
| 260 | G4DNAGenericIonsManager *instance; |
---|
| 261 | instance = G4DNAGenericIonsManager::Instance(); |
---|
| 262 | |
---|
| 263 | if ( |
---|
| 264 | particleDefinition != G4Proton::ProtonDefinition() |
---|
| 265 | && |
---|
| 266 | particleDefinition != instance->GetIon("alpha++") |
---|
| 267 | && |
---|
| 268 | particleDefinition != instance->GetIon("alpha+") |
---|
| 269 | ) |
---|
| 270 | |
---|
| 271 | return 0; |
---|
| 272 | |
---|
| 273 | G4double lowLim = 0; |
---|
| 274 | G4double highLim = 0; |
---|
| 275 | G4double crossSection = 0.; |
---|
| 276 | |
---|
| 277 | if (flagMaterialIsWater) |
---|
| 278 | { |
---|
| 279 | const G4String& particleName = particleDefinition->GetParticleName(); |
---|
| 280 | |
---|
| 281 | std::map< G4String,G4double,std::less<G4String> >::iterator pos1; |
---|
| 282 | pos1 = lowEnergyLimit.find(particleName); |
---|
| 283 | |
---|
| 284 | if (pos1 != lowEnergyLimit.end()) |
---|
| 285 | { |
---|
| 286 | lowLim = pos1->second; |
---|
| 287 | } |
---|
| 288 | |
---|
| 289 | std::map< G4String,G4double,std::less<G4String> >::iterator pos2; |
---|
| 290 | pos2 = highEnergyLimit.find(particleName); |
---|
| 291 | |
---|
| 292 | if (pos2 != highEnergyLimit.end()) |
---|
| 293 | { |
---|
| 294 | highLim = pos2->second; |
---|
| 295 | } |
---|
| 296 | |
---|
| 297 | if (k >= lowLim && k < highLim) |
---|
| 298 | { |
---|
| 299 | crossSection = Sum(k,particleDefinition); |
---|
| 300 | } |
---|
| 301 | |
---|
| 302 | if (verboseLevel > 3) |
---|
| 303 | { |
---|
| 304 | G4cout << "---> Kinetic energy(eV)=" << k/eV << G4endl; |
---|
| 305 | G4cout << " - Cross section per water molecule (cm^2)=" << crossSection/cm/cm << G4endl; |
---|
| 306 | G4cout << " - Cross section per water molecule (cm^-1)=" << crossSection*densityWater/(1./cm) << G4endl; |
---|
| 307 | } |
---|
| 308 | |
---|
| 309 | } // if (flagMaterialIsWater) |
---|
| 310 | |
---|
| 311 | return crossSection*densityWater; |
---|
| 312 | |
---|
| 313 | } |
---|
| 314 | |
---|
| 315 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 316 | |
---|
| 317 | void G4DNADingfelderChargeDecreaseModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect, |
---|
| 318 | const G4MaterialCutsCouple* /*couple*/, |
---|
| 319 | const G4DynamicParticle* aDynamicParticle, |
---|
| 320 | G4double, |
---|
| 321 | G4double) |
---|
| 322 | { |
---|
| 323 | if (verboseLevel > 3) |
---|
| 324 | G4cout << "Calling SampleSecondaries() of G4DNADingfelderChargeDecreaseModel" << G4endl; |
---|
| 325 | |
---|
| 326 | G4double inK = aDynamicParticle->GetKineticEnergy(); |
---|
| 327 | |
---|
| 328 | G4ParticleDefinition* definition = aDynamicParticle->GetDefinition(); |
---|
| 329 | |
---|
| 330 | G4int finalStateIndex = RandomSelect(inK,definition); |
---|
| 331 | |
---|
| 332 | G4int n = NumberOfFinalStates(definition, finalStateIndex); |
---|
| 333 | G4double waterBindingEnergy = WaterBindingEnergyConstant(definition, finalStateIndex); |
---|
| 334 | G4double outgoingParticleBindingEnergy = OutgoingParticleBindingEnergyConstant(definition, finalStateIndex); |
---|
| 335 | |
---|
| 336 | G4double outK = 0.; |
---|
| 337 | if (definition==G4Proton::Proton()) |
---|
| 338 | outK = inK - n*(inK*electron_mass_c2/proton_mass_c2) - waterBindingEnergy + outgoingParticleBindingEnergy; |
---|
| 339 | else |
---|
| 340 | outK = inK - n*(inK*electron_mass_c2/(3728*MeV)) - waterBindingEnergy + outgoingParticleBindingEnergy; |
---|
| 341 | |
---|
| 342 | if (outK<0) |
---|
| 343 | { |
---|
| 344 | G4String message; |
---|
| 345 | message="G4DNADingfelderChargeDecreaseModel::SampleSecondaries: final kinetic energy is below 0! Process "; |
---|
| 346 | G4Exception(message); |
---|
| 347 | } |
---|
| 348 | |
---|
| 349 | fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill); |
---|
| 350 | fParticleChangeForGamma->ProposeLocalEnergyDeposit(waterBindingEnergy); |
---|
| 351 | |
---|
| 352 | G4DynamicParticle* dp = new G4DynamicParticle (OutgoingParticleDefinition(definition, finalStateIndex), |
---|
| 353 | aDynamicParticle->GetMomentumDirection(), |
---|
| 354 | outK) ; |
---|
| 355 | fvect->push_back(dp); |
---|
| 356 | |
---|
| 357 | } |
---|
| 358 | |
---|
| 359 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 360 | |
---|
| 361 | G4int G4DNADingfelderChargeDecreaseModel::NumberOfFinalStates(G4ParticleDefinition* particleDefinition, |
---|
| 362 | G4int finalStateIndex ) |
---|
| 363 | |
---|
| 364 | { |
---|
| 365 | if (particleDefinition == G4Proton::Proton()) return 1; |
---|
| 366 | |
---|
| 367 | G4DNAGenericIonsManager*instance; |
---|
| 368 | instance = G4DNAGenericIonsManager::Instance(); |
---|
| 369 | |
---|
| 370 | if (particleDefinition == instance->GetIon("alpha++") ) |
---|
| 371 | { |
---|
| 372 | if (finalStateIndex==0) return 1; |
---|
| 373 | return 2; |
---|
| 374 | } |
---|
| 375 | |
---|
| 376 | if (particleDefinition == instance->GetIon("alpha+") ) return 1; |
---|
| 377 | |
---|
| 378 | return 0; |
---|
| 379 | } |
---|
| 380 | |
---|
| 381 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 382 | |
---|
| 383 | G4ParticleDefinition* G4DNADingfelderChargeDecreaseModel::OutgoingParticleDefinition (G4ParticleDefinition* particleDefinition, |
---|
| 384 | G4int finalStateIndex) |
---|
| 385 | { |
---|
| 386 | G4DNAGenericIonsManager * instance(G4DNAGenericIonsManager::Instance()); |
---|
| 387 | |
---|
| 388 | if (particleDefinition == G4Proton::Proton()) return instance->GetIon("hydrogen"); |
---|
| 389 | |
---|
| 390 | if (particleDefinition == instance->GetIon("alpha++") ) |
---|
| 391 | { |
---|
| 392 | if (finalStateIndex == 0) return instance->GetIon("alpha+"); |
---|
| 393 | return instance->GetIon("helium"); |
---|
| 394 | } |
---|
| 395 | |
---|
| 396 | if (particleDefinition == instance->GetIon("alpha+") ) return instance->GetIon("helium"); |
---|
| 397 | |
---|
| 398 | return 0; |
---|
| 399 | } |
---|
| 400 | |
---|
| 401 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 402 | |
---|
| 403 | G4double G4DNADingfelderChargeDecreaseModel::WaterBindingEnergyConstant(G4ParticleDefinition* particleDefinition, |
---|
| 404 | G4int finalStateIndex ) |
---|
| 405 | { |
---|
| 406 | // Ionization energy of first water shell |
---|
| 407 | // Rad. Phys. Chem. 59 p.267 by Dingf. et al. |
---|
| 408 | // W + 10.79 eV -> W+ + e- |
---|
| 409 | |
---|
| 410 | G4DNAGenericIonsManager * instance(G4DNAGenericIonsManager::Instance()); |
---|
| 411 | |
---|
| 412 | if(particleDefinition == G4Proton::Proton()) return 10.79*eV; |
---|
| 413 | |
---|
| 414 | if (particleDefinition == instance->GetIon("alpha++") ) |
---|
| 415 | { |
---|
| 416 | // Binding energy for W+ -> W++ + e- 10.79 eV |
---|
| 417 | // Binding energy for W -> W+ + e- 10.79 eV |
---|
| 418 | // |
---|
| 419 | // Ionization energy of first water shell |
---|
| 420 | // Rad. Phys. Chem. 59 p.267 by Dingf. et al. |
---|
| 421 | |
---|
| 422 | if (finalStateIndex == 0) return 10.79*eV; |
---|
| 423 | |
---|
| 424 | return 10.79*2*eV; |
---|
| 425 | } |
---|
| 426 | |
---|
| 427 | if (particleDefinition == instance->GetIon("alpha+") ) |
---|
| 428 | { |
---|
| 429 | // Binding energy for W+ -> W++ + e- 10.79 eV |
---|
| 430 | // Binding energy for W -> W+ + e- 10.79 eV |
---|
| 431 | // |
---|
| 432 | // Ionization energy of first water shell |
---|
| 433 | // Rad. Phys. Chem. 59 p.267 by Dingf. et al. |
---|
| 434 | |
---|
| 435 | return 10.79*eV; |
---|
| 436 | } |
---|
| 437 | |
---|
| 438 | return 0; |
---|
| 439 | } |
---|
| 440 | |
---|
| 441 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 442 | |
---|
| 443 | G4double G4DNADingfelderChargeDecreaseModel::OutgoingParticleBindingEnergyConstant(G4ParticleDefinition* particleDefinition, |
---|
| 444 | G4int finalStateIndex) |
---|
| 445 | { |
---|
| 446 | G4DNAGenericIonsManager * instance(G4DNAGenericIonsManager::Instance()); |
---|
| 447 | |
---|
| 448 | if(particleDefinition == G4Proton::Proton()) return 13.6*eV; |
---|
| 449 | |
---|
| 450 | if (particleDefinition == instance->GetIon("alpha++") ) |
---|
| 451 | { |
---|
| 452 | // Binding energy for He+ -> He++ + e- 54.509 eV |
---|
| 453 | // Binding energy for He -> He+ + e- 24.587 eV |
---|
| 454 | |
---|
| 455 | if (finalStateIndex==0) return 54.509*eV; |
---|
| 456 | |
---|
| 457 | return (54.509 + 24.587)*eV; |
---|
| 458 | } |
---|
| 459 | |
---|
| 460 | if (particleDefinition == instance->GetIon("alpha+") ) |
---|
| 461 | { |
---|
| 462 | // Binding energy for He+ -> He++ + e- 54.509 eV |
---|
| 463 | // Binding energy for He -> He+ + e- 24.587 eV |
---|
| 464 | |
---|
| 465 | return 24.587*eV; |
---|
| 466 | } |
---|
| 467 | |
---|
| 468 | return 0; |
---|
| 469 | } |
---|
| 470 | |
---|
| 471 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 472 | |
---|
| 473 | G4double G4DNADingfelderChargeDecreaseModel::PartialCrossSection(G4double k, G4int index, |
---|
| 474 | const G4ParticleDefinition* particleDefinition) |
---|
| 475 | { |
---|
| 476 | G4int particleTypeIndex = 0; |
---|
| 477 | G4DNAGenericIonsManager* instance; |
---|
| 478 | instance = G4DNAGenericIonsManager::Instance(); |
---|
| 479 | |
---|
| 480 | if (particleDefinition == G4Proton::ProtonDefinition()) particleTypeIndex=0; |
---|
| 481 | if (particleDefinition == instance->GetIon("alpha++")) particleTypeIndex=1; |
---|
| 482 | if (particleDefinition == instance->GetIon("alpha+")) particleTypeIndex=2; |
---|
| 483 | |
---|
| 484 | // |
---|
| 485 | // sigma(T) = f0 10 ^ y(log10(T/eV)) |
---|
| 486 | // |
---|
| 487 | // / a0 x + b0 x < x0 |
---|
| 488 | // | |
---|
| 489 | // y(x) = < a0 x + b0 - c0 (x - x0)^d0 x0 <= x < x1 |
---|
| 490 | // | |
---|
| 491 | // \ a1 x + b1 x >= x1 |
---|
| 492 | // |
---|
| 493 | // |
---|
| 494 | // f0, a0, a1, b0, b1, c0, d0, x0, x1 are parameters that change for protons and helium (0, +, ++) |
---|
| 495 | // |
---|
| 496 | // f0 has been added to the code in order to manage partial (shell-dependent) cross sections (if no shell dependence is present. f0=1. Sum of f0 over the considered shells should give 1) |
---|
| 497 | // |
---|
| 498 | // From Rad. Phys. and Chem. 59 (2000) 255-275, M. Dingfelder et al. |
---|
| 499 | // Inelastic-collision cross sections of liquid water for interactions of energetic proton |
---|
| 500 | // |
---|
| 501 | |
---|
| 502 | if (x1[index][particleTypeIndex]<x0[index][particleTypeIndex]) |
---|
| 503 | { |
---|
| 504 | // |
---|
| 505 | // if x1 < x0 means that x1 and b1 will be calculated with the following formula (this piece of code is run on all alphas and not on protons) |
---|
| 506 | // |
---|
| 507 | // x1 = x0 + ((a0 - a1)/(c0 * d0)) ^ (1 / (d0 - 1)) |
---|
| 508 | // |
---|
| 509 | // b1 = (a0 - a1) * x1 + b0 - c0 * (x1 - x0) ^ d0 |
---|
| 510 | // |
---|
| 511 | |
---|
| 512 | x1[index][particleTypeIndex]=x0[index][particleTypeIndex] + std::pow((a0[index][particleTypeIndex] - a1[index][particleTypeIndex]) |
---|
| 513 | / (c0[index][particleTypeIndex] * d0[index][particleTypeIndex]), 1. / (d0[index][particleTypeIndex] - 1.)); |
---|
| 514 | b1[index][particleTypeIndex]=(a0[index][particleTypeIndex] - a1[index][particleTypeIndex]) * x1[index][particleTypeIndex] |
---|
| 515 | + b0[index][particleTypeIndex] - c0[index][particleTypeIndex] * std::pow(x1[index][particleTypeIndex] |
---|
| 516 | - x0[index][particleTypeIndex], d0[index][particleTypeIndex]); |
---|
| 517 | } |
---|
| 518 | |
---|
| 519 | G4double x(std::log10(k/eV)); |
---|
| 520 | G4double y; |
---|
| 521 | |
---|
| 522 | if (x<x0[index][particleTypeIndex]) |
---|
| 523 | y=a0[index][particleTypeIndex] * x + b0[index][particleTypeIndex]; |
---|
| 524 | else if (x<x1[index][particleTypeIndex]) |
---|
| 525 | y=a0[index][particleTypeIndex] * x + b0[index][particleTypeIndex] - c0[index][particleTypeIndex] |
---|
| 526 | * std::pow(x - x0[index][particleTypeIndex], d0[index][particleTypeIndex]); |
---|
| 527 | else |
---|
| 528 | y=a1[index][particleTypeIndex] * x + b1[index][particleTypeIndex]; |
---|
| 529 | |
---|
| 530 | return f0[index][particleTypeIndex] * std::pow(10., y)*m*m; |
---|
| 531 | |
---|
| 532 | } |
---|
| 533 | |
---|
| 534 | G4int G4DNADingfelderChargeDecreaseModel::RandomSelect(G4double k, |
---|
| 535 | const G4ParticleDefinition* particleDefinition) |
---|
| 536 | { |
---|
| 537 | G4int particleTypeIndex = 0; |
---|
| 538 | G4DNAGenericIonsManager *instance; |
---|
| 539 | instance = G4DNAGenericIonsManager::Instance(); |
---|
| 540 | |
---|
| 541 | if (particleDefinition == G4Proton::ProtonDefinition()) particleTypeIndex = 0; |
---|
| 542 | if (particleDefinition == instance->GetIon("alpha++")) particleTypeIndex = 1; |
---|
| 543 | if (particleDefinition == instance->GetIon("alpha+")) particleTypeIndex = 2; |
---|
| 544 | |
---|
| 545 | const G4int n = numberOfPartialCrossSections[particleTypeIndex]; |
---|
| 546 | G4double* values(new G4double[n]); |
---|
| 547 | G4double value(0); |
---|
| 548 | G4int i = n; |
---|
| 549 | |
---|
| 550 | while (i>0) |
---|
| 551 | { |
---|
| 552 | i--; |
---|
| 553 | values[i]=PartialCrossSection(k, i, particleDefinition); |
---|
| 554 | value+=values[i]; |
---|
| 555 | } |
---|
| 556 | |
---|
| 557 | value*=G4UniformRand(); |
---|
| 558 | |
---|
| 559 | i=n; |
---|
| 560 | while (i>0) |
---|
| 561 | { |
---|
| 562 | i--; |
---|
| 563 | |
---|
| 564 | if (values[i]>value) |
---|
| 565 | break; |
---|
| 566 | |
---|
| 567 | value-=values[i]; |
---|
| 568 | } |
---|
| 569 | |
---|
| 570 | delete[] values; |
---|
| 571 | |
---|
| 572 | return i; |
---|
| 573 | } |
---|
| 574 | |
---|
| 575 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 576 | |
---|
| 577 | G4double G4DNADingfelderChargeDecreaseModel::Sum(G4double k, const G4ParticleDefinition* particleDefinition) |
---|
| 578 | { |
---|
| 579 | G4int particleTypeIndex = 0; |
---|
| 580 | G4DNAGenericIonsManager* instance; |
---|
| 581 | instance = G4DNAGenericIonsManager::Instance(); |
---|
| 582 | |
---|
| 583 | if (particleDefinition == G4Proton::ProtonDefinition()) particleTypeIndex=0; |
---|
| 584 | if (particleDefinition == instance->GetIon("alpha++")) particleTypeIndex=1; |
---|
| 585 | if (particleDefinition == instance->GetIon("alpha+")) particleTypeIndex=2; |
---|
| 586 | |
---|
| 587 | G4double totalCrossSection = 0.; |
---|
| 588 | |
---|
| 589 | for (G4int i=0; i<numberOfPartialCrossSections[particleTypeIndex]; i++) |
---|
| 590 | { |
---|
| 591 | totalCrossSection += PartialCrossSection(k,i,particleDefinition); |
---|
| 592 | } |
---|
| 593 | return totalCrossSection; |
---|
| 594 | } |
---|
| 595 | |
---|