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

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

update ti head

File size: 7.1 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//
[1340]26// $Id: G4BGGPionInelasticXS.cc,v 1.9 2010/10/12 06:03:49 dennis Exp $
27// GEANT4 tag $Name: hadr-cross-V09-03-12 $
[819]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"
[962]48#include "G4HadronNucleonXsc.hh"
[819]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;
[962]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;
[819]66}
67
68//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
69
70G4BGGPionInelasticXS::~G4BGGPionInelasticXS()
71{
72 delete fGlauber;
73 delete fPion;
[962]74 delete fHadron;
[819]75}
76
77//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
78
[1340]79G4double
80G4BGGPionInelasticXS::GetZandACrossSection(const G4DynamicParticle* dp,
81 G4int Z, G4int A, G4double)
[819]82{
83 G4double cross = 0.0;
84 G4double ekin = dp->GetKineticEnergy();
[1340]85// G4int iz = G4int(Z);
86 if(Z > 92) Z = 92;
[819]87
[962]88 if(ekin <= fLowEnergy) {
[1340]89 cross = theCoulombFac[Z];
[962]90 if(isPiplus) { cross *= CoulombFactor(ekin, A); }
91
[1340]92 } else if(Z == 1) {
[962]93
[1340]94 if( A < 2) {
[1228]95 fHadron->GetHadronNucleonXscNS(dp, G4Proton::Proton());
[962]96 cross = fHadron->GetInelasticHadronNucleonXsc();
97 } else {
[1228]98 fHadron->GetHadronNucleonXscNS(dp, G4Proton::Proton());
[962]99 cross = fHadron->GetInelasticHadronNucleonXsc();
[1228]100 fHadron->GetHadronNucleonXscNS(dp, G4Neutron::Neutron());
[962]101 cross += fHadron->GetInelasticHadronNucleonXsc();
102 }
103
104 } else if(ekin > fGlauberEnergy) {
[1340]105 cross = theGlauberFac[Z]*fGlauber->GetInelasticGlauberGribov(dp, Z, A);
[962]106
[819]107 } else {
[962]108 cross = fPion->GetIsoZACrossSection(dp, Z, A);
[819]109 }
110
111 if(verboseLevel > 1)
112 G4cout << "G4BGGPionInelasticXS::GetCrossSection for "
113 << dp->GetDefinition()->GetParticleName()
114 << " Ekin(GeV)= " << dp->GetKineticEnergy()
115 << " in nucleus Z= " << Z << " A= " << A
116 << " XS(b)= " << cross/barn
117 << G4endl;
118
119 return cross;
120}
121
122//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
123
[962]124void G4BGGPionInelasticXS::BuildPhysicsTable(const G4ParticleDefinition& p)
[819]125{
[962]126 if(&p == G4PionPlus::PionPlus() || &p == G4PionMinus::PionMinus()) {
127 particle = &p;
128 Initialise();
129 } else {
130 G4cout << "### G4BGGPionInelasticXS WARNING: is not applicable to "
131 << p.GetParticleName()
132 << G4endl;
133 }
[819]134}
135
136//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
137
138void G4BGGPionInelasticXS::DumpPhysicsTable(const G4ParticleDefinition&)
139{
140 G4cout << "G4BGGPionInelasticXS:"<<G4endl;
141}
142
143//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
144
145void G4BGGPionInelasticXS::Initialise()
146{
[962]147 if(isInitialized) return;
148 isInitialized = true;
149
150 fPion = new G4UPiNuclearCrossSection();
151 fGlauber = new G4GlauberGribovCrossSection();
152 fHadron = new G4HadronNucleonXsc();
153 fPion->BuildPhysicsTable(*particle);
154 fGlauber->BuildPhysicsTable(*particle);
155 if(particle == G4PionPlus::PionPlus()) isPiplus = true;
156
[819]157 G4ParticleDefinition* part = const_cast<G4ParticleDefinition*>(particle);
158 G4ThreeVector mom(0.0,0.0,1.0);
[962]159 G4DynamicParticle dp(part, mom, fGlauberEnergy);
[819]160
161 G4NistManager* nist = G4NistManager::Instance();
162
[1340]163 G4double csup, csdn;
164 G4int A;
[819]165
166 if(verboseLevel > 0) G4cout << "### G4BGGPionInelasticXS::Initialise for "
[962]167 << particle->GetParticleName()
168 << " isPiplus: " << isPiplus
169 << G4endl;
[819]170
171 for(G4int iz=2; iz<93; iz++) {
172
173 G4double Z = G4double(iz);
[1340]174 A = G4lrint(nist->GetAtomicMassAmu(iz));
[819]175
[1340]176 csup = fGlauber->GetInelasticGlauberGribov(&dp, iz, A);
177 csdn = fPion->GetInelasticCrossSection(&dp, iz, A);
[819]178
[962]179 theGlauberFac[iz] = csdn/csup;
[819]180 if(verboseLevel > 0) G4cout << "Z= " << Z << " A= " << A
[962]181 << " factor= " << theGlauberFac[iz] << G4endl;
[819]182 }
[962]183 dp.SetKineticEnergy(fLowEnergy);
184 fHadron->GetHadronNucleonXscNS(&dp, G4Proton::Proton());
185 theCoulombFac[1] = fHadron->GetInelasticHadronNucleonXsc();
[1340]186 if(isPiplus) { theCoulombFac[1] /= CoulombFactor(fLowEnergy,1); }
[962]187
188 for(G4int iz=2; iz<93; iz++) {
189
190 G4double Z = G4double(iz);
[1340]191 A = G4lrint(nist->GetAtomicMassAmu(iz));
[962]192
[1340]193 theCoulombFac[iz] = fPion->GetIsoZACrossSection(&dp, iz, A);
[962]194 if(isPiplus) { theCoulombFac[iz] /= CoulombFactor(fLowEnergy,A); }
195
196 if(verboseLevel > 0) G4cout << "Z= " << Z << " A= " << A
197 << " factor= " << theCoulombFac[iz] << G4endl;
198 }
[819]199}
200
201
[1340]202G4double G4BGGPionInelasticXS::CoulombFactor(G4double kinEnergy, G4int A)
[962]203{
204 G4double res= 0.0;
205 if(kinEnergy <= DBL_MIN) return res;
[1340]206 else if(A < 2) return kinEnergy*kinEnergy;
[962]207
208 G4double elog = std::log10(kinEnergy/GeV);
[1340]209 G4double aa = A;
[819]210
[962]211 // from G4ProtonInelasticCrossSection
[1340]212 G4double f1 = 8.0 - 8.0/aa - 0.008*aa;
213 G4double f2 = 2.34 - 5.4/aa - 0.0028*aa;
[962]214
215 res = 1.0/(1.0 + std::exp(-f1*(elog + f2)));
216
[1340]217 f1 = 5.6 - 0.016*aa;
218 f2 = 1.37 + 1.37/aa;
219 res *= ( 1.0 + (0.8 + 18./aa - 0.002*aa)/(1.0 + std::exp(f1*(elog + f2))));
[962]220 return res;
221}
222
223
Note: See TracBrowser for help on using the repository browser.