Ignore:
Timestamp:
Apr 6, 2009, 12:30:29 PM (15 years ago)
Author:
garnier
Message:

update processes

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/source/processes/hadronic/cross_sections/src/G4BGGPionInelasticXS.cc

    r819 r962  
    2424// ********************************************************************
    2525//
    26 // $Id: G4BGGPionInelasticXS.cc,v 1.1 2007/03/13 15:19:30 vnivanch Exp $
    27 // GEANT4 tag $Name: $
     26// $Id: G4BGGPionInelasticXS.cc,v 1.3 2008/12/01 16:50:23 vnivanch Exp $
     27// GEANT4 tag $Name: geant4-09-02-ref-02 $
    2828//
    2929// -------------------------------------------------------------------
     
    4646#include "G4GlauberGribovCrossSection.hh"
    4747#include "G4UPiNuclearCrossSection.hh"
     48#include "G4HadronNucleonXsc.hh"
    4849#include "G4PionPlus.hh"
    4950#include "G4PionMinus.hh"
     
    5556{
    5657  verboseLevel = 0;
    57   thEnergy     = 100.*GeV;
    58   if(p == G4PionPlus::PionPlus() || p == G4PionMinus::PionMinus()) {
    59     fPion = new G4UPiNuclearCrossSection();
    60     fGlauber = new G4GlauberGribovCrossSection();
    61     particle = p;
    62     Initialise();
    63   } else {
    64     fPion = 0;
    65     fGlauber = 0;
    66     particle = 0;
    67     if(p) G4cout << "### G4BGGPionInelasticXS WARNING: is not applicable to "
    68                  << p->GetParticleName()
    69                  << G4endl;
    70     else  G4cout << "### G4BGGPionInelasticXS WARNING: particle is not defined "
    71                  << G4endl;
    72   }
     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;
    7366}
    7467
     
    7972  delete fGlauber;
    8073  delete fPion;
     74  delete fHadron;
    8175}
    8276
     
    9084  G4double cross = 0.0;
    9185  G4double ekin = dp->GetKineticEnergy();
    92   G4int iz = G4int(Z + 0.5);
     86  G4int iz = G4int(Z);
    9387  if(iz > 92) iz = 92;
    9488
    95   if(ekin > thEnergy) {
    96     cross = theFac[iz]*fGlauber->GetInelasticGlauberGribov(dp, Z, A);
     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->GetHadronNucleonXscPDG(dp, G4Proton::Proton());
     97      //fHadron->GetHadronNucleonXscEL(dp, G4Proton::Proton());
     98      fHadron->GetHadronNucleonXscNS(dp, G4Proton::Proton());
     99      //fHadron->GetHadronNucleonXscVU(dp, G4Proton::Proton());
     100      //fHadron->GetHadronNucleonXscMK(dp, G4Proton::Proton());
     101      cross = fHadron->GetInelasticHadronNucleonXsc();
     102    } else {
     103      fHadron->GetHadronNucleonXscNS(dp, G4Proton::Proton());
     104      cross = fHadron->GetInelasticHadronNucleonXsc();
     105      fHadron->GetHadronNucleonXscNS(dp, G4Neutron::Neutron());
     106      cross += fHadron->GetInelasticHadronNucleonXsc();
     107    }
     108
     109  } else if(ekin > fGlauberEnergy) {
     110    cross = theGlauberFac[iz]*fGlauber->GetInelasticGlauberGribov(dp, Z, A);
     111
    97112  } else {
    98     cross = fPion->GetInelasticCrossSection(dp, Z, A);
     113    cross = fPion->GetIsoZACrossSection(dp, Z, A);
    99114  }
    100115
     
    112127//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
    113128
    114 void G4BGGPionInelasticXS::BuildPhysicsTable(const G4ParticleDefinition&)
    115 {
     129void G4BGGPionInelasticXS::BuildPhysicsTable(const G4ParticleDefinition& p)
     130{
     131  if(&p == G4PionPlus::PionPlus() || &p == G4PionMinus::PionMinus()) {
     132    particle = &p;
     133    Initialise();
     134  } else {
     135    G4cout << "### G4BGGPionInelasticXS WARNING: is not applicable to "
     136           << p.GetParticleName()
     137           << G4endl;
     138  }
    116139}
    117140
     
    127150void G4BGGPionInelasticXS::Initialise()
    128151{
     152  if(isInitialized) return;
     153  isInitialized = true;
     154
     155  fPion = new G4UPiNuclearCrossSection();
     156  fGlauber = new G4GlauberGribovCrossSection();
     157  fHadron  = new G4HadronNucleonXsc();
     158  fPion->BuildPhysicsTable(*particle);
     159  fGlauber->BuildPhysicsTable(*particle);
     160  if(particle == G4PionPlus::PionPlus()) isPiplus = true;
     161
    129162  G4ParticleDefinition* part = const_cast<G4ParticleDefinition*>(particle);
    130163  G4ThreeVector mom(0.0,0.0,1.0);
    131   G4DynamicParticle dp(part, mom, thEnergy);
     164  G4DynamicParticle dp(part, mom, fGlauberEnergy);
    132165
    133166  G4NistManager* nist = G4NistManager::Instance();
    134   G4double A = nist->GetAtomicMassAmu(2);
    135 
    136   G4double csup, csdn;
     167
     168  G4double A, csup, csdn;
    137169
    138170  if(verboseLevel > 0) G4cout << "### G4BGGPionInelasticXS::Initialise for "
    139                               << particle->GetParticleName() << G4endl;
     171                              << particle->GetParticleName()
     172                              << " isPiplus: " << isPiplus
     173                              << G4endl;
    140174
    141175  for(G4int iz=2; iz<93; iz++) {
     
    147181    csdn = fPion->GetInelasticCrossSection(&dp, Z, A);
    148182
    149     theFac[iz] = csdn/csup;
     183    theGlauberFac[iz] = csdn/csup;
    150184    if(verboseLevel > 0) G4cout << "Z= " << Z <<  "  A= " << A
    151                                 << " factor= " << theFac[iz] << G4endl;
    152   }
    153 }
    154 
    155 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
    156 
    157 
     185                                << " factor= " << theGlauberFac[iz] << G4endl;
     186  }
     187  dp.SetKineticEnergy(fLowEnergy);
     188  fHadron->GetHadronNucleonXscNS(&dp, G4Proton::Proton());
     189  theCoulombFac[1] = fHadron->GetInelasticHadronNucleonXsc();
     190  if(isPiplus) { theCoulombFac[1] /= CoulombFactor(fLowEnergy,1.0); }
     191     
     192  for(G4int iz=2; iz<93; iz++) {
     193
     194    G4double Z = G4double(iz);
     195    A = nist->GetAtomicMassAmu(iz);
     196
     197    theCoulombFac[iz] = fPion->GetIsoZACrossSection(&dp, Z, A);
     198    if(isPiplus) { theCoulombFac[iz] /= CoulombFactor(fLowEnergy,A); }
     199
     200    if(verboseLevel > 0) G4cout << "Z= " << Z <<  "  A= " << A
     201                                << " factor= " << theCoulombFac[iz] << G4endl;
     202  }
     203}
     204
     205//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
     206
     207G4double G4BGGPionInelasticXS::CoulombFactor(G4double kinEnergy, G4double A)
     208{
     209  G4double res= 0.0;
     210  if(kinEnergy <= DBL_MIN) return res;
     211  else if(A < 1.5) return kinEnergy*kinEnergy;
     212 
     213  G4double elog = std::log10(kinEnergy/GeV);
     214
     215  // from G4ProtonInelasticCrossSection
     216  G4double f1 = 8.0  - 8.0/A - 0.008*A;
     217  G4double f2 = 2.34 - 5.4/A - 0.0028*A;
     218
     219  res = 1.0/(1.0 + std::exp(-f1*(elog + f2)));
     220 
     221  f1 = 5.6 - 0.016*A;
     222  f2 = 1.37 + 1.37/A;
     223  res *= ( 1.0 + (0.8 + 18./A - 0.002*A)/(1.0 + std::exp(f1*(elog + f2))));
     224  return res; 
     225}
     226
     227//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
     228
     229
Note: See TracChangeset for help on using the changeset viewer.