Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4UrbanMscModel Class Reference

#include <G4UrbanMscModel.hh>

+ Inheritance diagram for G4UrbanMscModel:

Public Member Functions

 G4UrbanMscModel (const G4String &nam="UrbanMsc")
 
 ~G4UrbanMscModel () override
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
void StartTracking (G4Track *) override
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *particle, G4double KineticEnergy, G4double AtomicNumber, G4double AtomicWeight=0., G4double cut=0., G4double emax=DBL_MAX) override
 
G4ThreeVectorSampleScattering (const G4ThreeVector &, G4double safety) override
 
G4double ComputeTruePathLengthLimit (const G4Track &track, G4double &currentMinimalStep) override
 
G4double ComputeGeomPathLength (G4double truePathLength) override
 
G4double ComputeTrueStepLength (G4double geomStepLength) override
 
G4double ComputeTheta0 (G4double truePathLength, G4double KineticEnergy)
 
void SetDisplacementAlgorithm96 (const G4bool val)
 
void SetPositronCorrection (const G4bool val)
 
G4UrbanMscModeloperator= (const G4UrbanMscModel &right)=delete
 
 G4UrbanMscModel (const G4UrbanMscModel &)=delete
 
- Public Member Functions inherited from G4VMscModel
 G4VMscModel (const G4String &nam)
 
 ~G4VMscModel () override
 
void InitialiseParameters (const G4ParticleDefinition *)
 
void DumpParameters (std::ostream &out) const
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double tmax) override
 
void SetStepLimitType (G4MscStepLimitType)
 
void SetLateralDisplasmentFlag (G4bool val)
 
void SetRangeFactor (G4double)
 
void SetGeomFactor (G4double)
 
void SetSkin (G4double)
 
void SetLambdaLimit (G4double)
 
void SetSafetyFactor (G4double)
 
void SetSampleZ (G4bool)
 
G4VEnergyLossProcessGetIonisation () const
 
void SetIonisation (G4VEnergyLossProcess *, const G4ParticleDefinition *part)
 
G4double ComputeSafety (const G4ThreeVector &position, G4double limit=DBL_MAX)
 
G4double ComputeGeomLimit (const G4Track &, G4double &presafety, G4double limit)
 
G4double GetDEDX (const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
 
G4double GetDEDX (const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple, G4double logKineticEnergy)
 
G4double GetRange (const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
 
G4double GetRange (const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple, G4double logKineticEnergy)
 
G4double GetEnergy (const G4ParticleDefinition *part, G4double range, const G4MaterialCutsCouple *couple)
 
G4double GetTransportMeanFreePath (const G4ParticleDefinition *part, G4double kinEnergy)
 
G4double GetTransportMeanFreePath (const G4ParticleDefinition *part, G4double kinEnergy, G4double logKinEnergy)
 
G4VMscModeloperator= (const G4VMscModel &right)=delete
 
 G4VMscModel (const G4VMscModel &)=delete
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double &length, G4double &eloss)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
void SetLPMFlag (G4bool)
 
void SetMasterThread (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Additional Inherited Members

- Protected Member Functions inherited from G4VMscModel
G4ParticleChangeForMSCGetParticleChangeForMSC (const G4ParticleDefinition *p=nullptr)
 
G4double ConvertTrueToGeom (G4double &tLength, G4double &gLength)
 
void SetUseSplineForMSC (G4bool val)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 
- Protected Attributes inherited from G4VMscModel
G4double facrange = 0.04
 
G4double facgeom = 2.5
 
G4double facsafety = 0.6
 
G4double skin = 1.0
 
G4double dtrl = 0.05
 
G4double lambdalimit
 
G4double geomMin
 
G4double geomMax
 
G4ThreeVector fDisplacement
 
G4MscStepLimitType steppingAlgorithm
 
G4bool samplez = false
 
G4bool latDisplasment = true
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 
const G4MaterialpBaseMaterial = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4double inveplus
 
G4double pFactor = 1.0
 
std::size_t currentCoupleIndex = 0
 
std::size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Detailed Description

Definition at line 69 of file G4UrbanMscModel.hh.

Constructor & Destructor Documentation

◆ G4UrbanMscModel() [1/2]

G4UrbanMscModel::G4UrbanMscModel ( const G4String & nam = "UrbanMsc")
explicit

Definition at line 81 of file G4UrbanMscModel.cc.

82 : G4VMscModel(nam)
83{
84 masslimite = 0.6*CLHEP::MeV;
85 fr = 0.02;
86 taubig = 8.0;
87 tausmall = 1.e-16;
88 taulim = 1.e-6;
89 currentTau = taulim;
90 tlimitminfix = 0.01*CLHEP::nm;
91 tlimitminfix2 = 1.*CLHEP::nm;
92 stepmin = tlimitminfix;
93 smallstep = 1.e10;
94 currentRange = 0.;
95 rangeinit = 0.;
96 tlimit = 1.e10*CLHEP::mm;
97 tlimitmin = 10.*tlimitminfix;
98 tgeom = 1.e50*CLHEP::mm;
99 geombig = tgeom;
100 geommin = 1.e-3*CLHEP::mm;
101 geomlimit = geombig;
102 presafety = 0.;
103
104 positron = G4Positron::Positron();
105 rndmEngineMod = G4Random::getTheEngine();
106
107 drr = 0.35;
108 finalr = 10.*CLHEP::um;
109
110 tlow = 5.*CLHEP::keV;
111 invmev = 1.0/CLHEP::MeV;
112
113 skindepth = skin*stepmin;
114
115 mass = CLHEP::proton_mass_c2;
116 charge = chargeSquare = 1.0;
117 currentKinEnergy = currentRadLength = lambda0 = lambdaeff = tPathLength
118 = zPathLength = par1 = par2 = par3 = rndmarray[0] = rndmarray[1] = 0;
119 currentLogKinEnergy = LOG_EKIN_MIN;
120}
static G4Positron * Positron()
Definition G4Positron.cc:90
G4double skin
G4VMscModel(const G4String &nam)
#define LOG_EKIN_MIN
Definition templates.hh:98

Referenced by G4UrbanMscModel(), and operator=().

◆ ~G4UrbanMscModel()

G4UrbanMscModel::~G4UrbanMscModel ( )
override

Definition at line 124 of file G4UrbanMscModel.cc.

125{
126 if(isFirstInstance) {
127 for(auto const & ptr : msc) { delete ptr; }
128 msc.clear();
129 }
130}

◆ G4UrbanMscModel() [2/2]

G4UrbanMscModel::G4UrbanMscModel ( const G4UrbanMscModel & )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4UrbanMscModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition * particle,
G4double KineticEnergy,
G4double AtomicNumber,
G4double AtomicWeight = 0.,
G4double cut = 0.,
G4double emax = DBL_MAX )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 173 of file G4UrbanMscModel.cc.

178{
179 static const G4double epsmin = 1.e-4 , epsmax = 1.e10;
180
181 static const G4double Zdat[15] = { 4., 6., 13., 20., 26., 29., 32., 38.,47.,
182 50., 56., 64., 74., 79., 82. };
183
184 // corr. factors for e-/e+ lambda for T <= Tlim
185 static const G4double celectron[15][22] =
186 {{1.125,1.072,1.051,1.047,1.047,1.050,1.052,1.054,
187 1.054,1.057,1.062,1.069,1.075,1.090,1.105,1.111,
188 1.112,1.108,1.100,1.093,1.089,1.087 },
189 {1.408,1.246,1.143,1.096,1.077,1.059,1.053,1.051,
190 1.052,1.053,1.058,1.065,1.072,1.087,1.101,1.108,
191 1.109,1.105,1.097,1.090,1.086,1.082 },
192 {2.833,2.268,1.861,1.612,1.486,1.309,1.204,1.156,
193 1.136,1.114,1.106,1.106,1.109,1.119,1.129,1.132,
194 1.131,1.124,1.113,1.104,1.099,1.098 },
195 {3.879,3.016,2.380,2.007,1.818,1.535,1.340,1.236,
196 1.190,1.133,1.107,1.099,1.098,1.103,1.110,1.113,
197 1.112,1.105,1.096,1.089,1.085,1.098 },
198 {6.937,4.330,2.886,2.256,1.987,1.628,1.395,1.265,
199 1.203,1.122,1.080,1.065,1.061,1.063,1.070,1.073,
200 1.073,1.070,1.064,1.059,1.056,1.056 },
201 {9.616,5.708,3.424,2.551,2.204,1.762,1.485,1.330,
202 1.256,1.155,1.099,1.077,1.070,1.068,1.072,1.074,
203 1.074,1.070,1.063,1.059,1.056,1.052 },
204 {11.72,6.364,3.811,2.806,2.401,1.884,1.564,1.386,
205 1.300,1.180,1.112,1.082,1.073,1.066,1.068,1.069,
206 1.068,1.064,1.059,1.054,1.051,1.050 },
207 {18.08,8.601,4.569,3.183,2.662,2.025,1.646,1.439,
208 1.339,1.195,1.108,1.068,1.053,1.040,1.039,1.039,
209 1.039,1.037,1.034,1.031,1.030,1.036 },
210 {18.22,10.48,5.333,3.713,3.115,2.367,1.898,1.631,
211 1.498,1.301,1.171,1.105,1.077,1.048,1.036,1.033,
212 1.031,1.028,1.024,1.022,1.021,1.024 },
213 {14.14,10.65,5.710,3.929,3.266,2.453,1.951,1.669,
214 1.528,1.319,1.178,1.106,1.075,1.040,1.027,1.022,
215 1.020,1.017,1.015,1.013,1.013,1.020 },
216 {14.11,11.73,6.312,4.240,3.478,2.566,2.022,1.720,
217 1.569,1.342,1.186,1.102,1.065,1.022,1.003,0.997,
218 0.995,0.993,0.993,0.993,0.993,1.011 },
219 {22.76,20.01,8.835,5.287,4.144,2.901,2.219,1.855,
220 1.677,1.410,1.224,1.121,1.073,1.014,0.986,0.976,
221 0.974,0.972,0.973,0.974,0.975,0.987 },
222 {50.77,40.85,14.13,7.184,5.284,3.435,2.520,2.059,
223 1.837,1.512,1.283,1.153,1.091,1.010,0.969,0.954,
224 0.950,0.947,0.949,0.952,0.954,0.963 },
225 {65.87,59.06,15.87,7.570,5.567,3.650,2.682,2.182,
226 1.939,1.579,1.325,1.178,1.108,1.014,0.965,0.947,
227 0.941,0.938,0.940,0.944,0.946,0.954 },
228 {55.60,47.34,15.92,7.810,5.755,3.767,2.760,2.239,
229 1.985,1.609,1.343,1.188,1.113,1.013,0.960,0.939,
230 0.933,0.930,0.933,0.936,0.939,0.949 }};
231
232 static const G4double cpositron[15][22] = {
233 {2.589,2.044,1.658,1.446,1.347,1.217,1.144,1.110,
234 1.097,1.083,1.080,1.086,1.092,1.108,1.123,1.131,
235 1.131,1.126,1.117,1.108,1.103,1.100 },
236 {3.904,2.794,2.079,1.710,1.543,1.325,1.202,1.145,
237 1.122,1.096,1.089,1.092,1.098,1.114,1.130,1.137,
238 1.138,1.132,1.122,1.113,1.108,1.102 },
239 {7.970,6.080,4.442,3.398,2.872,2.127,1.672,1.451,
240 1.357,1.246,1.194,1.179,1.178,1.188,1.201,1.205,
241 1.203,1.190,1.173,1.159,1.151,1.145 },
242 {9.714,7.607,5.747,4.493,3.815,2.777,2.079,1.715,
243 1.553,1.353,1.253,1.219,1.211,1.214,1.225,1.228,
244 1.225,1.210,1.191,1.175,1.166,1.174 },
245 {17.97,12.95,8.628,6.065,4.849,3.222,2.275,1.820,
246 1.624,1.382,1.259,1.214,1.202,1.202,1.214,1.219,
247 1.217,1.203,1.184,1.169,1.160,1.151 },
248 {24.83,17.06,10.84,7.355,5.767,3.707,2.546,1.996,
249 1.759,1.465,1.311,1.252,1.234,1.228,1.238,1.241,
250 1.237,1.222,1.201,1.184,1.174,1.159 },
251 {23.26,17.15,11.52,8.049,6.375,4.114,2.792,2.155,
252 1.880,1.535,1.353,1.281,1.258,1.247,1.254,1.256,
253 1.252,1.234,1.212,1.194,1.183,1.170 },
254 {22.33,18.01,12.86,9.212,7.336,4.702,3.117,2.348,
255 2.015,1.602,1.385,1.297,1.268,1.251,1.256,1.258,
256 1.254,1.237,1.214,1.195,1.185,1.179 },
257 {33.91,24.13,15.71,10.80,8.507,5.467,3.692,2.808,
258 2.407,1.873,1.564,1.425,1.374,1.330,1.324,1.320,
259 1.312,1.288,1.258,1.235,1.221,1.205 },
260 {32.14,24.11,16.30,11.40,9.015,5.782,3.868,2.917,
261 2.490,1.925,1.596,1.447,1.391,1.342,1.332,1.327,
262 1.320,1.294,1.264,1.240,1.226,1.214 },
263 {29.51,24.07,17.19,12.28,9.766,6.238,4.112,3.066,
264 2.602,1.995,1.641,1.477,1.414,1.356,1.342,1.336,
265 1.328,1.302,1.270,1.245,1.231,1.233 },
266 {38.19,30.85,21.76,15.35,12.07,7.521,4.812,3.498,
267 2.926,2.188,1.763,1.563,1.484,1.405,1.382,1.371,
268 1.361,1.330,1.294,1.267,1.251,1.239 },
269 {49.71,39.80,27.96,19.63,15.36,9.407,5.863,4.155,
270 3.417,2.478,1.944,1.692,1.589,1.480,1.441,1.423,
271 1.409,1.372,1.330,1.298,1.280,1.258 },
272 {59.25,45.08,30.36,20.83,16.15,9.834,6.166,4.407,
273 3.641,2.648,2.064,1.779,1.661,1.531,1.482,1.459,
274 1.442,1.400,1.354,1.319,1.299,1.272 },
275 {56.38,44.29,30.50,21.18,16.51,10.11,6.354,4.542,
276 3.752,2.724,2.116,1.817,1.692,1.554,1.499,1.474,
277 1.456,1.412,1.364,1.328,1.307,1.282 }};
278
279 //data/corrections for T > Tlim
280
281 static const G4double hecorr[15] = {
282 120.70, 117.50, 105.00, 92.92, 79.23, 74.510, 68.29,
283 57.39, 41.97, 36.14, 24.53, 10.21, -7.855, -16.84,
284 -22.30};
285
286 G4double sigma;
287 SetParticle(part);
288
289 G4double Z23 = G4Pow::GetInstance()->Z23(G4lrint(atomicNumber));
290
291 // correction if particle .ne. e-/e+
292 // compute equivalent kinetic energy
293 // lambda depends on p*beta ....
294
295 G4double eKineticEnergy = kinEnergy;
296
297 if(mass > CLHEP::electron_mass_c2)
298 {
299 G4double TAU = kinEnergy/mass ;
300 G4double c = mass*TAU*(TAU+2.)/(CLHEP::electron_mass_c2*(TAU+1.)) ;
301 G4double w = c-2.;
302 G4double tau = 0.5*(w+std::sqrt(w*w+4.*c)) ;
303 eKineticEnergy = CLHEP::electron_mass_c2*tau ;
304 }
305
306 G4double eTotalEnergy = eKineticEnergy + CLHEP::electron_mass_c2 ;
307 G4double beta2 = eKineticEnergy*(eTotalEnergy+CLHEP::electron_mass_c2)
308 /(eTotalEnergy*eTotalEnergy);
309 G4double bg2 = eKineticEnergy*(eTotalEnergy+CLHEP::electron_mass_c2)
310 /(CLHEP::electron_mass_c2*CLHEP::electron_mass_c2);
311
312 static const G4double epsfactor = 2.*CLHEP::electron_mass_c2*
313 CLHEP::electron_mass_c2*CLHEP::Bohr_radius*CLHEP::Bohr_radius
314 /(CLHEP::hbarc*CLHEP::hbarc);
315 G4double eps = epsfactor*bg2/Z23;
316
317 if (eps<epsmin) sigma = 2.*eps*eps;
318 else if(eps<epsmax) sigma = G4Log(1.+2.*eps)-2.*eps/(1.+2.*eps);
319 else sigma = G4Log(2.*eps)-1.+1./eps;
320
321 sigma *= chargeSquare*atomicNumber*atomicNumber/(beta2*bg2);
322
323 // interpolate in AtomicNumber and beta2
324 G4double c1,c2,cc1;
325
326 // get bin number in Z
327 G4int iZ = 14;
328 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
329 while ((iZ>=0)&&(Zdat[iZ]>=atomicNumber)) { --iZ; }
330
331 iZ = std::min(std::max(iZ, 0), 13);
332
333 G4double ZZ1 = Zdat[iZ];
334 G4double ZZ2 = Zdat[iZ+1];
335 G4double ratZ = (atomicNumber-ZZ1)*(atomicNumber+ZZ1)/
336 ((ZZ2-ZZ1)*(ZZ2+ZZ1));
337
338 static const G4double Tlim = 10.*CLHEP::MeV;
339 static const G4double sigmafactor =
340 CLHEP::twopi*CLHEP::classic_electr_radius*CLHEP::classic_electr_radius;
341 static const G4double beta2lim = Tlim*(Tlim+2.*CLHEP::electron_mass_c2)/
342 ((Tlim+CLHEP::electron_mass_c2)*(Tlim+CLHEP::electron_mass_c2));
343 static const G4double bg2lim = Tlim*(Tlim+2.*CLHEP::electron_mass_c2)/
344 (CLHEP::electron_mass_c2*CLHEP::electron_mass_c2);
345
346 static const G4double sig0[15] = {
347 0.2672*CLHEP::barn, 0.5922*CLHEP::barn, 2.653*CLHEP::barn, 6.235*CLHEP::barn,
348 11.69*CLHEP::barn , 13.24*CLHEP::barn , 16.12*CLHEP::barn, 23.00*CLHEP::barn,
349 35.13*CLHEP::barn , 39.95*CLHEP::barn , 50.85*CLHEP::barn, 67.19*CLHEP::barn,
350 91.15*CLHEP::barn , 104.4*CLHEP::barn , 113.1*CLHEP::barn};
351
352 static const G4double Tdat[22] = {
353 100*CLHEP::eV, 200*CLHEP::eV, 400*CLHEP::eV, 700*CLHEP::eV,
354 1*CLHEP::keV, 2*CLHEP::keV, 4*CLHEP::keV, 7*CLHEP::keV,
355 10*CLHEP::keV, 20*CLHEP::keV, 40*CLHEP::keV, 70*CLHEP::keV,
356 100*CLHEP::keV, 200*CLHEP::keV, 400*CLHEP::keV, 700*CLHEP::keV,
357 1*CLHEP::MeV, 2*CLHEP::MeV, 4*CLHEP::MeV, 7*CLHEP::MeV,
358 10*CLHEP::MeV, 20*CLHEP::MeV};
359
360 if(eKineticEnergy <= Tlim)
361 {
362 // get bin number in T (beta2)
363 G4int iT = 21;
364 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
365 while ((iT>=0)&&(Tdat[iT]>=eKineticEnergy)) iT -= 1;
366
367 iT = std::min(std::max(iT, 0), 20);
368
369 // calculate betasquare values
370 G4double T = Tdat[iT];
371 G4double E = T + CLHEP::electron_mass_c2;
372 G4double b2small = T*(E+CLHEP::electron_mass_c2)/(E*E);
373
374 T = Tdat[iT+1];
375 E = T + CLHEP::electron_mass_c2;
376 G4double b2big = T*(E+CLHEP::electron_mass_c2)/(E*E);
377 G4double ratb2 = (beta2-b2small)/(b2big-b2small);
378
379 if (charge < 0.)
380 {
381 c1 = celectron[iZ][iT];
382 c2 = celectron[iZ+1][iT];
383 cc1 = c1+ratZ*(c2-c1);
384
385 c1 = celectron[iZ][iT+1];
386 c2 = celectron[iZ+1][iT+1];
387 }
388 else
389 {
390 c1 = cpositron[iZ][iT];
391 c2 = cpositron[iZ+1][iT];
392 cc1 = c1+ratZ*(c2-c1);
393
394 c1 = cpositron[iZ][iT+1];
395 c2 = cpositron[iZ+1][iT+1];
396 }
397 G4double cc2 = c1+ratZ*(c2-c1);
398 sigma *= sigmafactor/(cc1+ratb2*(cc2-cc1));
399 }
400 else
401 {
402 c1 = bg2lim*sig0[iZ]*(1.+hecorr[iZ]*(beta2-beta2lim))/bg2;
403 c2 = bg2lim*sig0[iZ+1]*(1.+hecorr[iZ+1]*(beta2-beta2lim))/bg2;
404 if((atomicNumber >= ZZ1) && (atomicNumber <= ZZ2))
405 sigma = c1+ratZ*(c2-c1) ;
406 else if(atomicNumber < ZZ1)
407 sigma = atomicNumber*atomicNumber*c1/(ZZ1*ZZ1);
408 else if(atomicNumber > ZZ2)
409 sigma = atomicNumber*atomicNumber*c2/(ZZ2*ZZ2);
410 }
411 // low energy correction based on theory
412 sigma *= (1.+0.30/(1.+std::sqrt(1000.*eKineticEnergy)));
413
414 return sigma;
415}
G4double G4Log(G4double x)
Definition G4Log.hh:227
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
static G4Pow * GetInstance()
Definition G4Pow.cc:41
G4double Z23(G4int Z) const
Definition G4Pow.hh:125
int G4lrint(double ad)
Definition templates.hh:134

◆ ComputeGeomPathLength()

G4double G4UrbanMscModel::ComputeGeomPathLength ( G4double truePathLength)
overridevirtual

Implements G4VMscModel.

Definition at line 671 of file G4UrbanMscModel.cc.

672{
673 lambdaeff = lambda0;
674 par1 = -1. ;
675 par2 = par3 = 0. ;
676
677 // this correction needed to run MSC with eIoni and eBrem inactivated
678 // and makes no harm for a normal run
679 tPathLength = std::min(tPathLength,currentRange);
680
681 // do the true -> geom transformation
682 zPathLength = tPathLength;
683
684 // z = t for very small tPathLength
685 if(tPathLength < tlimitminfix2) return zPathLength;
686
687 /*
688 G4cout << "ComputeGeomPathLength: tpl= " << tPathLength
689 << " R= " << currentRange << " L0= " << lambda0
690 << " E= " << currentKinEnergy << " "
691 << particle->GetParticleName() << G4endl;
692 */
693 G4double tau = tPathLength/lambda0 ;
694
695 if ((tau <= tausmall) || insideskin) {
696 zPathLength = std::min(tPathLength, lambda0);
697
698 } else if (tPathLength < currentRange*dtrl) {
699 zPathLength = (tau < taulim) ? tPathLength*(1.-0.5*tau)
700 : lambda0*(1.-G4Exp(-tau));
701
702 } else if(currentKinEnergy < mass || tPathLength == currentRange) {
703 par1 = 1./currentRange;
704 par2 = currentRange/lambda0;
705 par3 = 1.+par2;
706 if(tPathLength < currentRange) {
707 zPathLength =
708 (1.-G4Exp(par3*G4Log(1.-tPathLength/currentRange)))/(par1*par3);
709 } else {
710 zPathLength = 1./(par1*par3);
711 }
712
713 } else {
714 G4double rfin = std::max(currentRange-tPathLength, 0.01*currentRange);
715 G4double T1 = GetEnergy(particle,rfin,couple);
716 G4double lambda1 = GetTransportMeanFreePath(particle,T1);
717
718 par1 = (lambda0-lambda1)/(lambda0*tPathLength);
719 //G4cout << "par1= " << par1 << " L1= " << lambda1 << G4endl;
720 par2 = 1./(par1*lambda0);
721 par3 = 1.+par2;
722 zPathLength = (1.-G4Exp(par3*G4Log(lambda1/lambda0)))/(par1*par3);
723 }
724
725 zPathLength = std::min(zPathLength, lambda0);
726 //G4cout<< "zPathLength= "<< zPathLength<< " L0= " << lambda0 << G4endl;
727 return zPathLength;
728}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:180
G4double dtrl
G4double GetTransportMeanFreePath(const G4ParticleDefinition *part, G4double kinEnergy)
G4double GetEnergy(const G4ParticleDefinition *part, G4double range, const G4MaterialCutsCouple *couple)

◆ ComputeTheta0()

G4double G4UrbanMscModel::ComputeTheta0 ( G4double truePathLength,
G4double KineticEnergy )

Definition at line 972 of file G4UrbanMscModel.cc.

974{
975 // for all particles take the width of the central part
976 // from a parametrization similar to the Highland formula
977 // ( Highland formula: Particle Physics Booklet, July 2002, eq. 26.10)
978 G4double invbetacp = (kinEnergy+mass)/(kinEnergy*(kinEnergy+2.*mass));
979 if(currentKinEnergy != kinEnergy) {
980 invbetacp = std::sqrt(invbetacp*(currentKinEnergy+mass)/
981 (currentKinEnergy*(currentKinEnergy+2.*mass)));
982 }
983 G4double y = trueStepLength/currentRadLength;
984
985 if(fPosiCorrection && particle == positron)
986 {
987 static const G4double xl= 0.6;
988 static const G4double xh= 0.9;
989 static const G4double e = 113.0;
990 G4double corr;
991
992 G4double tau = std::sqrt(currentKinEnergy*kinEnergy)/mass;
993 G4double x = std::sqrt(tau*(tau+2.)/((tau+1.)*(tau+1.)));
994 G4double a = msc[idx]->posa;
995 G4double b = msc[idx]->posb;
996 G4double c = msc[idx]->posc;
997 G4double d = msc[idx]->posd;
998 if(x < xl) {
999 corr = a*(1.-G4Exp(-b*x));
1000 } else if(x > xh) {
1001 corr = c+d*G4Exp(e*(x-1.));
1002 } else {
1003 G4double yl = a*(1.-G4Exp(-b*xl));
1004 G4double yh = c+d*G4Exp(e*(xh-1.));
1005 G4double y0 = (yh-yl)/(xh-xl);
1006 G4double y1 = yl-y0*xl;
1007 corr = y0*x+y1;
1008 }
1009 //==================================================================
1010 y *= corr*msc[idx]->pose;
1011 }
1012
1013 static const G4double c_highland = 13.6*CLHEP::MeV;
1014 G4double theta0 = c_highland*std::abs(charge)*std::sqrt(y)*invbetacp;
1015
1016 // correction factor from e- scattering data
1017 theta0 *= (msc[idx]->coeffth1+msc[idx]->coeffth2*G4Log(y));
1018 return theta0;
1019}

◆ ComputeTruePathLengthLimit()

G4double G4UrbanMscModel::ComputeTruePathLengthLimit ( const G4Track & track,
G4double & currentMinimalStep )
overridevirtual

Implements G4VMscModel.

Definition at line 434 of file G4UrbanMscModel.cc.

437{
438 tPathLength = currentMinimalStep;
439 const G4DynamicParticle* dp = track.GetDynamicParticle();
440
441 G4StepPoint* sp = track.GetStep()->GetPreStepPoint();
442 G4StepStatus stepStatus = sp->GetStepStatus();
443 couple = track.GetMaterialCutsCouple();
444 SetCurrentCouple(couple);
445 idx = couple->GetIndex();
446 currentKinEnergy = dp->GetKineticEnergy();
447 currentLogKinEnergy = dp->GetLogKineticEnergy();
448 currentRange = GetRange(particle,currentKinEnergy,couple,currentLogKinEnergy);
449 lambda0 = GetTransportMeanFreePath(particle,currentKinEnergy,
450 currentLogKinEnergy);
451 tPathLength = std::min(tPathLength,currentRange);
452 /*
453 G4cout << "G4Urban::StepLimit tPathLength= " << tPathLength
454 << " range= " <<currentRange<< " lambda= "<<lambda0
455 <<G4endl;
456 */
457 // extreme small step
458 if(tPathLength < tlimitminfix) {
459 latDisplasment = false;
460 return ConvertTrueToGeom(tPathLength, currentMinimalStep);
461 }
462
463 presafety = (stepStatus == fGeomBoundary) ? sp->GetSafety()
464 : ComputeSafety(sp->GetPosition(), tPathLength);
465
466 // stop here if small step or range is less than safety
467 if((tPathLength == currentRange && tPathLength < presafety) ||
468 tPathLength < tlimitminfix) {
469 latDisplasment = false;
470 return ConvertTrueToGeom(tPathLength, currentMinimalStep);
471 }
472
473 // upper limit for the straight line distance the particle can travel
474 // for electrons and positrons
475 G4double distance = (mass < masslimite)
476 ? currentRange*msc[idx]->doverra
477 // for muons, hadrons
478 : currentRange*msc[idx]->doverrb;
479
480 /*
481 G4cout << "G4Urban::StepLimit tPathLength= "
482 <<tPathLength<<" safety= " << presafety
483 << " range= " <<currentRange<< " lambda= "<<lambda0
484 << " Alg: " << steppingAlgorithm <<G4endl;
485 */
486 // far from geometry boundary
487 if(distance < presafety)
488 {
489 latDisplasment = false;
490 return ConvertTrueToGeom(tPathLength, currentMinimalStep);
491 }
492
493 latDisplasment = latDisplasmentbackup;
494 // ----------------------------------------------------------------
495 // distance to boundary
497 {
498 //compute geomlimit and presafety
499 geomlimit = ComputeGeomLimit(track, presafety, currentRange);
500 /*
501 G4cout << "G4Urban::Distance to boundary geomlimit= "
502 <<geomlimit<<" safety= " << presafety<<G4endl;
503 */
504
505 smallstep += 1.;
506 insideskin = false;
507 tgeom = geombig;
508
509 // initialisation at first step and at the boundary
510 if(firstStep || (stepStatus == fGeomBoundary))
511 {
512 rangeinit = currentRange;
513 if(!firstStep) { smallstep = 1.; }
514
515 //stepmin ~ lambda_elastic
516 stepmin = ComputeStepmin();
517 skindepth = skin*stepmin;
518 tlimitmin = ComputeTlimitmin();
519 /*
520 G4cout << "rangeinit= " << rangeinit << " stepmin= " << stepmin
521 << " tlimitmin= " << tlimitmin << " geomlimit= "
522 << geomlimit <<G4endl;
523 */
524 }
525 // constraint from the geometry
526 if((geomlimit < geombig) && (geomlimit > geommin))
527 {
528 // geomlimit is a geometrical step length
529 // transform it to true path length (estimation)
530 if(lambda0 > geomlimit) {
531 geomlimit = -lambda0*G4Log(1.-geomlimit/lambda0)+tlimitmin;
532 }
533 tgeom = (stepStatus == fGeomBoundary) ? geomlimit/facgeom
534 : facrange*rangeinit + stepmin;
535 }
536
537 //step limit
538 tlimit = (currentRange > presafety) ?
539 std::max(facrange*rangeinit, facsafety*presafety) : currentRange;
540
541 //lower limit for tlimit
542 tlimit = std::min(std::max(tlimit,tlimitmin), tgeom);
543 /*
544 G4cout << "tgeom= " << tgeom << " geomlimit= " << geomlimit
545 << " tlimit= " << tlimit << " presafety= " << presafety << G4endl;
546 */
547 // shortcut
548 if((tPathLength < tlimit) && (tPathLength < presafety) &&
549 (smallstep > skin) && (tPathLength < geomlimit-0.999*skindepth))
550 {
551 return ConvertTrueToGeom(tPathLength, currentMinimalStep);
552 }
553
554 // step reduction near to boundary
555 if(smallstep <= skin)
556 {
557 tlimit = stepmin;
558 insideskin = true;
559 }
560 else if(geomlimit < geombig)
561 {
562 if(geomlimit > skindepth)
563 {
564 tlimit = std::min(tlimit, geomlimit-0.999*skindepth);
565 }
566 else
567 {
568 insideskin = true;
569 tlimit = std::min(tlimit, stepmin);
570 }
571 }
572
573 tlimit = std::max(tlimit, stepmin);
574
575 // randomise if not 'small' step and step determined by msc
576 tPathLength = (tlimit < tPathLength && smallstep > skin && !insideskin)
577 ? std::min(tPathLength, Randomizetlimit())
578 : std::min(tPathLength, tlimit);
579 }
580 // ----------------------------------------------------------------
581 // for simulation with or without magnetic field
582 // there no small step/single scattering at boundaries
583 else if(steppingAlgorithm == fUseSafety)
584 {
585 if(firstStep || (stepStatus == fGeomBoundary)) {
586 rangeinit = currentRange;
587 fr = facrange;
588 // stepping for e+/e- only (not for muons,hadrons)
589 if(mass < masslimite)
590 {
591 rangeinit = std::max(rangeinit, lambda0);
592 if(lambda0 > lambdalimit) {
593 fr *= (0.75+0.25*lambda0/lambdalimit);
594 }
595 }
596 //lower limit for tlimit
597 stepmin = ComputeStepmin();
598 tlimitmin = ComputeTlimitmin();
599 }
600
601 //step limit
602 tlimit = (currentRange > presafety) ?
603 std::max(fr*rangeinit, facsafety*presafety) : currentRange;
604
605 //lower limit for tlimit
606 tlimit = std::max(tlimit, tlimitmin);
607
608 // randomise if step determined by msc
609 tPathLength = (tlimit < tPathLength) ?
610 std::min(tPathLength, Randomizetlimit()) : tPathLength;
611 }
612 // ----------------------------------------------------------------
613 // for simulation with or without magnetic field
614 // there is small step/single scattering at boundaries
616 {
617 if(firstStep || (stepStatus == fGeomBoundary)) {
618 rangeinit = currentRange;
619 fr = facrange;
620 if(mass < masslimite)
621 {
622 if(lambda0 > lambdalimit) {
623 fr *= (0.84+0.16*lambda0/lambdalimit);
624 }
625 }
626 //lower limit for tlimit
627 stepmin = ComputeStepmin();
628 tlimitmin = ComputeTlimitmin();
629 }
630 //step limit
631 tlimit = (currentRange > presafety) ?
632 std::max(fr*rangeinit, facsafety*presafety) : currentRange;
633
634 //lower limit for tlimit
635 tlimit = std::max(tlimit, tlimitmin);
636
637 // condition for tPathLength from drr and finalr
638 if(currentRange > finalr) {
639 G4double tmax = drr*currentRange+
640 finalr*(1.-drr)*(2.-finalr/currentRange);
641 tPathLength = std::min(tPathLength,tmax);
642 }
643
644 // randomise if step determined by msc
645 tPathLength = (tlimit < tPathLength) ?
646 std::min(tPathLength, Randomizetlimit()) : tPathLength;
647 }
648
649 // ----------------------------------------------------------------
650 // simple step limitation
651 else
652 {
653 if (stepStatus == fGeomBoundary)
654 {
655 tlimit = (currentRange > lambda0)
656 ? facrange*currentRange : facrange*lambda0;
657 tlimit = std::max(tlimit, tlimitmin);
658 }
659 // randomise if step determined by msc
660 tPathLength = (tlimit < tPathLength) ?
661 std::min(tPathLength, Randomizetlimit()) : tPathLength;
662 }
663
664 // ----------------------------------------------------------------
665 firstStep = false;
666 return ConvertTrueToGeom(tPathLength, currentMinimalStep);
667}
@ fUseSafety
@ fUseSafetyPlus
@ fUseDistanceToBoundary
G4StepStatus
@ fGeomBoundary
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const
G4StepPoint * GetPreStepPoint() const
const G4DynamicParticle * GetDynamicParticle() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
const G4Step * GetStep() const
void SetCurrentCouple(const G4MaterialCutsCouple *)
G4double facrange
G4double ComputeGeomLimit(const G4Track &, G4double &presafety, G4double limit)
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4double lambdalimit
G4MscStepLimitType steppingAlgorithm
G4double ConvertTrueToGeom(G4double &tLength, G4double &gLength)
G4bool latDisplasment
G4double ComputeSafety(const G4ThreeVector &position, G4double limit=DBL_MAX)
G4double facsafety
G4double facgeom

◆ ComputeTrueStepLength()

G4double G4UrbanMscModel::ComputeTrueStepLength ( G4double geomStepLength)
overridevirtual

Implements G4VMscModel.

Definition at line 732 of file G4UrbanMscModel.cc.

733{
734 // step defined other than transportation
735 if(geomStepLength == zPathLength) {
736 //G4cout << "Urban::ComputeTrueLength: tPathLength= " << tPathLength
737 // << " step= " << geomStepLength << " *** " << G4endl;
738 return tPathLength;
739 }
740
741 zPathLength = geomStepLength;
742
743 // t = z for very small step
744 if(geomStepLength < tlimitminfix2) {
745 tPathLength = geomStepLength;
746
747 // recalculation
748 } else {
749
750 G4double tlength = geomStepLength;
751 if((geomStepLength > lambda0*tausmall) && !insideskin) {
752
753 if(par1 < 0.) {
754 tlength = -lambda0*G4Log(1.-geomStepLength/lambda0) ;
755 } else {
756 const G4double par4 = par1*par3;
757 if(par4*geomStepLength < 1.) {
758 tlength = (1.-G4Exp(G4Log(1.-par4*geomStepLength)/par3))/par1;
759 } else {
760 tlength = currentRange;
761 }
762 }
763
764 if(tlength < geomStepLength) { tlength = geomStepLength; }
765 else if(tlength > tPathLength) { tlength = tPathLength; }
766 }
767 tPathLength = tlength;
768 }
769 //G4cout << "Urban::ComputeTrueLength: tPathLength= " << tPathLength
770 // << " step= " << geomStepLength << " &&& " << G4endl;
771
772 return tPathLength;
773}

◆ Initialise()

void G4UrbanMscModel::Initialise ( const G4ParticleDefinition * p,
const G4DataVector &  )
overridevirtual

Implements G4VEmModel.

Definition at line 134 of file G4UrbanMscModel.cc.

136{
137 // set values of some data members
138 SetParticle(p);
139 fParticleChange = GetParticleChangeForMSC(p);
141
142 latDisplasmentbackup = latDisplasment;
143
144 // if model is locked parameters should be defined via Set methods
145 if(!IsLocked()) {
147 fPosiCorrection = G4EmParameters::Instance()->MscPositronCorrection();
148 }
149
150 // initialise cache only once
151 if(0 == msc.size()) {
152 G4AutoLock l(&theUrbanMutex);
153 if(0 == msc.size()) {
154 isFirstInstance = true;
155 msc.resize(1, nullptr);
156 }
157 l.unlock();
158 }
159 // initialise cache for each new run
160 if(isFirstInstance) { InitialiseModelCache(); }
161
162 /*
163 G4cout << "### G4UrbanMscModel::Initialise done for "
164 << p->GetParticleName() << " type= " << steppingAlgorithm << G4endl;
165 G4cout << " RangeFact= " << facrange << " GeomFact= " << facgeom
166 << " SafetyFact= " << facsafety << " LambdaLim= " << lambdalimit
167 << G4endl;
168 */
169}
G4TemplateAutoLock< G4Mutex > G4AutoLock
G4bool LateralDisplacementAlg96() const
static G4EmParameters * Instance()
G4bool MscPositronCorrection() const
G4bool IsLocked() const
G4ParticleChangeForMSC * GetParticleChangeForMSC(const G4ParticleDefinition *p=nullptr)
void InitialiseParameters(const G4ParticleDefinition *)

◆ operator=()

G4UrbanMscModel & G4UrbanMscModel::operator= ( const G4UrbanMscModel & right)
delete

◆ SampleScattering()

G4ThreeVector & G4UrbanMscModel::SampleScattering ( const G4ThreeVector & oldDirection,
G4double safety )
overridevirtual

Implements G4VMscModel.

Definition at line 778 of file G4UrbanMscModel.cc.

780{
781 fDisplacement.set(0.0,0.0,0.0);
782 if(tPathLength >= currentRange) { return fDisplacement; }
783
784 G4double kinEnergy = currentKinEnergy;
785 if (tPathLength > currentRange*dtrl) {
786 kinEnergy = GetEnergy(particle,currentRange-tPathLength,couple);
787 } else if(tPathLength > currentRange*0.01) {
788 kinEnergy -= tPathLength*GetDEDX(particle,currentKinEnergy,couple,
789 currentLogKinEnergy);
790 }
791
792 if((tPathLength <= tlimitminfix) || (tPathLength < tausmall*lambda0) ||
793 (kinEnergy <= CLHEP::eV)) { return fDisplacement; }
794
795 G4double cth = SampleCosineTheta(tPathLength,kinEnergy);
796
797 // protection against 'bad' cth values
798 if(std::abs(cth) >= 1.0) { return fDisplacement; }
799
800 G4double sth = std::sqrt((1.0 - cth)*(1.0 + cth));
801 G4double phi = CLHEP::twopi*rndmEngineMod->flat();
802 G4ThreeVector newDirection(sth*std::cos(phi),sth*std::sin(phi),cth);
803 newDirection.rotateUz(oldDirection);
804
805 fParticleChange->ProposeMomentumDirection(newDirection);
806 /*
807 G4cout << "G4UrbanMscModel::SampleSecondaries: e(MeV)= " << kineticEnergy
808 << " sinTheta= " << sth << " safety(mm)= " << safety
809 << " trueStep(mm)= " << tPathLength
810 << " geomStep(mm)= " << zPathLength
811 << G4endl;
812 */
813
814 if (latDisplasment && currentTau >= tausmall) {
815 if(dispAlg96) { SampleDisplacement(sth, phi); }
816 else { SampleDisplacementNew(cth, phi); }
817 fDisplacement.rotateUz(oldDirection);
818 }
819 return fDisplacement;
820}
CLHEP::Hep3Vector G4ThreeVector
G4double GetDEDX(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4ThreeVector fDisplacement

◆ SetDisplacementAlgorithm96()

void G4UrbanMscModel::SetDisplacementAlgorithm96 ( const G4bool val)
inline

Definition at line 268 of file G4UrbanMscModel.hh.

269{
270 dispAlg96 = val;
271}

◆ SetPositronCorrection()

void G4UrbanMscModel::SetPositronCorrection ( const G4bool val)
inline

Definition at line 275 of file G4UrbanMscModel.hh.

276{
277 fPosiCorrection = val;
278}

◆ StartTracking()

void G4UrbanMscModel::StartTracking ( G4Track * track)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 419 of file G4UrbanMscModel.cc.

420{
421 SetParticle(track->GetDynamicParticle()->GetDefinition());
422 firstStep = true;
423 insideskin = false;
424 fr = facrange;
425 tlimit = tgeom = rangeinit = geombig;
426 smallstep = 1.e10;
427 stepmin = tlimitminfix;
428 tlimitmin = 10.*tlimitminfix;
429 rndmEngineMod = G4Random::getTheEngine();
430}
G4ParticleDefinition * GetDefinition() const

The documentation for this class was generated from the following files: