source: trunk/source/processes/electromagnetic/lowenergy/src/G4CrossSectionExcitationMillerGreen.cc@ 830

Last change on this file since 830 was 819, checked in by garnier, 17 years ago

import all except CVS

File size: 7.1 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//
27// $Id: G4CrossSectionExcitationMillerGreen.cc,v 1.3 2007/11/09 20:11:04 pia Exp $
28// GEANT4 tag $Name: $
29//
30// Contact Author: Maria Grazia Pia (Maria.Grazia.Pia@cern.ch)
31//
32// Reference: TNS Geant4-DNA paper
33// Reference for implementation model: NIM. 155, pp. 145-156, 1978
34
35// History:
36// -----------
37// Date Name Modification
38// 28 Apr 2007 M.G. Pia Created in compliance with design described in TNS paper
39//
40// -------------------------------------------------------------------
41
42// Class description:
43// Geant4-DNA Cross total cross section for electron elastic scattering in water
44// Reference: TNS Geant4-DNA paper
45// S. Chauvie et al., Geant4 physics processes for microdosimetry simulation:
46// design foundation and implementation of the first set of models,
47// IEEE Trans. Nucl. Sci., vol. 54, no. 6, Dec. 2007.
48// Further documentation available from http://www.ge.infn.it/geant4/dna
49
50// -------------------------------------------------------------------
51
52
53#include "G4CrossSectionExcitationMillerGreen.hh"
54#include "G4Track.hh"
55#include "G4DynamicParticle.hh"
56#include "G4Proton.hh"
57#include "G4CrossSectionExcitationEmfietzoglouPartial.hh"
58#include "G4DNAGenericIonsManager.hh"
59
60G4CrossSectionExcitationMillerGreen::G4CrossSectionExcitationMillerGreen()
61{
62 // Default energy limits (defined for protection against anomalous behaviour only)
63 name = "ExcitationMillerGreen";
64 lowEnergyLimitDefault = 10 * eV;
65 highEnergyLimitDefault = 10 * MeV;
66
67 G4DNAGenericIonsManager *instance;
68 instance = G4DNAGenericIonsManager::Instance();
69 G4ParticleDefinition* protonDef = G4Proton::ProtonDefinition();
70 G4ParticleDefinition* alphaPlusPlusDef = instance->GetIon("alpha++");
71 G4ParticleDefinition* alphaPlusDef = instance->GetIon("alpha+");
72 G4ParticleDefinition* heliumDef = instance->GetIon("helium");
73
74 G4String proton;
75 G4String alphaPlusPlus;
76 G4String alphaPlus;
77 G4String helium;
78
79 if (protonDef != 0)
80 {
81 proton = protonDef->GetParticleName();
82 lowEnergyLimit[proton] = 10. * eV;
83 highEnergyLimit[proton] = 500. * keV;
84 }
85 else
86 {
87 G4Exception("G4CrossSectionExcitationMillerGreen Constructor: proton is not defined");
88 }
89
90 if (alphaPlusPlusDef != 0)
91 {
92 alphaPlusPlus = alphaPlusPlusDef->GetParticleName();
93 lowEnergyLimit[alphaPlusPlus] = 1. * keV;
94 highEnergyLimit[alphaPlusPlus] = 10. * MeV;
95 }
96 else
97 {
98 G4Exception("G4CrossSectionExcitationMillerGreen Constructor: alphaPlusPlus is not defined");
99 }
100
101 if (alphaPlusDef != 0)
102 {
103 alphaPlus = alphaPlusDef->GetParticleName();
104 lowEnergyLimit[alphaPlus] = 1. * keV;
105 highEnergyLimit[alphaPlus] = 10. * MeV;
106 }
107 else
108 {
109 G4Exception("G4CrossSectionExcitationMillerGreen Constructor: alphaPlus is not defined");
110 }
111
112 if (heliumDef != 0)
113 {
114 helium = heliumDef->GetParticleName();
115 lowEnergyLimit[helium] = 1. * keV;
116 highEnergyLimit[helium] = 10. * MeV;
117 }
118 else
119 {
120 G4Exception("G4CrossSectionExcitationMillerGreen Constructor: helium is not defined");
121 }
122
123
124}
125
126
127G4CrossSectionExcitationMillerGreen::~G4CrossSectionExcitationMillerGreen()
128{}
129
130
131G4double G4CrossSectionExcitationMillerGreen::CrossSection(const G4Track& track)
132{
133 G4DNAGenericIonsManager *instance;
134 instance = G4DNAGenericIonsManager::Instance();
135
136 if (
137 track.GetDefinition() != G4Proton::ProtonDefinition()
138 &&
139 track.GetDefinition() != instance->GetIon("alpha++")
140 &&
141 track.GetDefinition() != instance->GetIon("alpha+")
142 &&
143 track.GetDefinition() != instance->GetIon("helium")
144 )
145
146 G4Exception("G4CrossSectionMillerGreen: attempting to calculate cross section for wrong particle");
147
148 G4double lowLim = lowEnergyLimitDefault;
149 G4double highLim = highEnergyLimitDefault;
150
151 const G4DynamicParticle* particle = track.GetDynamicParticle();
152 G4double k = particle->GetKineticEnergy();
153
154 const G4String& particleName = particle->GetDefinition()->GetParticleName();
155
156 // Retrieve energy limits for the current particle type
157
158 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
159 pos1 = lowEnergyLimit.find(particleName);
160
161 // Lower limit
162 if (pos1 != lowEnergyLimit.end())
163 {
164 lowLim = pos1->second;
165 }
166
167 // Upper limit
168 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
169 pos2 = highEnergyLimit.find(particleName);
170
171 if (pos2 != highEnergyLimit.end())
172 {
173 highLim = pos2->second;
174 }
175
176 //
177 const G4ParticleDefinition* particleDefinition = track.GetDefinition();
178
179 G4double crossSection(0.);
180
181 if (k >= lowLim && k <= highLim)
182 {
183
184 crossSection = partialCrossSection.Sum(k,particleDefinition);
185 G4DNAGenericIonsManager *instance;
186 instance = G4DNAGenericIonsManager::Instance();
187
188 // add ONE or TWO electron-water excitation for alpha+ and helium
189
190 if ( particleDefinition == instance->GetIon("alpha+")
191 ||
192 particleDefinition == instance->GetIon("helium")
193 )
194 {
195 G4CrossSectionExcitationEmfietzoglouPartial * excitationXS =
196 new G4CrossSectionExcitationEmfietzoglouPartial();
197 G4double sigmaExcitation=0;
198 if (k*0.511/3728 > 7.4*eV && k*0.511/3728 < 10*keV) sigmaExcitation = excitationXS->Sum(k*0.511/3728);
199
200 if ( particleDefinition == instance->GetIon("alpha+") )
201 crossSection = crossSection + sigmaExcitation ;
202 if ( particleDefinition == instance->GetIon("helium") )
203 crossSection = crossSection + 2*sigmaExcitation ;
204 delete excitationXS;
205 }
206 }
207
208 return crossSection;
209}
210
Note: See TracBrowser for help on using the repository browser.