source: trunk/source/processes/hadronic/cross_sections/src/G4BGGPionInelasticXS.cc@ 1197

Last change on this file since 1197 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

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// $Id: G4BGGPionInelasticXS.cc,v 1.6 2009/11/19 12:24:41 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03-cand-01 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name: G4BGGPionInelasticXS
35//
36// Author: Vladimir Ivanchenko
37//
38// Creation date: 01.10.2003
39// Modifications:
40//
41//
42// -------------------------------------------------------------------
43//
44
45#include "G4BGGPionInelasticXS.hh"
46#include "G4GlauberGribovCrossSection.hh"
47#include "G4UPiNuclearCrossSection.hh"
48#include "G4HadronNucleonXsc.hh"
49#include "G4PionPlus.hh"
50#include "G4PionMinus.hh"
51#include "G4NistManager.hh"
52
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
54
55G4BGGPionInelasticXS::G4BGGPionInelasticXS(const G4ParticleDefinition* p)
56{
57 verboseLevel = 0;
58 fGlauberEnergy = 91.*GeV;
59 fLowEnergy = 20.*MeV;
60 fPion = 0;
61 fGlauber = 0;
62 fHadron = 0;
63 particle = p;
64 isPiplus = false;
65 isInitialized = false;
66}
67
68//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
69
70G4BGGPionInelasticXS::~G4BGGPionInelasticXS()
71{
72 delete fGlauber;
73 delete fPion;
74 delete fHadron;
75}
76
77//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
78
79G4double G4BGGPionInelasticXS::GetIsoZACrossSection(const G4DynamicParticle* dp,
80 G4double Z,
81 G4double A,
82 G4double)
83{
84 G4double cross = 0.0;
85 G4double ekin = dp->GetKineticEnergy();
86 G4int iz = G4int(Z);
87 if(iz > 92) iz = 92;
88
89 if(ekin <= fLowEnergy) {
90 cross = theCoulombFac[iz];
91 if(isPiplus) { cross *= CoulombFactor(ekin, A); }
92
93 } else if(iz == 1) {
94
95 if( A < 1.5) {
96 fHadron->GetHadronNucleonXscMK(dp, G4Proton::Proton());
97 cross = fHadron->GetInelasticHadronNucleonXsc();
98 } else {
99 fHadron->GetHadronNucleonXscMK(dp, G4Proton::Proton());
100 cross = fHadron->GetInelasticHadronNucleonXsc();
101 fHadron->GetHadronNucleonXscMK(dp, G4Neutron::Neutron());
102 cross += fHadron->GetInelasticHadronNucleonXsc();
103 }
104
105 } else if(ekin > fGlauberEnergy) {
106 cross = theGlauberFac[iz]*fGlauber->GetInelasticGlauberGribov(dp, Z, A);
107
108 } else {
109 cross = fPion->GetIsoZACrossSection(dp, Z, A);
110 }
111
112 if(verboseLevel > 1)
113 G4cout << "G4BGGPionInelasticXS::GetCrossSection for "
114 << dp->GetDefinition()->GetParticleName()
115 << " Ekin(GeV)= " << dp->GetKineticEnergy()
116 << " in nucleus Z= " << Z << " A= " << A
117 << " XS(b)= " << cross/barn
118 << G4endl;
119
120 return cross;
121}
122
123//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
124
125void G4BGGPionInelasticXS::BuildPhysicsTable(const G4ParticleDefinition& p)
126{
127 if(&p == G4PionPlus::PionPlus() || &p == G4PionMinus::PionMinus()) {
128 particle = &p;
129 Initialise();
130 } else {
131 G4cout << "### G4BGGPionInelasticXS WARNING: is not applicable to "
132 << p.GetParticleName()
133 << G4endl;
134 }
135}
136
137//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
138
139void G4BGGPionInelasticXS::DumpPhysicsTable(const G4ParticleDefinition&)
140{
141 G4cout << "G4BGGPionInelasticXS:"<<G4endl;
142}
143
144//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
145
146void G4BGGPionInelasticXS::Initialise()
147{
148 if(isInitialized) return;
149 isInitialized = true;
150
151 fPion = new G4UPiNuclearCrossSection();
152 fGlauber = new G4GlauberGribovCrossSection();
153 fHadron = new G4HadronNucleonXsc();
154 fPion->BuildPhysicsTable(*particle);
155 fGlauber->BuildPhysicsTable(*particle);
156 if(particle == G4PionPlus::PionPlus()) isPiplus = true;
157
158 G4ParticleDefinition* part = const_cast<G4ParticleDefinition*>(particle);
159 G4ThreeVector mom(0.0,0.0,1.0);
160 G4DynamicParticle dp(part, mom, fGlauberEnergy);
161
162 G4NistManager* nist = G4NistManager::Instance();
163
164 G4double A, csup, csdn;
165
166 if(verboseLevel > 0) G4cout << "### G4BGGPionInelasticXS::Initialise for "
167 << particle->GetParticleName()
168 << " isPiplus: " << isPiplus
169 << G4endl;
170
171 for(G4int iz=2; iz<93; iz++) {
172
173 G4double Z = G4double(iz);
174 A = nist->GetAtomicMassAmu(iz);
175
176 csup = fGlauber->GetInelasticGlauberGribov(&dp, Z, A);
177 csdn = fPion->GetInelasticCrossSection(&dp, Z, A);
178
179 theGlauberFac[iz] = csdn/csup;
180 if(verboseLevel > 0) G4cout << "Z= " << Z << " A= " << A
181 << " factor= " << theGlauberFac[iz] << G4endl;
182 }
183 dp.SetKineticEnergy(fLowEnergy);
184 fHadron->GetHadronNucleonXscNS(&dp, G4Proton::Proton());
185 theCoulombFac[1] = fHadron->GetInelasticHadronNucleonXsc();
186 if(isPiplus) { theCoulombFac[1] /= CoulombFactor(fLowEnergy,1.0); }
187
188 for(G4int iz=2; iz<93; iz++) {
189
190 G4double Z = G4double(iz);
191 A = nist->GetAtomicMassAmu(iz);
192
193 theCoulombFac[iz] = fPion->GetIsoZACrossSection(&dp, Z, A);
194 if(isPiplus) { theCoulombFac[iz] /= CoulombFactor(fLowEnergy,A); }
195
196 if(verboseLevel > 0) G4cout << "Z= " << Z << " A= " << A
197 << " factor= " << theCoulombFac[iz] << G4endl;
198 }
199}
200
201//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
202
203G4double G4BGGPionInelasticXS::CoulombFactor(G4double kinEnergy, G4double A)
204{
205 G4double res= 0.0;
206 if(kinEnergy <= DBL_MIN) return res;
207 else if(A < 1.5) return kinEnergy*kinEnergy;
208
209 G4double elog = std::log10(kinEnergy/GeV);
210
211 // from G4ProtonInelasticCrossSection
212 G4double f1 = 8.0 - 8.0/A - 0.008*A;
213 G4double f2 = 2.34 - 5.4/A - 0.0028*A;
214
215 res = 1.0/(1.0 + std::exp(-f1*(elog + f2)));
216
217 f1 = 5.6 - 0.016*A;
218 f2 = 1.37 + 1.37/A;
219 res *= ( 1.0 + (0.8 + 18./A - 0.002*A)/(1.0 + std::exp(f1*(elog + f2))));
220 return res;
221}
222
223//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
224
225
Note: See TracBrowser for help on using the repository browser.