source: trunk/source/processes/electromagnetic/lowenergy/src/G4DNARuddIonisationModel.cc

Last change on this file was 1347, checked in by garnier, 15 years ago

geant4 tag 9.4

File size: 30.2 KB
Line 
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// $Id: G4DNARuddIonisationModel.cc,v 1.21 2010/11/04 14:52:17 sincerti Exp $
27// GEANT4 tag $Name: geant4-09-04-ref-00 $
28//
29
30#include "G4DNARuddIonisationModel.hh"
31//#include "G4DynamicMolecule.hh"
32
33//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
34
35using namespace std;
36
37//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
38
39G4DNARuddIonisationModel::G4DNARuddIonisationModel(const G4ParticleDefinition*,
40 const G4String& nam)
41:G4VEmModel(nam),isInitialised(false)
42{
43 lowEnergyLimitForZ1 = 0 * eV;
44 lowEnergyLimitForZ2 = 0 * eV;
45 lowEnergyLimitOfModelForZ1 = 100 * eV;
46 lowEnergyLimitOfModelForZ2 = 1 * keV;
47 killBelowEnergyForZ1 = lowEnergyLimitOfModelForZ1;
48 killBelowEnergyForZ2 = lowEnergyLimitOfModelForZ2;
49
50 verboseLevel= 0;
51 // Verbosity scale:
52 // 0 = nothing
53 // 1 = warning for energy non-conservation
54 // 2 = details of energy budget
55 // 3 = calculation of cross sections, file openings, sampling of atoms
56 // 4 = entering in methods
57
58 if( verboseLevel>0 )
59 {
60 G4cout << "Rudd ionisation model is constructed " << G4endl;
61 }
62}
63
64//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
65
66G4DNARuddIonisationModel::~G4DNARuddIonisationModel()
67{
68 // Cross section
69
70 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
71 for (pos = tableData.begin(); pos != tableData.end(); ++pos)
72 {
73 G4DNACrossSectionDataSet* table = pos->second;
74 delete table;
75 }
76
77}
78
79//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
80
81void G4DNARuddIonisationModel::Initialise(const G4ParticleDefinition* particle,
82 const G4DataVector& /*cuts*/)
83{
84
85 if (verboseLevel > 3)
86 G4cout << "Calling G4DNARuddIonisationModel::Initialise()" << G4endl;
87
88 // Energy limits
89
90 G4String fileProton("dna/sigma_ionisation_p_rudd");
91 G4String fileHydrogen("dna/sigma_ionisation_h_rudd");
92 G4String fileAlphaPlusPlus("dna/sigma_ionisation_alphaplusplus_rudd");
93 G4String fileAlphaPlus("dna/sigma_ionisation_alphaplus_rudd");
94 G4String fileHelium("dna/sigma_ionisation_he_rudd");
95
96 G4DNAGenericIonsManager *instance;
97 instance = G4DNAGenericIonsManager::Instance();
98 G4ParticleDefinition* protonDef = G4Proton::ProtonDefinition();
99 G4ParticleDefinition* hydrogenDef = instance->GetIon("hydrogen");
100 G4ParticleDefinition* alphaPlusPlusDef = instance->GetIon("alpha++");
101 G4ParticleDefinition* alphaPlusDef = instance->GetIon("alpha+");
102 G4ParticleDefinition* heliumDef = instance->GetIon("helium");
103
104 G4String proton;
105 G4String hydrogen;
106 G4String alphaPlusPlus;
107 G4String alphaPlus;
108 G4String helium;
109
110 G4double scaleFactor = 1 * m*m;
111
112 if (protonDef != 0)
113 {
114 proton = protonDef->GetParticleName();
115 tableFile[proton] = fileProton;
116
117 lowEnergyLimit[proton] = lowEnergyLimitForZ1;
118 highEnergyLimit[proton] = 500. * keV;
119
120 // Cross section
121
122 G4DNACrossSectionDataSet* tableProton = new G4DNACrossSectionDataSet(new G4LogLogInterpolation,
123 eV,
124 scaleFactor );
125 tableProton->LoadData(fileProton);
126 tableData[proton] = tableProton;
127 }
128 else
129 {
130 G4Exception("G4DNARuddIonisationModel::Initialise: proton is not defined");
131 }
132
133 if (hydrogenDef != 0)
134 {
135 hydrogen = hydrogenDef->GetParticleName();
136 tableFile[hydrogen] = fileHydrogen;
137
138 lowEnergyLimit[hydrogen] = lowEnergyLimitForZ1;
139 highEnergyLimit[hydrogen] = 100. * MeV;
140
141 // Cross section
142
143 G4DNACrossSectionDataSet* tableHydrogen = new G4DNACrossSectionDataSet(new G4LogLogInterpolation,
144 eV,
145 scaleFactor );
146 tableHydrogen->LoadData(fileHydrogen);
147
148 tableData[hydrogen] = tableHydrogen;
149 }
150 else
151 {
152 G4Exception("G4DNARuddIonisationModel::Initialise: hydrogen is not defined");
153 }
154
155 if (alphaPlusPlusDef != 0)
156 {
157 alphaPlusPlus = alphaPlusPlusDef->GetParticleName();
158 tableFile[alphaPlusPlus] = fileAlphaPlusPlus;
159
160 lowEnergyLimit[alphaPlusPlus] = lowEnergyLimitForZ2;
161 highEnergyLimit[alphaPlusPlus] = 400. * MeV;
162
163 // Cross section
164
165 G4DNACrossSectionDataSet* tableAlphaPlusPlus = new G4DNACrossSectionDataSet(new G4LogLogInterpolation,
166 eV,
167 scaleFactor );
168 tableAlphaPlusPlus->LoadData(fileAlphaPlusPlus);
169
170 tableData[alphaPlusPlus] = tableAlphaPlusPlus;
171 }
172 else
173 {
174 G4Exception("G4DNARuddIonisationModel::Initialise: alphaPlusPlus is not defined");
175 }
176
177 if (alphaPlusDef != 0)
178 {
179 alphaPlus = alphaPlusDef->GetParticleName();
180 tableFile[alphaPlus] = fileAlphaPlus;
181
182 lowEnergyLimit[alphaPlus] = lowEnergyLimitForZ2;
183 highEnergyLimit[alphaPlus] = 400. * MeV;
184
185 // Cross section
186
187 G4DNACrossSectionDataSet* tableAlphaPlus = new G4DNACrossSectionDataSet(new G4LogLogInterpolation,
188 eV,
189 scaleFactor );
190 tableAlphaPlus->LoadData(fileAlphaPlus);
191 tableData[alphaPlus] = tableAlphaPlus;
192 }
193 else
194 {
195 G4Exception("G4DNARuddIonisationModel::Initialise: alphaPlus is not defined");
196 }
197
198 if (heliumDef != 0)
199 {
200 helium = heliumDef->GetParticleName();
201 tableFile[helium] = fileHelium;
202
203 lowEnergyLimit[helium] = lowEnergyLimitForZ2;
204 highEnergyLimit[helium] = 400. * MeV;
205
206 // Cross section
207
208 G4DNACrossSectionDataSet* tableHelium = new G4DNACrossSectionDataSet(new G4LogLogInterpolation,
209 eV,
210 scaleFactor );
211 tableHelium->LoadData(fileHelium);
212 tableData[helium] = tableHelium;
213 }
214 else
215 {
216 G4Exception("G4DNARuddIonisationModel::Initialise: helium is not defined");
217 }
218
219 if (particle==protonDef)
220 {
221 SetLowEnergyLimit(lowEnergyLimit[proton]);
222 SetHighEnergyLimit(highEnergyLimit[proton]);
223 }
224
225 if (particle==hydrogenDef)
226 {
227 SetLowEnergyLimit(lowEnergyLimit[hydrogen]);
228 SetHighEnergyLimit(highEnergyLimit[hydrogen]);
229 }
230
231 if (particle==heliumDef)
232 {
233 SetLowEnergyLimit(lowEnergyLimit[helium]);
234 SetHighEnergyLimit(highEnergyLimit[helium]);
235 }
236
237 if (particle==alphaPlusDef)
238 {
239 SetLowEnergyLimit(lowEnergyLimit[alphaPlus]);
240 SetHighEnergyLimit(highEnergyLimit[alphaPlus]);
241 }
242
243 if (particle==alphaPlusPlusDef)
244 {
245 SetLowEnergyLimit(lowEnergyLimit[alphaPlusPlus]);
246 SetHighEnergyLimit(highEnergyLimit[alphaPlusPlus]);
247 }
248
249 if( verboseLevel>0 )
250 {
251 G4cout << "Rudd ionisation model is initialized " << G4endl
252 << "Energy range: "
253 << LowEnergyLimit() / eV << " eV - "
254 << HighEnergyLimit() / keV << " keV for "
255 << particle->GetParticleName()
256 << G4endl;
257 }
258
259 //
260
261 if(!isInitialised)
262 {
263 isInitialised = true;
264
265 if(pParticleChange)
266 fParticleChangeForGamma = reinterpret_cast<G4ParticleChangeForGamma*>(pParticleChange);
267 else
268 fParticleChangeForGamma = new G4ParticleChangeForGamma();
269 }
270
271 // InitialiseElementSelectors(particle,cuts);
272
273}
274
275//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
276
277G4double G4DNARuddIonisationModel::CrossSectionPerVolume(const G4Material* material,
278 const G4ParticleDefinition* particleDefinition,
279 G4double k,
280 G4double,
281 G4double)
282{
283 if (verboseLevel > 3)
284 G4cout << "Calling CrossSectionPerVolume() of G4DNARuddIonisationModel" << G4endl;
285
286 // Calculate total cross section for model
287
288 G4DNAGenericIonsManager *instance;
289 instance = G4DNAGenericIonsManager::Instance();
290
291 if (
292 particleDefinition != G4Proton::ProtonDefinition()
293 &&
294 particleDefinition != instance->GetIon("hydrogen")
295 &&
296 particleDefinition != instance->GetIon("alpha++")
297 &&
298 particleDefinition != instance->GetIon("alpha+")
299 &&
300 particleDefinition != instance->GetIon("helium")
301 )
302
303 return 0;
304
305 G4double lowLim = 0;
306
307 if ( particleDefinition == G4Proton::ProtonDefinition()
308 || particleDefinition == instance->GetIon("hydrogen")
309 )
310
311 lowLim = lowEnergyLimitOfModelForZ1;
312
313 if ( particleDefinition == instance->GetIon("alpha++")
314 || particleDefinition == instance->GetIon("alpha+")
315 || particleDefinition == instance->GetIon("helium")
316 )
317
318 lowLim = lowEnergyLimitOfModelForZ2;
319
320 G4double highLim = 0;
321 G4double sigma=0;
322
323 if (material->GetName() == "G4_WATER")
324 {
325 const G4String& particleName = particleDefinition->GetParticleName();
326
327 // SI - the following is useless since lowLim is already defined
328 /*
329 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
330 pos1 = lowEnergyLimit.find(particleName);
331
332 if (pos1 != lowEnergyLimit.end())
333 {
334 lowLim = pos1->second;
335 }
336 */
337
338 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
339 pos2 = highEnergyLimit.find(particleName);
340
341 if (pos2 != highEnergyLimit.end())
342 {
343 highLim = pos2->second;
344 }
345
346 if (k < highLim)
347 {
348
349 //SI : XS must not be zero otherwise sampling of secondaries method ignored
350
351 if (k < lowLim) k = lowLim;
352
353 //
354
355 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
356 pos = tableData.find(particleName);
357
358 if (pos != tableData.end())
359 {
360 G4DNACrossSectionDataSet* table = pos->second;
361 if (table != 0)
362 {
363 sigma = table->FindValue(k);
364
365 }
366 }
367 else
368 {
369 G4Exception("G4DNARuddIonisationModel::CrossSectionPerVolume: attempting to calculate cross section for wrong particle");
370 }
371
372 } // if (k >= lowLim && k < highLim)
373
374 if (verboseLevel > 3)
375 {
376 G4cout << "---> Kinetic energy(eV)=" << k/eV << G4endl;
377 G4cout << " - Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
378 G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*
379 material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
380 }
381
382 } // if (waterMaterial)
383
384 return sigma*material->GetAtomicNumDensityVector()[1];
385
386}
387
388//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
389
390void G4DNARuddIonisationModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
391 const G4MaterialCutsCouple* /*couple*/,
392 const G4DynamicParticle* particle,
393 G4double,
394 G4double)
395{
396 if (verboseLevel > 3)
397 G4cout << "Calling SampleSecondaries() of G4DNARuddIonisationModel" << G4endl;
398
399 G4double lowLim = 0;
400 G4double highLim = 0;
401
402 G4DNAGenericIonsManager *instance;
403 instance = G4DNAGenericIonsManager::Instance();
404
405 if ( particle->GetDefinition() == G4Proton::ProtonDefinition()
406 || particle->GetDefinition() == instance->GetIon("hydrogen")
407 )
408
409 lowLim = killBelowEnergyForZ1;
410
411 if ( particle->GetDefinition() == instance->GetIon("alpha++")
412 || particle->GetDefinition() == instance->GetIon("alpha+")
413 || particle->GetDefinition() == instance->GetIon("helium")
414 )
415
416 lowLim = killBelowEnergyForZ2;
417
418 G4double k = particle->GetKineticEnergy();
419
420 const G4String& particleName = particle->GetDefinition()->GetParticleName();
421
422 // SI - the following is useless since lowLim is already defined
423 /*
424 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
425 pos1 = lowEnergyLimit.find(particleName);
426
427 if (pos1 != lowEnergyLimit.end())
428 {
429 lowLim = pos1->second;
430 }
431 */
432
433 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
434 pos2 = highEnergyLimit.find(particleName);
435
436 if (pos2 != highEnergyLimit.end())
437 {
438 highLim = pos2->second;
439 }
440
441 if (k >= lowLim && k < highLim)
442 {
443 G4ParticleDefinition* definition = particle->GetDefinition();
444 G4ParticleMomentum primaryDirection = particle->GetMomentumDirection();
445 /*
446 G4double particleMass = definition->GetPDGMass();
447 G4double totalEnergy = k + particleMass;
448 G4double pSquare = k*(totalEnergy+particleMass);
449 G4double totalMomentum = std::sqrt(pSquare);
450 */
451
452 G4int ionizationShell = RandomSelect(k,particleName);
453
454 G4double secondaryKinetic = RandomizeEjectedElectronEnergy(definition,k,ionizationShell);
455
456 G4double bindingEnergy = waterStructure.IonisationEnergy(ionizationShell);
457
458 G4double cosTheta = 0.;
459 G4double phi = 0.;
460 RandomizeEjectedElectronDirection(definition, k,secondaryKinetic, cosTheta, phi);
461
462 G4double sinTheta = std::sqrt(1.-cosTheta*cosTheta);
463 G4double dirX = sinTheta*std::cos(phi);
464 G4double dirY = sinTheta*std::sin(phi);
465 G4double dirZ = cosTheta;
466 G4ThreeVector deltaDirection(dirX,dirY,dirZ);
467 deltaDirection.rotateUz(primaryDirection);
468
469 // Ignored for ions on electrons
470 /*
471 G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
472
473 G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
474 G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
475 G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
476 G4double finalMomentum = std::sqrt(finalPx*finalPx+finalPy*finalPy+finalPz*finalPz);
477 finalPx /= finalMomentum;
478 finalPy /= finalMomentum;
479 finalPz /= finalMomentum;
480
481 G4ThreeVector direction;
482 direction.set(finalPx,finalPy,finalPz);
483
484 fParticleChangeForGamma->ProposeMomentumDirection(direction.unit()) ;
485 */
486 fParticleChangeForGamma->ProposeMomentumDirection(primaryDirection);
487
488 fParticleChangeForGamma->SetProposedKineticEnergy(k-bindingEnergy-secondaryKinetic);
489 fParticleChangeForGamma->ProposeLocalEnergyDeposit(bindingEnergy);
490
491 G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic) ;
492 fvect->push_back(dp);
493
494 }
495
496 // SI - not useful since low energy of model is 0 eV
497
498 if (k < lowLim)
499 {
500 fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill);
501 fParticleChangeForGamma->ProposeLocalEnergyDeposit(k);
502 }
503
504}
505
506//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
507
508G4double G4DNARuddIonisationModel::RandomizeEjectedElectronEnergy(G4ParticleDefinition* particleDefinition,
509 G4double k,
510 G4int shell)
511{
512 G4double maximumKineticEnergyTransfer = 0.;
513
514 G4DNAGenericIonsManager *instance;
515 instance = G4DNAGenericIonsManager::Instance();
516
517 if (particleDefinition == G4Proton::ProtonDefinition()
518 || particleDefinition == instance->GetIon("hydrogen"))
519 {
520 maximumKineticEnergyTransfer= 4.* (electron_mass_c2 / proton_mass_c2) * k;
521 }
522
523 if (particleDefinition == instance->GetIon("helium")
524 || particleDefinition == instance->GetIon("alpha+")
525 || particleDefinition == instance->GetIon("alpha++"))
526 {
527 maximumKineticEnergyTransfer= 4.* (0.511 / 3728) * k;
528 }
529
530 G4double crossSectionMaximum = 0.;
531
532 for(G4double value=waterStructure.IonisationEnergy(shell); value<=5.*waterStructure.IonisationEnergy(shell) && k>=value ; value+=0.1*eV)
533 {
534 G4double differentialCrossSection = DifferentialCrossSection(particleDefinition, k, value, shell);
535 if(differentialCrossSection >= crossSectionMaximum) crossSectionMaximum = differentialCrossSection;
536 }
537
538
539 G4double secElecKinetic = 0.;
540
541 do
542 {
543 secElecKinetic = G4UniformRand() * maximumKineticEnergyTransfer;
544 } while(G4UniformRand()*crossSectionMaximum > DifferentialCrossSection(particleDefinition,
545 k,
546 secElecKinetic+waterStructure.IonisationEnergy(shell),
547 shell));
548
549 return(secElecKinetic);
550}
551
552//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
553
554
555void G4DNARuddIonisationModel::RandomizeEjectedElectronDirection(G4ParticleDefinition* particleDefinition,
556 G4double k,
557 G4double secKinetic,
558 G4double & cosTheta,
559 G4double & phi )
560{
561 G4DNAGenericIonsManager *instance;
562 instance = G4DNAGenericIonsManager::Instance();
563
564 G4double maxSecKinetic = 0.;
565
566 if (particleDefinition == G4Proton::ProtonDefinition()
567 || particleDefinition == instance->GetIon("hydrogen"))
568 {
569 maxSecKinetic = 4.* (electron_mass_c2 / proton_mass_c2) * k;
570 }
571
572 if (particleDefinition == instance->GetIon("helium")
573 || particleDefinition == instance->GetIon("alpha+")
574 || particleDefinition == instance->GetIon("alpha++"))
575 {
576 maxSecKinetic = 4.* (0.511 / 3728) * k;
577 }
578
579 phi = twopi * G4UniformRand();
580
581 //cosTheta = std::sqrt(secKinetic / maxSecKinetic);
582
583 // Restriction below 100 eV from Emfietzoglou (2000)
584
585 if (secKinetic>100*eV) cosTheta = std::sqrt(secKinetic / maxSecKinetic);
586 else cosTheta = (2.*G4UniformRand())-1.;
587
588}
589
590//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
591
592
593G4double G4DNARuddIonisationModel::DifferentialCrossSection(G4ParticleDefinition* particleDefinition,
594 G4double k,
595 G4double energyTransfer,
596 G4int ionizationLevelIndex)
597{
598 // Shells ids are 0 1 2 3 4 (4 is k shell)
599 // !!Attention, "energyTransfer" here is the energy transfered to the electron which means
600 // that the secondary kinetic energy is w = energyTransfer - bindingEnergy
601 //
602 // ds S F1(nu) + w * F2(nu)
603 // ---- = G(k) * ---- -------------------------------------------
604 // dw Bj (1+w)^3 * [1 + exp{alpha * (w - wc) / nu}]
605 //
606 // w is the secondary electron kinetic Energy in eV
607 //
608 // All the other parameters can be found in Rudd's Papers
609 //
610 // M.Eugene Rudd, 1988, User-Friendly model for the energy distribution of
611 // electrons from protons or electron collisions. Nucl. Tracks Rad. Meas.Vol 16 N0 2/3 pp 219-218
612 //
613
614 const G4int j=ionizationLevelIndex;
615
616 G4double A1 ;
617 G4double B1 ;
618 G4double C1 ;
619 G4double D1 ;
620 G4double E1 ;
621 G4double A2 ;
622 G4double B2 ;
623 G4double C2 ;
624 G4double D2 ;
625 G4double alphaConst ;
626
627// const G4double Bj[5] = {12.61*eV, 14.73*eV, 18.55*eV, 32.20*eV, 539.7*eV};
628// The following values are provided by M. dingfelder (priv. comm)
629 const G4double Bj[5] = {12.60*eV, 14.70*eV, 18.40*eV, 32.20*eV, 540*eV};
630
631 if (j == 4)
632 {
633 //Data For Liquid Water K SHELL from Dingfelder (Protons in Water)
634 A1 = 1.25;
635 B1 = 0.5;
636 C1 = 1.00;
637 D1 = 1.00;
638 E1 = 3.00;
639 A2 = 1.10;
640 B2 = 1.30;
641 C2 = 1.00;
642 D2 = 0.00;
643 alphaConst = 0.66;
644 }
645 else
646 {
647 //Data For Liquid Water from Dingfelder (Protons in Water)
648 A1 = 1.02;
649 B1 = 82.0;
650 C1 = 0.45;
651 D1 = -0.80;
652 E1 = 0.38;
653 A2 = 1.07;
654 // Value provided by M. Dingfelder (priv. comm)
655 B2 = 11.6;
656 //
657 C2 = 0.60;
658 D2 = 0.04;
659 alphaConst = 0.64;
660 }
661
662 const G4double n = 2.;
663 const G4double Gj[5] = {0.99, 1.11, 1.11, 0.52, 1.};
664
665 G4DNAGenericIonsManager* instance;
666 instance = G4DNAGenericIonsManager::Instance();
667
668 G4double wBig = (energyTransfer - waterStructure.IonisationEnergy(ionizationLevelIndex));
669 if (wBig<0) return 0.;
670
671 G4double w = wBig / Bj[ionizationLevelIndex];
672 // Note that the following (j==4) cases are provided by M. Dingfelder (priv. comm)
673 if (j==4) w = wBig / waterStructure.IonisationEnergy(ionizationLevelIndex);
674
675 G4double Ry = 13.6*eV;
676
677 G4double tau = 0.;
678
679 if (particleDefinition == G4Proton::ProtonDefinition()
680 || particleDefinition == instance->GetIon("hydrogen"))
681 {
682 tau = (electron_mass_c2/proton_mass_c2) * k ;
683 }
684
685 if ( particleDefinition == instance->GetIon("helium")
686 || particleDefinition == instance->GetIon("alpha+")
687 || particleDefinition == instance->GetIon("alpha++"))
688 {
689 tau = (0.511/3728.) * k ;
690 }
691
692 G4double S = 4.*pi * Bohr_radius*Bohr_radius * n * std::pow((Ry/Bj[ionizationLevelIndex]),2);
693 if (j==4) S = 4.*pi * Bohr_radius*Bohr_radius * n * std::pow((Ry/waterStructure.IonisationEnergy(ionizationLevelIndex)),2);
694
695 G4double v2 = tau / Bj[ionizationLevelIndex];
696 if (j==4) v2 = tau / waterStructure.IonisationEnergy(ionizationLevelIndex);
697
698 G4double v = std::sqrt(v2);
699 G4double wc = 4.*v2 - 2.*v - (Ry/(4.*Bj[ionizationLevelIndex]));
700 if (j==4) wc = 4.*v2 - 2.*v - (Ry/(4.*waterStructure.IonisationEnergy(ionizationLevelIndex)));
701
702 G4double L1 = (C1* std::pow(v,(D1))) / (1.+ E1*std::pow(v, (D1+4.)));
703 G4double L2 = C2*std::pow(v,(D2));
704 G4double H1 = (A1*std::log(1.+v2)) / (v2+(B1/v2));
705 G4double H2 = (A2/v2) + (B2/(v2*v2));
706
707 G4double F1 = L1+H1;
708 G4double F2 = (L2*H2)/(L2+H2);
709
710 G4double sigma = CorrectionFactor(particleDefinition, k)
711 * Gj[j] * (S/Bj[ionizationLevelIndex])
712 * ( (F1+w*F2) / ( std::pow((1.+w),3) * ( 1.+std::exp(alphaConst*(w-wc)/v))) );
713
714 if (j==4) sigma = CorrectionFactor(particleDefinition, k)
715 * Gj[j] * (S/waterStructure.IonisationEnergy(ionizationLevelIndex))
716 * ( (F1+w*F2) / ( std::pow((1.+w),3) * ( 1.+std::exp(alphaConst*(w-wc)/v))) );
717
718 if ( (particleDefinition == instance->GetIon("hydrogen")) && (ionizationLevelIndex==4))
719
720// sigma = Gj[j] * (S/Bj[ionizationLevelIndex])
721 sigma = Gj[j] * (S/waterStructure.IonisationEnergy(ionizationLevelIndex))
722 * ( (F1+w*F2) / ( std::pow((1.+w),3) * ( 1.+std::exp(alphaConst*(w-wc)/v))) );
723
724 if ( particleDefinition == G4Proton::ProtonDefinition()
725 || particleDefinition == instance->GetIon("hydrogen")
726 )
727 {
728 return(sigma);
729 }
730
731 if (particleDefinition == instance->GetIon("alpha++") )
732 {
733 slaterEffectiveCharge[0]=0.;
734 slaterEffectiveCharge[1]=0.;
735 slaterEffectiveCharge[2]=0.;
736 sCoefficient[0]=0.;
737 sCoefficient[1]=0.;
738 sCoefficient[2]=0.;
739 }
740
741 if (particleDefinition == instance->GetIon("alpha+") )
742 {
743 slaterEffectiveCharge[0]=2.0;
744 // The following values are provided by M. Dingfelder (priv. comm)
745 slaterEffectiveCharge[1]=2.0;
746 slaterEffectiveCharge[2]=2.0;
747 //
748 sCoefficient[0]=0.7;
749 sCoefficient[1]=0.15;
750 sCoefficient[2]=0.15;
751 }
752
753 if (particleDefinition == instance->GetIon("helium") )
754 {
755 slaterEffectiveCharge[0]=1.7;
756 slaterEffectiveCharge[1]=1.15;
757 slaterEffectiveCharge[2]=1.15;
758 sCoefficient[0]=0.5;
759 sCoefficient[1]=0.25;
760 sCoefficient[2]=0.25;
761 }
762
763 if ( particleDefinition == instance->GetIon("helium")
764 || particleDefinition == instance->GetIon("alpha+")
765 || particleDefinition == instance->GetIon("alpha++")
766 )
767 {
768 sigma = Gj[j] * (S/Bj[ionizationLevelIndex]) * ( (F1+w*F2) / ( std::pow((1.+w),3) * ( 1.+std::exp(alphaConst*(w-wc)/v))) );
769
770 if (j==4) sigma = Gj[j] * (S/waterStructure.IonisationEnergy(ionizationLevelIndex))
771 * ( (F1+w*F2) / ( std::pow((1.+w),3) * ( 1.+std::exp(alphaConst*(w-wc)/v))) );
772
773 G4double zEff = particleDefinition->GetPDGCharge() / eplus + particleDefinition->GetLeptonNumber();
774
775 zEff -= ( sCoefficient[0] * S_1s(k, energyTransfer, slaterEffectiveCharge[0], 1.) +
776 sCoefficient[1] * S_2s(k, energyTransfer, slaterEffectiveCharge[1], 2.) +
777 sCoefficient[2] * S_2p(k, energyTransfer, slaterEffectiveCharge[2], 2.) );
778
779 return zEff * zEff * sigma ;
780 }
781
782 return 0;
783}
784
785//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
786
787G4double G4DNARuddIonisationModel::S_1s(G4double t,
788 G4double energyTransferred,
789 G4double slaterEffectiveChg,
790 G4double shellNumber)
791{
792 // 1 - e^(-2r) * ( 1 + 2 r + 2 r^2)
793 // Dingfelder, in Chattanooga 2005 proceedings, formula (7)
794
795 G4double r = R(t, energyTransferred, slaterEffectiveChg, shellNumber);
796 G4double value = 1. - std::exp(-2 * r) * ( ( 2. * r + 2. ) * r + 1. );
797
798 return value;
799}
800
801//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
802
803G4double G4DNARuddIonisationModel::S_2s(G4double t,
804 G4double energyTransferred,
805 G4double slaterEffectiveChg,
806 G4double shellNumber)
807{
808 // 1 - e^(-2 r) * ( 1 + 2 r + 2 r^2 + 2 r^4)
809 // Dingfelder, in Chattanooga 2005 proceedings, formula (8)
810
811 G4double r = R(t, energyTransferred, slaterEffectiveChg, shellNumber);
812 G4double value = 1. - std::exp(-2 * r) * (((2. * r * r + 2.) * r + 2.) * r + 1.);
813
814 return value;
815
816}
817
818//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
819
820G4double G4DNARuddIonisationModel::S_2p(G4double t,
821 G4double energyTransferred,
822 G4double slaterEffectiveChg,
823 G4double shellNumber)
824{
825 // 1 - e^(-2 r) * ( 1 + 2 r + 2 r^2 + 4/3 r^3 + 2/3 r^4)
826 // Dingfelder, in Chattanooga 2005 proceedings, formula (9)
827
828 G4double r = R(t, energyTransferred, slaterEffectiveChg, shellNumber);
829 G4double value = 1. - std::exp(-2 * r) * (((( 2./3. * r + 4./3.) * r + 2.) * r + 2.) * r + 1.);
830
831 return value;
832}
833
834//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
835
836G4double G4DNARuddIonisationModel::R(G4double t,
837 G4double energyTransferred,
838 G4double slaterEffectiveChg,
839 G4double shellNumber)
840{
841 // tElectron = m_electron / m_alpha * t
842 // Dingfelder, in Chattanooga 2005 proceedings, p 4
843
844 G4double tElectron = 0.511/3728. * t;
845 // The following values are provided by M. Dingfelder (priv. comm)
846 G4double H = 2.*13.60569172 * eV;
847 G4double value = std::sqrt ( 2. * tElectron / H ) / ( energyTransferred / H ) * (slaterEffectiveChg/shellNumber);
848
849 return value;
850}
851
852//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
853
854G4double G4DNARuddIonisationModel::CorrectionFactor(G4ParticleDefinition* particleDefinition, G4double k)
855{
856 G4DNAGenericIonsManager *instance;
857 instance = G4DNAGenericIonsManager::Instance();
858
859 if (particleDefinition == G4Proton::Proton())
860 {
861 return(1.);
862 }
863 else
864 if (particleDefinition == instance->GetIon("hydrogen"))
865 {
866 G4double value = (std::log10(k/eV)-4.2)/0.5;
867 // The following values are provided by M. Dingfelder (priv. comm)
868 return((0.6/(1+std::exp(value))) + 0.9);
869 }
870 else
871 {
872 return(1.);
873 }
874}
875
876//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
877
878G4int G4DNARuddIonisationModel::RandomSelect(G4double k, const G4String& particle )
879{
880
881 // BEGIN PART 1/2 OF ELECTRON CORRECTION
882
883 // add ONE or TWO electron-water ionisation for alpha+ and helium
884
885 G4DNAGenericIonsManager *instance;
886 instance = G4DNAGenericIonsManager::Instance();
887
888 G4int level = 0;
889
890 // Retrieve data table corresponding to the current particle type
891
892 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
893 pos = tableData.find(particle);
894
895 if (pos != tableData.end())
896 {
897 G4DNACrossSectionDataSet* table = pos->second;
898
899 if (table != 0)
900 {
901 G4double* valuesBuffer = new G4double[table->NumberOfComponents()];
902
903 const size_t n(table->NumberOfComponents());
904 size_t i(n);
905 G4double value = 0.;
906
907 while (i>0)
908 {
909 i--;
910 valuesBuffer[i] = table->GetComponent(i)->FindValue(k);
911 value += valuesBuffer[i];
912 }
913
914 value *= G4UniformRand();
915
916 i = n;
917
918 while (i > 0)
919 {
920 i--;
921
922
923 if (valuesBuffer[i] > value)
924 {
925 delete[] valuesBuffer;
926 return i;
927 }
928 value -= valuesBuffer[i];
929 }
930
931 if (valuesBuffer) delete[] valuesBuffer;
932
933 }
934 }
935 else
936 {
937 G4Exception("G4DNARuddIonisationModel::RandomSelect: attempting to calculate cross section for wrong particle");
938 }
939
940 return level;
941}
942
943//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
944
945G4double G4DNARuddIonisationModel::PartialCrossSection(const G4Track& track )
946{
947 G4double sigma = 0.;
948
949 const G4DynamicParticle* particle = track.GetDynamicParticle();
950 G4double k = particle->GetKineticEnergy();
951
952 G4double lowLim = 0;
953 G4double highLim = 0;
954
955 const G4String& particleName = particle->GetDefinition()->GetParticleName();
956
957 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
958 pos1 = lowEnergyLimit.find(particleName);
959
960 if (pos1 != lowEnergyLimit.end())
961 {
962 lowLim = pos1->second;
963 }
964
965 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
966 pos2 = highEnergyLimit.find(particleName);
967
968 if (pos2 != highEnergyLimit.end())
969 {
970 highLim = pos2->second;
971 }
972
973 if (k >= lowLim && k <= highLim)
974 {
975 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
976 pos = tableData.find(particleName);
977
978 if (pos != tableData.end())
979 {
980 G4DNACrossSectionDataSet* table = pos->second;
981 if (table != 0)
982 {
983 sigma = table->FindValue(k);
984 }
985 }
986 else
987 {
988 G4Exception("G4DNARuddIonisationModel::PartialCrossSection: attempting to calculate cross section for wrong particle");
989 }
990 }
991
992 return sigma;
993}
994
995//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
996
997G4double G4DNARuddIonisationModel::Sum(G4double /* energy */, const G4String& /* particle */)
998{
999 return 0;
1000}
1001
Note: See TracBrowser for help on using the repository browser.