Geant4 10.7.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)
 
G4UrbanMscModeloperator= (const G4UrbanMscModel &right)=delete
 
 G4UrbanMscModel (const G4UrbanMscModel &)=delete
 
- Public Member Functions inherited from G4VMscModel
 G4VMscModel (const G4String &nam)
 
 ~G4VMscModel () override
 
virtual G4double ComputeTruePathLengthLimit (const G4Track &track, G4double &stepLimit)=0
 
virtual G4double ComputeGeomPathLength (G4double truePathLength)=0
 
virtual G4double ComputeTrueStepLength (G4double geomPathLength)=0
 
virtual G4ThreeVectorSampleScattering (const G4ThreeVector &, G4double safety)=0
 
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 Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
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 StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
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 ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual 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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
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 LPMFlag () 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 SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (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)
 
- 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
 
G4double facgeom
 
G4double facsafety
 
G4double skin
 
G4double dtrl
 
G4double lambdalimit
 
G4double geomMin
 
G4double geomMax
 
G4ThreeVector fDisplacement
 
G4MscStepLimitType steppingAlgorithm
 
G4bool samplez
 
G4bool latDisplasment
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

Detailed Description

Definition at line 70 of file G4UrbanMscModel.hh.

Constructor & Destructor Documentation

◆ G4UrbanMscModel() [1/2]

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

Definition at line 82 of file G4UrbanMscModel.cc.

83 : G4VMscModel(nam)
84{
85 masslimite = 0.6*MeV;
86 fr = 0.02;
87 taubig = 8.0;
88 tausmall = 1.e-16;
89 taulim = 1.e-6;
90 currentTau = taulim;
91 tlimitminfix = 0.01*nm;
92 tlimitminfix2 = 1.*nm;
93 stepmin = tlimitminfix;
94 smallstep = 1.e10;
95 currentRange = 0. ;
96 rangeinit = 0.;
97 tlimit = 1.e10*mm;
98 tlimitmin = 10.*tlimitminfix;
99 tgeom = 1.e50*mm;
100 geombig = 1.e50*mm;
101 geommin = 1.e-3*mm;
102 geomlimit = geombig;
103 presafety = 0.*mm;
104
105 particle = nullptr;
106
107 positron = G4Positron::Positron();
108 theManager = G4LossTableManager::Instance();
109 rndmEngineMod = G4Random::getTheEngine();
110
111 firstStep = true;
112 insideskin = false;
113 latDisplasmentbackup = false;
114 dispAlg96 = true;
115
116 rangecut = geombig;
117 drr = 0.35;
118 finalr = 10.*um;
119
120 tlow = 5.*CLHEP::keV;
121 invmev = 1.0/CLHEP::MeV;
122
123 skindepth = skin*stepmin;
124
125 mass = proton_mass_c2;
126 charge = ChargeSquare = 1.0;
127 currentKinEnergy = currentRadLength = lambda0 = lambdaeff = tPathLength
128 = zPathLength = par1 = par2 = par3 = 0;
129 currentLogKinEnergy = LOG_EKIN_MIN;
130
131 idx = 0;
132 fParticleChange = nullptr;
133 couple = nullptr;
134}
static G4LossTableManager * Instance()
static G4Positron * Positron()
Definition: G4Positron.cc:93
G4double skin
Definition: G4VMscModel.hh:196
#define LOG_EKIN_MIN
Definition: templates.hh:98

◆ ~G4UrbanMscModel()

G4UrbanMscModel::~G4UrbanMscModel ( )
override

Definition at line 138 of file G4UrbanMscModel.cc.

139{
140 if(IsMaster()) {
141 for(auto & ptr : msc) { delete ptr; }
142 msc.clear();
143 }
144}
G4bool IsMaster() const
Definition: G4VEmModel.hh:736

◆ 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 170 of file G4UrbanMscModel.cc.

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

716{
717 lambdaeff = lambda0;
718 par1 = -1. ;
719 par2 = par3 = 0. ;
720
721 // this correction needed to run MSC with eIoni and eBrem inactivated
722 // and makes no harm for a normal run
723 tPathLength = std::min(tPathLength,currentRange);
724
725 // do the true -> geom transformation
726 zPathLength = tPathLength;
727
728 // z = t for very small tPathLength
729 if(tPathLength < tlimitminfix2) return zPathLength;
730
731 /*
732 G4cout << "ComputeGeomPathLength: tpl= " << tPathLength
733 << " R= " << currentRange << " L0= " << lambda0
734 << " E= " << currentKinEnergy << " "
735 << particle->GetParticleName() << G4endl;
736 */
737 G4double tau = tPathLength/lambda0 ;
738
739 if ((tau <= tausmall) || insideskin) {
740 zPathLength = std::min(tPathLength, lambda0);
741
742 } else if (tPathLength < currentRange*dtrl) {
743 if(tau < taulim) zPathLength = tPathLength*(1.-0.5*tau) ;
744 else zPathLength = lambda0*(1.-G4Exp(-tau));
745
746 } else if(currentKinEnergy < mass || tPathLength == currentRange) {
747 par1 = 1./currentRange ;
748 par2 = 1./(par1*lambda0) ;
749 par3 = 1.+par2 ;
750 if(tPathLength < currentRange) {
751 zPathLength =
752 (1.-G4Exp(par3*G4Log(1.-tPathLength/currentRange)))/(par1*par3);
753 } else {
754 zPathLength = 1./(par1*par3);
755 }
756
757 } else {
758 G4double rfin = std::max(currentRange-tPathLength, 0.01*currentRange);
759 G4double T1 = GetEnergy(particle,rfin,couple);
760 G4double lambda1 = GetTransportMeanFreePath(particle,T1);
761
762 par1 = (lambda0-lambda1)/(lambda0*tPathLength);
763 //G4cout << "par1= " << par1 << " L1= " << lambda1 << G4endl;
764 par2 = 1./(par1*lambda0);
765 par3 = 1.+par2 ;
766 zPathLength = (1.-G4Exp(par3*G4Log(lambda1/lambda0)))/(par1*par3);
767 }
768
769 zPathLength = std::min(zPathLength, lambda0);
770 //G4cout<< "zPathLength= "<< zPathLength<< " L0= " << lambda0 << G4endl;
771 return zPathLength;
772}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double dtrl
Definition: G4VMscModel.hh:197
G4double GetTransportMeanFreePath(const G4ParticleDefinition *part, G4double kinEnergy)
Definition: G4VMscModel.hh:405
G4double GetEnergy(const G4ParticleDefinition *part, G4double range, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.hh:368

◆ ComputeTheta0()

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

Definition at line 1056 of file G4UrbanMscModel.cc.

1058{
1059 // for all particles take the width of the central part
1060 // from a parametrization similar to the Highland formula
1061 // ( Highland formula: Particle Physics Booklet, July 2002, eq. 26.10)
1062 G4double invbetacp = std::sqrt((currentKinEnergy+mass)*(KineticEnergy+mass)/
1063 (currentKinEnergy*(currentKinEnergy+2.*mass)*
1064 KineticEnergy*(KineticEnergy+2.*mass)));
1065 G4double y = trueStepLength/currentRadLength;
1066
1067 if(particle == positron)
1068 {
1069 G4double Zeff = msc[idx]->Zeff;
1070 static const G4double xl= 0.6;
1071 static const G4double xh= 0.9;
1072 static const G4double e = 113.0;
1073 G4double corr;
1074
1075 G4double tau = std::sqrt(currentKinEnergy*KineticEnergy)/mass;
1076 G4double x = std::sqrt(tau*(tau+2.)/((tau+1.)*(tau+1.)));
1077 G4double a = 0.994-4.08e-3*Zeff;
1078 G4double b = 7.16+(52.6+365./Zeff)/Zeff;
1079 G4double c = 1.000-4.47e-3*Zeff;
1080 G4double d = 1.21e-3*Zeff;
1081 if(x < xl) {
1082 corr = a*(1.-G4Exp(-b*x));
1083 } else if(x > xh) {
1084 corr = c+d*G4Exp(e*(x-1.));
1085 } else {
1086 G4double yl = a*(1.-G4Exp(-b*xl));
1087 G4double yh = c+d*G4Exp(e*(xh-1.));
1088 G4double y0 = (yh-yl)/(xh-xl);
1089 G4double y1 = yl-y0*xl;
1090 corr = y0*x+y1;
1091 }
1092 //==================================================================
1093 y *= corr*(1.+Zeff*(1.84035e-4*Zeff-1.86427e-2)+0.41125);
1094 }
1095
1096 static const G4double c_highland = 13.6*CLHEP::MeV;
1097 G4double theta0 = c_highland*std::abs(charge)*std::sqrt(y)*invbetacp;
1098
1099 // correction factor from e- scattering data
1100 theta0 *= (msc[idx]->coeffth1+msc[idx]->coeffth2*G4Log(y));
1101 return theta0;
1102}

◆ ComputeTruePathLengthLimit()

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

Implements G4VMscModel.

Definition at line 438 of file G4UrbanMscModel.cc.

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

◆ ComputeTrueStepLength()

G4double G4UrbanMscModel::ComputeTrueStepLength ( G4double  geomStepLength)
overridevirtual

Implements G4VMscModel.

Definition at line 776 of file G4UrbanMscModel.cc.

777{
778 // step defined other than transportation
779 if(geomStepLength == zPathLength) {
780 //G4cout << "Urban::ComputeTrueLength: tPathLength= " << tPathLength
781 // << " step= " << geomStepLength << " *** " << G4endl;
782 return tPathLength;
783 }
784
785 zPathLength = geomStepLength;
786
787 // t = z for very small step
788 if(geomStepLength < tlimitminfix2) {
789 tPathLength = geomStepLength;
790
791 // recalculation
792 } else {
793
794 G4double tlength = geomStepLength;
795 if((geomStepLength > lambda0*tausmall) && !insideskin) {
796
797 if(par1 < 0.) {
798 tlength = -lambda0*G4Log(1.-geomStepLength/lambda0) ;
799 } else {
800 if(par1*par3*geomStepLength < 1.) {
801 tlength = (1.-G4Exp(G4Log(1.-par1*par3*geomStepLength)/par3))/par1 ;
802 } else {
803 tlength = currentRange;
804 }
805 }
806
807 if(tlength < geomStepLength) { tlength = geomStepLength; }
808 else if(tlength > tPathLength) { tlength = tPathLength; }
809 }
810 tPathLength = tlength;
811 }
812 //G4cout << "Urban::ComputeTrueLength: tPathLength= " << tPathLength
813 // << " step= " << geomStepLength << " &&& " << G4endl;
814
815 return tPathLength;
816}

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 148 of file G4UrbanMscModel.cc.

150{
151 // set values of some data members
152 SetParticle(p);
153 fParticleChange = GetParticleChangeForMSC(p);
155
156 latDisplasmentbackup = latDisplasment;
158 if(IsMaster() || msc.size() == 0) { InitialiseModelCache(); }
159 /*
160 G4cout << "### G4UrbanMscModel::Initialise done for "
161 << p->GetParticleName() << " type= " << steppingAlgorithm << G4endl;
162 G4cout << " RangeFact= " << facrange << " GeomFact= " << facgeom
163 << " SafetyFact= " << facsafety << " LambdaLim= " << lambdalimit
164 << G4endl;
165 */
166}
G4bool LateralDisplacementAlg96() const
static G4EmParameters * Instance()
G4ParticleChangeForMSC * GetParticleChangeForMSC(const G4ParticleDefinition *p=nullptr)
Definition: G4VMscModel.cc:91
void InitialiseParameters(const G4ParticleDefinition *)
Definition: G4VMscModel.cc:139

◆ operator=()

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

◆ SampleScattering()

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

Implements G4VMscModel.

Definition at line 821 of file G4UrbanMscModel.cc.

823{
824 fDisplacement.set(0.0,0.0,0.0);
825 G4double kineticEnergy = currentKinEnergy;
826 if (tPathLength > currentRange*dtrl) {
827 kineticEnergy = GetEnergy(particle,currentRange-tPathLength,couple);
828 } else {
829 kineticEnergy -= tPathLength*GetDEDX(particle,currentKinEnergy,couple,
830 currentLogKinEnergy);
831 }
832
833 if((kineticEnergy <= eV) || (tPathLength <= tlimitminfix) ||
834 (tPathLength < tausmall*lambda0)) { return fDisplacement; }
835
836 G4double cth = SampleCosineTheta(tPathLength,kineticEnergy);
837
838 // protection against 'bad' cth values
839 if(std::abs(cth) >= 1.0) { return fDisplacement; }
840
841 /*
842 if(cth < 1.0 - 1000*tPathLength/lambda0 && cth < 0.5 &&
843 kineticEnergy > 20*MeV) {
844 G4cout << "### G4UrbanMscModel::SampleScattering for "
845 << particle->GetParticleName()
846 << " E(MeV)= " << kineticEnergy/MeV
847 << " Step(mm)= " << tPathLength/mm
848 << " in " << CurrentCouple()->GetMaterial()->GetName()
849 << " CosTheta= " << cth << G4endl;
850 }
851 */
852 G4double sth = sqrt((1.0 - cth)*(1.0 + cth));
853 G4double phi = twopi*rndmEngineMod->flat();
854 G4double dirx = sth*cos(phi);
855 G4double diry = sth*sin(phi);
856
857 G4ThreeVector newDirection(dirx,diry,cth);
858 newDirection.rotateUz(oldDirection);
859
860 fParticleChange->ProposeMomentumDirection(newDirection);
861 /*
862 G4cout << "G4UrbanMscModel::SampleSecondaries: e(MeV)= " << kineticEnergy
863 << " sinTheta= " << sth << " safety(mm)= " << safety
864 << " trueStep(mm)= " << tPathLength
865 << " geomStep(mm)= " << zPathLength
866 << G4endl;
867 */
868
869 if (latDisplasment && currentTau >= tausmall) {
870 if(dispAlg96) { SampleDisplacement(sth, phi); }
871 else { SampleDisplacementNew(cth, phi); }
872 fDisplacement.rotateUz(oldDirection);
873 }
874 return fDisplacement;
875}
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
virtual double flat()=0
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
G4double GetDEDX(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.hh:300
G4ThreeVector fDisplacement
Definition: G4VMscModel.hh:202

◆ StartTracking()

void G4UrbanMscModel::StartTracking ( G4Track track)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 423 of file G4UrbanMscModel.cc.

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

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