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

#include <G4BetheBlochModel.hh>

Inheritance diagram for G4BetheBlochModel:
G4VEmModel G4BetheBlochIonGasModel G4BetheBlochNoDeltaModel

Public Member Functions

 G4BetheBlochModel (const G4ParticleDefinition *p=0, const G4String &nam="BetheBloch")
 
virtual ~G4BetheBlochModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *couple)
 
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 G4double GetChargeSquareRatio (const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double &eloss, G4double &, G4double length)
 
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 void StartTracking (G4Track *)
 
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)
 
G4double GetChargeSquareRatio () const
 
void SetChargeSquareRatio (G4double val)
 
- 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 72 of file G4BetheBlochModel.hh.

Constructor & Destructor Documentation

G4BetheBlochModel::G4BetheBlochModel ( const G4ParticleDefinition p = 0,
const G4String nam = "BetheBloch" 
)

Definition at line 75 of file G4BetheBlochModel.cc.

References G4Electron::Electron(), G4LossTableManager::EmCorrections(), G4NistManager::Instance(), G4LossTableManager::Instance(), python.hepunit::MeV, and G4VEmModel::SetLowEnergyLimit().

77  : G4VEmModel(nam),
78  particle(0),
79  tlimit(DBL_MAX),
80  twoln10(2.0*G4Log(10.0)),
81  bg2lim(0.0169),
82  taulim(8.4146e-3),
83  isIon(false),
84  isInitialised(false)
85 {
86  fParticleChange = 0;
87  theElectron = G4Electron::Electron();
88  if(p) {
89  SetGenericIon(p);
90  SetParticle(p);
91  } else {
92  SetParticle(theElectron);
93  }
95  nist = G4NistManager::Instance();
97 }
static G4LossTableManager * Instance()
static G4NistManager * 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
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
#define DBL_MAX
Definition: templates.hh:83
G4BetheBlochModel::~G4BetheBlochModel ( )
virtual

Definition at line 101 of file G4BetheBlochModel.cc.

102 {}

Member Function Documentation

G4double G4BetheBlochModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  Z,
G4double  A,
G4double  cutEnergy,
G4double  maxEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 226 of file G4BetheBlochModel.cc.

References ComputeCrossSectionPerElectron().

232 {
234  (p,kineticEnergy,cutEnergy,maxEnergy);
235  return cross;
236 }
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
double G4double
Definition: G4Types.hh:76
G4double G4BetheBlochModel::ComputeCrossSectionPerElectron ( const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
)
virtual

Definition at line 189 of file G4BetheBlochModel.cc.

References G4Log(), MaxSecondaryEnergy(), G4INCL::Math::min(), and python.hepunit::twopi_mc2_rcl2.

Referenced by ComputeCrossSectionPerAtom(), and CrossSectionPerVolume().

193 {
194  G4double cross = 0.0;
195  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
196  G4double maxEnergy = min(tmax,maxKinEnergy);
197  if(cutEnergy < maxEnergy) {
198 
199  G4double totEnergy = kineticEnergy + mass;
200  G4double energy2 = totEnergy*totEnergy;
201  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
202 
203  cross = 1.0/cutEnergy - 1.0/maxEnergy
204  - beta2*G4Log(maxEnergy/cutEnergy)/tmax;
205 
206  // +term for spin=1/2 particle
207  if( 0.5 == spin ) { cross += 0.5*(maxEnergy - cutEnergy)/energy2; }
208 
209  // High order correction different for hadrons and ions
210  // nevetheless they are applied to reduce high energy transfers
211  // if(!isIon)
212  //cross += corr->FiniteSizeCorrectionXS(p,currentMaterial,
213  // kineticEnergy,cutEnergy);
214 
215  cross *= twopi_mc2_rcl2*chargeSquare/beta2;
216  }
217 
218  // G4cout << "BB: e= " << kineticEnergy << " tmin= " << cutEnergy
219  // << " tmax= " << tmax << " cross= " << cross << G4endl;
220 
221  return cross;
222 }
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
G4double G4Log(G4double x)
Definition: G4Log.hh:227
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
double G4double
Definition: G4Types.hh:76
G4double G4BetheBlochModel::ComputeDEDXPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy 
)
virtual

Reimplemented from G4VEmModel.

Reimplemented in G4BetheBlochNoDeltaModel.

Definition at line 255 of file G4BetheBlochModel.cc.

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

Referenced by G4BetheBlochNoDeltaModel::ComputeDEDXPerVolume().

259 {
260  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
261  G4double cutEnergy = std::min(cut,tmax);
262 
263  G4double tau = kineticEnergy/mass;
264  G4double gam = tau + 1.0;
265  G4double bg2 = tau * (tau+2.0);
266  G4double beta2 = bg2/(gam*gam);
267 
268  G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
269  G4double eexc2 = eexc*eexc;
270 
271  G4double eDensity = material->GetElectronDensity();
272 
273  G4double dedx = G4Log(2.0*electron_mass_c2*bg2*cutEnergy/eexc2)
274  - (1.0 + cutEnergy/tmax)*beta2;
275 
276  if(0.5 == spin) {
277  G4double del = 0.5*cutEnergy/(kineticEnergy + mass);
278  dedx += del*del;
279  }
280 
281  // density correction
282  G4double x = G4Log(bg2)/twoln10;
283  dedx -= material->GetIonisation()->DensityCorrection(x);
284 
285  // shell correction
286  dedx -= 2.0*corr->ShellCorrection(p,material,kineticEnergy);
287 
288  // now compute the total ionization loss
289  dedx *= twopi_mc2_rcl2*chargeSquare*eDensity/beta2;
290 
291  //High order correction different for hadrons and ions
292  if(isIon) {
293  dedx += corr->IonBarkasCorrection(p,material,kineticEnergy);
294  } else {
295  dedx += corr->HighOrderCorrections(p,material,kineticEnergy,cutEnergy);
296  }
297 
298  if (dedx < 0.0) { dedx = 0.0; }
299 
300  //G4cout << "E(MeV)= " << kineticEnergy/MeV << " dedx= " << dedx
301  // << " " << material->GetName() << G4endl;
302 
303  return dedx;
304 }
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:224
G4double HighOrderCorrections(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy, G4double cutEnergy)
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
G4double GetElectronDensity() const
Definition: G4Material.hh:215
G4double IonBarkasCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
float electron_mass_c2
Definition: hepunit.py:274
G4double G4Log(G4double x)
Definition: G4Log.hh:227
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
void G4BetheBlochModel::CorrectionsAlongStep ( const G4MaterialCutsCouple couple,
const G4DynamicParticle dp,
G4double eloss,
G4double ,
G4double  length 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 308 of file G4BetheBlochModel.cc.

References G4EmCorrections::EffectiveChargeCorrection(), G4EmCorrections::EffectiveChargeSquareRatio(), G4DynamicParticle::GetDefinition(), G4DynamicParticle::GetKineticEnergy(), G4MaterialCutsCouple::GetMaterial(), G4VEmModel::GetModelOfFluctuations(), G4EmCorrections::IonHighOrderCorrections(), and G4VEmFluctuationModel::SetParticleAndCharge().

313 {
314  if(isIon) {
315  const G4ParticleDefinition* p = dp->GetDefinition();
316  const G4Material* mat = couple->GetMaterial();
317  G4double preKinEnergy = dp->GetKineticEnergy();
318  G4double e = preKinEnergy - eloss*0.5;
319  if(e < preKinEnergy*0.75) { e = preKinEnergy*0.75; }
320 
321  G4double q2 = corr->EffectiveChargeSquareRatio(p,mat,e);
323  G4double qfactor = q2*corr->EffectiveChargeCorrection(p,mat,e)/corrFactor;
324  G4double highOrder = length*corr->IonHighOrderCorrections(p,couple,e);
325  G4double elossnew = eloss*qfactor + highOrder;
326  if(elossnew > preKinEnergy) { elossnew = preKinEnergy; }
327  else if(elossnew < eloss*0.5) { elossnew = eloss*0.5; }
328  eloss = elossnew;
329  //G4cout << "G4BetheBlochModel::CorrectionsAlongStep: e= " << preKinEnergy
330  // << " qfactor= " << qfactor
331  // << " highOrder= " << highOrder << " ("
332  // << highOrder/eloss << ")" << G4endl;
333  }
334 }
G4double GetKineticEnergy() const
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
const char * p
Definition: xmltok.h:285
G4double EffectiveChargeCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4ParticleDefinition * GetDefinition() const
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:571
G4double IonHighOrderCorrections(const G4ParticleDefinition *, const G4MaterialCutsCouple *, G4double kineticEnergy)
double G4double
Definition: G4Types.hh:76
virtual void SetParticleAndCharge(const G4ParticleDefinition *, G4double q2)
const G4Material * GetMaterial() const
G4double G4BetheBlochModel::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
)
virtual

Reimplemented from G4VEmModel.

Reimplemented in G4BetheBlochNoDeltaModel.

Definition at line 240 of file G4BetheBlochModel.cc.

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

246 {
247  G4double eDensity = material->GetElectronDensity();
249  (p,kineticEnergy,cutEnergy,maxEnergy);
250  return cross;
251 }
G4double GetElectronDensity() const
Definition: G4Material.hh:215
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
double G4double
Definition: G4Types.hh:76
G4double G4BetheBlochModel::GetChargeSquareRatio ( const G4ParticleDefinition p,
const G4Material mat,
G4double  kineticEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 130 of file G4BetheBlochModel.cc.

References G4EmCorrections::EffectiveChargeCorrection(), and G4EmCorrections::EffectiveChargeSquareRatio().

133 {
134  // this method is called only for ions
135  G4double q2 = corr->EffectiveChargeSquareRatio(p,mat,kineticEnergy);
136  corrFactor = q2*corr->EffectiveChargeCorrection(p,mat,kineticEnergy);
137  return corrFactor;
138 }
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double EffectiveChargeCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
double G4double
Definition: G4Types.hh:76
G4double G4BetheBlochModel::GetChargeSquareRatio ( ) const
inlineprotected

Definition at line 196 of file G4BetheBlochModel.hh.

197 {
198  return chargeSquare;
199 }
G4double G4BetheBlochModel::GetParticleCharge ( const G4ParticleDefinition p,
const G4Material mat,
G4double  kineticEnergy 
)
virtual

Reimplemented from G4VEmModel.

Reimplemented in G4BetheBlochIonGasModel.

Definition at line 142 of file G4BetheBlochModel.cc.

References G4EmCorrections::GetParticleCharge().

145 {
146  //G4cout<<"G4BetheBlochModel::GetParticleCharge e= "<<kineticEnergy <<
147  // " q= " << corr->GetParticleCharge(p,mat,kineticEnergy) <<G4endl;
148  // this method is called only for ions, so no check if it is an ion
149  return corr->GetParticleCharge(p,mat,kineticEnergy);
150 }
G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
void G4BetheBlochModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 106 of file G4BetheBlochModel.cc.

References G4VEmModel::GetAngularDistribution(), G4VEmModel::GetParticleChangeForLoss(), G4VEmModel::SetAngularDistribution(), G4VEmModel::SetDeexcitationFlag(), and G4VEmModel::UseAngularGeneratorFlag().

108 {
109  SetGenericIon(p);
110  SetParticle(p);
111 
112  //G4cout << "G4BetheBlochModel::Initialise for " << p->GetParticleName()
113  // << " isIon= " << isIon
114  // << G4endl;
115 
116  // always false before the run
117  SetDeexcitationFlag(false);
118 
119  if(!isInitialised) {
120  isInitialised = true;
121  fParticleChange = GetParticleChangeForLoss();
124  }
125  }
126 }
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:107
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:578
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:655
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:585
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:739
G4double G4BetheBlochModel::MaxSecondaryEnergy ( const G4ParticleDefinition pd,
G4double  kinEnergy 
)
protectedvirtual

Reimplemented from G4VEmModel.

Definition at line 450 of file G4BetheBlochModel.cc.

References python.hepunit::electron_mass_c2, and G4INCL::Math::min().

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

452 {
453  // here particle type is checked for any method
454  SetParticle(pd);
455  G4double tau = kinEnergy/mass;
456  G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
457  (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
458  return std::min(tmax,tlimit);
459 }
float electron_mass_c2
Definition: hepunit.py:274
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
double G4double
Definition: G4Types.hh:76
G4double G4BetheBlochModel::MinEnergyCut ( const G4ParticleDefinition ,
const G4MaterialCutsCouple couple 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 180 of file G4BetheBlochModel.cc.

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

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

Implements G4VEmModel.

Definition at line 338 of file G4BetheBlochModel.cc.

References python.hepunit::electron_mass_c2, G4cout, G4endl, G4UniformRand, G4VEmModel::GetAngularDistribution(), G4DynamicParticle::GetDefinition(), G4DynamicParticle::GetKineticEnergy(), G4MaterialCutsCouple::GetMaterial(), G4DynamicParticle::GetMomentum(), G4DynamicParticle::GetMomentumDirection(), G4ParticleDefinition::GetParticleName(), G4DynamicParticle::GetTotalMomentum(), MaxSecondaryEnergy(), G4INCL::Math::min(), CLHEP::Hep3Vector::rotateUz(), G4VEmAngularDistribution::SampleDirection(), G4VEmModel::SelectRandomAtomNumber(), CLHEP::Hep3Vector::set(), G4ParticleChangeForLoss::SetProposedKineticEnergy(), G4ParticleChangeForLoss::SetProposedMomentumDirection(), python.hepunit::twopi, CLHEP::Hep3Vector::unit(), G4VEmModel::UseAngularGeneratorFlag(), and test::x.

343 {
344  G4double kineticEnergy = dp->GetKineticEnergy();
345  G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(),kineticEnergy);
346 
347  G4double maxKinEnergy = std::min(maxEnergy,tmax);
348  if(minKinEnergy >= maxKinEnergy) { return; }
349 
350  G4double totEnergy = kineticEnergy + mass;
351  G4double etot2 = totEnergy*totEnergy;
352  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
353 
354  G4double deltaKinEnergy, f;
355  G4double f1 = 0.0;
356  G4double fmax = 1.0;
357  if( 0.5 == spin ) { fmax += 0.5*maxKinEnergy*maxKinEnergy/etot2; }
358 
359  // sampling without nuclear size effect
360  do {
361  G4double q = G4UniformRand();
362  deltaKinEnergy = minKinEnergy*maxKinEnergy
363  /(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
364 
365  f = 1.0 - beta2*deltaKinEnergy/tmax;
366  if( 0.5 == spin ) {
367  f1 = 0.5*deltaKinEnergy*deltaKinEnergy/etot2;
368  f += f1;
369  }
370 
371  } while( fmax*G4UniformRand() > f);
372 
373  // projectile formfactor - suppresion of high energy
374  // delta-electron production at high energy
375 
376  G4double x = formfact*deltaKinEnergy;
377  if(x > 1.e-6) {
378 
379  G4double x1 = 1.0 + x;
380  G4double grej = 1.0/(x1*x1);
381  if( 0.5 == spin ) {
382  G4double x2 = 0.5*electron_mass_c2*deltaKinEnergy/(mass*mass);
383  grej *= (1.0 + magMoment2*(x2 - f1/f)/(1.0 + x2));
384  }
385  if(grej > 1.1) {
386  G4cout << "### G4BetheBlochModel WARNING: grej= " << grej
387  << " " << dp->GetDefinition()->GetParticleName()
388  << " Ekin(MeV)= " << kineticEnergy
389  << " delEkin(MeV)= " << deltaKinEnergy
390  << G4endl;
391  }
392  if(G4UniformRand() > grej) return;
393  }
394 
395  G4ThreeVector deltaDirection;
396 
398 
399  const G4Material* mat = couple->GetMaterial();
400  G4int Z = SelectRandomAtomNumber(mat);
401 
402  deltaDirection =
403  GetAngularDistribution()->SampleDirection(dp, deltaKinEnergy, Z, mat);
404 
405  } else {
406 
407  G4double deltaMomentum =
408  sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
409  G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
410  (deltaMomentum * dp->GetTotalMomentum());
411  if(cost > 1.0) { cost = 1.0; }
412  G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
413 
414  G4double phi = twopi * G4UniformRand() ;
415 
416  deltaDirection.set(sint*cos(phi),sint*sin(phi), cost) ;
417  deltaDirection.rotateUz(dp->GetMomentumDirection());
418  }
419 
420  /*
421  G4cout << "### G4BetheBlochModel "
422  << dp->GetDefinition()->GetParticleName()
423  << " Ekin(MeV)= " << kineticEnergy
424  << " delEkin(MeV)= " << deltaKinEnergy
425  << " tmin(MeV)= " << minKinEnergy
426  << " tmax(MeV)= " << maxKinEnergy
427  << " dir= " << dp->GetMomentumDirection()
428  << " dirDelta= " << deltaDirection
429  << G4endl;
430  }
431  */
432 
433  // create G4DynamicParticle object for delta ray
434  G4DynamicParticle* delta =
435  new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
436 
437  vdp->push_back(delta);
438 
439  // Change kinematics of primary particle
440  kineticEnergy -= deltaKinEnergy;
441  G4ThreeVector finalP = dp->GetMomentum() - delta->GetMomentum();
442  finalP = finalP.unit();
443 
444  fParticleChange->SetProposedKineticEnergy(kineticEnergy);
445  fParticleChange->SetProposedMomentumDirection(finalP);
446 }
void set(double x, double y, double z)
G4double GetKineticEnergy() const
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:578
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
G4double GetTotalMomentum() const
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:655
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)
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
G4ThreeVector GetMomentum() const
const G4Material * GetMaterial() const
G4int SelectRandomAtomNumber(const G4Material *)
Definition: G4VEmModel.hh:529
void G4BetheBlochModel::SetChargeSquareRatio ( G4double  val)
inlineprotected

Definition at line 203 of file G4BetheBlochModel.hh.

Referenced by G4BetheBlochIonGasModel::ChargeSquareRatio().

204 {
205  chargeSquare = val;
206 }

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