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

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

geant4 tag 9.4

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