Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions
G4MuBetheBlochModel Class Reference

#include <G4MuBetheBlochModel.hh>

Inheritance diagram for G4MuBetheBlochModel:
G4VEmModel

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 InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
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 *)
 
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)
 
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=0)
 
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 SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
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 ()
 
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
 

Detailed Description

Definition at line 69 of file G4MuBetheBlochModel.hh.

Constructor & Destructor Documentation

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

Definition at line 79 of file G4MuBetheBlochModel.cc.

References G4Electron::Electron(), G4LossTableManager::EmCorrections(), and G4LossTableManager::Instance().

81  : G4VEmModel(nam),
82  particle(0),
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 = 0;
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:65
G4EmCorrections * EmCorrections()
G4double G4Log(G4double x)
Definition: G4Log.hh:227
static G4Electron * Electron()
Definition: G4Electron.cc:94
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 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 185 of file G4MuBetheBlochModel.cc.

References ComputeCrossSectionPerElectron().

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
G4double G4MuBetheBlochModel::ComputeCrossSectionPerElectron ( const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
)
virtual

Definition at line 136 of file G4MuBetheBlochModel.cc.

References dcross(), python.hepunit::electron_mass_c2, G4Exp(), G4Log(), G4INCL::Math::max(), MaxSecondaryEnergy(), G4INCL::Math::min(), and python.hepunit::twopi_mc2_rcl2.

Referenced by ComputeCrossSectionPerAtom(), and CrossSectionPerVolume().

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 }
int G4int
Definition: G4Types.hh:78
float electron_mass_c2
Definition: hepunit.py:274
function dcross(V1, V2)
Definition: leptonew.f:3567
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
T max(const T t1, const T t2)
brief Return the largest of the two arguments
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
double G4double
Definition: G4Types.hh:76
G4double G4MuBetheBlochModel::ComputeDEDXPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 214 of file G4MuBetheBlochModel.cc.

References G4IonisParamMat::DensityCorrection(), python.hepunit::electron_mass_c2, G4Exp(), G4Log(), G4InuclParticleNames::gam, G4Material::GetElectronDensity(), G4Material::GetIonisation(), G4IonisParamMat::GetMeanExcitationEnergy(), G4EmCorrections::HighOrderCorrections(), MaxSecondaryEnergy(), G4INCL::Math::min(), G4EmCorrections::ShellCorrection(), python.hepunit::twopi_mc2_rcl2, and test::x.

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  //G4double cden = material->GetIonisation()->GetCdensity();
229  //G4double mden = material->GetIonisation()->GetMdensity();
230  //G4double aden = material->GetIonisation()->GetAdensity();
231  //G4double x0den = material->GetIonisation()->GetX0density();
232  //G4double x1den = material->GetIonisation()->GetX1density();
233 
234  G4double eDensity = material->GetElectronDensity();
235 
236  G4double dedx = G4Log(2.0*electron_mass_c2*bg2*cutEnergy/eexc2)
237  -(1.0 + cutEnergy/tmax)*beta2;
238 
239  G4double totEnergy = kineticEnergy + mass;
240  G4double del = 0.5*cutEnergy/totEnergy;
241  dedx += del*del;
242 
243  // density correction
244  G4double x = G4Log(bg2)/twoln10;
245  //if ( x >= x0den ) {
246  // dedx -= twoln10*x - cden ;
247  // if ( x < x1den ) dedx -= aden*pow((x1den-x),mden) ;
248  //}
249  dedx -= material->GetIonisation()->DensityCorrection(x);
250 
251  // shell correction
252  dedx -= 2.0*corr->ShellCorrection(p,material,kineticEnergy);
253 
254  // now compute the total ionization loss
255 
256  if (dedx < 0.0) dedx = 0.0 ;
257 
258  // radiative corrections of R. Kokoulin
259  if (cutEnergy > limitKinEnergy) {
260 
261  G4double logtmax = G4Log(cutEnergy);
262  G4double logstep = logtmax - logLimitKinEnergy;
263  G4double dloss = 0.0;
264  G4double ftot2= 0.5/(totEnergy*totEnergy);
265 
266  for (G4int ll=0; ll<8; ll++)
267  {
268  G4double ep = G4Exp(logLimitKinEnergy + xgi[ll]*logstep);
269  G4double a1 = G4Log(1.0 + 2.0*ep/electron_mass_c2);
270  G4double a3 = G4Log(4.0*totEnergy*(totEnergy - ep)/massSquare);
271  dloss += wgi[ll]*(1.0 - beta2*ep/tmax + ep*ep*ftot2)*a1*(a3 - a1);
272  }
273  dedx += dloss*logstep*alphaprime;
274  }
275 
276  dedx *= twopi_mc2_rcl2*eDensity/beta2;
277 
278  //High order corrections
279  dedx += corr->HighOrderCorrections(p,material,kineticEnergy,cutEnergy);
280 
281  return dedx;
282 }
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:224
G4double HighOrderCorrections(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy, G4double cutEnergy)
int G4int
Definition: G4Types.hh:78
G4double GetElectronDensity() const
Definition: G4Material.hh:215
float electron_mass_c2
Definition: hepunit.py:274
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
G4double ShellCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double DensityCorrection(G4double x)
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
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
G4double G4MuBetheBlochModel::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 199 of file G4MuBetheBlochModel.cc.

References ComputeCrossSectionPerElectron(), and G4Material::GetElectronDensity().

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:215
double G4double
Definition: G4Types.hh:76
void G4MuBetheBlochModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 127 of file G4MuBetheBlochModel.cc.

References G4VEmModel::GetParticleChangeForLoss().

129 {
130  if(p) { SetParticle(p); }
131  if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
132 }
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:107
G4double G4MuBetheBlochModel::MaxSecondaryEnergy ( const G4ParticleDefinition ,
G4double  kinEnergy 
)
protectedvirtual

Reimplemented from G4VEmModel.

Definition at line 116 of file G4MuBetheBlochModel.cc.

References python.hepunit::electron_mass_c2.

Referenced by ComputeCrossSectionPerElectron(), and ComputeDEDXPerVolume().

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
G4double G4MuBetheBlochModel::MinEnergyCut ( const G4ParticleDefinition ,
const G4MaterialCutsCouple couple 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 108 of file G4MuBetheBlochModel.cc.

References G4Material::GetIonisation(), G4MaterialCutsCouple::GetMaterial(), and G4IonisParamMat::GetMeanExcitationEnergy().

110 {
111  return couple->GetMaterial()->GetIonisation()->GetMeanExcitationEnergy();
112 }
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:224
G4double GetMeanExcitationEnergy() const
const G4Material * GetMaterial() const
void G4MuBetheBlochModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  vdp,
const G4MaterialCutsCouple ,
const G4DynamicParticle dp,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 286 of file G4MuBetheBlochModel.cc.

References python.hepunit::electron_mass_c2, G4cout, G4endl, G4Log(), G4UniformRand, G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), G4VEmModel::MaxSecondaryKinEnergy(), G4INCL::Math::min(), CLHEP::Hep3Vector::rotateUz(), G4ParticleChangeForLoss::SetProposedKineticEnergy(), G4ParticleChangeForLoss::SetProposedMomentumDirection(), python.hepunit::twopi, and CLHEP::Hep3Vector::unit().

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

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