[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 | // |
---|
[1055] | 26 | // $Id: G4UrbanMscModel.hh,v 1.35 2009/04/29 13:30:22 vnivanch Exp $ |
---|
| 27 | // GEANT4 tag $Name: geant4-09-03-beta-cand-01 $ |
---|
[819] | 28 | // |
---|
| 29 | // ------------------------------------------------------------------- |
---|
| 30 | // |
---|
| 31 | // |
---|
| 32 | // GEANT4 Class header file |
---|
| 33 | // |
---|
| 34 | // |
---|
| 35 | // File name: G4UrbanMscModel |
---|
| 36 | // |
---|
| 37 | // Author: Laszlo Urban |
---|
| 38 | // |
---|
[1055] | 39 | // Creation date: 06.03.2008 |
---|
[819] | 40 | // |
---|
| 41 | // Modifications: |
---|
| 42 | // |
---|
[1055] | 43 | // 28-04-09 move G4UrbanMscModel2 from the g49.2 to G4UrbanMscModel. |
---|
| 44 | // now it is frozen (V.Ivanchenk0) |
---|
[819] | 45 | // |
---|
| 46 | // Class Description: |
---|
| 47 | // |
---|
| 48 | // Implementation of the model of multiple scattering based on |
---|
| 49 | // H.W.Lewis Phys Rev 78 (1950) 526 and L.Urban model |
---|
| 50 | |
---|
| 51 | // ------------------------------------------------------------------- |
---|
| 52 | // |
---|
| 53 | |
---|
| 54 | #ifndef G4UrbanMscModel_h |
---|
| 55 | #define G4UrbanMscModel_h 1 |
---|
| 56 | |
---|
| 57 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 58 | |
---|
[961] | 59 | #include "G4VMscModel.hh" |
---|
[819] | 60 | #include "G4PhysicsTable.hh" |
---|
| 61 | #include "G4MscStepLimitType.hh" |
---|
| 62 | |
---|
| 63 | class G4ParticleChangeForMSC; |
---|
| 64 | class G4SafetyHelper; |
---|
| 65 | class G4LossTableManager; |
---|
| 66 | |
---|
| 67 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 68 | |
---|
[961] | 69 | class G4UrbanMscModel : public G4VMscModel |
---|
[819] | 70 | { |
---|
| 71 | |
---|
| 72 | public: |
---|
| 73 | |
---|
[961] | 74 | G4UrbanMscModel(const G4String& nam = "UrbanMscUni"); |
---|
[819] | 75 | |
---|
| 76 | virtual ~G4UrbanMscModel(); |
---|
| 77 | |
---|
[961] | 78 | void Initialise(const G4ParticleDefinition*, const G4DataVector&); |
---|
[1055] | 79 | |
---|
[819] | 80 | G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition* particle, |
---|
| 81 | G4double KineticEnergy, |
---|
| 82 | G4double AtomicNumber, |
---|
| 83 | G4double AtomicWeight=0., |
---|
| 84 | G4double cut =0., |
---|
| 85 | G4double emax=DBL_MAX); |
---|
| 86 | |
---|
| 87 | void SampleScattering(const G4DynamicParticle*, |
---|
| 88 | G4double safety); |
---|
| 89 | |
---|
| 90 | G4double ComputeTruePathLengthLimit(const G4Track& track, |
---|
| 91 | G4PhysicsTable* theLambdaTable, |
---|
| 92 | G4double currentMinimalStep); |
---|
| 93 | |
---|
| 94 | G4double ComputeGeomPathLength(G4double truePathLength); |
---|
| 95 | |
---|
| 96 | G4double ComputeTrueStepLength(G4double geomStepLength); |
---|
| 97 | |
---|
| 98 | G4double ComputeTheta0(G4double truePathLength, |
---|
| 99 | G4double KineticEnergy); |
---|
| 100 | |
---|
| 101 | private: |
---|
| 102 | |
---|
[1055] | 103 | G4double SimpleScattering(G4double xmeanth, G4double x2meanth); |
---|
| 104 | |
---|
[819] | 105 | G4double SampleCosineTheta(G4double trueStepLength, G4double KineticEnergy); |
---|
| 106 | |
---|
| 107 | G4double SampleDisplacement(); |
---|
| 108 | |
---|
| 109 | G4double LatCorrelation(); |
---|
| 110 | |
---|
[961] | 111 | inline G4double GetLambda(G4double kinEnergy); |
---|
[819] | 112 | |
---|
[961] | 113 | inline void SetParticle(const G4ParticleDefinition*); |
---|
| 114 | |
---|
[1055] | 115 | inline void UpdateCache(); |
---|
| 116 | |
---|
[819] | 117 | // hide assignment operator |
---|
| 118 | G4UrbanMscModel & operator=(const G4UrbanMscModel &right); |
---|
| 119 | G4UrbanMscModel(const G4UrbanMscModel&); |
---|
| 120 | |
---|
| 121 | const G4ParticleDefinition* particle; |
---|
| 122 | G4ParticleChangeForMSC* fParticleChange; |
---|
| 123 | |
---|
| 124 | G4PhysicsTable* theLambdaTable; |
---|
| 125 | const G4MaterialCutsCouple* couple; |
---|
| 126 | G4LossTableManager* theManager; |
---|
| 127 | |
---|
| 128 | G4double mass; |
---|
[1055] | 129 | G4double charge,ChargeSquare; |
---|
| 130 | G4double masslimite,lambdalimit,fr; |
---|
[819] | 131 | |
---|
| 132 | G4double taubig; |
---|
| 133 | G4double tausmall; |
---|
| 134 | G4double taulim; |
---|
| 135 | G4double currentTau; |
---|
| 136 | G4double tlimit; |
---|
| 137 | G4double tlimitmin; |
---|
| 138 | G4double tlimitminfix; |
---|
[1055] | 139 | G4double tgeom; |
---|
[819] | 140 | |
---|
| 141 | G4double geombig; |
---|
| 142 | G4double geommin; |
---|
| 143 | G4double geomlimit; |
---|
| 144 | G4double skindepth; |
---|
| 145 | G4double smallstep; |
---|
| 146 | |
---|
| 147 | G4double presafety; |
---|
| 148 | |
---|
| 149 | G4double lambda0; |
---|
| 150 | G4double lambdaeff; |
---|
| 151 | G4double tPathLength; |
---|
| 152 | G4double zPathLength; |
---|
[961] | 153 | G4double par1,par2,par3; |
---|
[819] | 154 | |
---|
[961] | 155 | G4double stepmin; |
---|
[819] | 156 | |
---|
| 157 | G4double currentKinEnergy; |
---|
| 158 | G4double currentRange; |
---|
[1055] | 159 | G4double rangeinit; |
---|
[819] | 160 | G4double currentRadLength; |
---|
| 161 | |
---|
[1055] | 162 | G4double theta0max,rellossmax; |
---|
| 163 | G4double third; |
---|
[819] | 164 | |
---|
| 165 | G4int currentMaterialIndex; |
---|
| 166 | |
---|
[1055] | 167 | G4double y; |
---|
| 168 | G4double Zold; |
---|
| 169 | G4double Zeff,Z2,Z23,lnZ; |
---|
| 170 | G4double coeffth1,coeffth2; |
---|
| 171 | G4double coeffc1,coeffc2; |
---|
| 172 | G4double scr1ini,scr2ini,scr1,scr2; |
---|
| 173 | |
---|
[819] | 174 | G4bool isInitialized; |
---|
| 175 | G4bool inside; |
---|
| 176 | G4bool insideskin; |
---|
| 177 | }; |
---|
| 178 | |
---|
| 179 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
[1055] | 180 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
[819] | 181 | |
---|
| 182 | inline |
---|
| 183 | G4double G4UrbanMscModel::GetLambda(G4double e) |
---|
| 184 | { |
---|
| 185 | G4double x; |
---|
| 186 | if(theLambdaTable) { |
---|
| 187 | G4bool b; |
---|
| 188 | x = ((*theLambdaTable)[currentMaterialIndex])->GetValue(e, b); |
---|
| 189 | } else { |
---|
| 190 | x = CrossSection(couple,particle,e); |
---|
| 191 | } |
---|
| 192 | if(x > DBL_MIN) x = 1./x; |
---|
| 193 | else x = DBL_MAX; |
---|
| 194 | return x; |
---|
| 195 | } |
---|
| 196 | |
---|
| 197 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 198 | |
---|
| 199 | inline |
---|
| 200 | void G4UrbanMscModel::SetParticle(const G4ParticleDefinition* p) |
---|
| 201 | { |
---|
| 202 | if (p != particle) { |
---|
| 203 | particle = p; |
---|
| 204 | mass = p->GetPDGMass(); |
---|
| 205 | charge = p->GetPDGCharge()/eplus; |
---|
[1055] | 206 | ChargeSquare = charge*charge; |
---|
[819] | 207 | } |
---|
| 208 | } |
---|
| 209 | |
---|
| 210 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 211 | |
---|
[1055] | 212 | inline |
---|
| 213 | void G4UrbanMscModel::UpdateCache() |
---|
| 214 | { |
---|
| 215 | lnZ = std::log(Zeff); |
---|
| 216 | coeffth1 = 0.885+lnZ*(0.104-0.0170*lnZ); |
---|
| 217 | coeffth2 = 0.028+lnZ*(0.012-0.00125*lnZ); |
---|
| 218 | coeffc1 = 2.134-lnZ*(0.1045-0.00602*lnZ); |
---|
| 219 | coeffc2 = 0.001126-lnZ*(0.0001089+0.0000247*lnZ); |
---|
| 220 | Z2 = Zeff*Zeff; |
---|
| 221 | Z23 = std::exp(2.*lnZ/3.); |
---|
| 222 | scr1 = scr1ini*Z23; |
---|
| 223 | scr2 = scr2ini*Z2*ChargeSquare; |
---|
| 224 | // lastMaterial = couple->GetMaterial(); |
---|
| 225 | Zold = Zeff; |
---|
| 226 | } |
---|
| 227 | |
---|
[819] | 228 | #endif |
---|
| 229 | |
---|