#include <G4DNARuddIonisationModel.hh>
Inheritance diagram for G4DNARuddIonisationModel:
Public Member Functions | |
G4DNARuddIonisationModel (const G4ParticleDefinition *p=0, const G4String &nam="DNARuddIonisationModel") | |
virtual | ~G4DNARuddIonisationModel () |
virtual void | Initialise (const G4ParticleDefinition *, const G4DataVector &) |
virtual G4double | CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) |
virtual void | SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) |
Protected Attributes | |
G4ParticleChangeForGamma * | fParticleChangeForGamma |
Definition at line 46 of file G4DNARuddIonisationModel.hh.
G4DNARuddIonisationModel::G4DNARuddIonisationModel | ( | const G4ParticleDefinition * | p = 0 , |
|
const G4String & | nam = "DNARuddIonisationModel" | |||
) |
Definition at line 44 of file G4DNARuddIonisationModel.cc.
References fParticleChangeForGamma, G4cout, G4endl, and G4VEmModel::SetDeexcitationFlag().
00046 :G4VEmModel(nam),isInitialised(false) 00047 { 00048 // nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER"); 00049 fpWaterDensity = 0; 00050 00051 slaterEffectiveCharge[0]=0.; 00052 slaterEffectiveCharge[1]=0.; 00053 slaterEffectiveCharge[2]=0.; 00054 sCoefficient[0]=0.; 00055 sCoefficient[1]=0.; 00056 sCoefficient[2]=0.; 00057 00058 lowEnergyLimitForZ1 = 0 * eV; 00059 lowEnergyLimitForZ2 = 0 * eV; 00060 lowEnergyLimitOfModelForZ1 = 100 * eV; 00061 lowEnergyLimitOfModelForZ2 = 1 * keV; 00062 killBelowEnergyForZ1 = lowEnergyLimitOfModelForZ1; 00063 killBelowEnergyForZ2 = lowEnergyLimitOfModelForZ2; 00064 00065 verboseLevel= 0; 00066 // Verbosity scale: 00067 // 0 = nothing 00068 // 1 = warning for energy non-conservation 00069 // 2 = details of energy budget 00070 // 3 = calculation of cross sections, file openings, sampling of atoms 00071 // 4 = entering in methods 00072 00073 if( verboseLevel>0 ) 00074 { 00075 G4cout << "Rudd ionisation model is constructed " << G4endl; 00076 } 00077 00078 //Mark this model as "applicable" for atomic deexcitation 00079 SetDeexcitationFlag(true); 00080 fAtomDeexcitation = 0; 00081 fParticleChangeForGamma = 0; 00082 }
G4DNARuddIonisationModel::~G4DNARuddIonisationModel | ( | ) | [virtual] |
Definition at line 86 of file G4DNARuddIonisationModel.cc.
00087 { 00088 // Cross section 00089 00090 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos; 00091 for (pos = tableData.begin(); pos != tableData.end(); ++pos) 00092 { 00093 G4DNACrossSectionDataSet* table = pos->second; 00094 delete table; 00095 } 00096 00097 00098 // The following removal is forbidden since G4VEnergyLossmodel takes care of deletion 00099 // Coverity however will signal this as an error 00100 //if (fAtomDeexcitation) {delete fAtomDeexcitation;} 00101 00102 }
G4double G4DNARuddIonisationModel::CrossSectionPerVolume | ( | const G4Material * | material, | |
const G4ParticleDefinition * | p, | |||
G4double | ekin, | |||
G4double | emin, | |||
G4double | emax | |||
) | [virtual] |
Reimplemented from G4VEmModel.
Definition at line 278 of file G4DNARuddIonisationModel.cc.
References FatalException, G4cout, G4endl, G4Exception(), G4Material::GetIndex(), G4DNAGenericIonsManager::GetIon(), G4ParticleDefinition::GetParticleName(), G4DNAGenericIonsManager::Instance(), and G4Proton::ProtonDefinition().
00283 { 00284 if (verboseLevel > 3) 00285 G4cout << "Calling CrossSectionPerVolume() of G4DNARuddIonisationModel" << G4endl; 00286 00287 // Calculate total cross section for model 00288 00289 G4DNAGenericIonsManager *instance; 00290 instance = G4DNAGenericIonsManager::Instance(); 00291 00292 if ( 00293 particleDefinition != G4Proton::ProtonDefinition() 00294 && 00295 particleDefinition != instance->GetIon("hydrogen") 00296 && 00297 particleDefinition != instance->GetIon("alpha++") 00298 && 00299 particleDefinition != instance->GetIon("alpha+") 00300 && 00301 particleDefinition != instance->GetIon("helium") 00302 ) 00303 00304 return 0; 00305 00306 G4double lowLim = 0; 00307 00308 if ( particleDefinition == G4Proton::ProtonDefinition() 00309 || particleDefinition == instance->GetIon("hydrogen") 00310 ) 00311 00312 lowLim = lowEnergyLimitOfModelForZ1; 00313 00314 if ( particleDefinition == instance->GetIon("alpha++") 00315 || particleDefinition == instance->GetIon("alpha+") 00316 || particleDefinition == instance->GetIon("helium") 00317 ) 00318 00319 lowLim = lowEnergyLimitOfModelForZ2; 00320 00321 G4double highLim = 0; 00322 G4double sigma=0; 00323 00324 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()]; 00325 00326 if(waterDensity!= 0.0) 00327 // if (material == nistwater || material->GetBaseMaterial() == nistwater) 00328 { 00329 const G4String& particleName = particleDefinition->GetParticleName(); 00330 00331 // SI - the following is useless since lowLim is already defined 00332 /* 00333 std::map< G4String,G4double,std::less<G4String> >::iterator pos1; 00334 pos1 = lowEnergyLimit.find(particleName); 00335 00336 if (pos1 != lowEnergyLimit.end()) 00337 { 00338 lowLim = pos1->second; 00339 } 00340 */ 00341 00342 std::map< G4String,G4double,std::less<G4String> >::iterator pos2; 00343 pos2 = highEnergyLimit.find(particleName); 00344 00345 if (pos2 != highEnergyLimit.end()) 00346 { 00347 highLim = pos2->second; 00348 } 00349 00350 if (k <= highLim) 00351 { 00352 //SI : XS must not be zero otherwise sampling of secondaries method ignored 00353 00354 if (k < lowLim) k = lowLim; 00355 00356 // 00357 00358 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos; 00359 pos = tableData.find(particleName); 00360 00361 if (pos != tableData.end()) 00362 { 00363 G4DNACrossSectionDataSet* table = pos->second; 00364 if (table != 0) 00365 { 00366 sigma = table->FindValue(k); 00367 } 00368 } 00369 else 00370 { 00371 G4Exception("G4DNARuddIonisationModel::CrossSectionPerVolume","em0002", 00372 FatalException,"Model not applicable to particle type."); 00373 } 00374 00375 } // if (k >= lowLim && k < highLim) 00376 00377 if (verboseLevel > 2) 00378 { 00379 G4cout << "__________________________________" << G4endl; 00380 G4cout << "°°° G4DNARuddIonisationModel - XS INFO START" << G4endl; 00381 G4cout << "°°° Kinetic energy(eV)=" << k/eV << " particle : " << particleDefinition->GetParticleName() << G4endl; 00382 G4cout << "°°° Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl; 00383 G4cout << "°°° Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl; 00384 // G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl; 00385 G4cout << "°°° G4DNARuddIonisationModel - XS INFO END" << G4endl; 00386 } 00387 00388 } // if (waterMaterial) 00389 else 00390 { 00391 if (verboseLevel > 2) 00392 { 00393 G4cout << "Warning : RuddIonisationModel: WATER DENSITY IS NULL" << G4endl; 00394 } 00395 } 00396 00397 return sigma*waterDensity; 00398 // return sigma*material->GetAtomicNumDensityVector()[1]; 00399 00400 }
void G4DNARuddIonisationModel::Initialise | ( | const G4ParticleDefinition * | , | |
const G4DataVector & | ||||
) | [virtual] |
Implements G4VEmModel.
Definition at line 106 of file G4DNARuddIonisationModel.cc.
References G4LossTableManager::AtomDeexcitation(), fParticleChangeForGamma, G4cout, G4endl, G4DNAGenericIonsManager::GetIon(), G4Material::GetMaterial(), G4DNAMolecularMaterial::GetNumMolPerVolTableFor(), G4VEmModel::GetParticleChangeForGamma(), G4ParticleDefinition::GetParticleName(), G4VEmModel::HighEnergyLimit(), G4LossTableManager::Instance(), G4DNAMolecularMaterial::Instance(), G4DNAGenericIonsManager::Instance(), G4DNACrossSectionDataSet::LoadData(), G4VEmModel::LowEnergyLimit(), G4InuclParticleNames::proton, G4Proton::ProtonDefinition(), G4VEmModel::SetHighEnergyLimit(), and G4VEmModel::SetLowEnergyLimit().
00108 { 00109 00110 if (verboseLevel > 3) 00111 G4cout << "Calling G4DNARuddIonisationModel::Initialise()" << G4endl; 00112 00113 // Energy limits 00114 00115 G4String fileProton("dna/sigma_ionisation_p_rudd"); 00116 G4String fileHydrogen("dna/sigma_ionisation_h_rudd"); 00117 G4String fileAlphaPlusPlus("dna/sigma_ionisation_alphaplusplus_rudd"); 00118 G4String fileAlphaPlus("dna/sigma_ionisation_alphaplus_rudd"); 00119 G4String fileHelium("dna/sigma_ionisation_he_rudd"); 00120 00121 G4DNAGenericIonsManager *instance; 00122 instance = G4DNAGenericIonsManager::Instance(); 00123 G4ParticleDefinition* protonDef = G4Proton::ProtonDefinition(); 00124 G4ParticleDefinition* hydrogenDef = instance->GetIon("hydrogen"); 00125 G4ParticleDefinition* alphaPlusPlusDef = instance->GetIon("alpha++"); 00126 G4ParticleDefinition* alphaPlusDef = instance->GetIon("alpha+"); 00127 G4ParticleDefinition* heliumDef = instance->GetIon("helium"); 00128 00129 G4String proton; 00130 G4String hydrogen; 00131 G4String alphaPlusPlus; 00132 G4String alphaPlus; 00133 G4String helium; 00134 00135 G4double scaleFactor = 1 * m*m; 00136 00137 // LIMITS AND DATA 00138 00139 // ******************************************************** 00140 00141 proton = protonDef->GetParticleName(); 00142 tableFile[proton] = fileProton; 00143 00144 lowEnergyLimit[proton] = lowEnergyLimitForZ1; 00145 highEnergyLimit[proton] = 500. * keV; 00146 00147 // Cross section 00148 00149 G4DNACrossSectionDataSet* tableProton = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, 00150 eV, 00151 scaleFactor ); 00152 tableProton->LoadData(fileProton); 00153 tableData[proton] = tableProton; 00154 00155 // ******************************************************** 00156 00157 hydrogen = hydrogenDef->GetParticleName(); 00158 tableFile[hydrogen] = fileHydrogen; 00159 00160 lowEnergyLimit[hydrogen] = lowEnergyLimitForZ1; 00161 highEnergyLimit[hydrogen] = 100. * MeV; 00162 00163 // Cross section 00164 00165 G4DNACrossSectionDataSet* tableHydrogen = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, 00166 eV, 00167 scaleFactor ); 00168 tableHydrogen->LoadData(fileHydrogen); 00169 00170 tableData[hydrogen] = tableHydrogen; 00171 00172 // ******************************************************** 00173 00174 alphaPlusPlus = alphaPlusPlusDef->GetParticleName(); 00175 tableFile[alphaPlusPlus] = fileAlphaPlusPlus; 00176 00177 lowEnergyLimit[alphaPlusPlus] = lowEnergyLimitForZ2; 00178 highEnergyLimit[alphaPlusPlus] = 400. * MeV; 00179 00180 // Cross section 00181 00182 G4DNACrossSectionDataSet* tableAlphaPlusPlus = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, 00183 eV, 00184 scaleFactor ); 00185 tableAlphaPlusPlus->LoadData(fileAlphaPlusPlus); 00186 00187 tableData[alphaPlusPlus] = tableAlphaPlusPlus; 00188 00189 // ******************************************************** 00190 00191 alphaPlus = alphaPlusDef->GetParticleName(); 00192 tableFile[alphaPlus] = fileAlphaPlus; 00193 00194 lowEnergyLimit[alphaPlus] = lowEnergyLimitForZ2; 00195 highEnergyLimit[alphaPlus] = 400. * MeV; 00196 00197 // Cross section 00198 00199 G4DNACrossSectionDataSet* tableAlphaPlus = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, 00200 eV, 00201 scaleFactor ); 00202 tableAlphaPlus->LoadData(fileAlphaPlus); 00203 tableData[alphaPlus] = tableAlphaPlus; 00204 00205 // ******************************************************** 00206 00207 helium = heliumDef->GetParticleName(); 00208 tableFile[helium] = fileHelium; 00209 00210 lowEnergyLimit[helium] = lowEnergyLimitForZ2; 00211 highEnergyLimit[helium] = 400. * MeV; 00212 00213 // Cross section 00214 00215 G4DNACrossSectionDataSet* tableHelium = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, 00216 eV, 00217 scaleFactor ); 00218 tableHelium->LoadData(fileHelium); 00219 tableData[helium] = tableHelium; 00220 00221 // 00222 00223 if (particle==protonDef) 00224 { 00225 SetLowEnergyLimit(lowEnergyLimit[proton]); 00226 SetHighEnergyLimit(highEnergyLimit[proton]); 00227 } 00228 00229 if (particle==hydrogenDef) 00230 { 00231 SetLowEnergyLimit(lowEnergyLimit[hydrogen]); 00232 SetHighEnergyLimit(highEnergyLimit[hydrogen]); 00233 } 00234 00235 if (particle==heliumDef) 00236 { 00237 SetLowEnergyLimit(lowEnergyLimit[helium]); 00238 SetHighEnergyLimit(highEnergyLimit[helium]); 00239 } 00240 00241 if (particle==alphaPlusDef) 00242 { 00243 SetLowEnergyLimit(lowEnergyLimit[alphaPlus]); 00244 SetHighEnergyLimit(highEnergyLimit[alphaPlus]); 00245 } 00246 00247 if (particle==alphaPlusPlusDef) 00248 { 00249 SetLowEnergyLimit(lowEnergyLimit[alphaPlusPlus]); 00250 SetHighEnergyLimit(highEnergyLimit[alphaPlusPlus]); 00251 } 00252 00253 if( verboseLevel>0 ) 00254 { 00255 G4cout << "Rudd ionisation model is initialized " << G4endl 00256 << "Energy range: " 00257 << LowEnergyLimit() / eV << " eV - " 00258 << HighEnergyLimit() / keV << " keV for " 00259 << particle->GetParticleName() 00260 << G4endl; 00261 } 00262 00263 // Initialize water density pointer 00264 fpWaterDensity = G4DNAMolecularMaterial::Instance()->GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER")); 00265 00266 // 00267 00268 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation(); 00269 00270 if (isInitialised) { return; } 00271 fParticleChangeForGamma = GetParticleChangeForGamma(); 00272 isInitialised = true; 00273 00274 }
void G4DNARuddIonisationModel::SampleSecondaries | ( | std::vector< G4DynamicParticle * > * | , | |
const G4MaterialCutsCouple * | , | |||
const G4DynamicParticle * | , | |||
G4double | tmin, | |||
G4double | maxEnergy | |||
) | [virtual] |
Implements G4VEmModel.
Definition at line 404 of file G4DNARuddIonisationModel.cc.
References G4InuclSpecialFunctions::bindingEnergy(), G4DNAChemistryManager::CreateWaterMolecule(), eIonizedMolecule, G4Electron::Electron(), fKShell, fParticleChangeForGamma, fStopAndKill, G4cout, G4endl, G4VAtomDeexcitation::GenerateParticles(), G4VAtomDeexcitation::GetAtomicShell(), G4ParticleChangeForGamma::GetCurrentTrack(), G4DynamicParticle::GetDefinition(), G4DNAGenericIonsManager::GetIon(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), G4ParticleDefinition::GetParticleName(), G4DNAChemistryManager::Instance(), G4DNAGenericIonsManager::Instance(), G4DNAWaterIonisationStructure::IonisationEnergy(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChangeForGamma::ProposeMomentumDirection(), G4VParticleChange::ProposeTrackStatus(), G4Proton::ProtonDefinition(), and G4ParticleChangeForGamma::SetProposedKineticEnergy().
00409 { 00410 if (verboseLevel > 3) 00411 G4cout << "Calling SampleSecondaries() of G4DNARuddIonisationModel" << G4endl; 00412 00413 G4double lowLim = 0; 00414 G4double highLim = 0; 00415 00416 G4DNAGenericIonsManager *instance; 00417 instance = G4DNAGenericIonsManager::Instance(); 00418 00419 if ( particle->GetDefinition() == G4Proton::ProtonDefinition() 00420 || particle->GetDefinition() == instance->GetIon("hydrogen") 00421 ) 00422 00423 lowLim = killBelowEnergyForZ1; 00424 00425 if ( particle->GetDefinition() == instance->GetIon("alpha++") 00426 || particle->GetDefinition() == instance->GetIon("alpha+") 00427 || particle->GetDefinition() == instance->GetIon("helium") 00428 ) 00429 00430 lowLim = killBelowEnergyForZ2; 00431 00432 G4double k = particle->GetKineticEnergy(); 00433 00434 const G4String& particleName = particle->GetDefinition()->GetParticleName(); 00435 00436 // SI - the following is useless since lowLim is already defined 00437 /* 00438 std::map< G4String,G4double,std::less<G4String> >::iterator pos1; 00439 pos1 = lowEnergyLimit.find(particleName); 00440 00441 if (pos1 != lowEnergyLimit.end()) 00442 { 00443 lowLim = pos1->second; 00444 } 00445 */ 00446 00447 std::map< G4String,G4double,std::less<G4String> >::iterator pos2; 00448 pos2 = highEnergyLimit.find(particleName); 00449 00450 if (pos2 != highEnergyLimit.end()) 00451 { 00452 highLim = pos2->second; 00453 } 00454 00455 if (k >= lowLim && k <= highLim) 00456 { 00457 G4ParticleDefinition* definition = particle->GetDefinition(); 00458 G4ParticleMomentum primaryDirection = particle->GetMomentumDirection(); 00459 /* 00460 G4double particleMass = definition->GetPDGMass(); 00461 G4double totalEnergy = k + particleMass; 00462 G4double pSquare = k*(totalEnergy+particleMass); 00463 G4double totalMomentum = std::sqrt(pSquare); 00464 */ 00465 00466 G4int ionizationShell = RandomSelect(k,particleName); 00467 00468 // sample deexcitation 00469 // here we assume that H_{2}O electronic levels are the same of Oxigen. 00470 // this can be considered true with a rough 10% error in energy on K-shell, 00471 00472 G4int secNumberInit = 0; // need to know at a certain point the enrgy of secondaries 00473 G4int secNumberFinal = 0; // So I'll make the diference and then sum the energies 00474 00475 G4double bindingEnergy = 0; 00476 bindingEnergy = waterStructure.IonisationEnergy(ionizationShell); 00477 00478 if(fAtomDeexcitation) { 00479 G4int Z = 8; 00480 G4AtomicShellEnumerator as = fKShell; 00481 00482 if (ionizationShell <5 && ionizationShell >1) 00483 { 00484 as = G4AtomicShellEnumerator(4-ionizationShell); 00485 } 00486 else if (ionizationShell <2) 00487 { 00488 as = G4AtomicShellEnumerator(3); 00489 } 00490 00491 00492 00493 // DEBUG 00494 // if (ionizationShell == 4) { 00495 // 00496 // G4cout << "Z: " << Z << " as: " << as 00497 // << " ionizationShell: " << ionizationShell << " bindingEnergy: "<< bindingEnergy/eV << G4endl; 00498 // G4cout << "Press <Enter> key to continue..." << G4endl; 00499 // G4cin.ignore(); 00500 // } 00501 00502 const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as); 00503 secNumberInit = fvect->size(); 00504 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, 0, 0); 00505 secNumberFinal = fvect->size(); 00506 } 00507 00508 00509 G4double secondaryKinetic = RandomizeEjectedElectronEnergy(definition,k,ionizationShell); 00510 00511 G4double cosTheta = 0.; 00512 G4double phi = 0.; 00513 RandomizeEjectedElectronDirection(definition, k,secondaryKinetic, cosTheta, phi); 00514 00515 G4double sinTheta = std::sqrt(1.-cosTheta*cosTheta); 00516 G4double dirX = sinTheta*std::cos(phi); 00517 G4double dirY = sinTheta*std::sin(phi); 00518 G4double dirZ = cosTheta; 00519 G4ThreeVector deltaDirection(dirX,dirY,dirZ); 00520 deltaDirection.rotateUz(primaryDirection); 00521 00522 // Ignored for ions on electrons 00523 /* 00524 G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 )); 00525 00526 G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x(); 00527 G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y(); 00528 G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z(); 00529 G4double finalMomentum = std::sqrt(finalPx*finalPx+finalPy*finalPy+finalPz*finalPz); 00530 finalPx /= finalMomentum; 00531 finalPy /= finalMomentum; 00532 finalPz /= finalMomentum; 00533 00534 G4ThreeVector direction; 00535 direction.set(finalPx,finalPy,finalPz); 00536 00537 fParticleChangeForGamma->ProposeMomentumDirection(direction.unit()) ; 00538 */ 00539 fParticleChangeForGamma->ProposeMomentumDirection(primaryDirection); 00540 00541 G4double scatteredEnergy = k-bindingEnergy-secondaryKinetic; 00542 G4double deexSecEnergy = 0; 00543 for (G4int j=secNumberInit; j < secNumberFinal; j++) { 00544 00545 deexSecEnergy = deexSecEnergy + (*fvect)[j]->GetKineticEnergy(); 00546 00547 } 00548 00549 fParticleChangeForGamma->SetProposedKineticEnergy(scatteredEnergy); 00550 fParticleChangeForGamma->ProposeLocalEnergyDeposit(k-scatteredEnergy-secondaryKinetic-deexSecEnergy); 00551 00552 // debug 00553 // k-scatteredEnergy-secondaryKinetic-deexSecEnergy = k-(k-bindingEnergy-secondaryKinetic)-secondaryKinetic-deexSecEnergy = 00554 // = k-k+bindingEnergy+secondaryKinetic-secondaryKinetic-deexSecEnergy= 00555 // = bindingEnergy-deexSecEnergy 00556 // SO deexSecEnergy=0 => LocalEnergyDeposit = bindingEnergy 00557 00558 G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic) ; 00559 fvect->push_back(dp); 00560 00561 const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack(); 00562 G4DNAChemistryManager::Instance()->CreateWaterMolecule(eIonizedMolecule, 00563 ionizationShell, 00564 theIncomingTrack); 00565 } 00566 00567 // SI - not useful since low energy of model is 0 eV 00568 00569 if (k < lowLim) 00570 { 00571 fParticleChangeForGamma->SetProposedKineticEnergy(0.); 00572 fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill); 00573 fParticleChangeForGamma->ProposeLocalEnergyDeposit(k); 00574 } 00575 00576 }
Definition at line 72 of file G4DNARuddIonisationModel.hh.
Referenced by G4DNARuddIonisationModel(), Initialise(), and SampleSecondaries().