[1316] | 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 | // |
---|
[1337] | 26 | // $Id: G4WentzelOKandVIxSection.hh,v 1.6 2010/06/25 09:41:38 gunter Exp $ |
---|
| 27 | // GEANT4 tag $Name: geant4-09-04-beta-01 $ |
---|
[1316] | 28 | // |
---|
| 29 | // ------------------------------------------------------------------- |
---|
| 30 | // |
---|
| 31 | // |
---|
| 32 | // GEANT4 Class header file |
---|
| 33 | // |
---|
| 34 | // |
---|
| 35 | // File name: G4WentzelOKandVIxSection |
---|
| 36 | // |
---|
| 37 | // Authors: V.Ivanchenko and O.Kadri |
---|
| 38 | // |
---|
| 39 | // Creation date: 21.05.2010 |
---|
| 40 | // |
---|
| 41 | // Modifications: |
---|
| 42 | // |
---|
| 43 | // |
---|
| 44 | // Class Description: |
---|
| 45 | // |
---|
| 46 | // Implementation of the computation of total and transport cross sections, |
---|
| 47 | // sample scattering angle for the single scattering case. |
---|
| 48 | // to be used by single and multiple scattering models. References: |
---|
| 49 | // 1) G.Wentzel, Z. Phys. 40 (1927) 590. |
---|
| 50 | // 2) J.M. Fernandez-Varea et al., NIM B73 (1993) 447. |
---|
| 51 | // |
---|
| 52 | // ------------------------------------------------------------------- |
---|
| 53 | // |
---|
| 54 | |
---|
| 55 | #ifndef G4WentzelOKandVIxSection_h |
---|
| 56 | #define G4WentzelOKandVIxSection_h 1 |
---|
| 57 | |
---|
| 58 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 59 | |
---|
| 60 | #include "globals.hh" |
---|
| 61 | #include "G4Material.hh" |
---|
| 62 | #include "G4Element.hh" |
---|
| 63 | #include "G4ElementVector.hh" |
---|
| 64 | #include "G4NistManager.hh" |
---|
| 65 | #include "G4ThreeVector.hh" |
---|
| 66 | #include "G4Pow.hh" |
---|
| 67 | |
---|
| 68 | class G4ParticleDefinition; |
---|
| 69 | |
---|
| 70 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 71 | |
---|
| 72 | class G4WentzelOKandVIxSection |
---|
| 73 | { |
---|
| 74 | |
---|
| 75 | public: |
---|
| 76 | |
---|
| 77 | G4WentzelOKandVIxSection(); |
---|
| 78 | |
---|
| 79 | virtual ~G4WentzelOKandVIxSection(); |
---|
| 80 | |
---|
| 81 | void Initialise(const G4ParticleDefinition*, G4double CosThetaLim); |
---|
| 82 | |
---|
| 83 | void SetupParticle(const G4ParticleDefinition*); |
---|
| 84 | |
---|
| 85 | // return cos(ThetaMax) for msc and cos(thetaMin) for single scattering |
---|
| 86 | // cut = DBL_MAX means no scattering off electrons |
---|
| 87 | G4double SetupTarget(G4int Z, G4double cut = DBL_MAX); |
---|
| 88 | |
---|
| 89 | G4double ComputeTransportCrossSectionPerAtom(G4double CosThetaMax); |
---|
| 90 | |
---|
| 91 | G4ThreeVector SampleSingleScattering(G4double CosThetaMin, |
---|
| 92 | G4double CosThetaMax, |
---|
| 93 | G4double elecRatio = 0.0); |
---|
| 94 | |
---|
| 95 | inline G4double ComputeNuclearCrossSection(G4double CosThetaMin, |
---|
| 96 | G4double CosThetaMax); |
---|
| 97 | |
---|
| 98 | inline G4double ComputeElectronCrossSection(G4double CosThetaMin, |
---|
| 99 | G4double CosThetaMax); |
---|
| 100 | |
---|
| 101 | inline G4double SetupKinematic(G4double kinEnergy, const G4Material* mat); |
---|
| 102 | |
---|
| 103 | inline void SetTargetMass(G4double value); |
---|
| 104 | |
---|
| 105 | inline void SetRelativisticMass(G4double value); |
---|
| 106 | |
---|
| 107 | inline G4double GetMomentumSquare() const; |
---|
| 108 | |
---|
| 109 | inline G4double GetCosThetaNuc() const; |
---|
| 110 | |
---|
| 111 | inline G4double GetCosThetaElec() const; |
---|
| 112 | |
---|
| 113 | private: |
---|
| 114 | |
---|
| 115 | void ComputeMaxElectronScattering(G4double cut); |
---|
| 116 | |
---|
| 117 | // hide assignment operator |
---|
| 118 | G4WentzelOKandVIxSection & operator=(const G4WentzelOKandVIxSection &right); |
---|
| 119 | G4WentzelOKandVIxSection(const G4WentzelOKandVIxSection&); |
---|
| 120 | |
---|
| 121 | const G4ParticleDefinition* theProton; |
---|
| 122 | const G4ParticleDefinition* theElectron; |
---|
| 123 | const G4ParticleDefinition* thePositron; |
---|
| 124 | const G4Material* currentMaterial; |
---|
| 125 | |
---|
| 126 | G4NistManager* fNistManager; |
---|
| 127 | G4Pow* fG4pow; |
---|
| 128 | |
---|
| 129 | G4double numlimit; |
---|
| 130 | |
---|
| 131 | G4double elecXSRatio; |
---|
| 132 | |
---|
| 133 | // integer parameters |
---|
| 134 | G4int nwarnings; |
---|
| 135 | G4int nwarnlimit; |
---|
| 136 | |
---|
| 137 | // single scattering parameters |
---|
| 138 | G4double coeff; |
---|
| 139 | G4double cosTetMaxElec; |
---|
| 140 | G4double cosTetMaxNuc; |
---|
| 141 | G4double cosThetaMax; |
---|
| 142 | G4double alpha2; |
---|
| 143 | |
---|
| 144 | // projectile |
---|
| 145 | const G4ParticleDefinition* particle; |
---|
| 146 | |
---|
| 147 | G4double chargeSquare; |
---|
| 148 | G4double charge3; |
---|
| 149 | G4double spin; |
---|
| 150 | G4double mass; |
---|
| 151 | G4double tkin; |
---|
| 152 | G4double mom2; |
---|
| 153 | G4double invbeta2; |
---|
| 154 | G4double kinFactor; |
---|
| 155 | G4double etag; |
---|
| 156 | G4double ecut; |
---|
| 157 | G4double lowEnergyLimit; |
---|
| 158 | |
---|
| 159 | // target |
---|
| 160 | G4int targetZ; |
---|
| 161 | G4double targetMass; |
---|
| 162 | G4double screenZ; |
---|
| 163 | G4double formfactA; |
---|
| 164 | G4double factorA2; |
---|
| 165 | G4double factB; |
---|
| 166 | G4double factD; |
---|
| 167 | |
---|
| 168 | static G4double ScreenRSquare[100]; |
---|
| 169 | static G4double FormFactor[100]; |
---|
| 170 | }; |
---|
| 171 | |
---|
| 172 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 173 | |
---|
| 174 | inline G4double |
---|
| 175 | G4WentzelOKandVIxSection::SetupKinematic(G4double ekin, const G4Material* mat) |
---|
| 176 | { |
---|
| 177 | if(ekin != tkin || mat != currentMaterial) { |
---|
| 178 | currentMaterial = mat; |
---|
| 179 | tkin = ekin; |
---|
| 180 | mom2 = tkin*(tkin + 2.0*mass); |
---|
| 181 | invbeta2 = 1.0 + mass*mass/mom2; |
---|
| 182 | factB = spin/invbeta2; |
---|
| 183 | cosTetMaxNuc = cosThetaMax; |
---|
| 184 | if(std::fabs(cosThetaMax) < 1.0) { |
---|
| 185 | cosTetMaxNuc = |
---|
| 186 | std::max(cosThetaMax,1.-factorA2*mat->GetIonisation()->GetInvA23()/mom2); |
---|
| 187 | } |
---|
| 188 | } |
---|
| 189 | return cosTetMaxNuc; |
---|
| 190 | } |
---|
| 191 | |
---|
| 192 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 193 | |
---|
| 194 | inline void G4WentzelOKandVIxSection::SetTargetMass(G4double value) |
---|
| 195 | { |
---|
| 196 | targetMass = value; |
---|
[1337] | 197 | factD = std::sqrt(mom2)/value; |
---|
[1316] | 198 | } |
---|
| 199 | |
---|
| 200 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 201 | |
---|
| 202 | inline void G4WentzelOKandVIxSection::SetRelativisticMass(G4double value) |
---|
| 203 | { |
---|
| 204 | mass = value; |
---|
| 205 | } |
---|
| 206 | |
---|
| 207 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 208 | |
---|
| 209 | inline G4double G4WentzelOKandVIxSection::GetMomentumSquare() const |
---|
| 210 | { |
---|
| 211 | return mom2; |
---|
| 212 | } |
---|
| 213 | |
---|
| 214 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 215 | |
---|
| 216 | inline G4double G4WentzelOKandVIxSection::GetCosThetaNuc() const |
---|
| 217 | { |
---|
| 218 | return cosTetMaxNuc; |
---|
| 219 | } |
---|
| 220 | |
---|
| 221 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 222 | |
---|
| 223 | inline G4double G4WentzelOKandVIxSection::GetCosThetaElec() const |
---|
| 224 | { |
---|
| 225 | return cosTetMaxElec; |
---|
| 226 | } |
---|
| 227 | |
---|
| 228 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 229 | |
---|
| 230 | inline G4double |
---|
| 231 | G4WentzelOKandVIxSection::ComputeNuclearCrossSection(G4double cosTMin, |
---|
| 232 | G4double cosTMax) |
---|
| 233 | { |
---|
| 234 | G4double xsec = 0.0; |
---|
| 235 | if(cosTMax < cosTMin) { |
---|
| 236 | xsec = targetZ*kinFactor*(cosTMin - cosTMax)/ |
---|
| 237 | ((1.0 - cosTMin + screenZ)*(1.0 - cosTMax + screenZ)); |
---|
| 238 | } |
---|
| 239 | return xsec; |
---|
| 240 | } |
---|
| 241 | |
---|
| 242 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 243 | |
---|
| 244 | inline G4double |
---|
| 245 | G4WentzelOKandVIxSection::ComputeElectronCrossSection(G4double cosTMin, |
---|
| 246 | G4double cosTMax) |
---|
| 247 | { |
---|
| 248 | G4double xsec = 0.0; |
---|
| 249 | G4double cost1 = std::max(cosTMin,cosTetMaxElec); |
---|
| 250 | G4double cost2 = std::max(cosTMax,cosTetMaxElec); |
---|
| 251 | if(cost1 > cost2) { |
---|
| 252 | xsec = kinFactor*(cost1 - cost2)/((1.0 - cost1 + screenZ)*(1.0 - cost2 + screenZ)); |
---|
| 253 | } |
---|
| 254 | return xsec; |
---|
| 255 | } |
---|
| 256 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 257 | |
---|
| 258 | #endif |
---|
| 259 | |
---|