Geant4 10.7.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=nullptr, const G4String &nam="MuBetheBloch")
 
 ~G4MuBetheBlochModel ()=default
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *) override
 
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) override
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
- 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
 

Protected Member Functions

virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kinEnergy) override
 
- 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
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 68 of file G4MuBetheBlochModel.hh.

Constructor & Destructor Documentation

◆ G4MuBetheBlochModel()

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

Definition at line 78 of file G4MuBetheBlochModel.cc.

80 : G4VEmModel(nam),
81 particle(nullptr),
82 limitKinEnergy(100.*keV),
83 logLimitKinEnergy(G4Log(limitKinEnergy)),
84 twoln10(2.0*G4Log(10.0)),
85 //bg2lim(0.0169),
86 //taulim(8.4146e-3),
87 alphaprime(fine_structure_const/twopi)
88{
89 theElectron = G4Electron::Electron();
91 fParticleChange = nullptr;
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}
G4double G4Log(G4double x)
Definition: G4Log.hh:226
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()

◆ ~G4MuBetheBlochModel()

G4MuBetheBlochModel::~G4MuBetheBlochModel ( )
default

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 179 of file G4MuBetheBlochModel.cc.

185{
187 (p,kineticEnergy,cutEnergy,maxEnergy);
188 return cross;
189}
double G4double
Definition: G4Types.hh:83
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 130 of file G4MuBetheBlochModel.cc.

135{
136 G4double cross = 0.0;
137 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
138 G4double maxEnergy = min(tmax,maxKinEnergy);
139 if(cutEnergy < maxEnergy) {
140
141 G4double totEnergy = kineticEnergy + mass;
142 G4double energy2 = totEnergy*totEnergy;
143 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
144
145 cross = 1.0/cutEnergy - 1.0/maxEnergy - beta2*G4Log(maxEnergy/cutEnergy)/tmax
146 + 0.5*(maxEnergy - cutEnergy)/energy2;
147
148 // radiative corrections of R. Kokoulin
149 if (maxEnergy > limitKinEnergy) {
150
151 G4double logtmax = G4Log(maxEnergy);
152 G4double logtmin = G4Log(max(cutEnergy,limitKinEnergy));
153 G4double logstep = logtmax - logtmin;
154 G4double dcross = 0.0;
155
156 for (G4int ll=0; ll<8; ll++)
157 {
158 G4double ep = G4Exp(logtmin + xgi[ll]*logstep);
159 G4double a1 = G4Log(1.0 + 2.0*ep/electron_mass_c2);
160 G4double a3 = G4Log(4.0*totEnergy*(totEnergy - ep)/massSquare);
161 dcross += wgi[ll]*(1.0/ep - beta2/tmax + 0.5*ep/energy2)*a1*(a3 - a1);
162 }
163
164 cross += dcross*logstep*alphaprime;
165 }
166
167 cross *= twopi_mc2_rcl2/beta2;
168
169 }
170
171 // G4cout << "tmin= " << cutEnergy << " tmax= " << tmax
172 // << " cross= " << cross << G4endl;
173
174 return cross;
175}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
int G4int
Definition: G4Types.hh:85
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments

Referenced by ComputeCrossSectionPerAtom(), and CrossSectionPerVolume().

◆ ComputeDEDXPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 208 of file G4MuBetheBlochModel.cc.

212{
213 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
214 G4double tau = kineticEnergy/mass;
215 G4double cutEnergy = min(cut,tmax);
216 G4double gam = tau + 1.0;
217 G4double bg2 = tau * (tau+2.0);
218 G4double beta2 = bg2/(gam*gam);
219
220 G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
221 G4double eexc2 = eexc*eexc;
222
223 G4double eDensity = material->GetElectronDensity();
224
225 G4double dedx = G4Log(2.0*electron_mass_c2*bg2*cutEnergy/eexc2)
226 -(1.0 + cutEnergy/tmax)*beta2;
227
228 G4double totEnergy = kineticEnergy + mass;
229 G4double del = 0.5*cutEnergy/totEnergy;
230 dedx += del*del;
231
232 // density correction
233 G4double x = G4Log(bg2)/twoln10;
234 dedx -= material->GetIonisation()->DensityCorrection(x);
235
236 // shell correction
237 dedx -= 2.0*corr->ShellCorrection(p,material,kineticEnergy);
238
239 // now compute the total ionization loss
240
241 if (dedx < 0.0) dedx = 0.0 ;
242
243 // radiative corrections of R. Kokoulin
244 if (cutEnergy > limitKinEnergy) {
245
246 G4double logtmax = G4Log(cutEnergy);
247 G4double logstep = logtmax - logLimitKinEnergy;
248 G4double dloss = 0.0;
249 G4double ftot2= 0.5/(totEnergy*totEnergy);
250
251 for (G4int ll=0; ll<8; ll++)
252 {
253 G4double ep = G4Exp(logLimitKinEnergy + xgi[ll]*logstep);
254 G4double a1 = G4Log(1.0 + 2.0*ep/electron_mass_c2);
255 G4double a3 = G4Log(4.0*totEnergy*(totEnergy - ep)/massSquare);
256 dloss += wgi[ll]*(1.0 - beta2*ep/tmax + ep*ep*ftot2)*a1*(a3 - a1);
257 }
258 dedx += dloss*logstep*alphaprime;
259 }
260
261 dedx *= twopi_mc2_rcl2*eDensity/beta2;
262
263 //High order corrections
264 dedx += corr->HighOrderCorrections(p,material,kineticEnergy,cutEnergy);
265
266 return dedx;
267}
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:224
G4double GetElectronDensity() const
Definition: G4Material.hh:215

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 193 of file G4MuBetheBlochModel.cc.

199{
200 G4double eDensity = material->GetElectronDensity();
202 (p,kineticEnergy,cutEnergy,maxEnergy);
203 return cross;
204}

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 121 of file G4MuBetheBlochModel.cc.

123{
124 if(p) { SetParticle(p); }
125 if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
126}
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:118

◆ MaxSecondaryEnergy()

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

Reimplemented from G4VEmModel.

Definition at line 110 of file G4MuBetheBlochModel.cc.

112{
113 G4double tau = kinEnergy/mass;
114 G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
115 (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
116 return tmax;
117}

Referenced by ComputeCrossSectionPerElectron(), and ComputeDEDXPerVolume().

◆ MinEnergyCut()

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

Reimplemented from G4VEmModel.

Definition at line 102 of file G4MuBetheBlochModel.cc.

104{
106}
const G4Material * GetMaterial() const

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 271 of file G4MuBetheBlochModel.cc.

276{
278 G4double maxKinEnergy = min(maxEnergy,tmax);
279 if(minKinEnergy >= maxKinEnergy) { return; }
280
281 G4double kineticEnergy = dp->GetKineticEnergy();
282 G4double totEnergy = kineticEnergy + mass;
283 G4double etot2 = totEnergy*totEnergy;
284 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
285
286 G4double grej = 1.;
287 if(tmax > limitKinEnergy) {
288 G4double a0 = G4Log(2.*totEnergy/mass);
289 grej += alphaprime*a0*a0;
290 }
291
292 G4double deltaKinEnergy, f;
293
294 // sampling follows ...
295 do {
297 deltaKinEnergy = minKinEnergy*maxKinEnergy
298 /(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
299
300
301 f = 1.0 - beta2*deltaKinEnergy/tmax
302 + 0.5*deltaKinEnergy*deltaKinEnergy/etot2;
303
304 if(deltaKinEnergy > limitKinEnergy) {
305 G4double a1 = G4Log(1.0 + 2.0*deltaKinEnergy/electron_mass_c2);
306 G4double a3 = G4Log(4.0*totEnergy*(totEnergy - deltaKinEnergy)/massSquare);
307 f *= (1. + alphaprime*a1*(a3 - a1));
308 }
309
310 if(f > grej) {
311 G4cout << "G4MuBetheBlochModel::SampleSecondary Warning! "
312 << "Majorant " << grej << " < "
313 << f << " for edelta= " << deltaKinEnergy
314 << " tmin= " << minKinEnergy << " max= " << maxKinEnergy
315 << G4endl;
316 }
317 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
318 } while( grej*G4UniformRand() > f );
319
320 G4double deltaMomentum =
321 sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
322 G4double totalMomentum = totEnergy*sqrt(beta2);
323 G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
324 (deltaMomentum * totalMomentum);
325
326 G4double sint = sqrt(1.0 - cost*cost);
327
328 G4double phi = twopi * G4UniformRand() ;
329
330 G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost) ;
331 G4ThreeVector direction = dp->GetMomentumDirection();
332 deltaDirection.rotateUz(direction);
333
334 // primary change
335 kineticEnergy -= deltaKinEnergy;
336 G4ThreeVector dir = totalMomentum*direction - deltaMomentum*deltaDirection;
337 direction = dir.unit();
338 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
339 fParticleChange->SetProposedMomentumDirection(direction);
340
341 // create G4DynamicParticle object for delta ray
342 G4DynamicParticle* delta = new G4DynamicParticle(theElectron,
343 deltaDirection,deltaKinEnergy);
344 vdp->push_back(delta);
345}
const G4double a0
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:510

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