Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BetheBlochModel Class Reference

#include <G4BetheBlochModel.hh>

+ Inheritance diagram for G4BetheBlochModel:

Public Member Functions

 G4BetheBlochModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="BetheBloch")
 
 ~G4BetheBlochModel () override
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
 
virtual 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
 
G4double GetChargeSquareRatio (const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
 
G4double GetParticleCharge (const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
 
void CorrectionsAlongStep (const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, const G4double &length, G4double &eloss) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
G4BetheBlochModeloperator= (const G4BetheBlochModel &right)=delete
 
 G4BetheBlochModel (const G4BetheBlochModel &)=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
 
G4double GetChargeSquareRatio () const
 
void SetChargeSquareRatio (G4double val)
 
- 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 60 of file G4BetheBlochModel.hh.

Constructor & Destructor Documentation

◆ G4BetheBlochModel() [1/2]

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

Definition at line 76 of file G4BetheBlochModel.cc.

78 : G4VEmModel(nam),
79 twoln10(2.0*G4Log(10.0)),
80 fAlphaTlimit(1*CLHEP::GeV),
81 fProtonTlimit(10*CLHEP::GeV)
82{
83 theElectron = G4Electron::Electron();
86 SetLowEnergyLimit(2.0*CLHEP::MeV);
87}
G4double G4Log(G4double x)
Definition: G4Log.hh:227
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
static G4NistManager * Instance()
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753

◆ ~G4BetheBlochModel()

G4BetheBlochModel::~G4BetheBlochModel ( )
overridedefault

◆ G4BetheBlochModel() [2/2]

G4BetheBlochModel::G4BetheBlochModel ( const G4BetheBlochModel )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 217 of file G4BetheBlochModel.cc.

223{
224 return Z*ComputeCrossSectionPerElectron(p,kinEnergy,cutEnergy,maxEnergy);
225}
const G4int Z[17]
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)

◆ ComputeCrossSectionPerElectron()

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

Definition at line 186 of file G4BetheBlochModel.cc.

190{
191 G4double cross = 0.0;
192 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
193 G4double maxEnergy = std::min(tmax, maxKinEnergy);
194 if(cutEnergy < maxEnergy) {
195
196 G4double totEnergy = kineticEnergy + mass;
197 G4double energy2 = totEnergy*totEnergy;
198 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
199
200 cross = (maxEnergy - cutEnergy)/(cutEnergy*maxEnergy)
201 - beta2*G4Log(maxEnergy/cutEnergy)/tmax;
202
203 // +term for spin=1/2 particle
204 if( 0.0 < spin ) { cross += 0.5*(maxEnergy - cutEnergy)/energy2; }
205
206 cross *= CLHEP::twopi_mc2_rcl2*chargeSquare/beta2;
207 }
208
209 // G4cout << "BB: e= " << kineticEnergy << " tmin= " << cutEnergy
210 // << " tmax= " << tmax << " cross= " << cross << G4endl;
211
212 return cross;
213}
double G4double
Definition: G4Types.hh:83
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override

Referenced by ComputeCrossSectionPerAtom(), G4LindhardSorensenIonModel::ComputeCrossSectionPerElectron(), and CrossSectionPerVolume().

◆ ComputeDEDXPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 246 of file G4BetheBlochModel.cc.

250{
251 const G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
252 // projectile formfactor limit energy loss
253 const G4double cutEnergy = std::min(std::min(cut,tmax), tlimit);
254
255 G4double tau = kineticEnergy/mass;
256 G4double gam = tau + 1.0;
257 G4double bg2 = tau * (tau+2.0);
258 G4double beta2 = bg2/(gam*gam);
259 G4double xc = cutEnergy/tmax;
260
261 G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
262 G4double eexc2 = eexc*eexc;
263
264 G4double eDensity = material->GetElectronDensity();
265
266 // added ICRU90 stopping data for limited list of materials
267 /*
268 G4cout << "### DEDX ICRI90:" << (nullptr != fICRU90)
269 << " Ekin=" << kineticEnergy
270 << " " << p->GetParticleName()
271 << " q2=" << chargeSquare
272 << " inside " << material->GetName() << G4endl;
273 */
274 if(nullptr != fICRU90 && kineticEnergy < fProtonTlimit) {
275 if(material != currentMaterial) {
276 currentMaterial = material;
277 baseMaterial = material->GetBaseMaterial()
278 ? material->GetBaseMaterial() : material;
279 iICRU90 = fICRU90->GetIndex(baseMaterial);
280 }
281 if(iICRU90 >= 0) {
282 G4double dedx = 0.0;
283 // only for alpha
284 if(isAlpha) {
285 if(kineticEnergy <= fAlphaTlimit) {
286 dedx = fICRU90->GetElectronicDEDXforAlpha(iICRU90, kineticEnergy);
287 } else {
288 const G4double e = kineticEnergy*CLHEP::proton_mass_c2/mass;
289 dedx = fICRU90->GetElectronicDEDXforProton(iICRU90, e)*chargeSquare;
290 }
291 } else {
292 dedx = fICRU90->GetElectronicDEDXforProton(iICRU90, kineticEnergy)
293 *chargeSquare;
294 }
295 dedx *= material->GetDensity();
296 if(cutEnergy < tmax) {
297 dedx += (G4Log(xc) + (1.0 - xc)*beta2)*CLHEP::twopi_mc2_rcl2
298 *(eDensity*chargeSquare/beta2);
299 }
300 //G4cout << " iICRU90=" << iICRU90 << " dedx=" << dedx << G4endl;
301 if(dedx > 0.0) { return dedx; }
302 }
303 }
304 // general Bethe-Bloch formula
305 G4double dedx = G4Log(2.0*CLHEP::electron_mass_c2*bg2*cutEnergy/eexc2)
306 - (1.0 + xc)*beta2;
307
308 if(0.0 < spin) {
309 G4double del = 0.5*cutEnergy/(kineticEnergy + mass);
310 dedx += del*del;
311 }
312
313 // density correction
314 G4double x = G4Log(bg2)/twoln10;
315 dedx -= material->GetIonisation()->DensityCorrection(x);
316
317 // shell correction
318 dedx -= 2.0*corr->ShellCorrection(p,material,kineticEnergy);
319
320 // now compute the total ionization loss
321 dedx *= CLHEP::twopi_mc2_rcl2*chargeSquare*eDensity/beta2;
322
323 //High order correction different for hadrons and ions
324 if(isIon) {
325 dedx += corr->IonBarkasCorrection(p,material,kineticEnergy);
326 } else {
327 dedx += corr->HighOrderCorrections(p,material,kineticEnergy,cutEnergy);
328 }
329
330 dedx = std::max(dedx, 0.0);
331 /*
332 G4cout << "E(MeV)= " << kineticEnergy/CLHEP::MeV << " dedx= " << dedx
333 << " " << material->GetName() << G4endl;
334 */
335 return dedx;
336}
G4double IonBarkasCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double HighOrderCorrections(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy, G4double cutEnergy)
G4double ShellCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double GetElectronicDEDXforProton(const G4Material *, G4double kinEnergy) const
G4int GetIndex(const G4Material *) const
G4double GetElectronicDEDXforAlpha(const G4Material *, G4double scaledKinEnergy) const
G4double DensityCorrection(G4double x) const
G4double GetMeanExcitationEnergy() const
G4double GetDensity() const
Definition: G4Material.hh:175
const G4Material * GetBaseMaterial() const
Definition: G4Material.hh:228
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:221
G4double GetElectronDensity() const
Definition: G4Material.hh:212

Referenced by G4BetheBlochNoDeltaModel::ComputeDEDXPerVolume(), and G4LindhardSorensenIonModel::ComputeDEDXPerVolume().

◆ CorrectionsAlongStep()

void G4BetheBlochModel::CorrectionsAlongStep ( const G4MaterialCutsCouple couple,
const G4DynamicParticle dp,
const G4double length,
G4double eloss 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 340 of file G4BetheBlochModel.cc.

344{
345 // no correction for alpha
346 if(isAlpha) { return; }
347
348 // no correction at the last step or at small step
349 const G4double preKinEnergy = dp->GetKineticEnergy();
350 if(eloss >= preKinEnergy || eloss < preKinEnergy*0.05) { return; }
351
352 // corrections for all charged particles with Q > 1
353 const G4ParticleDefinition* p = dp->GetDefinition();
354 if(p != particle) { SetupParameters(p); }
355 if(!isIon) { return; }
356
357 // effective energy and charge at a step
358 const G4double e = std::max(preKinEnergy - eloss*0.5, preKinEnergy*0.5);
359 const G4Material* mat = couple->GetMaterial();
360 const G4double q20 = corr->EffectiveChargeSquareRatio(p, mat, preKinEnergy);
361 const G4double q2 = corr->EffectiveChargeSquareRatio(p, mat, e);
362 const G4double qfactor = q2/q20;
363
364 /*
365 G4cout << "G4BetheBlochModel::CorrectionsAlongStep: Epre(MeV)="
366 << preKinEnergy << " Eeff(MeV)=" << e
367 << " eloss=" << eloss << " elossnew=" << eloss*qfactor
368 << " qfactor=" << qfactor << " Qpre=" << q20
369 << p->GetParticleName() <<G4endl;
370 */
371 eloss *= qfactor;
372}
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
const G4Material * GetMaterial() const

◆ CrossSectionPerVolume()

G4double G4BetheBlochModel::CrossSectionPerVolume ( const G4Material mat,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 229 of file G4BetheBlochModel.cc.

235{
236 G4double sigma = mat->GetElectronDensity()
237 *ComputeCrossSectionPerElectron(p,kinEnergy,cutEnergy,maxEnergy);
238 if(isAlpha) {
239 sigma *= corr->EffectiveChargeSquareRatio(p,mat,kinEnergy)/chargeSquare;
240 }
241 return sigma;
242}

◆ GetChargeSquareRatio() [1/2]

G4double G4BetheBlochModel::GetChargeSquareRatio ( ) const
inlineprotected

Definition at line 162 of file G4BetheBlochModel.hh.

163{
164 return chargeSquare;
165}

◆ GetChargeSquareRatio() [2/2]

G4double G4BetheBlochModel::GetChargeSquareRatio ( const G4ParticleDefinition p,
const G4Material mat,
G4double  kineticEnergy 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 128 of file G4BetheBlochModel.cc.

131{
132 // this method is called only for ions, so no check if it is an ion
133 return
134 (!isAlpha) ? corr->EffectiveChargeSquareRatio(p,mat,kineticEnergy) : 1.0;
135}

◆ GetParticleCharge()

G4double G4BetheBlochModel::GetParticleCharge ( const G4ParticleDefinition p,
const G4Material mat,
G4double  kineticEnergy 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 139 of file G4BetheBlochModel.cc.

142{
143 // this method is called only for ions, so no check if it is an ion
144 return corr->GetParticleCharge(p,mat,kineticEnergy);
145}
G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 95 of file G4BetheBlochModel.cc.

97{
98 if(p != particle) { SetupParameters(p); }
99
100 //G4cout << "G4BetheBlochModel::Initialise for " << p->GetParticleName()
101 // << " isIon= " << isIon
102 // << G4endl;
103
104 // always false before the run
105 SetDeexcitationFlag(false);
106
107 // initialisation once
108 if(nullptr == fParticleChange) {
109 const G4String& pname = particle->GetParticleName();
110 if(IsMaster() && G4EmParameters::Instance()->UseICRU90Data() &&
111 (pname == "proton" || pname == "GenericIon" || pname == "alpha")) {
112 fICRU90 = nist->GetICRU90StoppingData();
113 fICRU90->Initialise();
114 }
115 if(particle->GetPDGCharge() > CLHEP::eplus ||
116 pname == "GenericIon") { isIon = true; }
117 if(pname == "alpha") { isAlpha = true; }
118
119 fParticleChange = GetParticleChangeForLoss();
120 if(UseAngularGeneratorFlag() && nullptr == GetAngularDistribution()) {
122 }
123 }
124}
static G4EmParameters * Instance()
G4ICRU90StoppingData * GetICRU90StoppingData()
G4double GetPDGCharge() const
const G4String & GetParticleName() const
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:600
G4bool IsMaster() const
Definition: G4VEmModel.hh:725
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:802
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:607
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:697
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:109

Referenced by G4LindhardSorensenIonModel::Initialise().

◆ MaxSecondaryEnergy()

G4double G4BetheBlochModel::MaxSecondaryEnergy ( const G4ParticleDefinition pd,
G4double  kinEnergy 
)
overrideprotectedvirtual

Reimplemented from G4VEmModel.

Definition at line 487 of file G4BetheBlochModel.cc.

489{
490 // here particle type is checked for the case,
491 // when this model is shared between particles
492 if(pd != particle) { SetupParameters(pd); }
493 G4double tau = kinEnergy/mass;
494 return 2.0*CLHEP::electron_mass_c2*tau*(tau + 2.) /
495 (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
496}

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

◆ MinEnergyCut()

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

Reimplemented from G4VEmModel.

Definition at line 177 of file G4BetheBlochModel.cc.

179{
181}

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 376 of file G4BetheBlochModel.cc.

381{
382 G4double kineticEnergy = dp->GetKineticEnergy();
383 const G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(),kineticEnergy);
384 const G4double maxKinEnergy = std::min(maxEnergy,tmax);
385 if(minKinEnergy >= maxKinEnergy) { return; }
386
387 //G4cout << "G4BetheBlochModel::SampleSecondaries Emin= " << minKinEnergy
388 // << " Emax= " << maxKinEnergy << G4endl;
389
390 const G4double totEnergy = kineticEnergy + mass;
391 const G4double etot2 = totEnergy*totEnergy;
392 const G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
393
394 G4double deltaKinEnergy, f;
395 G4double f1 = 0.0;
396 G4double fmax = 1.0;
397 if( 0.0 < spin ) { fmax += 0.5*maxKinEnergy*maxKinEnergy/etot2; }
398
399 CLHEP::HepRandomEngine* rndmEngineMod = G4Random::getTheEngine();
400 G4double rndm[2];
401
402 // sampling without nuclear size effect
403 do {
404 rndmEngineMod->flatArray(2, rndm);
405 deltaKinEnergy = minKinEnergy*maxKinEnergy
406 /(minKinEnergy*(1.0 - rndm[0]) + maxKinEnergy*rndm[0]);
407
408 f = 1.0 - beta2*deltaKinEnergy/tmax;
409 if( 0.0 < spin ) {
410 f1 = 0.5*deltaKinEnergy*deltaKinEnergy/etot2;
411 f += f1;
412 }
413
414 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
415 } while( fmax*rndm[1] > f);
416
417 // projectile formfactor - suppresion of high energy
418 // delta-electron production at high energy
419
420 G4double x = formfact*deltaKinEnergy;
421 if(x > 1.e-6) {
422
423 G4double x1 = 1.0 + x;
424 G4double grej = 1.0/(x1*x1);
425 if( 0.0 < spin ) {
426 G4double x2 = 0.5*electron_mass_c2*deltaKinEnergy/(mass*mass);
427 grej *= (1.0 + magMoment2*(x2 - f1/f)/(1.0 + x2));
428 }
429 if(grej > 1.1) {
430 G4cout << "### G4BetheBlochModel WARNING: grej= " << grej
431 << " " << dp->GetDefinition()->GetParticleName()
432 << " Ekin(MeV)= " << kineticEnergy
433 << " delEkin(MeV)= " << deltaKinEnergy
434 << G4endl;
435 }
436 if(rndmEngineMod->flat() > grej) { return; }
437 }
438
439 G4ThreeVector deltaDirection;
440
442 const G4Material* mat = couple->GetMaterial();
443 deltaDirection =
444 GetAngularDistribution()->SampleDirection(dp, deltaKinEnergy,
446 mat);
447 } else {
448
449 G4double deltaMomentum =
450 sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
451 G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
452 (deltaMomentum * dp->GetTotalMomentum());
453 cost = std::min(cost, 1.0);
454 const G4double sint = std::sqrt((1.0 - cost)*(1.0 + cost));
455 const G4double phi = twopi*rndmEngineMod->flat();
456
457 deltaDirection.set(sint*std::cos(phi),sint*std::sin(phi), cost) ;
458 deltaDirection.rotateUz(dp->GetMomentumDirection());
459 }
460 /*
461 G4cout << "### G4BetheBlochModel "
462 << dp->GetDefinition()->GetParticleName()
463 << " Ekin(MeV)= " << kineticEnergy
464 << " delEkin(MeV)= " << deltaKinEnergy
465 << " tmin(MeV)= " << minKinEnergy
466 << " tmax(MeV)= " << maxKinEnergy
467 << " dir= " << dp->GetMomentumDirection()
468 << " dirDelta= " << deltaDirection
469 << G4endl;
470 */
471 // create G4DynamicParticle object for delta ray
472 auto delta = new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
473
474 vdp->push_back(delta);
475
476 // Change kinematics of primary particle
477 kineticEnergy -= deltaKinEnergy;
478 G4ThreeVector finalP = dp->GetMomentum() - delta->GetMomentum();
479 finalP = finalP.unit();
480
481 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
482 fParticleChange->SetProposedMomentumDirection(finalP);
483}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Hep3Vector unit() const
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
virtual double flat()=0
virtual void flatArray(const int size, double *vect)=0
const G4ThreeVector & GetMomentumDirection() const
G4ThreeVector GetMomentum() const
G4double GetTotalMomentum() 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
Definition: G4VEmModel.cc:253

◆ SetChargeSquareRatio()

void G4BetheBlochModel::SetChargeSquareRatio ( G4double  val)
inlineprotected

Definition at line 169 of file G4BetheBlochModel.hh.

170{
171 chargeSquare = val;
172}

Referenced by G4BetheBlochIonGasModel::ChargeSquareRatio().


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