#include <G4DNABornIonisationModel.hh>
Inheritance diagram for G4DNABornIonisationModel:
Public Member Functions | |
G4DNABornIonisationModel (const G4ParticleDefinition *p=0, const G4String &nam="DNABornIonisationModel") | |
virtual | ~G4DNABornIonisationModel () |
virtual void | Initialise (const G4ParticleDefinition *, const G4DataVector &=*(new 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) |
double | DifferentialCrossSection (G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell) |
Protected Attributes | |
G4ParticleChangeForGamma * | fParticleChangeForGamma |
Definition at line 48 of file G4DNABornIonisationModel.hh.
G4DNABornIonisationModel::G4DNABornIonisationModel | ( | const G4ParticleDefinition * | p = 0 , |
|
const G4String & | nam = "DNABornIonisationModel" | |||
) |
Definition at line 43 of file G4DNABornIonisationModel.cc.
References fParticleChangeForGamma, G4cout, G4endl, and G4VEmModel::SetDeexcitationFlag().
00045 :G4VEmModel(nam),isInitialised(false) 00046 { 00047 // nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER"); 00048 00049 verboseLevel= 0; 00050 // Verbosity scale: 00051 // 0 = nothing 00052 // 1 = warning for energy non-conservation 00053 // 2 = details of energy budget 00054 // 3 = calculation of cross sections, file openings, sampling of atoms 00055 // 4 = entering in methods 00056 00057 if( verboseLevel>0 ) 00058 { 00059 G4cout << "Born ionisation model is constructed " << G4endl; 00060 } 00061 00062 //Mark this model as "applicable" for atomic deexcitation 00063 SetDeexcitationFlag(true); 00064 fAtomDeexcitation = 0; 00065 fParticleChangeForGamma = 0; 00066 fpMolWaterDensity = 0; 00067 }
G4DNABornIonisationModel::~G4DNABornIonisationModel | ( | ) | [virtual] |
Definition at line 71 of file G4DNABornIonisationModel.cc.
00072 { 00073 // Cross section 00074 00075 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos; 00076 for (pos = tableData.begin(); pos != tableData.end(); ++pos) 00077 { 00078 G4DNACrossSectionDataSet* table = pos->second; 00079 delete table; 00080 } 00081 00082 // Final state 00083 00084 eVecm.clear(); 00085 pVecm.clear(); 00086 00087 }
G4double G4DNABornIonisationModel::CrossSectionPerVolume | ( | const G4Material * | material, | |
const G4ParticleDefinition * | p, | |||
G4double | ekin, | |||
G4double | emin, | |||
G4double | emax | |||
) | [virtual] |
Reimplemented from G4VEmModel.
Definition at line 243 of file G4DNABornIonisationModel.cc.
References G4Electron::ElectronDefinition(), FatalException, G4cout, G4endl, G4Exception(), G4Material::GetIndex(), G4ParticleDefinition::GetParticleName(), and G4Proton::ProtonDefinition().
00248 { 00249 if (verboseLevel > 3) 00250 G4cout << "Calling CrossSectionPerVolume() of G4DNABornIonisationModel" << G4endl; 00251 00252 if ( 00253 particleDefinition != G4Proton::ProtonDefinition() 00254 && 00255 particleDefinition != G4Electron::ElectronDefinition() 00256 ) 00257 00258 return 0; 00259 00260 // Calculate total cross section for model 00261 00262 G4double lowLim = 0; 00263 G4double highLim = 0; 00264 G4double sigma=0; 00265 00266 G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()]; 00267 00268 if(waterDensity!= 0.0) 00269 // if (material == nistwater || material->GetBaseMaterial() == nistwater) 00270 { 00271 const G4String& particleName = particleDefinition->GetParticleName(); 00272 00273 std::map< G4String,G4double,std::less<G4String> >::iterator pos1; 00274 pos1 = lowEnergyLimit.find(particleName); 00275 if (pos1 != lowEnergyLimit.end()) 00276 { 00277 lowLim = pos1->second; 00278 } 00279 00280 std::map< G4String,G4double,std::less<G4String> >::iterator pos2; 00281 pos2 = highEnergyLimit.find(particleName); 00282 if (pos2 != highEnergyLimit.end()) 00283 { 00284 highLim = pos2->second; 00285 } 00286 00287 if (ekin >= lowLim && ekin < highLim) 00288 { 00289 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos; 00290 pos = tableData.find(particleName); 00291 00292 if (pos != tableData.end()) 00293 { 00294 G4DNACrossSectionDataSet* table = pos->second; 00295 if (table != 0) 00296 { 00297 sigma = table->FindValue(ekin); 00298 } 00299 } 00300 else 00301 { 00302 G4Exception("G4DNABornIonisationModel::CrossSectionPerVolume","em0002", 00303 FatalException,"Model not applicable to particle type."); 00304 } 00305 } 00306 00307 if (verboseLevel > 2) 00308 { 00309 G4cout << "__________________________________" << G4endl; 00310 G4cout << "°°° G4DNABornIonisationModel - XS INFO START" << G4endl; 00311 G4cout << "°°° Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl; 00312 G4cout << "°°° Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl; 00313 G4cout << "°°° Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl; 00314 // G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl; 00315 G4cout << "°°° G4DNABornIonisationModel - XS INFO END" << G4endl; 00316 } 00317 00318 } // if (waterMaterial) 00319 00320 // return sigma*material->GetAtomicNumDensityVector()[1]; 00321 return sigma*waterDensity; 00322 }
double G4DNABornIonisationModel::DifferentialCrossSection | ( | G4ParticleDefinition * | aParticleDefinition, | |
G4double | k, | |||
G4double | energyTransfer, | |||
G4int | shell | |||
) |
Definition at line 584 of file G4DNABornIonisationModel.cc.
References G4Electron::ElectronDefinition(), G4DNAWaterIonisationStructure::IonisationEnergy(), and G4Proton::ProtonDefinition().
00588 { 00589 G4double sigma = 0.; 00590 00591 if (energyTransfer >= waterStructure.IonisationEnergy(ionizationLevelIndex)) 00592 { 00593 G4double valueT1 = 0; 00594 G4double valueT2 = 0; 00595 G4double valueE21 = 0; 00596 G4double valueE22 = 0; 00597 G4double valueE12 = 0; 00598 G4double valueE11 = 0; 00599 00600 G4double xs11 = 0; 00601 G4double xs12 = 0; 00602 G4double xs21 = 0; 00603 G4double xs22 = 0; 00604 00605 if (particleDefinition == G4Electron::ElectronDefinition()) 00606 { 00607 // k should be in eV and energy transfer eV also 00608 00609 std::vector<double>::iterator t2 = std::upper_bound(eTdummyVec.begin(),eTdummyVec.end(), k); 00610 00611 std::vector<double>::iterator t1 = t2-1; 00612 00613 // SI : the following condition avoids situations where energyTransfer >last vector element 00614 if (energyTransfer <= eVecm[(*t1)].back() && energyTransfer <= eVecm[(*t2)].back() ) 00615 { 00616 std::vector<double>::iterator e12 = std::upper_bound(eVecm[(*t1)].begin(),eVecm[(*t1)].end(), energyTransfer); 00617 std::vector<double>::iterator e11 = e12-1; 00618 00619 std::vector<double>::iterator e22 = std::upper_bound(eVecm[(*t2)].begin(),eVecm[(*t2)].end(), energyTransfer); 00620 std::vector<double>::iterator e21 = e22-1; 00621 00622 valueT1 =*t1; 00623 valueT2 =*t2; 00624 valueE21 =*e21; 00625 valueE22 =*e22; 00626 valueE12 =*e12; 00627 valueE11 =*e11; 00628 00629 xs11 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE11]; 00630 xs12 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE12]; 00631 xs21 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE21]; 00632 xs22 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE22]; 00633 } 00634 00635 } 00636 00637 if (particleDefinition == G4Proton::ProtonDefinition()) 00638 { 00639 // k should be in eV and energy transfer eV also 00640 std::vector<double>::iterator t2 = std::upper_bound(pTdummyVec.begin(),pTdummyVec.end(), k); 00641 std::vector<double>::iterator t1 = t2-1; 00642 00643 std::vector<double>::iterator e12 = std::upper_bound(pVecm[(*t1)].begin(),pVecm[(*t1)].end(), energyTransfer); 00644 std::vector<double>::iterator e11 = e12-1; 00645 00646 std::vector<double>::iterator e22 = std::upper_bound(pVecm[(*t2)].begin(),pVecm[(*t2)].end(), energyTransfer); 00647 std::vector<double>::iterator e21 = e22-1; 00648 00649 valueT1 =*t1; 00650 valueT2 =*t2; 00651 valueE21 =*e21; 00652 valueE22 =*e22; 00653 valueE12 =*e12; 00654 valueE11 =*e11; 00655 00656 xs11 = pDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE11]; 00657 xs12 = pDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE12]; 00658 xs21 = pDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE21]; 00659 xs22 = pDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE22]; 00660 00661 } 00662 00663 G4double xsProduct = xs11 * xs12 * xs21 * xs22; 00664 if (xsProduct != 0.) 00665 { 00666 sigma = QuadInterpolator( valueE11, valueE12, 00667 valueE21, valueE22, 00668 xs11, xs12, 00669 xs21, xs22, 00670 valueT1, valueT2, 00671 k, energyTransfer); 00672 } 00673 00674 } 00675 00676 return sigma; 00677 }
void G4DNABornIonisationModel::Initialise | ( | const G4ParticleDefinition * | , | |
const G4DataVector & | = *(new G4DataVector()) | |||
) | [virtual] |
Implements G4VEmModel.
Definition at line 91 of file G4DNABornIonisationModel.cc.
References G4LossTableManager::AtomDeexcitation(), G4Electron::ElectronDefinition(), FatalException, fParticleChangeForGamma, G4cout, G4endl, G4Exception(), G4Material::GetMaterial(), G4DNAMolecularMaterial::GetNumMolPerVolTableFor(), G4VEmModel::GetParticleChangeForGamma(), G4ParticleDefinition::GetParticleName(), G4VEmModel::HighEnergyLimit(), G4LossTableManager::Instance(), G4DNAMolecularMaterial::Instance(), G4DNACrossSectionDataSet::LoadData(), G4VEmModel::LowEnergyLimit(), G4InuclParticleNames::proton, G4Proton::ProtonDefinition(), G4VEmModel::SetHighEnergyLimit(), and G4VEmModel::SetLowEnergyLimit().
00093 { 00094 00095 if (verboseLevel > 3) 00096 G4cout << "Calling G4DNABornIonisationModel::Initialise()" << G4endl; 00097 00098 // Energy limits 00099 00100 G4String fileElectron("dna/sigma_ionisation_e_born"); 00101 G4String fileProton("dna/sigma_ionisation_p_born"); 00102 00103 G4ParticleDefinition* electronDef = G4Electron::ElectronDefinition(); 00104 G4ParticleDefinition* protonDef = G4Proton::ProtonDefinition(); 00105 00106 G4String electron; 00107 G4String proton; 00108 00109 G4double scaleFactor = (1.e-22 / 3.343) * m*m; 00110 00111 char *path = getenv("G4LEDATA"); 00112 00113 // *** ELECTRON 00114 00115 electron = electronDef->GetParticleName(); 00116 00117 tableFile[electron] = fileElectron; 00118 00119 lowEnergyLimit[electron] = 11. * eV; 00120 highEnergyLimit[electron] = 1. * MeV; 00121 00122 // Cross section 00123 00124 G4DNACrossSectionDataSet* tableE = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV,scaleFactor ); 00125 tableE->LoadData(fileElectron); 00126 00127 tableData[electron] = tableE; 00128 00129 // Final state 00130 00131 std::ostringstream eFullFileName; 00132 eFullFileName << path << "/dna/sigmadiff_ionisation_e_born.dat"; 00133 std::ifstream eDiffCrossSection(eFullFileName.str().c_str()); 00134 00135 if (!eDiffCrossSection) 00136 { 00137 G4Exception("G4DNABornIonisationModel::Initialise","em0003", 00138 FatalException,"Missing data file:/dna/sigmadiff_ionisation_e_born.dat"); 00139 } 00140 00141 eTdummyVec.push_back(0.); 00142 while(!eDiffCrossSection.eof()) 00143 { 00144 double tDummy; 00145 double eDummy; 00146 eDiffCrossSection>>tDummy>>eDummy; 00147 if (tDummy != eTdummyVec.back()) eTdummyVec.push_back(tDummy); 00148 for (int j=0; j<5; j++) 00149 { 00150 eDiffCrossSection>>eDiffCrossSectionData[j][tDummy][eDummy]; 00151 00152 // SI - only if eof is not reached ! 00153 if (!eDiffCrossSection.eof()) eDiffCrossSectionData[j][tDummy][eDummy]*=scaleFactor; 00154 00155 eVecm[tDummy].push_back(eDummy); 00156 00157 } 00158 } 00159 00160 // *** PROTON 00161 00162 proton = protonDef->GetParticleName(); 00163 00164 tableFile[proton] = fileProton; 00165 00166 lowEnergyLimit[proton] = 500. * keV; 00167 highEnergyLimit[proton] = 100. * MeV; 00168 00169 // Cross section 00170 00171 G4DNACrossSectionDataSet* tableP = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV,scaleFactor ); 00172 tableP->LoadData(fileProton); 00173 00174 tableData[proton] = tableP; 00175 00176 // Final state 00177 00178 std::ostringstream pFullFileName; 00179 pFullFileName << path << "/dna/sigmadiff_ionisation_p_born.dat"; 00180 std::ifstream pDiffCrossSection(pFullFileName.str().c_str()); 00181 00182 if (!pDiffCrossSection) 00183 { 00184 G4Exception("G4DNABornIonisationModel::Initialise","em0003", 00185 FatalException,"Missing data file:/dna/sigmadiff_ionisation_p_born.dat"); 00186 } 00187 00188 pTdummyVec.push_back(0.); 00189 while(!pDiffCrossSection.eof()) 00190 { 00191 double tDummy; 00192 double eDummy; 00193 pDiffCrossSection>>tDummy>>eDummy; 00194 if (tDummy != pTdummyVec.back()) pTdummyVec.push_back(tDummy); 00195 for (int j=0; j<5; j++) 00196 { 00197 pDiffCrossSection>>pDiffCrossSectionData[j][tDummy][eDummy]; 00198 00199 // SI - only if eof is not reached ! 00200 if (!pDiffCrossSection.eof()) pDiffCrossSectionData[j][tDummy][eDummy]*=scaleFactor; 00201 00202 pVecm[tDummy].push_back(eDummy); 00203 } 00204 } 00205 00206 // 00207 00208 if (particle==electronDef) 00209 { 00210 SetLowEnergyLimit(lowEnergyLimit[electron]); 00211 SetHighEnergyLimit(highEnergyLimit[electron]); 00212 } 00213 00214 if (particle==protonDef) 00215 { 00216 SetLowEnergyLimit(lowEnergyLimit[proton]); 00217 SetHighEnergyLimit(highEnergyLimit[proton]); 00218 } 00219 00220 if( verboseLevel>0 ) 00221 { 00222 G4cout << "Born ionisation model is initialized " << G4endl 00223 << "Energy range: " 00224 << LowEnergyLimit() / eV << " eV - " 00225 << HighEnergyLimit() / keV << " keV for " 00226 << particle->GetParticleName() 00227 << G4endl; 00228 } 00229 00230 // Initialize water density pointer 00231 fpMolWaterDensity = G4DNAMolecularMaterial::Instance()->GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER")); 00232 00233 // 00234 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation(); 00235 00236 if (isInitialised) { return; } 00237 fParticleChangeForGamma = GetParticleChangeForGamma(); 00238 isInitialised = true; 00239 }
void G4DNABornIonisationModel::SampleSecondaries | ( | std::vector< G4DynamicParticle * > * | , | |
const G4MaterialCutsCouple * | , | |||
const G4DynamicParticle * | , | |||
G4double | tmin, | |||
G4double | maxEnergy | |||
) | [virtual] |
Implements G4VEmModel.
Definition at line 326 of file G4DNABornIonisationModel.cc.
References G4InuclSpecialFunctions::bindingEnergy(), G4DNAChemistryManager::CreateWaterMolecule(), eIonizedMolecule, G4Electron::Electron(), G4Electron::ElectronDefinition(), fKShell, fParticleChangeForGamma, G4cout, G4endl, G4VAtomDeexcitation::GenerateParticles(), G4VAtomDeexcitation::GetAtomicShell(), G4ParticleChangeForGamma::GetCurrentTrack(), G4DynamicParticle::GetDefinition(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetPDGMass(), G4DNAChemistryManager::Instance(), G4DNAWaterIonisationStructure::IonisationEnergy(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChangeForGamma::ProposeMomentumDirection(), and G4ParticleChangeForGamma::SetProposedKineticEnergy().
00331 { 00332 00333 if (verboseLevel > 3) 00334 G4cout << "Calling SampleSecondaries() of G4DNABornIonisationModel" << G4endl; 00335 00336 G4double lowLim = 0; 00337 G4double highLim = 0; 00338 00339 G4double k = particle->GetKineticEnergy(); 00340 00341 const G4String& particleName = particle->GetDefinition()->GetParticleName(); 00342 00343 std::map< G4String,G4double,std::less<G4String> >::iterator pos1; 00344 pos1 = lowEnergyLimit.find(particleName); 00345 00346 if (pos1 != lowEnergyLimit.end()) 00347 { 00348 lowLim = pos1->second; 00349 } 00350 00351 std::map< G4String,G4double,std::less<G4String> >::iterator pos2; 00352 pos2 = highEnergyLimit.find(particleName); 00353 00354 if (pos2 != highEnergyLimit.end()) 00355 { 00356 highLim = pos2->second; 00357 } 00358 00359 if (k >= lowLim && k < highLim) 00360 { 00361 G4ParticleMomentum primaryDirection = particle->GetMomentumDirection(); 00362 G4double particleMass = particle->GetDefinition()->GetPDGMass(); 00363 G4double totalEnergy = k + particleMass; 00364 G4double pSquare = k * (totalEnergy + particleMass); 00365 G4double totalMomentum = std::sqrt(pSquare); 00366 00367 G4int ionizationShell = RandomSelect(k,particleName); 00368 00369 // sample deexcitation 00370 // here we assume that H_{2}O electronic levels are the same of Oxigen. 00371 // this can be considered true with a rough 10% error in energy on K-shell, 00372 00373 G4int secNumberInit = 0; // need to know at a certain point the enrgy of secondaries 00374 G4int secNumberFinal = 0; // So I'll make the diference and then sum the energies 00375 00376 G4double bindingEnergy = 0; 00377 bindingEnergy = waterStructure.IonisationEnergy(ionizationShell); 00378 00379 if(fAtomDeexcitation) { 00380 G4int Z = 8; 00381 G4AtomicShellEnumerator as = fKShell; 00382 00383 if (ionizationShell <5 && ionizationShell >1) 00384 { 00385 as = G4AtomicShellEnumerator(4-ionizationShell); 00386 } 00387 else if (ionizationShell <2) 00388 { 00389 as = G4AtomicShellEnumerator(3); 00390 } 00391 00392 // FOR DEBUG ONLY 00393 // if (ionizationShell == 4) { 00394 // 00395 // G4cout << "Z: " << Z << " as: " << as 00396 // << " ionizationShell: " << ionizationShell << " bindingEnergy: "<< bindingEnergy/eV << G4endl; 00397 // G4cout << "Press <Enter> key to continue..." << G4endl; 00398 // G4cin.ignore(); 00399 // } 00400 00401 const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as); 00402 secNumberInit = fvect->size(); 00403 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, 0, 0); 00404 secNumberFinal = fvect->size(); 00405 } 00406 00407 G4double secondaryKinetic = RandomizeEjectedElectronEnergy(particle->GetDefinition(),k,ionizationShell); 00408 00409 G4double cosTheta = 0.; 00410 G4double phi = 0.; 00411 RandomizeEjectedElectronDirection(particle->GetDefinition(), k,secondaryKinetic, cosTheta, phi); 00412 00413 G4double sinTheta = std::sqrt(1.-cosTheta*cosTheta); 00414 G4double dirX = sinTheta*std::cos(phi); 00415 G4double dirY = sinTheta*std::sin(phi); 00416 G4double dirZ = cosTheta; 00417 G4ThreeVector deltaDirection(dirX,dirY,dirZ); 00418 deltaDirection.rotateUz(primaryDirection); 00419 00420 if (particle->GetDefinition() == G4Electron::ElectronDefinition()) 00421 { 00422 G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 )); 00423 00424 G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x(); 00425 G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y(); 00426 G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z(); 00427 G4double finalMomentum = std::sqrt(finalPx*finalPx + finalPy*finalPy + finalPz*finalPz); 00428 finalPx /= finalMomentum; 00429 finalPy /= finalMomentum; 00430 finalPz /= finalMomentum; 00431 00432 G4ThreeVector direction; 00433 direction.set(finalPx,finalPy,finalPz); 00434 00435 fParticleChangeForGamma->ProposeMomentumDirection(direction.unit()) ; 00436 } 00437 00438 else fParticleChangeForGamma->ProposeMomentumDirection(primaryDirection) ; 00439 00440 // note thta secondaryKinetic is the nergy of the delta ray, not of all secondaries. 00441 G4double scatteredEnergy = k-bindingEnergy-secondaryKinetic; 00442 G4double deexSecEnergy = 0; 00443 for (G4int j=secNumberInit; j < secNumberFinal; j++) { 00444 00445 deexSecEnergy = deexSecEnergy + (*fvect)[j]->GetKineticEnergy(); 00446 00447 } 00448 00449 fParticleChangeForGamma->SetProposedKineticEnergy(scatteredEnergy); 00450 fParticleChangeForGamma->ProposeLocalEnergyDeposit(k-scatteredEnergy-secondaryKinetic-deexSecEnergy); 00451 00452 G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic) ; 00453 fvect->push_back(dp); 00454 00455 00456 const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack(); 00457 G4DNAChemistryManager::Instance()->CreateWaterMolecule(eIonizedMolecule, 00458 ionizationShell, 00459 theIncomingTrack); 00460 } 00461 00462 }
Definition at line 76 of file G4DNABornIonisationModel.hh.
Referenced by G4DNABornIonisationModel(), Initialise(), and SampleSecondaries().