Geant4 11.2.2
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 InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void 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 void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
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)
 
void SetLPMFlag (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 ()
 
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
 
std::size_t currentCoupleIndex = 0
 
std::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 77 of file G4MuBetheBlochModel.cc.

79 : G4VEmModel(nam),
80 limitRadCorrection(250.*CLHEP::MeV),
81 limitKinEnergy(100.*CLHEP::keV),
82 logLimitKinEnergy(G4Log(limitKinEnergy)),
83 twoln10(2.0*G4Log(10.0)),
84 alphaprime(CLHEP::fine_structure_const/CLHEP::twopi)
85{
86 theElectron = G4Electron::Electron();
88 if(nullptr != p) { SetParticle(p); }
89}
G4double G4Log(G4double x)
Definition G4Log.hh:227
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67

◆ ~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 182 of file G4MuBetheBlochModel.cc.

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

143{
144 G4double cross = 0.0;
145 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
146 G4double maxEnergy = std::min(tmax, maxKinEnergy);
147 if(cutEnergy < maxEnergy) {
148
149 G4double totEnergy = kineticEnergy + mass;
150 G4double energy2 = totEnergy*totEnergy;
151 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
152
153 cross = 1.0/cutEnergy - 1.0/maxEnergy -
154 beta2*G4Log(maxEnergy/cutEnergy)/tmax +
155 0.5*(maxEnergy - cutEnergy)/energy2;
156
157 // radiative corrections of R. Kokoulin
158 if (maxEnergy > limitKinEnergy && kineticEnergy > limitRadCorrection) {
159
160 G4double logtmax = G4Log(maxEnergy);
161 G4double logtmin = G4Log(std::max(cutEnergy,limitKinEnergy));
162 G4double logstep = logtmax - logtmin;
163 G4double dcross = 0.0;
164
165 for (G4int ll=0; ll<8; ++ll) {
166 G4double ep = G4Exp(logtmin + xgi[ll]*logstep);
167 G4double a1 = G4Log(1.0 + 2.0*ep/CLHEP::electron_mass_c2);
168 G4double a3 = G4Log(4.0*totEnergy*(totEnergy - ep)/massSquare);
169 dcross += wgi[ll]*(1.0/ep - beta2/tmax + 0.5*ep/energy2)*a1*(a3 - a1);
170 }
171 cross += dcross*logstep*alphaprime;
172 }
173 cross *= CLHEP::twopi_mc2_rcl2/beta2;
174 }
175 // G4cout << "tmin= " << cutEnergy << " tmax= " << tmax
176 // << " cross= " << cross << G4endl;
177 return cross;
178}
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 211 of file G4MuBetheBlochModel.cc.

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

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 196 of file G4MuBetheBlochModel.cc.

202{
203 G4double eDensity = material->GetElectronDensity();
205 (p,kineticEnergy,cutEnergy,maxEnergy);
206 return cross;
207}

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 124 of file G4MuBetheBlochModel.cc.

126{
127 SetParticle(p);
128 if(nullptr == fParticleChange) {
129 fParticleChange = GetParticleChangeForLoss();
130 if(UseAngularGeneratorFlag() && nullptr == GetAngularDistribution()) {
132 }
133 }
134}
G4VEmAngularDistribution * GetAngularDistribution()
void SetAngularDistribution(G4VEmAngularDistribution *)
G4bool UseAngularGeneratorFlag() const
G4ParticleChangeForLoss * GetParticleChangeForLoss()

◆ MaxSecondaryEnergy()

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

Reimplemented from G4VEmModel.

Definition at line 101 of file G4MuBetheBlochModel.cc.

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

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

◆ MinEnergyCut()

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

Reimplemented from G4VEmModel.

Definition at line 93 of file G4MuBetheBlochModel.cc.

95{
97}
const G4Material * GetMaterial() const

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 268 of file G4MuBetheBlochModel.cc.

274{
275 G4double kineticEnergy = dp->GetKineticEnergy();
276 G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(), kineticEnergy);
277 G4double maxKinEnergy = std::min(maxEnergy, tmax);
278 if(minKinEnergy >= maxKinEnergy) { return; }
279
280 G4double totEnergy = kineticEnergy + mass;
281 G4double etot2 = totEnergy*totEnergy;
282 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
283
284 G4double grej = 1.;
285 G4bool radC = (tmax > limitKinEnergy && kineticEnergy > limitRadCorrection);
286 if(radC) {
287 G4double a0 = G4Log(2.*totEnergy/mass);
288 grej += alphaprime*a0*a0;
289 }
290
291 G4double tkin, f;
292
293 // sampling follows ...
294 do {
296 tkin = minKinEnergy*maxKinEnergy/(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
297 f = 1.0 - beta2*tkin/tmax + 0.5*tkin*tkin/etot2;
298
299 if(radC && tkin > limitKinEnergy) {
300 G4double a1 = G4Log(1.0 + 2.0*tkin/CLHEP::electron_mass_c2);
301 G4double a3 = G4Log(4.0*totEnergy*(totEnergy - tkin)/massSquare);
302 f *= (1. + alphaprime*a1*(a3 - a1));
303 }
304
305 if(f > grej) {
306 G4cout << "G4MuBetheBlochModel::SampleSecondary Warning! "
307 << "Majorant " << grej << " < "
308 << f << " for edelta= " << tkin
309 << " tmin= " << minKinEnergy << " max= " << maxKinEnergy
310 << G4endl;
311 }
312 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
313 } while( grej*G4UniformRand() > f );
314
315 G4ThreeVector deltaDirection;
316
318 const G4Material* mat = couple->GetMaterial();
319 deltaDirection = GetAngularDistribution()->SampleDirection(dp, tkin,
320 SelectRandomAtomNumber(mat), mat);
321 } else {
322
323 G4double deltaMom = std::sqrt(tkin * (tkin + 2.0*CLHEP::electron_mass_c2));
324 G4double totalMom = totEnergy*std::sqrt(beta2);
325 G4double cost = tkin * (totEnergy + CLHEP::electron_mass_c2) /
326 (deltaMom * totalMom);
327 cost = std::min(cost, 1.0);
328 const G4double sint = std::sqrt((1.0 - cost)*(1.0 + cost));
329 const G4double phi = twopi*G4UniformRand();
330
331 deltaDirection.set(sint*std::cos(phi),sint*std::sin(phi), cost) ;
332 deltaDirection.rotateUz(dp->GetMomentumDirection());
333 }
334 // create G4DynamicParticle object for delta ray
335 auto delta = new G4DynamicParticle(theElectron, deltaDirection, tkin);
336 vdp->push_back(delta);
337
338 // primary change
339 kineticEnergy -= tkin;
340 G4ThreeVector dir = dp->GetMomentum() - delta->GetMomentum();
341 dir = dir.unit();
342 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
343 fParticleChange->SetProposedMomentumDirection(dir);
344}
const G4double a0
bool G4bool
Definition G4Types.hh:86
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector unit() const
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4int SelectRandomAtomNumber(const G4Material *) const

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