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

Last change on this file since 1058 was 1055, checked in by garnier, 17 years ago

maj sur la beta de geant 4.9.3

File size: 6.8 KB
RevLine 
[819]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//
[1055]26// $Id: G4CrossSectionExcitationMillerGreen.cc,v 1.5 2009/05/02 15:07:47 sincerti Exp $
27// GEANT4 tag $Name: geant4-09-03-beta-cand-01 $
[819]28
[961]29#include "G4CrossSectionExcitationMillerGreen.hh"
[819]30
[961]31//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
[819]32
33G4CrossSectionExcitationMillerGreen::G4CrossSectionExcitationMillerGreen()
34{
35 lowEnergyLimitDefault = 10 * eV;
36 highEnergyLimitDefault = 10 * MeV;
37
38 G4DNAGenericIonsManager *instance;
39 instance = G4DNAGenericIonsManager::Instance();
40 G4ParticleDefinition* protonDef = G4Proton::ProtonDefinition();
41 G4ParticleDefinition* alphaPlusPlusDef = instance->GetIon("alpha++");
42 G4ParticleDefinition* alphaPlusDef = instance->GetIon("alpha+");
43 G4ParticleDefinition* heliumDef = instance->GetIon("helium");
44
45 G4String proton;
46 G4String alphaPlusPlus;
47 G4String alphaPlus;
48 G4String helium;
49
50 if (protonDef != 0)
[961]51 {
52 proton = protonDef->GetParticleName();
53 lowEnergyLimit[proton] = 10. * eV;
54 highEnergyLimit[proton] = 500. * keV;
55 }
[819]56 else
[961]57 {
58 G4Exception("G4CrossSectionExcitationMillerGreen Constructor: proton is not defined");
59 }
[819]60
61 if (alphaPlusPlusDef != 0)
[961]62 {
63 alphaPlusPlus = alphaPlusPlusDef->GetParticleName();
64 lowEnergyLimit[alphaPlusPlus] = 1. * keV;
65 highEnergyLimit[alphaPlusPlus] = 10. * MeV;
66 }
[819]67 else
[961]68 {
[819]69 G4Exception("G4CrossSectionExcitationMillerGreen Constructor: alphaPlusPlus is not defined");
[961]70 }
[819]71
72 if (alphaPlusDef != 0)
[961]73 {
74 alphaPlus = alphaPlusDef->GetParticleName();
75 lowEnergyLimit[alphaPlus] = 1. * keV;
76 highEnergyLimit[alphaPlus] = 10. * MeV;
77 }
[819]78 else
[961]79 {
80 G4Exception("G4CrossSectionExcitationMillerGreen Constructor: alphaPlus is not defined");
81 }
[819]82
83 if (heliumDef != 0)
[961]84 {
85 helium = heliumDef->GetParticleName();
86 lowEnergyLimit[helium] = 1. * keV;
87 highEnergyLimit[helium] = 10. * MeV;
88 }
[819]89 else
[961]90 {
91 G4Exception("G4CrossSectionExcitationMillerGreen Constructor: helium is not defined");
92 }
93
[1055]94
95 G4cout << G4endl;
96 G4cout << "*******************************************************************************" << G4endl;
97 G4cout << "*******************************************************************************" << G4endl;
98 G4cout << " The class G4CrossSectionExcitationMillerGreen is NOT SUPPORTED ANYMORE. " << G4endl;
99 G4cout << " It will be REMOVED with the next major release of Geant4. " << G4endl;
100 G4cout << " Please consult: https://twiki.cern.ch/twiki/bin/view/Geant4/LoweProcesses" << G4endl;
101 G4cout << "*******************************************************************************" << G4endl;
102 G4cout << "*******************************************************************************" << G4endl;
103 G4cout << G4endl;
[819]104}
105
[961]106//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
[819]107
108G4CrossSectionExcitationMillerGreen::~G4CrossSectionExcitationMillerGreen()
109{}
110
[961]111//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
[819]112
113G4double G4CrossSectionExcitationMillerGreen::CrossSection(const G4Track& track)
114{
115 G4DNAGenericIonsManager *instance;
116 instance = G4DNAGenericIonsManager::Instance();
117
118 if (
119 track.GetDefinition() != G4Proton::ProtonDefinition()
120 &&
121 track.GetDefinition() != instance->GetIon("alpha++")
122 &&
123 track.GetDefinition() != instance->GetIon("alpha+")
124 &&
125 track.GetDefinition() != instance->GetIon("helium")
126 )
127
128 G4Exception("G4CrossSectionMillerGreen: attempting to calculate cross section for wrong particle");
129
130 G4double lowLim = lowEnergyLimitDefault;
131 G4double highLim = highEnergyLimitDefault;
132
133 const G4DynamicParticle* particle = track.GetDynamicParticle();
134 G4double k = particle->GetKineticEnergy();
135
136 const G4String& particleName = particle->GetDefinition()->GetParticleName();
137
138 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
139 pos1 = lowEnergyLimit.find(particleName);
140
141 if (pos1 != lowEnergyLimit.end())
[961]142 {
143 lowLim = pos1->second;
144 }
[819]145
146 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
147 pos2 = highEnergyLimit.find(particleName);
148
149 if (pos2 != highEnergyLimit.end())
[961]150 {
151 highLim = pos2->second;
152 }
[819]153
154 const G4ParticleDefinition* particleDefinition = track.GetDefinition();
155
156 G4double crossSection(0.);
157
158 if (k >= lowLim && k <= highLim)
[961]159 {
160 crossSection = partialCrossSection.Sum(k,particleDefinition);
161 G4DNAGenericIonsManager *instance;
162 instance = G4DNAGenericIonsManager::Instance();
[819]163
[961]164 // add ONE or TWO electron-water excitation for alpha+ and helium
[819]165
[961]166 if ( particleDefinition == instance->GetIon("alpha+")
167 ||
168 particleDefinition == instance->GetIon("helium")
169 )
170 {
[819]171 G4CrossSectionExcitationEmfietzoglouPartial * excitationXS =
172 new G4CrossSectionExcitationEmfietzoglouPartial();
[961]173
[819]174 G4double sigmaExcitation=0;
175 if (k*0.511/3728 > 7.4*eV && k*0.511/3728 < 10*keV) sigmaExcitation = excitationXS->Sum(k*0.511/3728);
176
177 if ( particleDefinition == instance->GetIon("alpha+") )
178 crossSection = crossSection + sigmaExcitation ;
[961]179
[819]180 if ( particleDefinition == instance->GetIon("helium") )
181 crossSection = crossSection + 2*sigmaExcitation ;
[961]182
[819]183 delete excitationXS;
[961]184 }
185 }
[819]186
187 return crossSection;
188}
189
Note: See TracBrowser for help on using the repository browser.