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

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

import all except CVS

File size: 6.9 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: G4CrossSectionIonisationBorn.cc,v 1.2 2007/11/08 18:51:34 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 "G4CrossSectionIonisationBorn.hh"
54#include "G4ParticleDefinition.hh"
55#include "G4Electron.hh"
56#include "G4Proton.hh"
57#include "G4Track.hh"
58#include "G4LogLogInterpolation.hh"
59#include "G4SystemOfUnits.hh"
60
61
62G4CrossSectionIonisationBorn::G4CrossSectionIonisationBorn()
63{
64
65 name = "IonisationBorn";
66
67 // Default energy limits (defined for protection against anomalous behaviour only)
68 lowEnergyLimitDefault = 25 * eV;
69 highEnergyLimitDefault = 30 * keV;
70
71 G4String fileElectron("dna/sigma_ionisation_e_born");
72 G4String fileProton("dna/sigma_ionisation_p_born");
73
74 G4ParticleDefinition* electronDef = G4Electron::ElectronDefinition();
75 G4ParticleDefinition* protonDef = G4Proton::ProtonDefinition();
76
77 G4String electron;
78 G4String proton;
79
80 // Factor to scale microscopic/macroscopic cross section data in water
81 // ---- MGP ---- Hardcoded (taken from prototype code); to be replaced with proper calculation
82 G4double scaleFactor = (1.e-22 / 3.343) * m*m;
83
84
85 // Data members for electrons
86
87 if (electronDef != 0)
88 {
89 electron = electronDef->GetParticleName();
90 tableFile[electron] = fileElectron;
91
92 // Energy limits
93 lowEnergyLimit[electron] = 25. * eV;
94 highEnergyLimit[electron] = 30. * keV;
95
96 // Create data set with electron cross section data and load values stored in file
97 G4DNACrossSectionDataSet* tableE = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV,scaleFactor );
98 tableE->LoadData(fileElectron);
99
100 // Insert key-table pair in map
101 tableData[electron] = tableE;
102 }
103 else
104 {
105 G4Exception("G4CrossSectionIonisationBorn Constructor: electron is not defined");
106 }
107
108 // Data members for protons
109
110 if (protonDef != 0)
111 {
112 proton = protonDef->GetParticleName();
113 tableFile[proton] = fileProton;
114
115 // Energy limits
116 lowEnergyLimit[proton] = 500. * keV;
117 highEnergyLimit[proton] = 10. * MeV;
118
119 // Create data set with proton cross section data and load values stored in file
120 G4DNACrossSectionDataSet* tableP = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV,scaleFactor );
121 tableP->LoadData(fileProton);
122
123 // Insert key-table pair in map
124 tableData[proton] = tableP;
125 }
126 else
127 {
128 G4Exception("G4CrossSectionIonisationBorn Constructor: proton is not defined");
129 }
130}
131
132
133G4CrossSectionIonisationBorn::~G4CrossSectionIonisationBorn()
134{
135 // Destroy the content of the map
136 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
137 for (pos = tableData.begin(); pos != tableData.end(); ++pos)
138 {
139 G4DNACrossSectionDataSet* table = pos->second;
140 delete table;
141 }
142}
143
144
145
146G4double G4CrossSectionIonisationBorn::CrossSection(const G4Track& track )
147{
148 const G4DynamicParticle* particle = track.GetDynamicParticle();
149 G4double k = particle->GetKineticEnergy();
150
151 // Cross section = 0 outside the energy validity limits set in the constructor
152 G4double sigma = 0.;
153
154 // ---- MGP ---- Better handling of these limits to be set in a following design iteration
155 G4double lowLim = lowEnergyLimitDefault;
156 G4double highLim = highEnergyLimitDefault;
157
158 const G4String& particleName = particle->GetDefinition()->GetParticleName();
159
160 // Retrieve energy limits for the current particle type
161
162 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
163 pos1 = lowEnergyLimit.find(particleName);
164
165 // Lower limit
166 if (pos1 != lowEnergyLimit.end())
167 {
168 lowLim = pos1->second;
169 }
170
171 // Upper limit
172 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
173 pos2 = highEnergyLimit.find(particleName);
174
175 if (pos2 != highEnergyLimit.end())
176 {
177 highLim = pos2->second;
178 }
179
180 // Verify that the current track is within the energy limits of validity of the cross section model
181 if (k > lowLim && k < highLim)
182 {
183 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
184 pos = tableData.find(particleName);
185
186 if (pos != tableData.end())
187 {
188 G4DNACrossSectionDataSet* table = pos->second;
189 if (table != 0)
190 {
191 // Cross section
192 sigma = table->FindValue(k);
193 }
194 }
195 else
196 {
197 // The track does not corresponds to a particle pertinent to this model
198 G4Exception("G4CrossSectionIonisationBorn: attempting to calculate cross section for wrong particle");
199 }
200 }
201
202 return sigma;
203}
204
Note: See TracBrowser for help on using the repository browser.