| 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 | // $Id: G4BGGPionElasticXS.cc,v 1.8 2010/10/12 06:03:37 dennis Exp $
|
|---|
| 27 | // GEANT4 tag $Name: hadr-cross-V09-03-12 $
|
|---|
| 28 | //
|
|---|
| 29 | // -------------------------------------------------------------------
|
|---|
| 30 | //
|
|---|
| 31 | // GEANT4 Class file
|
|---|
| 32 | //
|
|---|
| 33 | //
|
|---|
| 34 | // File name: G4BGGPionElasticXS
|
|---|
| 35 | //
|
|---|
| 36 | // Author: Vladimir Ivanchenko
|
|---|
| 37 | //
|
|---|
| 38 | // Creation date: 01.10.2003
|
|---|
| 39 | // Modifications:
|
|---|
| 40 | //
|
|---|
| 41 | //
|
|---|
| 42 | // -------------------------------------------------------------------
|
|---|
| 43 | //
|
|---|
| 44 |
|
|---|
| 45 | #include "G4BGGPionElasticXS.hh"
|
|---|
| 46 | #include "G4GlauberGribovCrossSection.hh"
|
|---|
| 47 | #include "G4UPiNuclearCrossSection.hh"
|
|---|
| 48 | #include "G4HadronNucleonXsc.hh"
|
|---|
| 49 | #include "G4PionPlus.hh"
|
|---|
| 50 | #include "G4PionMinus.hh"
|
|---|
| 51 | #include "G4NistManager.hh"
|
|---|
| 52 |
|
|---|
| 53 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 54 |
|
|---|
| 55 | G4BGGPionElasticXS::G4BGGPionElasticXS(const G4ParticleDefinition*)
|
|---|
| 56 | {
|
|---|
| 57 | verboseLevel = 0;
|
|---|
| 58 | fGlauberEnergy = 91.*GeV;
|
|---|
| 59 | fLowEnergy = 20.*MeV;
|
|---|
| 60 | fPion = 0;
|
|---|
| 61 | fGlauber = 0;
|
|---|
| 62 | fHadron = 0;
|
|---|
| 63 | particle = 0;
|
|---|
| 64 | isPiplus = false;
|
|---|
| 65 | isInitialized = false;
|
|---|
| 66 | }
|
|---|
| 67 |
|
|---|
| 68 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 69 |
|
|---|
| 70 | G4BGGPionElasticXS::~G4BGGPionElasticXS()
|
|---|
| 71 | {
|
|---|
| 72 | delete fGlauber;
|
|---|
| 73 | delete fPion;
|
|---|
| 74 | delete fHadron;
|
|---|
| 75 | }
|
|---|
| 76 |
|
|---|
| 77 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 78 |
|
|---|
| 79 | G4double
|
|---|
| 80 | G4BGGPionElasticXS::GetZandACrossSection(const G4DynamicParticle* dp,
|
|---|
| 81 | G4int Z, G4int A, G4double)
|
|---|
| 82 | {
|
|---|
| 83 | G4double cross = 0.0;
|
|---|
| 84 | G4double ekin = dp->GetKineticEnergy();
|
|---|
| 85 | // G4int iz = G4int(Z);
|
|---|
| 86 | if(Z > 92) Z = 92;
|
|---|
| 87 |
|
|---|
| 88 | if(ekin <= fLowEnergy) {
|
|---|
| 89 | cross = theCoulombFac[Z];
|
|---|
| 90 | if(isPiplus) { cross *= CoulombFactor(ekin, A); }
|
|---|
| 91 |
|
|---|
| 92 | } else if(Z == 1) {
|
|---|
| 93 | if( A < 2) {
|
|---|
| 94 | fHadron->GetHadronNucleonXscNS(dp, G4Proton::Proton());
|
|---|
| 95 | cross = fHadron->GetElasticHadronNucleonXsc();
|
|---|
| 96 | } else {
|
|---|
| 97 | fHadron->GetHadronNucleonXscNS(dp, G4Proton::Proton());
|
|---|
| 98 | cross = fHadron->GetElasticHadronNucleonXsc();
|
|---|
| 99 | fHadron->GetHadronNucleonXscNS(dp, G4Neutron::Neutron());
|
|---|
| 100 | cross += fHadron->GetElasticHadronNucleonXsc();
|
|---|
| 101 | }
|
|---|
| 102 | } else if(ekin > fGlauberEnergy) {
|
|---|
| 103 | cross = theGlauberFac[Z]*fGlauber->GetElasticGlauberGribov(dp, Z, A);
|
|---|
| 104 | } else {
|
|---|
| 105 | cross = fPion->GetElasticCrossSection(dp, Z, A);
|
|---|
| 106 | }
|
|---|
| 107 |
|
|---|
| 108 | if(verboseLevel > 1)
|
|---|
| 109 | G4cout << "G4BGGPionElasticXS::GetCrossSection for "
|
|---|
| 110 | << dp->GetDefinition()->GetParticleName()
|
|---|
| 111 | << " Ekin(GeV)= " << dp->GetKineticEnergy()
|
|---|
| 112 | << " in nucleus Z= " << Z << " A= " << A
|
|---|
| 113 | << " XS(b)= " << cross/barn
|
|---|
| 114 | << G4endl;
|
|---|
| 115 |
|
|---|
| 116 | return cross;
|
|---|
| 117 | }
|
|---|
| 118 |
|
|---|
| 119 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 120 |
|
|---|
| 121 | void G4BGGPionElasticXS::BuildPhysicsTable(const G4ParticleDefinition& p)
|
|---|
| 122 | {
|
|---|
| 123 | if(&p == G4PionPlus::PionPlus() || &p == G4PionMinus::PionMinus()) {
|
|---|
| 124 | particle = &p;
|
|---|
| 125 | Initialise();
|
|---|
| 126 | } else {
|
|---|
| 127 | G4cout << "### G4BGGPionElasticXS WARNING: is not applicable to "
|
|---|
| 128 | << p.GetParticleName()
|
|---|
| 129 | << G4endl;
|
|---|
| 130 | }
|
|---|
| 131 | }
|
|---|
| 132 |
|
|---|
| 133 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 134 |
|
|---|
| 135 | void G4BGGPionElasticXS::DumpPhysicsTable(const G4ParticleDefinition&)
|
|---|
| 136 | {
|
|---|
| 137 | G4cout << "G4BGGPionElasticXS:"<<G4endl;
|
|---|
| 138 | }
|
|---|
| 139 |
|
|---|
| 140 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 141 |
|
|---|
| 142 | void G4BGGPionElasticXS::Initialise()
|
|---|
| 143 | {
|
|---|
| 144 | if(isInitialized) return;
|
|---|
| 145 | isInitialized = true;
|
|---|
| 146 |
|
|---|
| 147 | fPion = new G4UPiNuclearCrossSection();
|
|---|
| 148 | fGlauber = new G4GlauberGribovCrossSection();
|
|---|
| 149 | fHadron = new G4HadronNucleonXsc();
|
|---|
| 150 | fPion->BuildPhysicsTable(*particle);
|
|---|
| 151 | fGlauber->BuildPhysicsTable(*particle);
|
|---|
| 152 | if(particle == G4PionPlus::PionPlus()) isPiplus = true;
|
|---|
| 153 |
|
|---|
| 154 | G4ParticleDefinition* part = const_cast<G4ParticleDefinition*>(particle);
|
|---|
| 155 | G4ThreeVector mom(0.0,0.0,1.0);
|
|---|
| 156 | G4DynamicParticle dp(part, mom, fGlauberEnergy);
|
|---|
| 157 |
|
|---|
| 158 | G4NistManager* nist = G4NistManager::Instance();
|
|---|
| 159 |
|
|---|
| 160 | G4double csup, csdn;
|
|---|
| 161 | G4int A;
|
|---|
| 162 |
|
|---|
| 163 | if(verboseLevel > 0) G4cout << "### G4BGGPionElasticXS::Initialise for "
|
|---|
| 164 | << particle->GetParticleName() << G4endl;
|
|---|
| 165 |
|
|---|
| 166 | for(G4int iz=2; iz<93; iz++) {
|
|---|
| 167 |
|
|---|
| 168 | G4double Z = G4double(iz);
|
|---|
| 169 | A = G4lrint(nist->GetAtomicMassAmu(iz));
|
|---|
| 170 |
|
|---|
| 171 | csup = fGlauber->GetElasticGlauberGribov(&dp, iz, A);
|
|---|
| 172 | csdn = fPion->GetElasticCrossSection(&dp, iz, A);
|
|---|
| 173 |
|
|---|
| 174 | theGlauberFac[iz] = csdn/csup;
|
|---|
| 175 | if(verboseLevel > 0) G4cout << "Z= " << Z << " A= " << A
|
|---|
| 176 | << " factor= " << theGlauberFac[iz] << G4endl;
|
|---|
| 177 | }
|
|---|
| 178 | dp.SetKineticEnergy(fLowEnergy);
|
|---|
| 179 | fHadron->GetHadronNucleonXscNS(&dp, G4Proton::Proton());
|
|---|
| 180 | theCoulombFac[1] = fHadron->GetElasticHadronNucleonXsc();
|
|---|
| 181 | if(isPiplus) { theCoulombFac[1] /= CoulombFactor(fLowEnergy,1); }
|
|---|
| 182 |
|
|---|
| 183 | for(G4int iz=2; iz<93; iz++) {
|
|---|
| 184 |
|
|---|
| 185 | G4double Z = G4double(iz);
|
|---|
| 186 | A = G4lrint(nist->GetAtomicMassAmu(iz));
|
|---|
| 187 |
|
|---|
| 188 | theCoulombFac[iz] = fPion->GetElasticCrossSection(&dp, iz, A);
|
|---|
| 189 | if(isPiplus) { theCoulombFac[iz] /= CoulombFactor(fLowEnergy,A); }
|
|---|
| 190 | if(verboseLevel > 0) G4cout << "Z= " << Z << " A= " << A
|
|---|
| 191 | << " factor= " << theCoulombFac[iz] << G4endl;
|
|---|
| 192 | }
|
|---|
| 193 | }
|
|---|
| 194 |
|
|---|
| 195 |
|
|---|
| 196 | G4double G4BGGPionElasticXS::CoulombFactor(G4double kinEnergy, G4int A)
|
|---|
| 197 | {
|
|---|
| 198 | G4double res= 0.0;
|
|---|
| 199 | if(kinEnergy <= DBL_MIN) return res;
|
|---|
| 200 | else if(A < 2) return kinEnergy*kinEnergy;
|
|---|
| 201 |
|
|---|
| 202 | G4double elog = std::log10(kinEnergy/GeV);
|
|---|
| 203 | G4double aa = A;
|
|---|
| 204 |
|
|---|
| 205 | // from G4ProtonInelasticCrossSection
|
|---|
| 206 | G4double f1 = 8.0 - 8.0/aa - 0.008*aa;
|
|---|
| 207 | G4double f2 = 2.34 - 5.4/aa - 0.0028*aa;
|
|---|
| 208 |
|
|---|
| 209 | res = 1.0/(1.0 + std::exp(-f1*(elog + f2)));
|
|---|
| 210 |
|
|---|
| 211 | f1 = 5.6 - 0.016*aa;
|
|---|
| 212 | f2 = 1.37 + 1.37/aa;
|
|---|
| 213 | res *= ( 1.0 + (0.8 + 18./aa - 0.002*aa)/(1.0 + std::exp(f1*(elog + f2))));
|
|---|
| 214 | return res;
|
|---|
| 215 | }
|
|---|
| 216 |
|
|---|