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

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

update ti head

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