[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 | // $Id: G4UrbanMscModel90.hh,v 1.1 2007/12/07 17:35:52 vnivanch Exp $ |
---|
| 27 | // GEANT4 tag $Name: geant4-09-01-patch-02 $ |
---|
| 28 | // |
---|
| 29 | // ------------------------------------------------------------------- |
---|
| 30 | // |
---|
| 31 | // |
---|
| 32 | // GEANT4 Class header file |
---|
| 33 | // |
---|
| 34 | // |
---|
| 35 | // File name: G4UrbanMscModel90 |
---|
| 36 | // |
---|
| 37 | // Author: V.Ivanchenko clone Laszlo Urban model |
---|
| 38 | // |
---|
| 39 | // Creation date: 07.12.2007 |
---|
| 40 | // |
---|
| 41 | // Modifications: |
---|
| 42 | // |
---|
| 43 | // |
---|
| 44 | // Class Description: |
---|
| 45 | // |
---|
| 46 | // Implementation of the model of multiple scattering based on |
---|
| 47 | // H.W.Lewis Phys Rev 78 (1950) 526 and L.Urban model |
---|
| 48 | |
---|
| 49 | // ------------------------------------------------------------------- |
---|
| 50 | // |
---|
| 51 | |
---|
| 52 | #ifndef G4UrbanMscModel90_h |
---|
| 53 | #define G4UrbanMscModel90_h 1 |
---|
| 54 | |
---|
| 55 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 56 | |
---|
| 57 | #include "G4VEmModel.hh" |
---|
| 58 | #include "G4PhysicsTable.hh" |
---|
| 59 | #include "G4MscStepLimitType.hh" |
---|
| 60 | |
---|
| 61 | class G4ParticleChangeForMSC; |
---|
| 62 | class G4SafetyHelper; |
---|
| 63 | class G4LossTableManager; |
---|
| 64 | |
---|
| 65 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 66 | |
---|
| 67 | class G4UrbanMscModel90 : public G4VEmModel |
---|
| 68 | { |
---|
| 69 | |
---|
| 70 | public: |
---|
| 71 | |
---|
| 72 | G4UrbanMscModel90(G4double facrange, G4double dtrl, G4double lambdalimit, |
---|
| 73 | G4double facgeom,G4double skin, |
---|
| 74 | G4bool samplez, G4MscStepLimitType stepAlg, |
---|
| 75 | const G4String& nam = "UrbanMscUni"); |
---|
| 76 | |
---|
| 77 | virtual ~G4UrbanMscModel90(); |
---|
| 78 | |
---|
| 79 | virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&); |
---|
| 80 | |
---|
| 81 | G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition* particle, |
---|
| 82 | G4double KineticEnergy, |
---|
| 83 | G4double AtomicNumber, |
---|
| 84 | G4double AtomicWeight=0., |
---|
| 85 | G4double cut =0., |
---|
| 86 | G4double emax=DBL_MAX); |
---|
| 87 | |
---|
| 88 | void SampleSecondaries(std::vector<G4DynamicParticle*>*, |
---|
| 89 | const G4MaterialCutsCouple*, |
---|
| 90 | const G4DynamicParticle*, |
---|
| 91 | G4double, |
---|
| 92 | G4double); |
---|
| 93 | |
---|
| 94 | void SampleScattering(const G4DynamicParticle*, |
---|
| 95 | G4double safety); |
---|
| 96 | |
---|
| 97 | G4double ComputeTruePathLengthLimit(const G4Track& track, |
---|
| 98 | G4PhysicsTable* theLambdaTable, |
---|
| 99 | G4double currentMinimalStep); |
---|
| 100 | |
---|
| 101 | G4double ComputeGeomPathLength(G4double truePathLength); |
---|
| 102 | |
---|
| 103 | G4double ComputeTrueStepLength(G4double geomStepLength); |
---|
| 104 | |
---|
| 105 | G4double ComputeTheta0(G4double truePathLength, |
---|
| 106 | G4double KineticEnergy); |
---|
| 107 | |
---|
| 108 | void SetStepLimitType(G4MscStepLimitType); |
---|
| 109 | |
---|
| 110 | void SetLateralDisplasmentFlag(G4bool val); |
---|
| 111 | |
---|
| 112 | void SetRangeFactor(G4double); |
---|
| 113 | |
---|
| 114 | void SetGeomFactor(G4double); |
---|
| 115 | |
---|
| 116 | void SetSkin(G4double); |
---|
| 117 | |
---|
| 118 | private: |
---|
| 119 | |
---|
| 120 | G4double SampleCosineTheta(G4double trueStepLength, G4double KineticEnergy); |
---|
| 121 | |
---|
| 122 | G4double SampleDisplacement(); |
---|
| 123 | |
---|
| 124 | G4double LatCorrelation(); |
---|
| 125 | |
---|
| 126 | G4double GetLambda(G4double kinEnergy); |
---|
| 127 | |
---|
| 128 | void GeomLimit(const G4Track& track); |
---|
| 129 | |
---|
| 130 | void SetParticle(const G4ParticleDefinition* p); |
---|
| 131 | |
---|
| 132 | // hide assignment operator |
---|
| 133 | G4UrbanMscModel90 & operator=(const G4UrbanMscModel90 &right); |
---|
| 134 | G4UrbanMscModel90(const G4UrbanMscModel90&); |
---|
| 135 | |
---|
| 136 | const G4ParticleDefinition* particle; |
---|
| 137 | G4ParticleChangeForMSC* fParticleChange; |
---|
| 138 | |
---|
| 139 | G4SafetyHelper* safetyHelper; |
---|
| 140 | G4PhysicsTable* theLambdaTable; |
---|
| 141 | const G4MaterialCutsCouple* couple; |
---|
| 142 | G4LossTableManager* theManager; |
---|
| 143 | |
---|
| 144 | |
---|
| 145 | G4double mass; |
---|
| 146 | G4double charge; |
---|
| 147 | |
---|
| 148 | G4double masslimite,masslimitmu; |
---|
| 149 | |
---|
| 150 | G4double taubig; |
---|
| 151 | G4double tausmall; |
---|
| 152 | G4double taulim; |
---|
| 153 | G4double currentTau; |
---|
| 154 | G4double dtrl; |
---|
| 155 | |
---|
| 156 | G4double lambdalimit; |
---|
| 157 | G4double facrange; |
---|
| 158 | G4double frscaling1,frscaling2; |
---|
| 159 | G4double tlimit; |
---|
| 160 | G4double tlimitmin; |
---|
| 161 | G4double tlimitminfix; |
---|
| 162 | |
---|
| 163 | G4double nstepmax; |
---|
| 164 | G4double geombig; |
---|
| 165 | G4double geommin; |
---|
| 166 | G4double geomlimit; |
---|
| 167 | G4double facgeom; |
---|
| 168 | G4double skin; |
---|
| 169 | G4double skindepth; |
---|
| 170 | G4double smallstep; |
---|
| 171 | |
---|
| 172 | G4double presafety; |
---|
| 173 | G4double facsafety; |
---|
| 174 | |
---|
| 175 | G4double lambda0; |
---|
| 176 | G4double lambdaeff; |
---|
| 177 | G4double tPathLength; |
---|
| 178 | G4double zPathLength; |
---|
| 179 | G4double par1,par2,par3 ; |
---|
| 180 | |
---|
| 181 | G4double stepmin ; |
---|
| 182 | |
---|
| 183 | G4double currentKinEnergy; |
---|
| 184 | G4double currentRange; |
---|
| 185 | G4double currentRadLength; |
---|
| 186 | |
---|
| 187 | G4double Zeff; |
---|
| 188 | |
---|
| 189 | G4int currentMaterialIndex; |
---|
| 190 | |
---|
| 191 | G4MscStepLimitType steppingAlgorithm; |
---|
| 192 | |
---|
| 193 | G4bool samplez; |
---|
| 194 | G4bool latDisplasment; |
---|
| 195 | G4bool isInitialized; |
---|
| 196 | |
---|
| 197 | G4bool inside; |
---|
| 198 | G4bool insideskin; |
---|
| 199 | |
---|
| 200 | }; |
---|
| 201 | |
---|
| 202 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 203 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 204 | |
---|
| 205 | inline |
---|
| 206 | void G4UrbanMscModel90::SetLateralDisplasmentFlag(G4bool val) |
---|
| 207 | { |
---|
| 208 | latDisplasment = val; |
---|
| 209 | } |
---|
| 210 | |
---|
| 211 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 212 | |
---|
| 213 | inline |
---|
| 214 | void G4UrbanMscModel90::SetSkin(G4double val) |
---|
| 215 | { |
---|
| 216 | skin = val; |
---|
| 217 | stepmin = tlimitminfix; |
---|
| 218 | skindepth = skin*stepmin; |
---|
| 219 | } |
---|
| 220 | |
---|
| 221 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 222 | |
---|
| 223 | inline |
---|
| 224 | void G4UrbanMscModel90::SetRangeFactor(G4double val) |
---|
| 225 | { |
---|
| 226 | facrange = val; |
---|
| 227 | } |
---|
| 228 | |
---|
| 229 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 230 | |
---|
| 231 | inline |
---|
| 232 | void G4UrbanMscModel90::SetGeomFactor(G4double val) |
---|
| 233 | { |
---|
| 234 | facgeom = val; |
---|
| 235 | } |
---|
| 236 | |
---|
| 237 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 238 | |
---|
| 239 | inline |
---|
| 240 | void G4UrbanMscModel90::SetStepLimitType(G4MscStepLimitType val) |
---|
| 241 | { |
---|
| 242 | steppingAlgorithm = val; |
---|
| 243 | } |
---|
| 244 | |
---|
| 245 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 246 | |
---|
| 247 | inline |
---|
| 248 | G4double G4UrbanMscModel90::GetLambda(G4double e) |
---|
| 249 | { |
---|
| 250 | G4double x; |
---|
| 251 | if(theLambdaTable) { |
---|
| 252 | G4bool b; |
---|
| 253 | x = ((*theLambdaTable)[currentMaterialIndex])->GetValue(e, b); |
---|
| 254 | } else { |
---|
| 255 | x = CrossSection(couple,particle,e); |
---|
| 256 | } |
---|
| 257 | if(x > DBL_MIN) x = 1./x; |
---|
| 258 | else x = DBL_MAX; |
---|
| 259 | return x; |
---|
| 260 | } |
---|
| 261 | |
---|
| 262 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 263 | |
---|
| 264 | inline |
---|
| 265 | void G4UrbanMscModel90::SetParticle(const G4ParticleDefinition* p) |
---|
| 266 | { |
---|
| 267 | if (p != particle) { |
---|
| 268 | particle = p; |
---|
| 269 | mass = p->GetPDGMass(); |
---|
| 270 | charge = p->GetPDGCharge()/eplus; |
---|
| 271 | } |
---|
| 272 | } |
---|
| 273 | |
---|
| 274 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 275 | |
---|
| 276 | #endif |
---|
| 277 | |
---|