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