source: trunk/source/processes/hadronic/cross_sections/src/G4UInelasticCrossSection.cc@ 1350

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

geant4 tag 9.4

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