| [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: G4DNAMillerGreenExcitationModel.cc,v 1.6 2009/08/13 11:32:47 sincerti Exp $
|
|---|
| [1196] | 27 | // GEANT4 tag $Name: geant4-09-03-cand-01 $
|
|---|
| [1058] | 28 | //
|
|---|
| 29 |
|
|---|
| 30 | #include "G4DNAMillerGreenExcitationModel.hh"
|
|---|
| 31 |
|
|---|
| 32 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 33 |
|
|---|
| 34 | using namespace std;
|
|---|
| 35 |
|
|---|
| 36 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 37 |
|
|---|
| 38 | G4DNAMillerGreenExcitationModel::G4DNAMillerGreenExcitationModel(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 << "Miller & Green excitation model is constructed " << G4endl;
|
|---|
| 54 | }
|
|---|
| [1058] | 55 | }
|
|---|
| 56 |
|
|---|
| 57 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 58 |
|
|---|
| 59 | G4DNAMillerGreenExcitationModel::~G4DNAMillerGreenExcitationModel()
|
|---|
| 60 | {}
|
|---|
| 61 |
|
|---|
| 62 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 63 |
|
|---|
| 64 | void G4DNAMillerGreenExcitationModel::Initialise(const G4ParticleDefinition* particle,
|
|---|
| 65 | const G4DataVector& /*cuts*/)
|
|---|
| 66 | {
|
|---|
| 67 |
|
|---|
| 68 | if (verboseLevel > 3)
|
|---|
| 69 | G4cout << "Calling G4DNAMillerGreenExcitationModel::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 | G4ParticleDefinition* heliumDef = instance->GetIon("helium");
|
|---|
| 79 |
|
|---|
| 80 | G4String proton;
|
|---|
| 81 | G4String alphaPlusPlus;
|
|---|
| 82 | G4String alphaPlus;
|
|---|
| 83 | G4String helium;
|
|---|
| 84 |
|
|---|
| 85 | if (protonDef != 0)
|
|---|
| 86 | {
|
|---|
| 87 | proton = protonDef->GetParticleName();
|
|---|
| 88 | lowEnergyLimit[proton] = 10. * eV;
|
|---|
| 89 | highEnergyLimit[proton] = 500. * keV;
|
|---|
| 90 |
|
|---|
| 91 | kineticEnergyCorrection[0] = 1.;
|
|---|
| 92 | slaterEffectiveCharge[0][0] = 0.;
|
|---|
| 93 | slaterEffectiveCharge[1][0] = 0.;
|
|---|
| 94 | slaterEffectiveCharge[2][0] = 0.;
|
|---|
| 95 | sCoefficient[0][0] = 0.;
|
|---|
| 96 | sCoefficient[1][0] = 0.;
|
|---|
| 97 | sCoefficient[2][0] = 0.;
|
|---|
| 98 | }
|
|---|
| 99 | else
|
|---|
| 100 | {
|
|---|
| 101 | G4Exception("G4DNAMillerGreenExcitationModel::Initialise: proton is not defined");
|
|---|
| 102 | }
|
|---|
| 103 |
|
|---|
| 104 | if (alphaPlusPlusDef != 0)
|
|---|
| 105 | {
|
|---|
| 106 | alphaPlusPlus = alphaPlusPlusDef->GetParticleName();
|
|---|
| 107 | lowEnergyLimit[alphaPlusPlus] = 1. * keV;
|
|---|
| 108 | highEnergyLimit[alphaPlusPlus] = 10. * MeV;
|
|---|
| 109 |
|
|---|
| 110 | kineticEnergyCorrection[1] = 0.9382723/3.727417;
|
|---|
| 111 | slaterEffectiveCharge[0][1]=0.;
|
|---|
| 112 | slaterEffectiveCharge[1][1]=0.;
|
|---|
| 113 | slaterEffectiveCharge[2][1]=0.;
|
|---|
| 114 | sCoefficient[0][1]=0.;
|
|---|
| 115 | sCoefficient[1][1]=0.;
|
|---|
| 116 | sCoefficient[2][1]=0.;
|
|---|
| 117 | }
|
|---|
| 118 | else
|
|---|
| 119 | {
|
|---|
| 120 | G4Exception("G4DNAMillerGreenExcitationModel::Initialise: alphaPlusPlus is not defined");
|
|---|
| 121 | }
|
|---|
| 122 |
|
|---|
| 123 | if (alphaPlusDef != 0)
|
|---|
| 124 | {
|
|---|
| 125 | alphaPlus = alphaPlusDef->GetParticleName();
|
|---|
| 126 | lowEnergyLimit[alphaPlus] = 1. * keV;
|
|---|
| 127 | highEnergyLimit[alphaPlus] = 10. * MeV;
|
|---|
| 128 |
|
|---|
| 129 | kineticEnergyCorrection[2] = 0.9382723/3.727417;
|
|---|
| 130 | slaterEffectiveCharge[0][2]=2.0;
|
|---|
| 131 | slaterEffectiveCharge[1][2]=1.15;
|
|---|
| 132 | slaterEffectiveCharge[2][2]=1.15;
|
|---|
| 133 | sCoefficient[0][2]=0.7;
|
|---|
| 134 | sCoefficient[1][2]=0.15;
|
|---|
| 135 | sCoefficient[2][2]=0.15;
|
|---|
| 136 | }
|
|---|
| 137 | else
|
|---|
| 138 | {
|
|---|
| 139 | G4Exception("G4DNAMillerGreenExcitationModel::Initialise: alphaPlus is not defined");
|
|---|
| 140 | }
|
|---|
| 141 |
|
|---|
| 142 | if (heliumDef != 0)
|
|---|
| 143 | {
|
|---|
| 144 | helium = heliumDef->GetParticleName();
|
|---|
| 145 | lowEnergyLimit[helium] = 1. * keV;
|
|---|
| 146 | highEnergyLimit[helium] = 10. * MeV;
|
|---|
| 147 |
|
|---|
| 148 | kineticEnergyCorrection[3] = 0.9382723/3.727417;
|
|---|
| 149 | slaterEffectiveCharge[0][3]=1.7;
|
|---|
| 150 | slaterEffectiveCharge[1][3]=1.15;
|
|---|
| 151 | slaterEffectiveCharge[2][3]=1.15;
|
|---|
| 152 | sCoefficient[0][3]=0.5;
|
|---|
| 153 | sCoefficient[1][3]=0.25;
|
|---|
| 154 | sCoefficient[2][3]=0.25;
|
|---|
| 155 | }
|
|---|
| 156 | else
|
|---|
| 157 | {
|
|---|
| 158 | G4Exception("G4DNAMillerGreenExcitationModel::Initialise: helium is not defined");
|
|---|
| 159 | }
|
|---|
| 160 |
|
|---|
| 161 | if (particle==protonDef)
|
|---|
| 162 | {
|
|---|
| 163 | SetLowEnergyLimit(lowEnergyLimit[proton]);
|
|---|
| 164 | SetHighEnergyLimit(highEnergyLimit[proton]);
|
|---|
| 165 | }
|
|---|
| 166 |
|
|---|
| 167 | if (particle==alphaPlusPlusDef)
|
|---|
| 168 | {
|
|---|
| 169 | SetLowEnergyLimit(lowEnergyLimit[alphaPlusPlus]);
|
|---|
| 170 | SetHighEnergyLimit(highEnergyLimit[alphaPlusPlus]);
|
|---|
| 171 | }
|
|---|
| 172 |
|
|---|
| 173 | if (particle==alphaPlusDef)
|
|---|
| 174 | {
|
|---|
| 175 | SetLowEnergyLimit(lowEnergyLimit[alphaPlus]);
|
|---|
| 176 | SetHighEnergyLimit(highEnergyLimit[alphaPlus]);
|
|---|
| 177 | }
|
|---|
| 178 |
|
|---|
| 179 | if (particle==heliumDef)
|
|---|
| 180 | {
|
|---|
| 181 | SetLowEnergyLimit(lowEnergyLimit[helium]);
|
|---|
| 182 | SetHighEnergyLimit(highEnergyLimit[helium]);
|
|---|
| 183 | }
|
|---|
| 184 |
|
|---|
| 185 | //
|
|---|
| 186 |
|
|---|
| 187 | nLevels = waterExcitation.NumberOfLevels();
|
|---|
| 188 |
|
|---|
| 189 | //
|
|---|
| [1192] | 190 | if( verboseLevel>0 )
|
|---|
| 191 | {
|
|---|
| 192 | G4cout << "Miller & Green excitation model is initialized " << G4endl
|
|---|
| 193 | << "Energy range: "
|
|---|
| 194 | << LowEnergyLimit() / eV << " eV - "
|
|---|
| 195 | << HighEnergyLimit() / keV << " keV for "
|
|---|
| 196 | << particle->GetParticleName()
|
|---|
| 197 | << G4endl;
|
|---|
| 198 | }
|
|---|
| [1058] | 199 |
|
|---|
| 200 | if(!isInitialised)
|
|---|
| 201 | {
|
|---|
| 202 | isInitialised = true;
|
|---|
| 203 |
|
|---|
| 204 | if(pParticleChange)
|
|---|
| 205 | fParticleChangeForGamma = reinterpret_cast<G4ParticleChangeForGamma*>(pParticleChange);
|
|---|
| 206 | else
|
|---|
| 207 | fParticleChangeForGamma = new G4ParticleChangeForGamma();
|
|---|
| 208 | }
|
|---|
| 209 |
|
|---|
| 210 | // InitialiseElementSelectors(particle,cuts);
|
|---|
| 211 |
|
|---|
| 212 | // Test if water material
|
|---|
| 213 |
|
|---|
| 214 | flagMaterialIsWater= false;
|
|---|
| 215 | densityWater = 0;
|
|---|
| 216 |
|
|---|
| 217 | const G4ProductionCutsTable* theCoupleTable = G4ProductionCutsTable::GetProductionCutsTable();
|
|---|
| 218 |
|
|---|
| 219 | if(theCoupleTable)
|
|---|
| 220 | {
|
|---|
| 221 | G4int numOfCouples = theCoupleTable->GetTableSize();
|
|---|
| 222 |
|
|---|
| 223 | if(numOfCouples>0)
|
|---|
| 224 | {
|
|---|
| 225 | for (G4int i=0; i<numOfCouples; i++)
|
|---|
| 226 | {
|
|---|
| 227 | const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(i);
|
|---|
| 228 | const G4Material* material = couple->GetMaterial();
|
|---|
| 229 |
|
|---|
| [1192] | 230 | if (material->GetName() == "G4_WATER")
|
|---|
| [1058] | 231 | {
|
|---|
| [1192] | 232 | G4double density = material->GetAtomicNumDensityVector()[1];
|
|---|
| 233 | flagMaterialIsWater = true;
|
|---|
| 234 | densityWater = density;
|
|---|
| 235 |
|
|---|
| 236 | if (verboseLevel > 3)
|
|---|
| 237 | G4cout << "****** Water material is found with density(cm^-3)=" << density/(cm*cm*cm) << G4endl;
|
|---|
| [1058] | 238 | }
|
|---|
| 239 |
|
|---|
| 240 | }
|
|---|
| 241 |
|
|---|
| [1192] | 242 | } // if(numOfCouples>0)
|
|---|
| 243 |
|
|---|
| [1058] | 244 | } // if (theCoupleTable)
|
|---|
| 245 |
|
|---|
| 246 | }
|
|---|
| 247 |
|
|---|
| 248 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 249 |
|
|---|
| 250 | G4double G4DNAMillerGreenExcitationModel::CrossSectionPerVolume(const G4Material* material,
|
|---|
| 251 | const G4ParticleDefinition* particleDefinition,
|
|---|
| 252 | G4double k,
|
|---|
| 253 | G4double,
|
|---|
| 254 | G4double)
|
|---|
| 255 | {
|
|---|
| 256 | if (verboseLevel > 3)
|
|---|
| 257 | G4cout << "Calling CrossSectionPerVolume() of G4DNAMillerGreenExcitationModel" << G4endl;
|
|---|
| 258 |
|
|---|
| 259 | // Calculate total cross section for model
|
|---|
| 260 |
|
|---|
| 261 | G4DNAGenericIonsManager *instance;
|
|---|
| 262 | instance = G4DNAGenericIonsManager::Instance();
|
|---|
| 263 |
|
|---|
| 264 | if (
|
|---|
| 265 | particleDefinition != G4Proton::ProtonDefinition()
|
|---|
| 266 | &&
|
|---|
| 267 | particleDefinition != instance->GetIon("alpha++")
|
|---|
| 268 | &&
|
|---|
| 269 | particleDefinition != instance->GetIon("alpha+")
|
|---|
| 270 | &&
|
|---|
| 271 | particleDefinition != instance->GetIon("helium")
|
|---|
| 272 | )
|
|---|
| 273 |
|
|---|
| 274 | return 0;
|
|---|
| 275 |
|
|---|
| 276 | G4double lowLim = 0;
|
|---|
| 277 | G4double highLim = 0;
|
|---|
| 278 | G4double crossSection = 0.;
|
|---|
| 279 |
|
|---|
| 280 | if (flagMaterialIsWater)
|
|---|
| 281 | {
|
|---|
| 282 | const G4String& particleName = particleDefinition->GetParticleName();
|
|---|
| 283 |
|
|---|
| 284 | std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
|
|---|
| 285 | pos1 = lowEnergyLimit.find(particleName);
|
|---|
| 286 |
|
|---|
| 287 | if (pos1 != lowEnergyLimit.end())
|
|---|
| 288 | {
|
|---|
| 289 | lowLim = pos1->second;
|
|---|
| 290 | }
|
|---|
| 291 |
|
|---|
| 292 | std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
|
|---|
| 293 | pos2 = highEnergyLimit.find(particleName);
|
|---|
| 294 |
|
|---|
| 295 | if (pos2 != highEnergyLimit.end())
|
|---|
| 296 | {
|
|---|
| 297 | highLim = pos2->second;
|
|---|
| 298 | }
|
|---|
| 299 |
|
|---|
| 300 | if (k >= lowLim && k < highLim)
|
|---|
| 301 | {
|
|---|
| 302 | crossSection = Sum(k,particleDefinition);
|
|---|
| 303 |
|
|---|
| 304 | G4DNAGenericIonsManager *instance;
|
|---|
| 305 | instance = G4DNAGenericIonsManager::Instance();
|
|---|
| 306 |
|
|---|
| 307 | // add ONE or TWO electron-water excitation for alpha+ and helium
|
|---|
| 308 |
|
|---|
| 309 | if ( particleDefinition == instance->GetIon("alpha+")
|
|---|
| 310 | ||
|
|---|
| 311 | particleDefinition == instance->GetIon("helium")
|
|---|
| 312 | )
|
|---|
| 313 | {
|
|---|
| 314 | G4DNAEmfietzoglouExcitationModel * excitationXS = new G4DNAEmfietzoglouExcitationModel();
|
|---|
| 315 |
|
|---|
| 316 | G4double sigmaExcitation=0;
|
|---|
| 317 | G4double tmp =0.;
|
|---|
| 318 |
|
|---|
| 319 | if (k*0.511/3728 > 7.4*eV && k*0.511/3728 < 10*keV) sigmaExcitation =
|
|---|
| 320 | excitationXS->CrossSectionPerVolume(material,particleDefinition,k*0.511/3728,tmp,tmp)/densityWater;
|
|---|
| 321 |
|
|---|
| 322 | if ( particleDefinition == instance->GetIon("alpha+") )
|
|---|
| 323 | crossSection = crossSection + sigmaExcitation ;
|
|---|
| 324 |
|
|---|
| 325 | if ( particleDefinition == instance->GetIon("helium") )
|
|---|
| 326 | crossSection = crossSection + 2*sigmaExcitation ;
|
|---|
| 327 |
|
|---|
| 328 | delete excitationXS;
|
|---|
| 329 | }
|
|---|
| 330 |
|
|---|
| 331 | }
|
|---|
| 332 |
|
|---|
| 333 | if (verboseLevel > 3)
|
|---|
| 334 | {
|
|---|
| 335 | G4cout << "---> Kinetic energy(eV)=" << k/eV << G4endl;
|
|---|
| 336 | G4cout << " - Cross section per water molecule (cm^2)=" << crossSection/cm/cm << G4endl;
|
|---|
| 337 | G4cout << " - Cross section per water molecule (cm^-1)=" << crossSection*densityWater/(1./cm) << G4endl;
|
|---|
| 338 | }
|
|---|
| 339 |
|
|---|
| 340 | } // if (flagMaterialIsWater)
|
|---|
| 341 |
|
|---|
| 342 | return crossSection*densityWater;
|
|---|
| 343 |
|
|---|
| 344 | }
|
|---|
| 345 |
|
|---|
| 346 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 347 |
|
|---|
| 348 | void G4DNAMillerGreenExcitationModel::SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/,
|
|---|
| 349 | const G4MaterialCutsCouple* /*couple*/,
|
|---|
| 350 | const G4DynamicParticle* aDynamicParticle,
|
|---|
| 351 | G4double,
|
|---|
| 352 | G4double)
|
|---|
| 353 | {
|
|---|
| 354 |
|
|---|
| 355 | if (verboseLevel > 3)
|
|---|
| 356 | G4cout << "Calling SampleSecondaries() of G4DNAMillerGreenExcitationModel" << G4endl;
|
|---|
| 357 |
|
|---|
| 358 | G4double particleEnergy0 = aDynamicParticle->GetKineticEnergy();
|
|---|
| 359 |
|
|---|
| 360 | G4int level = RandomSelect(particleEnergy0,aDynamicParticle->GetDefinition());
|
|---|
| 361 |
|
|---|
| 362 | G4double excitationEnergy = waterExcitation.ExcitationEnergy(level);
|
|---|
| 363 | G4double newEnergy = particleEnergy0 - excitationEnergy;
|
|---|
| 364 |
|
|---|
| 365 | if (newEnergy>0)
|
|---|
| 366 | {
|
|---|
| 367 | fParticleChangeForGamma->ProposeMomentumDirection(aDynamicParticle->GetMomentumDirection());
|
|---|
| 368 | fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy);
|
|---|
| 369 | fParticleChangeForGamma->ProposeLocalEnergyDeposit(excitationEnergy);
|
|---|
| 370 | }
|
|---|
| 371 |
|
|---|
| 372 | }
|
|---|
| 373 |
|
|---|
| 374 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 375 |
|
|---|
| 376 | G4double G4DNAMillerGreenExcitationModel::PartialCrossSection(G4double k, G4int excitationLevel,
|
|---|
| 377 | const G4ParticleDefinition* particleDefinition)
|
|---|
| 378 | {
|
|---|
| 379 | // ( ( z * aj ) ^ omegaj ) * ( t - ej ) ^ nu
|
|---|
| 380 | // sigma(t) = zEff^2 * sigma0 * --------------------------------------------
|
|---|
| 381 | // jj ^ ( omegaj + nu ) + t ^ ( omegaj + nu )
|
|---|
| 382 | //
|
|---|
| 383 | // where t is the kinetic energy corrected by Helium mass over proton mass for Helium ions
|
|---|
| 384 | //
|
|---|
| 385 | // zEff is:
|
|---|
| 386 | // 1 for protons
|
|---|
| 387 | // 2 for alpha++
|
|---|
| 388 | // and 2 - c1 S_1s - c2 S_2s - c3 S_2p for alpha+ and He
|
|---|
| 389 | //
|
|---|
| 390 | // Dingfelder et al., RPC 59, 255-275, 2000 from Miller and Green (1973)
|
|---|
| 391 | // Formula (34) and Table 2
|
|---|
| 392 |
|
|---|
| 393 | const G4double sigma0(1.E+8 * barn);
|
|---|
| 394 | const G4double nu(1.);
|
|---|
| 395 | const G4double aj[]={876.*eV, 2084.* eV, 1373.*eV, 692.*eV, 900.*eV};
|
|---|
| 396 | const G4double jj[]={19820.*eV, 23490.*eV, 27770.*eV, 30830.*eV, 33080.*eV};
|
|---|
| 397 | const G4double omegaj[]={0.85, 0.88, 0.88, 0.78, 0.78};
|
|---|
| 398 |
|
|---|
| 399 | G4int particleTypeIndex = 0;
|
|---|
| 400 | G4DNAGenericIonsManager* instance;
|
|---|
| 401 | instance = G4DNAGenericIonsManager::Instance();
|
|---|
| 402 |
|
|---|
| 403 | if (particleDefinition == G4Proton::ProtonDefinition()) particleTypeIndex=0;
|
|---|
| 404 | if (particleDefinition == instance->GetIon("alpha++")) particleTypeIndex=1;
|
|---|
| 405 | if (particleDefinition == instance->GetIon("alpha+")) particleTypeIndex=2;
|
|---|
| 406 | if (particleDefinition == instance->GetIon("helium")) particleTypeIndex=3;
|
|---|
| 407 |
|
|---|
| 408 | G4double tCorrected;
|
|---|
| 409 | tCorrected = k * kineticEnergyCorrection[particleTypeIndex];
|
|---|
| 410 |
|
|---|
| 411 | // SI - added protection
|
|---|
| 412 | if (tCorrected < waterExcitation.ExcitationEnergy(excitationLevel)) return 0;
|
|---|
| 413 | //
|
|---|
| 414 |
|
|---|
| 415 | G4int z = 10;
|
|---|
| 416 |
|
|---|
| 417 | G4double numerator;
|
|---|
| 418 | numerator = std::pow(z * aj[excitationLevel], omegaj[excitationLevel]) *
|
|---|
| 419 | std::pow(tCorrected - waterExcitation.ExcitationEnergy(excitationLevel), nu);
|
|---|
| 420 |
|
|---|
| 421 | G4double power;
|
|---|
| 422 | power = omegaj[excitationLevel] + nu;
|
|---|
| 423 |
|
|---|
| 424 | G4double denominator;
|
|---|
| 425 | denominator = std::pow(jj[excitationLevel], power) + std::pow(tCorrected, power);
|
|---|
| 426 |
|
|---|
| 427 | G4double zEff = particleDefinition->GetPDGCharge() / eplus + particleDefinition->GetLeptonNumber();
|
|---|
| 428 |
|
|---|
| 429 | zEff -= ( sCoefficient[0][particleTypeIndex] * S_1s(k, waterExcitation.ExcitationEnergy(excitationLevel), slaterEffectiveCharge[0][particleTypeIndex], 1.) +
|
|---|
| 430 | sCoefficient[1][particleTypeIndex] * S_2s(k, waterExcitation.ExcitationEnergy(excitationLevel), slaterEffectiveCharge[1][particleTypeIndex], 2.) +
|
|---|
| 431 | sCoefficient[2][particleTypeIndex] * S_2p(k, waterExcitation.ExcitationEnergy(excitationLevel), slaterEffectiveCharge[2][particleTypeIndex], 2.) );
|
|---|
| 432 |
|
|---|
| 433 | G4double cross = sigma0 * zEff * zEff * numerator / denominator;
|
|---|
| 434 |
|
|---|
| 435 | return cross;
|
|---|
| 436 | }
|
|---|
| 437 |
|
|---|
| 438 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 439 |
|
|---|
| 440 | G4int G4DNAMillerGreenExcitationModel::RandomSelect(G4double k,const G4ParticleDefinition* particle)
|
|---|
| 441 | {
|
|---|
| 442 | G4int i = nLevels;
|
|---|
| 443 | G4double value = 0.;
|
|---|
| 444 | std::deque<double> values;
|
|---|
| 445 |
|
|---|
| 446 | G4DNAGenericIonsManager *instance;
|
|---|
| 447 | instance = G4DNAGenericIonsManager::Instance();
|
|---|
| 448 |
|
|---|
| 449 | if ( particle == instance->GetIon("alpha++") ||
|
|---|
| 450 | particle == G4Proton::ProtonDefinition() )
|
|---|
| 451 | {
|
|---|
| 452 | while (i > 0)
|
|---|
| 453 | {
|
|---|
| 454 | i--;
|
|---|
| 455 | G4double partial = PartialCrossSection(k,i,particle);
|
|---|
| 456 | values.push_front(partial);
|
|---|
| 457 | value += partial;
|
|---|
| 458 | }
|
|---|
| 459 |
|
|---|
| 460 | value *= G4UniformRand();
|
|---|
| 461 |
|
|---|
| 462 | i = nLevels;
|
|---|
| 463 |
|
|---|
| 464 | while (i > 0)
|
|---|
| 465 | {
|
|---|
| 466 | i--;
|
|---|
| 467 | if (values[i] > value) return i;
|
|---|
| 468 | value -= values[i];
|
|---|
| 469 | }
|
|---|
| 470 | }
|
|---|
| 471 |
|
|---|
| 472 | // add ONE or TWO electron-water excitation for alpha+ and helium
|
|---|
| 473 |
|
|---|
| 474 | if ( particle == instance->GetIon("alpha+")
|
|---|
| 475 | ||
|
|---|
| 476 | particle == instance->GetIon("helium")
|
|---|
| 477 | )
|
|---|
| 478 | {
|
|---|
| 479 | while (i>0)
|
|---|
| 480 | {
|
|---|
| 481 | i--;
|
|---|
| 482 |
|
|---|
| 483 | G4DNAEmfietzoglouExcitationModel * excitationXS = new G4DNAEmfietzoglouExcitationModel();
|
|---|
| 484 |
|
|---|
| 485 | G4double sigmaExcitation=0;
|
|---|
| 486 |
|
|---|
| 487 | if (k*0.511/3728 > 7.4*eV && k*0.511/3728 < 10*keV) sigmaExcitation = excitationXS->PartialCrossSection(k*0.511/3728,i);
|
|---|
| 488 |
|
|---|
| 489 | G4double partial = PartialCrossSection(k,i,particle);
|
|---|
| 490 | if (particle == instance->GetIon("alpha+")) partial = PartialCrossSection(k,i,particle) + sigmaExcitation;
|
|---|
| 491 | if (particle == instance->GetIon("helium")) partial = PartialCrossSection(k,i,particle) + 2*sigmaExcitation;
|
|---|
| 492 | values.push_front(partial);
|
|---|
| 493 | value += partial;
|
|---|
| 494 | delete excitationXS;
|
|---|
| 495 | }
|
|---|
| 496 |
|
|---|
| 497 | value*=G4UniformRand();
|
|---|
| 498 |
|
|---|
| 499 | i=5;
|
|---|
| 500 | while (i>0)
|
|---|
| 501 | {
|
|---|
| 502 | i--;
|
|---|
| 503 |
|
|---|
| 504 | if (values[i]>value) return i;
|
|---|
| 505 |
|
|---|
| 506 | value-=values[i];
|
|---|
| 507 | }
|
|---|
| 508 | }
|
|---|
| 509 |
|
|---|
| 510 | return 0;
|
|---|
| 511 | }
|
|---|
| 512 |
|
|---|
| 513 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 514 |
|
|---|
| 515 | G4double G4DNAMillerGreenExcitationModel::Sum(G4double k, const G4ParticleDefinition* particle)
|
|---|
| 516 | {
|
|---|
| 517 | G4double totalCrossSection = 0.;
|
|---|
| 518 |
|
|---|
| 519 | for (G4int i=0; i<nLevels; i++)
|
|---|
| 520 | {
|
|---|
| 521 | totalCrossSection += PartialCrossSection(k,i,particle);
|
|---|
| 522 | }
|
|---|
| 523 | return totalCrossSection;
|
|---|
| 524 | }
|
|---|
| 525 |
|
|---|
| 526 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 527 |
|
|---|
| 528 | G4double G4DNAMillerGreenExcitationModel::S_1s(G4double t,
|
|---|
| 529 | G4double energyTransferred,
|
|---|
| 530 | G4double slaterEffectiveCharge,
|
|---|
| 531 | G4double shellNumber)
|
|---|
| 532 | {
|
|---|
| 533 | // 1 - e^(-2r) * ( 1 + 2 r + 2 r^2)
|
|---|
| 534 | // Dingfelder, in Chattanooga 2005 proceedings, formula (7)
|
|---|
| 535 |
|
|---|
| 536 | G4double r = R(t, energyTransferred, slaterEffectiveCharge, shellNumber);
|
|---|
| 537 | G4double value = 1. - std::exp(-2 * r) * ( ( 2. * r + 2. ) * r + 1. );
|
|---|
| 538 |
|
|---|
| 539 | return value;
|
|---|
| 540 | }
|
|---|
| 541 |
|
|---|
| 542 |
|
|---|
| 543 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 544 |
|
|---|
| 545 | G4double G4DNAMillerGreenExcitationModel::S_2s(G4double t,
|
|---|
| 546 | G4double energyTransferred,
|
|---|
| 547 | G4double slaterEffectiveCharge,
|
|---|
| 548 | G4double shellNumber)
|
|---|
| 549 | {
|
|---|
| 550 | // 1 - e^(-2 r) * ( 1 + 2 r + 2 r^2 + 2 r^4)
|
|---|
| 551 | // Dingfelder, in Chattanooga 2005 proceedings, formula (8)
|
|---|
| 552 |
|
|---|
| 553 | G4double r = R(t, energyTransferred, slaterEffectiveCharge, shellNumber);
|
|---|
| 554 | G4double value = 1. - std::exp(-2 * r) * (((2. * r * r + 2.) * r + 2.) * r + 1.);
|
|---|
| 555 |
|
|---|
| 556 | return value;
|
|---|
| 557 |
|
|---|
| 558 | }
|
|---|
| 559 |
|
|---|
| 560 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 561 |
|
|---|
| 562 | G4double G4DNAMillerGreenExcitationModel::S_2p(G4double t,
|
|---|
| 563 | G4double energyTransferred,
|
|---|
| 564 | G4double slaterEffectiveCharge,
|
|---|
| 565 | G4double shellNumber)
|
|---|
| 566 | {
|
|---|
| 567 | // 1 - e^(-2 r) * ( 1 + 2 r + 2 r^2 + 4/3 r^3 + 2/3 r^4)
|
|---|
| 568 | // Dingfelder, in Chattanooga 2005 proceedings, formula (9)
|
|---|
| 569 |
|
|---|
| 570 | G4double r = R(t, energyTransferred, slaterEffectiveCharge, shellNumber);
|
|---|
| 571 | G4double value = 1. - std::exp(-2 * r) * (((( 2./3. * r + 4./3.) * r + 2.) * r + 2.) * r + 1.);
|
|---|
| 572 |
|
|---|
| 573 | return value;
|
|---|
| 574 | }
|
|---|
| 575 |
|
|---|
| 576 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 577 |
|
|---|
| 578 | G4double G4DNAMillerGreenExcitationModel::R(G4double t,
|
|---|
| 579 | G4double energyTransferred,
|
|---|
| 580 | G4double slaterEffectiveCharge,
|
|---|
| 581 | G4double shellNumber)
|
|---|
| 582 | {
|
|---|
| 583 | // tElectron = m_electron / m_alpha * t
|
|---|
| 584 | // Dingfelder, in Chattanooga 2005 proceedings, p 4
|
|---|
| 585 |
|
|---|
| 586 | G4double tElectron = 0.511/3728. * t;
|
|---|
| 587 | G4double value = 2. * tElectron * slaterEffectiveCharge / (energyTransferred * shellNumber);
|
|---|
| 588 |
|
|---|
| 589 | return value;
|
|---|
| 590 | }
|
|---|
| 591 |
|
|---|