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

#include <G4MuBetheBlochModel.hh>

+ Inheritance diagram for G4MuBetheBlochModel:

Public Member Functions

 G4MuBetheBlochModel (const G4ParticleDefinition *p=0, const G4String &nam="MuBetheBloch")
 
virtual ~G4MuBetheBlochModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual G4double ComputeCrossSectionPerElectron (const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
- 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 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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., 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 *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
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)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
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
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void ForceBuildTable (G4bool val)
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Member Functions

virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kinEnergy)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Additional Inherited Members

- Protected Attributes inherited from G4VEmModel
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

Detailed Description

Definition at line 69 of file G4MuBetheBlochModel.hh.

Constructor & Destructor Documentation

◆ G4MuBetheBlochModel()

G4MuBetheBlochModel::G4MuBetheBlochModel ( const G4ParticleDefinition p = 0,
const G4String nam = "MuBetheBloch" 
)

Definition at line 78 of file G4MuBetheBlochModel.cc.

80 : G4VEmModel(nam),
81 particle(0),
82 limitKinEnergy(100.*keV),
83 logLimitKinEnergy(log(limitKinEnergy)),
84 twoln10(2.0*log(10.0)),
85 bg2lim(0.0169),
86 taulim(8.4146e-3),
87 alphaprime(fine_structure_const/twopi)
88{
89 theElectron = G4Electron::Electron();
91 fParticleChange = 0;
92
93 // initial initialisation of memeber should be overwritten
94 // by SetParticle
95 mass = massSquare = ratio = 1.0;
96
97 if(p) { SetParticle(p); }
98}
static G4Electron * Electron()
Definition: G4Electron.cc:94
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()

◆ ~G4MuBetheBlochModel()

G4MuBetheBlochModel::~G4MuBetheBlochModel ( )
virtual

Definition at line 102 of file G4MuBetheBlochModel.cc.

103{}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4MuBetheBlochModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  Z,
G4double  A,
G4double  cutEnergy,
G4double  maxEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 184 of file G4MuBetheBlochModel.cc.

190{
192 (p,kineticEnergy,cutEnergy,maxEnergy);
193 return cross;
194}
double G4double
Definition: G4Types.hh:64
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)

◆ ComputeCrossSectionPerElectron()

G4double G4MuBetheBlochModel::ComputeCrossSectionPerElectron ( const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
)
virtual

Definition at line 135 of file G4MuBetheBlochModel.cc.

140{
141 G4double cross = 0.0;
142 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
143 G4double maxEnergy = min(tmax,maxKinEnergy);
144 if(cutEnergy < maxEnergy) {
145
146 G4double totEnergy = kineticEnergy + mass;
147 G4double energy2 = totEnergy*totEnergy;
148 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
149
150 cross = 1.0/cutEnergy - 1.0/maxEnergy - beta2*log(maxEnergy/cutEnergy)/tmax
151 + 0.5*(maxEnergy - cutEnergy)/energy2;
152
153 // radiative corrections of R. Kokoulin
154 if (maxEnergy > limitKinEnergy) {
155
156 G4double logtmax = log(maxEnergy);
157 G4double logtmin = log(max(cutEnergy,limitKinEnergy));
158 G4double logstep = logtmax - logtmin;
159 G4double dcross = 0.0;
160
161 for (G4int ll=0; ll<8; ll++)
162 {
163 G4double ep = exp(logtmin + xgi[ll]*logstep);
164 G4double a1 = log(1.0 + 2.0*ep/electron_mass_c2);
165 G4double a3 = log(4.0*totEnergy*(totEnergy - ep)/massSquare);
166 dcross += wgi[ll]*(1.0/ep - beta2/tmax + 0.5*ep/energy2)*a1*(a3 - a1);
167 }
168
169 cross += dcross*logstep*alphaprime;
170 }
171
172 cross *= twopi_mc2_rcl2/beta2;
173
174 }
175
176 // G4cout << "tmin= " << cutEnergy << " tmax= " << tmax
177 // << " cross= " << cross << G4endl;
178
179 return cross;
180}
int G4int
Definition: G4Types.hh:66
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)

Referenced by ComputeCrossSectionPerAtom(), and CrossSectionPerVolume().

◆ ComputeDEDXPerVolume()

G4double G4MuBetheBlochModel::ComputeDEDXPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 213 of file G4MuBetheBlochModel.cc.

217{
218 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
219 G4double tau = kineticEnergy/mass;
220 G4double cutEnergy = min(cut,tmax);
221 G4double gam = tau + 1.0;
222 G4double bg2 = tau * (tau+2.0);
223 G4double beta2 = bg2/(gam*gam);
224
225 G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
226 G4double eexc2 = eexc*eexc;
227 //G4double cden = material->GetIonisation()->GetCdensity();
228 //G4double mden = material->GetIonisation()->GetMdensity();
229 //G4double aden = material->GetIonisation()->GetAdensity();
230 //G4double x0den = material->GetIonisation()->GetX0density();
231 //G4double x1den = material->GetIonisation()->GetX1density();
232
233 G4double eDensity = material->GetElectronDensity();
234
235 G4double dedx = log(2.0*electron_mass_c2*bg2*cutEnergy/eexc2)
236 -(1.0 + cutEnergy/tmax)*beta2;
237
238 G4double totEnergy = kineticEnergy + mass;
239 G4double del = 0.5*cutEnergy/totEnergy;
240 dedx += del*del;
241
242 // density correction
243 G4double x = log(bg2)/twoln10;
244 //if ( x >= x0den ) {
245 // dedx -= twoln10*x - cden ;
246 // if ( x < x1den ) dedx -= aden*pow((x1den-x),mden) ;
247 //}
248 dedx -= material->GetIonisation()->DensityCorrection(x);
249
250 // shell correction
251 dedx -= 2.0*corr->ShellCorrection(p,material,kineticEnergy);
252
253 // now compute the total ionization loss
254
255 if (dedx < 0.0) dedx = 0.0 ;
256
257 // radiative corrections of R. Kokoulin
258 if (cutEnergy > limitKinEnergy) {
259
260 G4double logtmax = log(cutEnergy);
261 G4double logstep = logtmax - logLimitKinEnergy;
262 G4double dloss = 0.0;
263 G4double ftot2= 0.5/(totEnergy*totEnergy);
264
265 for (G4int ll=0; ll<8; ll++)
266 {
267 G4double ep = exp(logLimitKinEnergy + xgi[ll]*logstep);
268 G4double a1 = log(1.0 + 2.0*ep/electron_mass_c2);
269 G4double a3 = log(4.0*totEnergy*(totEnergy - ep)/massSquare);
270 dloss += wgi[ll]*(1.0 - beta2*ep/tmax + ep*ep*ftot2)*a1*(a3 - a1);
271 }
272 dedx += dloss*logstep*alphaprime;
273 }
274
275 dedx *= twopi_mc2_rcl2*eDensity/beta2;
276
277 //High order corrections
278 dedx += corr->HighOrderCorrections(p,material,kineticEnergy,cutEnergy);
279
280 return dedx;
281}
G4double HighOrderCorrections(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy, G4double cutEnergy)
G4double ShellCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double GetMeanExcitationEnergy() const
G4double DensityCorrection(G4double x)
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:225
G4double GetElectronDensity() const
Definition: G4Material.hh:216

◆ CrossSectionPerVolume()

G4double G4MuBetheBlochModel::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 198 of file G4MuBetheBlochModel.cc.

204{
205 G4double eDensity = material->GetElectronDensity();
207 (p,kineticEnergy,cutEnergy,maxEnergy);
208 return cross;
209}

◆ Initialise()

void G4MuBetheBlochModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 126 of file G4MuBetheBlochModel.cc.

128{
129 if(p) { SetParticle(p); }
130 if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
131}
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:95

◆ MaxSecondaryEnergy()

G4double G4MuBetheBlochModel::MaxSecondaryEnergy ( const G4ParticleDefinition ,
G4double  kinEnergy 
)
protectedvirtual

Reimplemented from G4VEmModel.

Definition at line 115 of file G4MuBetheBlochModel.cc.

117{
118 G4double tau = kinEnergy/mass;
119 G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
120 (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
121 return tmax;
122}

Referenced by ComputeCrossSectionPerElectron(), and ComputeDEDXPerVolume().

◆ MinEnergyCut()

G4double G4MuBetheBlochModel::MinEnergyCut ( const G4ParticleDefinition ,
const G4MaterialCutsCouple couple 
)
virtual

Definition at line 107 of file G4MuBetheBlochModel.cc.

109{
111}
const G4Material * GetMaterial() const

◆ SampleSecondaries()

void G4MuBetheBlochModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  vdp,
const G4MaterialCutsCouple ,
const G4DynamicParticle dp,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 285 of file G4MuBetheBlochModel.cc.

290{
292 G4double maxKinEnergy = min(maxEnergy,tmax);
293 if(minKinEnergy >= maxKinEnergy) { return; }
294
295 G4double kineticEnergy = dp->GetKineticEnergy();
296 G4double totEnergy = kineticEnergy + mass;
297 G4double etot2 = totEnergy*totEnergy;
298 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
299
300 G4double grej = 1.;
301 if(tmax > limitKinEnergy) {
302 G4double a0 = log(2.*totEnergy/mass);
303 grej += alphaprime*a0*a0;
304 }
305
306 G4double deltaKinEnergy, f;
307
308 // sampling follows ...
309 do {
311 deltaKinEnergy = minKinEnergy*maxKinEnergy
312 /(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
313
314
315 f = 1.0 - beta2*deltaKinEnergy/tmax
316 + 0.5*deltaKinEnergy*deltaKinEnergy/etot2;
317
318 if(deltaKinEnergy > limitKinEnergy) {
319 G4double a1 = log(1.0 + 2.0*deltaKinEnergy/electron_mass_c2);
320 G4double a3 = log(4.0*totEnergy*(totEnergy - deltaKinEnergy)/massSquare);
321 f *= (1. + alphaprime*a1*(a3 - a1));
322 }
323
324 if(f > grej) {
325 G4cout << "G4MuBetheBlochModel::SampleSecondary Warning! "
326 << "Majorant " << grej << " < "
327 << f << " for edelta= " << deltaKinEnergy
328 << " tmin= " << minKinEnergy << " max= " << maxKinEnergy
329 << G4endl;
330 }
331
332
333 } while( grej*G4UniformRand() > f );
334
335 G4double deltaMomentum =
336 sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
337 G4double totalMomentum = totEnergy*sqrt(beta2);
338 G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
339 (deltaMomentum * totalMomentum);
340
341 G4double sint = sqrt(1.0 - cost*cost);
342
343 G4double phi = twopi * G4UniformRand() ;
344
345 G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost) ;
346 G4ThreeVector direction = dp->GetMomentumDirection();
347 deltaDirection.rotateUz(direction);
348
349 // primary change
350 kineticEnergy -= deltaKinEnergy;
351 G4ThreeVector dir = totalMomentum*direction - deltaMomentum*deltaDirection;
352 direction = dir.unit();
353 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
354 fParticleChange->SetProposedMomentumDirection(direction);
355
356 // create G4DynamicParticle object for delta ray
357 G4DynamicParticle* delta = new G4DynamicParticle(theElectron,
358 deltaDirection,deltaKinEnergy);
359 vdp->push_back(delta);
360}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector unit() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:399

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