| [966] | 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 | //
|
|---|
| [1007] | 26 | // $Id: G4eBremsstrahlungRelModel.hh,v 1.10 2008/11/14 09:25:19 vnivanch Exp $
|
|---|
| 27 | // GEANT4 tag $Name: geant4-09-02 $
|
|---|
| [966] | 28 | //
|
|---|
| 29 | // -------------------------------------------------------------------
|
|---|
| 30 | //
|
|---|
| 31 | // GEANT4 Class header file
|
|---|
| 32 | //
|
|---|
| 33 | //
|
|---|
| 34 | // File name: G4eBremsstrahlungRelModel
|
|---|
| 35 | // extention of standard G4eBremsstrahlungModel
|
|---|
| 36 | //
|
|---|
| 37 | // Author: Andreas Schaelicke
|
|---|
| 38 | //
|
|---|
| 39 | // Creation date: 28.03.2008
|
|---|
| 40 | //
|
|---|
| 41 | // Modifications:
|
|---|
| 42 | //
|
|---|
| 43 | //
|
|---|
| 44 | // Class Description:
|
|---|
| 45 | //
|
|---|
| 46 | // Implementation of energy loss for gamma emission by electrons and
|
|---|
| 47 | // positrons including an improved version of the LPM effect
|
|---|
| 48 |
|
|---|
| 49 | // -------------------------------------------------------------------
|
|---|
| 50 | //
|
|---|
| 51 |
|
|---|
| 52 | #ifndef G4eBremsstrahlungRelModel_h
|
|---|
| 53 | #define G4eBremsstrahlungRelModel_h 1
|
|---|
| 54 |
|
|---|
| 55 | #include "G4VEmModel.hh"
|
|---|
| 56 | #include "G4NistManager.hh"
|
|---|
| 57 |
|
|---|
| 58 | class G4ParticleChangeForLoss;
|
|---|
| 59 | class G4PhysicsVector;
|
|---|
| 60 |
|
|---|
| 61 | class G4eBremsstrahlungRelModel : public G4VEmModel
|
|---|
| 62 | {
|
|---|
| 63 |
|
|---|
| 64 | public:
|
|---|
| 65 |
|
|---|
| 66 | G4eBremsstrahlungRelModel(const G4ParticleDefinition* p = 0,
|
|---|
| 67 | const G4String& nam = "eBremRel");
|
|---|
| 68 |
|
|---|
| 69 | virtual ~G4eBremsstrahlungRelModel();
|
|---|
| 70 |
|
|---|
| 71 | virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
|
|---|
| 72 |
|
|---|
| [1007] | 73 | G4double MinEnergyCut(const G4ParticleDefinition*,
|
|---|
| 74 | const G4MaterialCutsCouple*);
|
|---|
| [966] | 75 |
|
|---|
| 76 | virtual G4double ComputeDEDXPerVolume(const G4Material*,
|
|---|
| 77 | const G4ParticleDefinition*,
|
|---|
| 78 | G4double kineticEnergy,
|
|---|
| 79 | G4double cutEnergy);
|
|---|
| 80 |
|
|---|
| 81 | virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
|
|---|
| 82 | G4double tkin,
|
|---|
| 83 | G4double Z, G4double,
|
|---|
| 84 | G4double cutEnergy,
|
|---|
| 85 | G4double maxEnergy = DBL_MAX);
|
|---|
| 86 |
|
|---|
| 87 | virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
|
|---|
| 88 | const G4MaterialCutsCouple*,
|
|---|
| 89 | const G4DynamicParticle*,
|
|---|
| 90 | G4double cutEnergy,
|
|---|
| 91 | G4double maxEnergy);
|
|---|
| 92 |
|
|---|
| 93 | virtual void SetupForMaterial(const G4ParticleDefinition*,
|
|---|
| 94 | const G4Material*,G4double);
|
|---|
| 95 |
|
|---|
| 96 | inline void SetLPMconstant(G4double val);
|
|---|
| 97 | inline G4double LPMconstant() const;
|
|---|
| 98 |
|
|---|
| [1007] | 99 | protected:
|
|---|
| 100 |
|
|---|
| 101 | inline G4double MaxSecondaryEnergy(const G4ParticleDefinition*,
|
|---|
| 102 | G4double kineticEnergy);
|
|---|
| 103 |
|
|---|
| [966] | 104 | private:
|
|---|
| 105 |
|
|---|
| 106 | void InitialiseConstants();
|
|---|
| 107 |
|
|---|
| 108 | void CalcLPMFunctions(G4double gammaEnergy);
|
|---|
| 109 |
|
|---|
| 110 | G4double ComputeBremLoss(G4double cutEnergy);
|
|---|
| 111 |
|
|---|
| 112 | G4double ComputeXSectionPerAtom(G4double cutEnergy);
|
|---|
| 113 |
|
|---|
| 114 | G4double ComputeDXSectionPerAtom(G4double gammaEnergy);
|
|---|
| 115 |
|
|---|
| 116 | G4double ComputeRelDXSectionPerAtom(G4double gammaEnergy);
|
|---|
| 117 |
|
|---|
| 118 | void SetParticle(const G4ParticleDefinition* p);
|
|---|
| 119 |
|
|---|
| 120 | // * fast inline functions *
|
|---|
| 121 | inline void SetCurrentElement(const G4double);
|
|---|
| 122 |
|
|---|
| 123 | inline G4double Phi1(G4double,G4double);
|
|---|
| 124 | inline G4double Phi1M2(G4double,G4double);
|
|---|
| 125 | inline G4double Psi1(G4double,G4double);
|
|---|
| 126 | inline G4double Psi1M2(G4double,G4double);
|
|---|
| 127 |
|
|---|
| 128 | // hide assignment operator
|
|---|
| 129 | G4eBremsstrahlungRelModel & operator=(const G4eBremsstrahlungRelModel &right);
|
|---|
| 130 | G4eBremsstrahlungRelModel(const G4eBremsstrahlungRelModel&);
|
|---|
| 131 |
|
|---|
| 132 | protected:
|
|---|
| 133 |
|
|---|
| 134 | G4NistManager* nist;
|
|---|
| 135 | const G4ParticleDefinition* particle;
|
|---|
| 136 | G4ParticleDefinition* theGamma;
|
|---|
| 137 | G4ParticleChangeForLoss* fParticleChange;
|
|---|
| 138 |
|
|---|
| 139 | static const G4double xgi[8], wgi[8];
|
|---|
| 140 | static const G4double Fel_light[5];
|
|---|
| 141 | static const G4double Finel_light[5];
|
|---|
| 142 |
|
|---|
| 143 | G4double minThreshold;
|
|---|
| 144 |
|
|---|
| 145 | // cash
|
|---|
| 146 | G4double particleMass;
|
|---|
| 147 | G4double kinEnergy;
|
|---|
| 148 | G4double totalEnergy;
|
|---|
| 149 | G4double currentZ;
|
|---|
| 150 | G4double z13, z23, lnZ;
|
|---|
| 151 | G4double Fel, Finel, fCoulomb, fMax;
|
|---|
| 152 | G4double densityFactor;
|
|---|
| 153 | G4double densityCorr;
|
|---|
| 154 |
|
|---|
| 155 | // LPM effect
|
|---|
| 156 | G4double lpmEnergy;
|
|---|
| 157 | G4PhysicsVector *fXiLPM, *fPhiLPM, *fGLPM;
|
|---|
| 158 | G4double xiLPM, phiLPM, gLPM;
|
|---|
| 159 |
|
|---|
| 160 | // critical gamma energies
|
|---|
| 161 | G4double klpm, kp;
|
|---|
| 162 |
|
|---|
| 163 | G4bool isElectron;
|
|---|
| 164 |
|
|---|
| 165 | private:
|
|---|
| 166 | // consts
|
|---|
| 167 | G4double highKinEnergy;
|
|---|
| 168 | G4double lowKinEnergy;
|
|---|
| 169 | G4double fMigdalConstant;
|
|---|
| 170 | G4double fLPMconstant;
|
|---|
| 171 | G4double bremFactor;
|
|---|
| 172 | G4double energyThresholdLPM;
|
|---|
| 173 | G4double facFel, facFinel;
|
|---|
| 174 | G4double preS1,logTwo;
|
|---|
| 175 |
|
|---|
| 176 | G4bool use_completescreening;
|
|---|
| 177 |
|
|---|
| 178 | G4bool isInitialised;
|
|---|
| 179 | };
|
|---|
| 180 |
|
|---|
| 181 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 182 |
|
|---|
| 183 | inline void G4eBremsstrahlungRelModel::SetCurrentElement(const G4double Z)
|
|---|
| 184 | {
|
|---|
| 185 | if(Z != currentZ) {
|
|---|
| 186 | currentZ = Z;
|
|---|
| 187 |
|
|---|
| 188 | G4int iz = G4int(Z);
|
|---|
| 189 | z13 = nist->GetZ13(iz);
|
|---|
| 190 | z23 = z13*z13;
|
|---|
| 191 | lnZ = nist->GetLOGZ(iz);
|
|---|
| 192 |
|
|---|
| 193 | if (iz <= 4) {
|
|---|
| 194 | Fel = Fel_light[iz];
|
|---|
| 195 | Finel = Finel_light[iz] ;
|
|---|
| 196 | }
|
|---|
| 197 | else {
|
|---|
| 198 | Fel = facFel - lnZ/3. ;
|
|---|
| 199 | Finel = facFinel - 2.*lnZ/3. ;
|
|---|
| 200 | }
|
|---|
| 201 |
|
|---|
| 202 | fCoulomb=GetCurrentElement()->GetfCoulomb();
|
|---|
| 203 | fMax = Fel-fCoulomb + Finel/currentZ + (1.+1./currentZ)/12.;
|
|---|
| 204 | }
|
|---|
| 205 | }
|
|---|
| 206 |
|
|---|
| 207 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 208 |
|
|---|
| [1007] | 209 | inline G4double
|
|---|
| 210 | G4eBremsstrahlungRelModel::MaxSecondaryEnergy(const G4ParticleDefinition*,
|
|---|
| 211 | G4double kineticEnergy)
|
|---|
| 212 | {
|
|---|
| 213 | return kineticEnergy;
|
|---|
| 214 | }
|
|---|
| [966] | 215 |
|
|---|
| [1007] | 216 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 217 |
|
|---|
| 218 |
|
|---|
| [966] | 219 | inline G4double G4eBremsstrahlungRelModel::Phi1(G4double gg, G4double)
|
|---|
| 220 | {
|
|---|
| 221 | // Thomas-Fermi FF from Tsai, eq.(3.38) for Z>=5
|
|---|
| 222 | return 20.863 - 2.*std::log(1. + sqr(0.55846*gg) )
|
|---|
| 223 | - 4.*( 1. - 0.6*std::exp(-0.9*gg) - 0.4*std::exp(-1.5*gg) );
|
|---|
| 224 | }
|
|---|
| 225 |
|
|---|
| 226 | inline G4double G4eBremsstrahlungRelModel::Phi1M2(G4double gg, G4double)
|
|---|
| 227 | {
|
|---|
| 228 | // Thomas-Fermi FF from Tsai, eq. (3.39) for Z>=5
|
|---|
| 229 | // return Phi1(gg,Z) -
|
|---|
| 230 | return 2./(3.*(1. + 6.5*gg +6.*gg*gg) );
|
|---|
| 231 | }
|
|---|
| 232 |
|
|---|
| 233 | inline G4double G4eBremsstrahlungRelModel::Psi1(G4double eps, G4double)
|
|---|
| 234 | {
|
|---|
| 235 | // Thomas-Fermi FF from Tsai, eq.(3.40) for Z>=5
|
|---|
| 236 | return 28.340 - 2.*std::log(1. + sqr(3.621*eps) )
|
|---|
| 237 | - 4.*( 1. - 0.7*std::exp(-8*eps) - 0.3*std::exp(-29.*eps) );
|
|---|
| 238 | }
|
|---|
| 239 |
|
|---|
| 240 | inline G4double G4eBremsstrahlungRelModel::Psi1M2(G4double eps, G4double)
|
|---|
| 241 | {
|
|---|
| 242 | // Thomas-Fermi FF from Tsai, eq. (3.41) for Z>=5
|
|---|
| 243 | return 2./(3.*(1. + 40.*eps +400.*eps*eps) );
|
|---|
| 244 | }
|
|---|
| 245 |
|
|---|
| 246 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 247 |
|
|---|
| 248 | inline
|
|---|
| 249 | void G4eBremsstrahlungRelModel::SetLPMconstant(G4double val)
|
|---|
| 250 | {
|
|---|
| 251 | fLPMconstant = val;
|
|---|
| 252 | }
|
|---|
| 253 |
|
|---|
| 254 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 255 |
|
|---|
| 256 | inline
|
|---|
| 257 | G4double G4eBremsstrahlungRelModel::LPMconstant() const
|
|---|
| 258 | {
|
|---|
| 259 | return fLPMconstant;
|
|---|
| 260 | }
|
|---|
| 261 |
|
|---|
| 262 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 263 |
|
|---|
| 264 |
|
|---|
| 265 | #endif
|
|---|