Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4MuBetheBlochModel Class Reference

#include <G4MuBetheBlochModel.hh>

Inheritance diagram for G4MuBetheBlochModel:
Collaboration diagram for G4MuBetheBlochModel:

Public Member Functions

 G4MuBetheBlochModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="MuBetheBloch")
 
virtual ~G4MuBetheBlochModel ()
 
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 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 *, 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 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)
 
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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
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
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
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 SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () 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)
 

Protected Member Functions

virtual 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
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 69 of file G4MuBetheBlochModel.hh.

Constructor & Destructor Documentation

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

Definition at line 79 of file G4MuBetheBlochModel.cc.

81  : G4VEmModel(nam),
82  particle(nullptr),
83  limitKinEnergy(100.*keV),
84  logLimitKinEnergy(G4Log(limitKinEnergy)),
85  twoln10(2.0*G4Log(10.0)),
86  //bg2lim(0.0169),
87  //taulim(8.4146e-3),
88  alphaprime(fine_structure_const/twopi)
89 {
90  theElectron = G4Electron::Electron();
92  fParticleChange = nullptr;
93 
94  // initial initialisation of memeber should be overwritten
95  // by SetParticle
96  mass = massSquare = ratio = 1.0;
97 
98  if(p) { SetParticle(p); }
99 }
static G4LossTableManager * Instance()
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
static constexpr double twopi
Definition: G4SIunits.hh:76
G4EmCorrections * EmCorrections()
G4double G4Log(G4double x)
Definition: G4Log.hh:230
static G4Electron * Electron()
Definition: G4Electron.cc:94
static constexpr double keV
Definition: G4SIunits.hh:216

Here is the call graph for this function:

G4MuBetheBlochModel::~G4MuBetheBlochModel ( )
virtual

Definition at line 103 of file G4MuBetheBlochModel.cc.

104 {}

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 185 of file G4MuBetheBlochModel.cc.

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

Here is the call graph for this function:

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

Definition at line 136 of file G4MuBetheBlochModel.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 214 of file G4MuBetheBlochModel.cc.

218 {
219  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
220  G4double tau = kineticEnergy/mass;
221  G4double cutEnergy = min(cut,tmax);
222  G4double gam = tau + 1.0;
223  G4double bg2 = tau * (tau+2.0);
224  G4double beta2 = bg2/(gam*gam);
225 
226  G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
227  G4double eexc2 = eexc*eexc;
228 
229  G4double eDensity = material->GetElectronDensity();
230 
231  G4double dedx = G4Log(2.0*electron_mass_c2*bg2*cutEnergy/eexc2)
232  -(1.0 + cutEnergy/tmax)*beta2;
233 
234  G4double totEnergy = kineticEnergy + mass;
235  G4double del = 0.5*cutEnergy/totEnergy;
236  dedx += del*del;
237 
238  // density correction
239  G4double x = G4Log(bg2)/twoln10;
240  dedx -= material->GetIonisation()->DensityCorrection(x);
241 
242  // shell correction
243  dedx -= 2.0*corr->ShellCorrection(p,material,kineticEnergy);
244 
245  // now compute the total ionization loss
246 
247  if (dedx < 0.0) dedx = 0.0 ;
248 
249  // radiative corrections of R. Kokoulin
250  if (cutEnergy > limitKinEnergy) {
251 
252  G4double logtmax = G4Log(cutEnergy);
253  G4double logstep = logtmax - logLimitKinEnergy;
254  G4double dloss = 0.0;
255  G4double ftot2= 0.5/(totEnergy*totEnergy);
256 
257  for (G4int ll=0; ll<8; ll++)
258  {
259  G4double ep = G4Exp(logLimitKinEnergy + xgi[ll]*logstep);
260  G4double a1 = G4Log(1.0 + 2.0*ep/electron_mass_c2);
261  G4double a3 = G4Log(4.0*totEnergy*(totEnergy - ep)/massSquare);
262  dloss += wgi[ll]*(1.0 - beta2*ep/tmax + ep*ep*ftot2)*a1*(a3 - a1);
263  }
264  dedx += dloss*logstep*alphaprime;
265  }
266 
267  dedx *= twopi_mc2_rcl2*eDensity/beta2;
268 
269  //High order corrections
270  dedx += corr->HighOrderCorrections(p,material,kineticEnergy,cutEnergy);
271 
272  return dedx;
273 }
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:226
G4double HighOrderCorrections(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy, G4double cutEnergy)
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override
tuple x
Definition: test.py:50
int G4int
Definition: G4Types.hh:78
G4double GetElectronDensity() const
Definition: G4Material.hh:217
float electron_mass_c2
Definition: hepunit.py:274
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4double ShellCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double DensityCorrection(G4double x)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double GetMeanExcitationEnergy() const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 199 of file G4MuBetheBlochModel.cc.

205 {
206  G4double eDensity = material->GetElectronDensity();
208  (p,kineticEnergy,cutEnergy,maxEnergy);
209  return cross;
210 }
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4double GetElectronDensity() const
Definition: G4Material.hh:217
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 127 of file G4MuBetheBlochModel.cc.

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

Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 116 of file G4MuBetheBlochModel.cc.

118 {
119  G4double tau = kinEnergy/mass;
120  G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
121  (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
122  return tmax;
123 }
float electron_mass_c2
Definition: hepunit.py:274
double G4double
Definition: G4Types.hh:76

Here is the caller graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 108 of file G4MuBetheBlochModel.cc.

110 {
111  return couple->GetMaterial()->GetIonisation()->GetMeanExcitationEnergy();
112 }
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:226
G4double GetMeanExcitationEnergy() const
const G4Material * GetMaterial() const

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 277 of file G4MuBetheBlochModel.cc.

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

Here is the call graph for this function:


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