source: trunk/source/processes/electromagnetic/lowenergy/src/G4CrossSectionChargeIncrease.cc@ 1005

Last change on this file since 1005 was 991, checked in by garnier, 17 years ago

update

  • Property svn:executable set to *
File size: 5.3 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: G4CrossSectionChargeIncrease.cc,v 1.4 2008/07/14 20:47:34 sincerti Exp $
27// GEANT4 tag $Name: geant4-09-02 $
28
29#include "G4CrossSectionChargeIncrease.hh"
30
31//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32
33G4CrossSectionChargeIncrease::G4CrossSectionChargeIncrease()
34{
35 lowEnergyLimitDefault = 1 * keV;
36 highEnergyLimitDefault = 10 * MeV;
37
38 G4DNAGenericIonsManager *instance;
39 instance = G4DNAGenericIonsManager::Instance();
40 G4ParticleDefinition* hydrogenDef = instance->GetIon("hydrogen");
41 G4ParticleDefinition* alphaPlusDef = instance->GetIon("alpha+");
42 G4ParticleDefinition* heliumDef = instance->GetIon("helium");
43
44 G4String hydrogen;
45 G4String alphaPlus;
46 G4String helium;
47
48 if (hydrogenDef != 0)
49 {
50 hydrogen = hydrogenDef->GetParticleName();
51 lowEnergyLimit[hydrogen] = 1. * keV;
52 highEnergyLimit[hydrogen] = 10. * MeV;
53 }
54 else
55 {
56 G4Exception("G4CrossSectionChargeIncrease Constructor: hydrogen is not defined");
57 }
58
59 if (alphaPlusDef != 0)
60 {
61 alphaPlus = alphaPlusDef->GetParticleName();
62 lowEnergyLimit[alphaPlus] = 1. * keV;
63 highEnergyLimit[alphaPlus] = 10. * MeV;
64 }
65 else
66 {
67 G4Exception("G4CrossSectionChargeIncrease Constructor: alphaPlus is not defined");
68 }
69
70 if (heliumDef != 0)
71 {
72 helium = heliumDef->GetParticleName();
73 lowEnergyLimit[helium] = 1. * keV;
74 highEnergyLimit[helium] = 10. * MeV;
75 }
76 else
77 {
78 G4Exception("G4CrossSectionChargeIncrease Constructor: helium is not defined");
79 }
80
81}
82
83//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
84
85G4CrossSectionChargeIncrease::~G4CrossSectionChargeIncrease()
86{}
87
88//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
89
90G4double G4CrossSectionChargeIncrease::CrossSection(const G4Track& track)
91{
92 G4double lowLim = lowEnergyLimitDefault;
93 G4double highLim = highEnergyLimitDefault;
94
95 const G4DynamicParticle* particle = track.GetDynamicParticle();
96 G4double k = particle->GetKineticEnergy();
97
98 const G4String& particleName = particle->GetDefinition()->GetParticleName();
99
100 G4DNAGenericIonsManager *instance;
101 instance = G4DNAGenericIonsManager::Instance();
102
103 const G4ParticleDefinition* particleDefinition = track.GetDefinition();
104
105 if (
106 particleDefinition != instance->GetIon("hydrogen")
107 &&
108 particleDefinition != instance->GetIon("alpha+")
109 &&
110 particleDefinition != instance->GetIon("helium")
111 )
112
113 G4Exception("G4CrossSectionChargeIncrease: attempting to calculate cross section for wrong particle");
114
115 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
116 pos1 = lowEnergyLimit.find(particleName);
117
118 if (pos1 != lowEnergyLimit.end())
119 {
120 lowLim = pos1->second;
121 }
122
123 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
124 pos2 = highEnergyLimit.find(particleName);
125
126 if (pos2 != highEnergyLimit.end())
127 {
128 highLim = pos2->second;
129 }
130
131 G4double totalCrossSection = 0.;
132
133 if (k >= lowLim && k <= highLim)
134 {
135 //HYDROGEN
136 if (particleDefinition == instance->GetIon("hydrogen"))
137 {
138 const G4double aa = 2.835;
139 const G4double bb = 0.310;
140 const G4double cc = 2.100;
141 const G4double dd = 0.760;
142 const G4double fac = 1.0e-18;
143 const G4double rr = 13.606 * eV;
144
145 G4double t = k / (proton_mass_c2/electron_mass_c2);
146 G4double x = t / rr;
147 G4double temp = 4.0 * pi * Bohr_radius/nm * Bohr_radius/nm * fac;
148 G4double sigmal = temp * cc * (std::pow(x,dd));
149 G4double sigmah = temp * (aa * std::log(1.0 + x) + bb) / x;
150 totalCrossSection = 1.0/(1.0/sigmal + 1.0/sigmah) *m*m;
151 }
152 else
153 {
154 totalCrossSection = partialCrossSection.Sum(k,particleDefinition);
155 }
156 }
157
158 return totalCrossSection;
159}
160
Note: See TracBrowser for help on using the repository browser.