source: trunk/source/processes/hadronic/cross_sections/src/G4UElasticCrossSection.cc@ 1340

Last change on this file since 1340 was 1340, checked in by garnier, 15 years ago

update ti head

File size: 7.2 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// 12.08.06 V.Ivanchenko - first implementation
28// 22.01.07 V.Ivanchenko - add GetIsoZACrossSection
29// 05.03.07 V.Ivanchenko - use G4NucleonNuclearCrossSection
30// 06.03.07 V.Ivanchenko - add Initialise function
31//
32//
33
34
35#include "G4UElasticCrossSection.hh"
36
37#include "G4ParticleTable.hh"
38#include "G4GlauberGribovCrossSection.hh"
39#include "G4NucleonNuclearCrossSection.hh"
40#include "G4UPiNuclearCrossSection.hh"
41#include "G4HadronCrossSections.hh"
42#include "G4ParticleDefinition.hh"
43#include "G4Element.hh"
44#include "G4Proton.hh"
45#include "G4Neutron.hh"
46#include "G4PionPlus.hh"
47#include "G4PionMinus.hh"
48#include "G4NistManager.hh"
49
50//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
51
52G4UElasticCrossSection::G4UElasticCrossSection(const G4ParticleDefinition*)
53{
54 verboseLevel = 0;
55 hasGlauber = false;
56 thEnergy = 90.*GeV;
57 fGlauber = new G4GlauberGribovCrossSection();
58 fGheisha = G4HadronCrossSections::Instance();
59 fNucleon = 0;
60 fUPi = 0;
61}
62
63//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
64
65G4UElasticCrossSection::~G4UElasticCrossSection()
66{
67 delete fGlauber;
68 if(fNucleon) delete fNucleon;
69 if(fUPi) delete fUPi;
70}
71
72//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
73
74G4bool G4UElasticCrossSection::IsApplicable(const G4DynamicParticle* dp,
75 const G4Element* elm)
76{
77 return IsZAApplicable(dp, elm->GetZ(), elm->GetN());
78}
79
80//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
81
82// G4bool G4UElasticCrossSection::IsZAApplicable(const G4DynamicParticle* dp,
83// G4double Z, G4double A)
84G4bool G4UElasticCrossSection::IsIsoApplicable(const G4DynamicParticle* dp,
85 G4int Z, G4int A)
86{
87 G4bool res = false;
88 if(fNucleon || fUPi) res = true;
89 else res = fGheisha->IsApplicable(dp, Z, A);
90 if(verboseLevel > 1)
91 G4cout << "G4UElasticCrossSection::IsApplicable for "
92 << dp->GetDefinition()->GetParticleName()
93 << " Ekin(GeV)= " << dp->GetKineticEnergy()
94 << G4endl;
95 return res;
96}
97
98
99G4double
100G4UElasticCrossSection::GetCrossSection(const G4DynamicParticle* dp,
101 const G4Element* elm, G4double temp)
102{
103 G4int Z = G4lrint(elm->GetZ());
104 G4int N = G4lrint(elm->GetN());
105
106 return GetZandACrossSection(dp, Z, N, temp);
107}
108
109
110G4double
111G4UElasticCrossSection::GetZandACrossSection(const G4DynamicParticle* dp,
112 G4int Z, G4int A, G4double)
113{
114 G4double cross = 0.0;
115 G4double ekin = dp->GetKineticEnergy();
116// G4int iz = G4int(Z);
117 if(Z > 92) Z = 92;
118
119 // proton and neutron
120 if(fNucleon) {
121 if(Z == 1) cross = fGheisha->GetElasticCrossSection(dp, Z, A);
122 else if(ekin > thEnergy) {
123 cross = theFac[Z]*fGlauber->GetElasticGlauberGribov(dp, Z, A);
124 } else {
125 cross = fNucleon->GetElasticCrossSection(dp, Z, A);
126 }
127
128 // pions
129 } else if(fUPi) {
130 if(Z == 1) cross = fGheisha->GetElasticCrossSection(dp, Z, A);
131 else if(ekin > thEnergy) {
132 cross = theFac[Z]*fGlauber->GetElasticGlauberGribov(dp, Z, A);
133 } else {
134 cross = fUPi->GetElasticCrossSection(dp, Z, A);
135 }
136
137 //others
138 } else {
139 if(hasGlauber && ekin > thEnergy) {
140 cross = theFac[Z]*fGlauber->GetElasticGlauberGribov(dp, Z, A);
141 } else if(fGheisha->IsApplicable(dp, Z, A)) {
142 cross = fGheisha->GetElasticCrossSection(dp, Z, A);
143 }
144 }
145
146 if(verboseLevel > 1)
147 G4cout << "G4UElasticCrossSection::GetCrossSection for "
148 << dp->GetDefinition()->GetParticleName()
149 << " Ekin(GeV)= " << dp->GetKineticEnergy()
150 << " in nucleus Z= " << Z << " A= " << A
151 << " XS(b)= " << cross/barn
152 << G4endl;
153
154 return cross;
155}
156
157//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
158
159void G4UElasticCrossSection::BuildPhysicsTable(const G4ParticleDefinition& p)
160{
161 if(&p == G4Proton::Proton() || &p == G4Neutron::Neutron()) {
162 fNucleon = new G4NucleonNuclearCrossSection();
163 } else if(&p == G4PionPlus::PionPlus() || &p == G4PionMinus::PionMinus()) {
164 fUPi = new G4UPiNuclearCrossSection();
165 }
166 Initialise(&p);
167}
168
169//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
170
171void G4UElasticCrossSection::DumpPhysicsTable(const G4ParticleDefinition&)
172{
173 G4cout << "G4UElasticCrossSection:"<<G4endl;
174}
175
176
177//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
178
179void G4UElasticCrossSection::Initialise(const G4ParticleDefinition* p)
180{
181 G4ParticleDefinition* part = const_cast<G4ParticleDefinition*>(p);
182 G4ThreeVector mom(0.0,0.0,1.0);
183 G4DynamicParticle dp(part, mom, thEnergy);
184
185 G4NistManager* nist = G4NistManager::Instance();
186 G4int A = G4lrint(nist->GetAtomicMassAmu(2));
187
188 if(fGlauber->IsZAApplicable(&dp, 2.0, A)) {
189 hasGlauber = true;
190
191 if(verboseLevel > 0) G4cout << "### G4UElasticCrossSection::Initialise for "
192 << p->GetParticleName() << G4endl;
193 G4double csup, csdn;
194 for(G4int iz=2; iz<93; iz++) {
195
196 G4double Z = G4double(iz);
197 A = G4lrint(nist->GetAtomicMassAmu(iz));
198 csup = fGlauber->GetElasticGlauberGribov(&dp, iz, A);
199 // proton and neutron
200 if(fNucleon) {
201 csdn = fNucleon->GetElasticCrossSection(&dp, iz, A);
202
203 // pions
204 } else if(fUPi) {
205 csdn = fUPi->GetElasticCrossSection(&dp, iz, A);
206
207 // other
208 } else if(fGheisha->IsApplicable(&dp, iz, A)) {
209 csdn = fGheisha->GetElasticCrossSection(&dp, iz, A);
210
211 } else {
212 csdn = csup;
213 }
214 theFac[iz] = csdn/csup;
215 if(verboseLevel > 0) G4cout << "Z= " << Z << " A= " << A
216 << " factor= " << theFac[iz] << G4endl;
217 }
218 }
219}
220
221
Note: See TracBrowser for help on using the repository browser.