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

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

import all except CVS

  • Property svn:executable set to *
File size: 5.8 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: G4CrossSectionChargeDecrease.cc,v 1.2 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 "G4CrossSectionChargeDecrease.hh"
54#include "G4Track.hh"
55#include "G4DynamicParticle.hh"
56#include "G4Proton.hh"
57#include "G4CrossSectionExcitationEmfietzoglouPartial.hh"
58#include "G4DNAGenericIonsManager.hh"
59
60G4CrossSectionChargeDecrease::G4CrossSectionChargeDecrease()
61{
62 // Default energy limits (defined for protection against anomalous behaviour only)
63 name = "ChargeDecrease";
64 lowEnergyLimitDefault = 1 * keV;
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
73 G4String proton;
74 G4String alphaPlusPlus;
75 G4String alphaPlus;
76
77 if (protonDef != 0)
78 {
79 proton = protonDef->GetParticleName();
80 lowEnergyLimit[proton] = 1. * keV;
81 highEnergyLimit[proton] = 10. * keV;
82 }
83 else
84 {
85 G4Exception("G4CrossSectionChargeDecrease Constructor: proton is not defined");
86 }
87
88 if (alphaPlusPlusDef != 0)
89 {
90 alphaPlusPlus = alphaPlusPlusDef->GetParticleName();
91 lowEnergyLimit[alphaPlusPlus] = 1. * keV;
92 highEnergyLimit[alphaPlusPlus] = 10. * MeV;
93 }
94 else
95 {
96 G4Exception("G4CrossSectionChargeDecrease Constructor: alphaPlusPlus is not defined");
97 }
98
99 if (alphaPlusDef != 0)
100 {
101 alphaPlus = alphaPlusDef->GetParticleName();
102 lowEnergyLimit[alphaPlus] = 1. * keV;
103 highEnergyLimit[alphaPlus] = 10. * MeV;
104 }
105 else
106 {
107 G4Exception("G4CrossSectionChargeDecrease Constructor: alphaPlus is not defined");
108 }
109
110}
111
112
113G4CrossSectionChargeDecrease::~G4CrossSectionChargeDecrease()
114{}
115
116
117G4double G4CrossSectionChargeDecrease::CrossSection(const G4Track& track)
118{
119 G4double lowLim = lowEnergyLimitDefault;
120 G4double highLim = highEnergyLimitDefault;
121
122 const G4DynamicParticle* particle = track.GetDynamicParticle();
123 G4double k = particle->GetKineticEnergy();
124
125 const G4ParticleDefinition* particleDefinition = track.GetDefinition();
126
127 const G4String& particleName = particleDefinition->GetParticleName();
128
129 G4DNAGenericIonsManager *instance;
130 instance = G4DNAGenericIonsManager::Instance();
131
132 if (
133 particleDefinition != G4Proton::ProtonDefinition()
134 &&
135 particleDefinition != instance->GetIon("alpha++")
136 &&
137 particleDefinition != instance->GetIon("alpha+")
138 )
139
140 G4Exception("G4CrossSectionChargeDecrease: attempting to calculate cross section for wrong particle");
141
142 // Retrieve energy limits for the current particle type
143
144 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
145 pos1 = lowEnergyLimit.find(particleName);
146
147 // Lower limit
148 if (pos1 != lowEnergyLimit.end())
149 {
150 lowLim = pos1->second;
151 }
152
153 // Upper limit
154 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
155 pos2 = highEnergyLimit.find(particleName);
156
157 if (pos2 != highEnergyLimit.end())
158 {
159 highLim = pos2->second;
160 }
161
162 //
163
164 G4double crossSection(0.);
165 if (k >= lowLim && k <= highLim)
166 {
167 crossSection = partialCrossSection.Sum(k,particleDefinition);
168 }
169
170 return crossSection;
171}
172
Note: See TracBrowser for help on using the repository browser.