Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/cross_sections/src/G4BGGNucleonElasticXS.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

Diff markup

Differences between /processes/hadronic/cross_sections/src/G4BGGNucleonElasticXS.cc (Version 11.3.0) and /processes/hadronic/cross_sections/src/G4BGGNucleonElasticXS.cc (Version 10.0.p2)


  1 //                                                  1 //
  2 // *******************************************      2 // ********************************************************************
  3 // * License and Disclaimer                         3 // * License and Disclaimer                                           *
  4 // *                                                4 // *                                                                  *
  5 // * The  Geant4 software  is  copyright of th      5 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
  6 // * the Geant4 Collaboration.  It is provided      6 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
  7 // * conditions of the Geant4 Software License      7 // * conditions of the Geant4 Software License,  included in the file *
  8 // * LICENSE and available at  http://cern.ch/      8 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
  9 // * include a list of copyright holders.           9 // * include a list of copyright holders.                             *
 10 // *                                               10 // *                                                                  *
 11 // * Neither the authors of this software syst     11 // * Neither the authors of this software system, nor their employing *
 12 // * institutes,nor the agencies providing fin     12 // * institutes,nor the agencies providing financial support for this *
 13 // * work  make  any representation or  warran     13 // * work  make  any representation or  warranty, express or implied, *
 14 // * regarding  this  software system or assum     14 // * regarding  this  software system or assume any liability for its *
 15 // * use.  Please see the license in the file      15 // * use.  Please see the license in the file  LICENSE  and URL above *
 16 // * for the full disclaimer and the limitatio     16 // * for the full disclaimer and the limitation of liability.         *
 17 // *                                               17 // *                                                                  *
 18 // * This  code  implementation is the result      18 // * This  code  implementation is the result of  the  scientific and *
 19 // * technical work of the GEANT4 collaboratio     19 // * technical work of the GEANT4 collaboration.                      *
 20 // * By using,  copying,  modifying or  distri     20 // * By using,  copying,  modifying or  distributing the software (or *
 21 // * any work based  on the software)  you  ag     21 // * any work based  on the software)  you  agree  to acknowledge its *
 22 // * use  in  resulting  scientific  publicati     22 // * use  in  resulting  scientific  publications,  and indicate your *
 23 // * acceptance of all terms of the Geant4 Sof     23 // * acceptance of all terms of the Geant4 Software license.          *
 24 // *******************************************     24 // ********************************************************************
 25 //                                                 25 //
                                                   >>  26 // $Id: G4BGGNucleonElasticXS.cc 78189 2013-12-04 16:34:08Z gcosmo $
                                                   >>  27 //
 26 // -------------------------------------------     28 // -------------------------------------------------------------------
 27 //                                                 29 //
 28 // GEANT4 Class file                               30 // GEANT4 Class file
 29 //                                                 31 //
 30 //                                                 32 //
 31 // File name:     G4BGGNucleonElasticXS            33 // File name:     G4BGGNucleonElasticXS
 32 //                                                 34 //
 33 // Author:        Vladimir Ivanchenko              35 // Author:        Vladimir Ivanchenko
 34 //                                                 36 //
 35 // Creation date: 13.03.2007                       37 // Creation date: 13.03.2007
                                                   >>  38 // Modifications:
                                                   >>  39 //
 36 //                                                 40 //
 37 // -------------------------------------------     41 // -------------------------------------------------------------------
 38 //                                                 42 //
 39                                                    43 
 40 #include "G4BGGNucleonElasticXS.hh"                44 #include "G4BGGNucleonElasticXS.hh"
 41 #include "G4SystemOfUnits.hh"                      45 #include "G4SystemOfUnits.hh"
 42 #include "G4ComponentGGHadronNucleusXsc.hh"    <<  46 #include "G4GlauberGribovCrossSection.hh"
 43 #include "G4NucleonNuclearCrossSection.hh"         47 #include "G4NucleonNuclearCrossSection.hh"
 44 #include "G4HadronNucleonXsc.hh"                   48 #include "G4HadronNucleonXsc.hh"
 45 #include "G4NuclearRadii.hh"                   <<  49 #include "G4ComponentSAIDTotalXS.hh"
 46 #include "G4Proton.hh"                             50 #include "G4Proton.hh"
 47 #include "G4Neutron.hh"                            51 #include "G4Neutron.hh"
 48 #include "G4NistManager.hh"                        52 #include "G4NistManager.hh"
 49 #include "G4NuclearRadii.hh"                   <<  53 #include "G4Log.hh"
                                                   >>  54 #include "G4Exp.hh"
 50                                                    55 
 51 #include "G4CrossSectionDataSetRegistry.hh"        56 #include "G4CrossSectionDataSetRegistry.hh"
 52                                                    57 
 53 G4double G4BGGNucleonElasticXS::theGlauberFacP <<  58 const G4double llog10 = G4Log(10.);
 54 G4double G4BGGNucleonElasticXS::theCoulombFacP << 
 55 G4double G4BGGNucleonElasticXS::theGlauberFacN << 
 56 G4double G4BGGNucleonElasticXS::theCoulombFacN << 
 57 G4int G4BGGNucleonElasticXS::theA[93] = {0};   << 
 58                                                    59 
 59 G4BGGNucleonElasticXS::G4BGGNucleonElasticXS(c     60 G4BGGNucleonElasticXS::G4BGGNucleonElasticXS(const G4ParticleDefinition* p)
 60  : G4VCrossSectionDataSet("BarashenkovGlauberG <<  61  : G4VCrossSectionDataSet("Barashenkov-Glauber") 
 61 {                                                  62 {
 62   verboseLevel = 0;                                63   verboseLevel = 0;
 63   fGlauberEnergy = 91.*GeV;                        64   fGlauberEnergy = 91.*GeV;
 64   fLowEnergy = 14.0*MeV;                       <<  65   fPDGEnergy = 5*GeV;
 65   fNucleon = new G4NucleonNuclearCrossSection( <<  66   fLowEnergy = 14.*MeV;
 66   fGlauber = new G4ComponentGGHadronNucleusXsc <<  67   fSAIDLowEnergyLimit = 1*MeV;
 67   fHadron = new G4HadronNucleonXsc();          <<  68   fSAIDHighEnergyLimit = 1.3*GeV;
 68                                                <<  69   fLowestXSection = millibarn;
 69   theProton = G4Proton::Proton();              <<  70   for (G4int i = 0; i < 93; ++i) {
 70   isProton = (theProton == p);                 <<  71     theGlauberFac[i] = 0.0;
 71   SetForAllAtomsAndEnergies(true);             <<  72     theCoulombFac[i] = 0.0;
 72                                                <<  73     theA[i] = 1;
 73   if (0 == theA[0]) { Initialise(); }          <<  74   }
                                                   >>  75   fNucleon = 0;
                                                   >>  76   fGlauber = 0;
                                                   >>  77   fHadron  = 0;
                                                   >>  78   fSAID    = 0;
                                                   >>  79   particle = p;
                                                   >>  80   theProton= G4Proton::Proton();
                                                   >>  81   isProton = false;
                                                   >>  82   isInitialized = false;
 74 }                                                  83 }
 75                                                    84 
 76 //....oooOO0OOooo........oooOO0OOooo........oo     85 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 77                                                    86 
 78 G4BGGNucleonElasticXS::~G4BGGNucleonElasticXS(     87 G4BGGNucleonElasticXS::~G4BGGNucleonElasticXS()
 79 {                                                  88 {
 80   delete fHadron;                                  89   delete fHadron;
                                                   >>  90   delete fSAID;
 81 }                                                  91 }
 82                                                    92 
 83 //....oooOO0OOooo........oooOO0OOooo........oo     93 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 84                                                    94 
 85 G4bool                                             95 G4bool 
 86 G4BGGNucleonElasticXS::IsElementApplicable(con     96 G4BGGNucleonElasticXS::IsElementApplicable(const G4DynamicParticle*, G4int,
 87                                            con <<  97              const G4Material*)
 88 {                                                  98 {
 89   return true;                                     99   return true;
 90 }                                                 100 }
 91                                                   101 
 92 //....oooOO0OOooo........oooOO0OOooo........oo    102 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 93                                                   103 
 94 G4bool G4BGGNucleonElasticXS::IsIsoApplicable(    104 G4bool G4BGGNucleonElasticXS::IsIsoApplicable(const G4DynamicParticle*, 
 95                                                << 105                 G4int Z, G4int,  
 96                                                << 106                 const G4Element*,
 97                                                << 107                 const G4Material*)
 98 {                                                 108 {
 99   return (1 == Z);                                109   return (1 == Z);
100 }                                                 110 }
101                                                   111 
102 //....oooOO0OOooo........oooOO0OOooo........oo    112 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
103                                                   113 
104 G4double                                          114 G4double
105 G4BGGNucleonElasticXS::GetElementCrossSection(    115 G4BGGNucleonElasticXS::GetElementCrossSection(const G4DynamicParticle* dp,
106                                                << 116                 G4int ZZ, const G4Material*)
107 {                                                 117 {
108   // this method should be called only for Z >    118   // this method should be called only for Z > 1
109                                                   119 
110   G4double cross = 0.0;                           120   G4double cross = 0.0;
111   G4int Z = std::min(ZZ, 92);                  << 
112   G4double ekin = dp->GetKineticEnergy();         121   G4double ekin = dp->GetKineticEnergy();
                                                   >> 122   G4int Z = ZZ;
113   if(1 == Z) {                                    123   if(1 == Z) {
114     cross = 1.0115*GetIsoCrossSection(dp,1,1);    124     cross = 1.0115*GetIsoCrossSection(dp,1,1);
115   } else {                                        125   } else {
                                                   >> 126     if(Z > 92) { Z = 92; }
                                                   >> 127 
116     if(ekin <= fLowEnergy) {                      128     if(ekin <= fLowEnergy) {
117       cross = (isProton) ? theCoulombFacP[Z] : << 129       cross = theCoulombFac[Z]*CoulombFactor(ekin, Z);
118       cross *= CoulombFactor(ekin, Z);         << 130 
119     } else if(ekin > fGlauberEnergy) {            131     } else if(ekin > fGlauberEnergy) {
120       cross = (isProton) ? theGlauberFacP[Z] : << 132       cross = theGlauberFac[Z]*fGlauber->GetElasticGlauberGribov(dp, Z, theA[Z]);
121       cross *= fGlauber->GetElasticGlauberGrib << 
122     } else {                                      133     } else {
123       cross = fNucleon->GetElasticCrossSection    134       cross = fNucleon->GetElasticCrossSection(dp, Z);
124     }                                             135     }
125   }                                               136   }
126 #ifdef G4VERBOSE                               << 137 
127   if (verboseLevel > 1) {                      << 138   if(verboseLevel > 1) {
128     G4cout << "G4BGGNucleonElasticXS::GetEleme    139     G4cout << "G4BGGNucleonElasticXS::GetElementCrossSection  for "
129            << dp->GetDefinition()->GetParticle << 140      << dp->GetDefinition()->GetParticleName()
130            << "  Ekin(GeV)= " << dp->GetKineti << 141      << "  Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
131            << " in nucleus Z= " << Z << "  A=  << 142      << " in nucleus Z= " << Z << "  A= " << theA[Z]
132            << " XS(b)= " << cross/barn         << 143      << " XS(b)= " << cross/barn 
133            << G4endl;                          << 144      << G4endl;
134   }                                               145   }
135 #endif                                         << 146   //AR-04Dec2013  if(cross <= fLowestXSection) { cross = 0.0; }
136   return cross;                                   147   return cross;
137 }                                                 148 }
138                                                   149 
139 //....oooOO0OOooo........oooOO0OOooo........oo    150 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
140                                                   151 
141 G4double                                          152 G4double
142 G4BGGNucleonElasticXS::GetIsoCrossSection(cons    153 G4BGGNucleonElasticXS::GetIsoCrossSection(const G4DynamicParticle* dp, 
143                                           G4in << 154             G4int Z, G4int A, 
144                                           cons << 155             const G4Isotope*,
145                                           cons << 156             const G4Element*,
146                                           cons << 157             const G4Material*)
147 {                                                 158 {
148   // this method should be called only for Z =    159   // this method should be called only for Z = 1
149   fHadron->HadronNucleonXscNS(dp->GetDefinitio << 
150                               dp->GetKineticEn << 
151   G4double cross = A*fHadron->GetElasticHadron << 
152                                                   160 
153 #ifdef G4VERBOSE                               << 161   G4double cross = 0.0;
154   if (verboseLevel > 1) {                      << 162   G4double ekin = dp->GetKineticEnergy();
                                                   >> 163 
                                                   >> 164   if(ekin <= fSAIDLowEnergyLimit) {
                                                   >> 165     cross = theCoulombFac[0]*CoulombFactor(ekin, 1);
                                                   >> 166   } else if(ekin <= fSAIDHighEnergyLimit) {
                                                   >> 167     cross = fSAID->GetElasticIsotopeCrossSection(particle, ekin, 1, 1);
                                                   >> 168   } else if(ekin <= fPDGEnergy) {
                                                   >> 169     G4double cross1 = 
                                                   >> 170       fSAID->GetElasticIsotopeCrossSection(particle, fPDGEnergy, 1, 1);
                                                   >> 171     G4double cross2 = theCoulombFac[1]; 
                                                   >> 172     cross = cross1 + (cross2 - cross1)*(ekin - fSAIDHighEnergyLimit)
                                                   >> 173       /(fPDGEnergy - fSAIDHighEnergyLimit);
                                                   >> 174   } else {
                                                   >> 175     fHadron->GetHadronNucleonXscPDG(dp, theProton);
                                                   >> 176     cross = fHadron->GetElasticHadronNucleonXsc();
                                                   >> 177   } 
                                                   >> 178   cross *= A;
                                                   >> 179 
                                                   >> 180   if(verboseLevel > 1) {
155     G4cout << "G4BGGNucleonElasticXS::GetIsoCr    181     G4cout << "G4BGGNucleonElasticXS::GetIsoCrossSection  for "
156            << dp->GetDefinition()->GetParticle << 182      << dp->GetDefinition()->GetParticleName()
157            << "  Ekin(GeV)= " << dp->GetKineti << 183      << "  Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
158            << " in nucleus  Z=1  A=" << A      << 184      << " in nucleus Z= " << Z << "  A= " << A
159            << " XS(b)= " << cross/barn         << 185      << " XS(b)= " << cross/barn 
160            << G4endl;                          << 186      << G4endl;
161   }                                               187   }
162 #endif                                         << 188   //AR-04Dec2013  if(cross <= fLowestXSection) { cross = 0.0; }
163   return cross;                                   189   return cross;
164 }                                                 190 }
165                                                   191 
166 //....oooOO0OOooo........oooOO0OOooo........oo    192 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
167                                                   193 
168 void G4BGGNucleonElasticXS::BuildPhysicsTable(    194 void G4BGGNucleonElasticXS::BuildPhysicsTable(const G4ParticleDefinition& p)
169 {                                                 195 {
170   if(&p == theProton || &p == G4Neutron::Neutr    196   if(&p == theProton || &p == G4Neutron::Neutron()) {
171     isProton = (theProton == &p);              << 197     particle = &p;
172                                                   198 
173   } else {                                        199   } else {
174     G4ExceptionDescription ed;                 << 200     G4cout << "### G4BGGNucleonElasticXS WARNING: is not applicable to " 
175     ed << "This BGG cross section is applicabl << 201      << p.GetParticleName()
176        << p.GetParticleName() << G4endl;       << 202      << G4endl;
177     G4Exception("G4BGGNucleonElasticXS::BuildP << 203     throw G4HadronicException(__FILE__, __LINE__,
178                 FatalException, ed);           << 204     "G4BGGNucleonElasticXS::BuildPhysicsTable is used for wrong particle");
                                                   >> 205     return;
179   }                                               206   }
180 }                                              << 
181                                                   207 
182 //....oooOO0OOooo........oooOO0OOooo........oo << 208   if(isInitialized) { return; }
                                                   >> 209   isInitialized = true;
183                                                   210 
184 void G4BGGNucleonElasticXS::Initialise()       << 211   fNucleon = (G4NucleonNuclearCrossSection*)G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4NucleonNuclearCrossSection::Default_Name());
185 {                                              << 212   fGlauber = (G4GlauberGribovCrossSection*)G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4GlauberGribovCrossSection::Default_Name());
186   theA[0] = theA[1] = 1;                       << 213     
                                                   >> 214   fHadron  = new G4HadronNucleonXsc();
                                                   >> 215   fSAID    = new G4ComponentSAIDTotalXS();
                                                   >> 216   fNucleon->BuildPhysicsTable(*particle);
                                                   >> 217   fGlauber->BuildPhysicsTable(*particle);
                                                   >> 218   if(particle == theProton) { 
                                                   >> 219     isProton = true; 
                                                   >> 220     fSAIDHighEnergyLimit = 3*GeV;
                                                   >> 221   }
                                                   >> 222 
                                                   >> 223   G4ParticleDefinition* part = const_cast<G4ParticleDefinition*>(particle);
187   G4ThreeVector mom(0.0,0.0,1.0);                 224   G4ThreeVector mom(0.0,0.0,1.0);
188   G4DynamicParticle dp(theProton, mom, fGlaube << 225   G4DynamicParticle dp(part, mom, fGlauberEnergy);
189                                                   226 
190   G4NistManager* nist = G4NistManager::Instanc    227   G4NistManager* nist = G4NistManager::Instance();
                                                   >> 228 
191   G4double csup, csdn;                            229   G4double csup, csdn;
                                                   >> 230   G4int A;
                                                   >> 231 
                                                   >> 232   if(verboseLevel > 0) {
                                                   >> 233     G4cout << "### G4BGGNucleonElasticXS::Initialise for "
                                                   >> 234      << particle->GetParticleName() << G4endl;
                                                   >> 235   }
                                                   >> 236 
                                                   >> 237   for(G4int iz=2; iz<93; iz++) {
192                                                   238 
193   for (G4int iz=2; iz<93; ++iz) {              << 239     A = G4lrint(nist->GetAtomicMassAmu(iz));
194     G4int A = G4lrint(nist->GetAtomicMassAmu(i << 
195     theA[iz] = A;                                 240     theA[iz] = A;
196                                                   241 
197     csup = fGlauber->GetElasticGlauberGribov(&    242     csup = fGlauber->GetElasticGlauberGribov(&dp, iz, A);
198     csdn = fNucleon->GetElasticCrossSection(&d    243     csdn = fNucleon->GetElasticCrossSection(&dp, iz);
199     theGlauberFacP[iz] = csdn/csup;            << 
200   }                                            << 
201                                                << 
202   dp.SetDefinition(G4Neutron::Neutron());      << 
203   for (G4int iz=2; iz<93; ++iz) {              << 
204     csup = fGlauber->GetElasticGlauberGribov(& << 
205     csdn = fNucleon->GetElasticCrossSection(&d << 
206     theGlauberFacN[iz] = csdn/csup;            << 
207                                                   244 
208     if (verboseLevel > 1) {                    << 245     theGlauberFac[iz] = csdn/csup;
209         G4cout << "G4BGGNucleonElasticXS::Init << 246     if(verboseLevel > 0) { 
210                << " GFactorP=" << theGlauberFa << 247       G4cout << "Z= " << iz <<  "  A= " << A 
211                << " GFactorN=" << theGlauberFa << 248        << " factor= " << theGlauberFac[iz] << G4endl;
212     }                                             249     } 
213   }                                               250   }
214                                                   251 
215   theCoulombFacP[0] = theCoulombFacP[1] =      << 252   theCoulombFac[0] = 
216     theCoulombFacN[0] = theCoulombFacN[1] = 1. << 253     fSAID->GetElasticIsotopeCrossSection(particle,fSAIDLowEnergyLimit,1,1)
217   dp.SetDefinition(theProton);                 << 254     /CoulombFactor(fSAIDLowEnergyLimit, 1);
                                                   >> 255 
                                                   >> 256   dp.SetKineticEnergy(fPDGEnergy);
                                                   >> 257   fHadron->GetHadronNucleonXscPDG(&dp, theProton);
                                                   >> 258   theCoulombFac[1] = fHadron->GetElasticHadronNucleonXsc();
                                                   >> 259 
                                                   >> 260   if(verboseLevel > 0) {
                                                   >> 261     G4cout << "Z=1  A=1 " << " factor0= " << theCoulombFac[0] 
                                                   >> 262      << " factor1= " << theCoulombFac[1] 
                                                   >> 263      << G4endl; 
                                                   >> 264   }
                                                   >> 265 
218   dp.SetKineticEnergy(fLowEnergy);                266   dp.SetKineticEnergy(fLowEnergy);
219   for (G4int iz=2; iz<93; ++iz) {              << 267   for(G4int iz=2; iz<93; iz++) {
220     theCoulombFacP[iz] = fNucleon->GetElasticC << 268     theCoulombFac[iz] = 
221       /CoulombFactor(fLowEnergy, iz);          << 269       fNucleon->GetElasticCrossSection(&dp, iz)/CoulombFactor(fLowEnergy, iz);
222   }                                            << 270     if(verboseLevel > 0) {
223   dp.SetDefinition(G4Neutron::Neutron());      << 271       G4cout << "Z= " << iz <<  "  A= " << theA[iz]
224   for(G4int iz=2; iz<93; ++iz) {               << 272        << " factor= " << theCoulombFac[iz] << G4endl; 
225     theCoulombFacN[iz] = fNucleon->GetElasticC << 
226       /CoulombFactor(fLowEnergy, iz);          << 
227                                                << 
228     if (verboseLevel > 1) {                    << 
229       G4cout << "G4BGGNucleonElasticXS::Initia << 
230        << " CFactorP=" << theCoulombFacP[iz]   << 
231        << " CFactorN=" << theCoulombFacN[iz] < << 
232     }                                             273     }
233   }                                               274   }
234 }                                                 275 }
235                                                   276 
236 //....oooOO0OOooo........oooOO0OOooo........oo    277 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
237                                                   278 
238 G4double G4BGGNucleonElasticXS::CoulombFactor(    279 G4double G4BGGNucleonElasticXS::CoulombFactor(G4double kinEnergy, G4int Z)
239 {                                                 280 {
240   return (isProton) ?                          << 281   G4double res= 1.0;
241     G4NuclearRadii::CoulombFactor(Z, theA[Z],  << 282   
                                                   >> 283   // from G4ProtonInelasticCrossSection
                                                   >> 284   if(isProton) {
                                                   >> 285 
                                                   >> 286     if (Z <= 1) { return kinEnergy*kinEnergy; }
                                                   >> 287 
                                                   >> 288     G4double elog = G4Log(kinEnergy/GeV)/llog10;
                                                   >> 289     G4double aa = theA[Z];
                                                   >> 290 
                                                   >> 291     G4double ff1 = 5.6  - 0.016*aa;    // slope of the drop at medium energies.
                                                   >> 292     G4double ff2 = 1.37 + 1.37/aa;     // start of the slope.
                                                   >> 293     G4double ff3 = 0.8  + 18./aa - 0.002*aa;   // stephight
                                                   >> 294     res = 1.0 + ff3*(1.0 - (1.0/(1+G4Exp(-ff1*(elog + ff2)))));
                                                   >> 295 
                                                   >> 296     ff1 = 8.   - 8./aa  - 0.008*aa; // slope of the rise
                                                   >> 297     ff2 = 2.34 - 5.4/aa - 0.0028*aa; // start of the rise
                                                   >> 298     res /= (1.0 + G4Exp(-ff1*(elog + ff2)));
                                                   >> 299 
                                                   >> 300   } 
                                                   >> 301   return res;  
242 }                                                 302 }
243                                                   303 
244 //....oooOO0OOooo........oooOO0OOooo........oo    304 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
245                                                   305 
246 void G4BGGNucleonElasticXS::CrossSectionDescri    306 void G4BGGNucleonElasticXS::CrossSectionDescription(std::ostream& outFile) const
247 {                                                 307 {
248   outFile << "The Barashenkov-Glauber-Gribov c    308   outFile << "The Barashenkov-Glauber-Gribov cross section handles elastic\n"
249           << "scattering of protons and neutro    309           << "scattering of protons and neutrons from nuclei using the\n"
250           << "Barashenkov parameterization bel    310           << "Barashenkov parameterization below 91 GeV and the Glauber-Gribov\n"
251           << "parameterization above 91 GeV. n    311           << "parameterization above 91 GeV. n";
252 }                                                 312 }
253                                                   313 
254 //....oooOO0OOooo........oooOO0OOooo........oo    314 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
255                                                   315