Geant4 11.1.1
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
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *) override
 
G4double ComputeCrossSectionPerElectron (const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
 
G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
 
G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
G4MuBetheBlochModeloperator= (const G4MuBetheBlochModel &right)=delete
 
 G4MuBetheBlochModel (const G4MuBetheBlochModel &)=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 *, 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 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 *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Member Functions

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 = 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
 
size_t currentCoupleIndex = 0
 
size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Detailed Description

Definition at line 66 of file G4MuBetheBlochModel.hh.

Constructor & Destructor Documentation

◆ G4MuBetheBlochModel() [1/2]

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

Definition at line 76 of file G4MuBetheBlochModel.cc.

78 : G4VEmModel(nam),
79 limitKinEnergy(100.*CLHEP::keV),
80 logLimitKinEnergy(G4Log(limitKinEnergy)),
81 twoln10(2.0*G4Log(10.0)),
82 alphaprime(CLHEP::fine_structure_const/CLHEP::twopi)
83{
84 theElectron = G4Electron::Electron();
86 if(nullptr != p) { SetParticle(p); }
87}
G4double G4Log(G4double x)
Definition: G4Log.hh:227
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()

◆ ~G4MuBetheBlochModel()

G4MuBetheBlochModel::~G4MuBetheBlochModel ( )
default

◆ G4MuBetheBlochModel() [2/2]

G4MuBetheBlochModel::G4MuBetheBlochModel ( const G4MuBetheBlochModel )
delete

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 177 of file G4MuBetheBlochModel.cc.

183{
185 (p,kineticEnergy,cutEnergy,maxEnergy);
186 return cross;
187}
double G4double
Definition: G4Types.hh:83
const G4int Z[17]
G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)

◆ ComputeCrossSectionPerElectron()

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

Definition at line 133 of file G4MuBetheBlochModel.cc.

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

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 206 of file G4MuBetheBlochModel.cc.

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

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 191 of file G4MuBetheBlochModel.cc.

197{
198 G4double eDensity = material->GetElectronDensity();
200 (p,kineticEnergy,cutEnergy,maxEnergy);
201 return cross;
202}

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 122 of file G4MuBetheBlochModel.cc.

124{
125 SetParticle(p);
126 if(nullptr == fParticleChange) {
127 fParticleChange = GetParticleChangeForLoss();
128 }
129}
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:109

◆ MaxSecondaryEnergy()

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

Reimplemented from G4VEmModel.

Definition at line 99 of file G4MuBetheBlochModel.cc.

101{
102 G4double tau = kinEnergy/mass;
103 G4double tmax = 2.0*CLHEP::electron_mass_c2*tau*(tau + 2.) /
104 (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
105 return tmax;
106}

Referenced by ComputeCrossSectionPerElectron(), ComputeDEDXPerVolume(), and SampleSecondaries().

◆ MinEnergyCut()

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

Reimplemented from G4VEmModel.

Definition at line 91 of file G4MuBetheBlochModel.cc.

93{
95}
const G4Material * GetMaterial() const

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 265 of file G4MuBetheBlochModel.cc.

271{
272 G4double kineticEnergy = dp->GetKineticEnergy();
273 G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(), kineticEnergy);
274 G4double maxKinEnergy = std::min(maxEnergy, tmax);
275 if(minKinEnergy >= maxKinEnergy) { return; }
276
277 G4double totEnergy = kineticEnergy + mass;
278 G4double etot2 = totEnergy*totEnergy;
279 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
280
281 G4double grej = 1.;
282 if(tmax > limitKinEnergy) {
283 G4double a0 = G4Log(2.*totEnergy/mass);
284 grej += alphaprime*a0*a0;
285 }
286
287 G4double tkin, f;
288
289 // sampling follows ...
290 do {
292 tkin = minKinEnergy*maxKinEnergy/(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
293 f = 1.0 - beta2*tkin/tmax + 0.5*tkin*tkin/etot2;
294
295 if(tkin > limitKinEnergy) {
296 G4double a1 = G4Log(1.0 + 2.0*tkin/CLHEP::electron_mass_c2);
297 G4double a3 = G4Log(4.0*totEnergy*(totEnergy - tkin)/massSquare);
298 f *= (1. + alphaprime*a1*(a3 - a1));
299 }
300
301 if(f > grej) {
302 G4cout << "G4MuBetheBlochModel::SampleSecondary Warning! "
303 << "Majorant " << grej << " < "
304 << f << " for edelta= " << tkin
305 << " tmin= " << minKinEnergy << " max= " << maxKinEnergy
306 << G4endl;
307 }
308 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
309 } while( grej*G4UniformRand() > f );
310
311 G4double deltaMomentum =
312 std::sqrt(tkin * (tkin + 2.0*CLHEP::electron_mass_c2));
313 G4double totalMomentum = totEnergy*std::sqrt(beta2);
314 G4double cost = tkin * (totEnergy + CLHEP::electron_mass_c2) /
315 (deltaMomentum * totalMomentum);
316
317 G4double sint = std::sqrt(1.0 - cost*cost);
318 G4double phi = CLHEP::twopi * G4UniformRand();
319 G4ThreeVector deltaDirection(sint*std::cos(phi), sint*std::sin(phi), cost);
320 G4ThreeVector direction = dp->GetMomentumDirection();
321 deltaDirection.rotateUz(direction);
322
323 // primary change
324 kineticEnergy -= tkin;
325 G4ThreeVector dir = totalMomentum*direction - deltaMomentum*deltaDirection;
326 direction = dir.unit();
327 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
328 fParticleChange->SetProposedMomentumDirection(direction);
329
330 // create G4DynamicParticle object for delta ray
331 G4DynamicParticle* delta =
332 new G4DynamicParticle(theElectron, deltaDirection, tkin);
333 vdp->push_back(delta);
334}
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
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)

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