Geant4 11.3.0
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 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 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 *)
 
virtual void FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
void SetLPMFlag (G4bool)
 
void SetMasterThread (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 ()
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Additional Inherited Members

- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 
const G4MaterialpBaseMaterial = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4double inveplus
 
G4double pFactor = 1.0
 
std::size_t currentCoupleIndex = 0
 
std::size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Detailed Description

Definition at line 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 74 of file G4BetheBlochModel.cc.

76 : G4VEmModel(nam),
77 twoln10(2.0*G4Log(10.0)),
78 fAlphaTlimit(1*CLHEP::GeV),
79 fProtonTlimit(10*CLHEP::GeV)
80{
81 theElectron = G4Electron::Electron();
84 SetLowEnergyLimit(2.0*CLHEP::MeV);
85}
G4double G4Log(G4double x)
Definition G4Log.hh:227
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
static G4NistManager * Instance()
void SetLowEnergyLimit(G4double)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67

Referenced by G4BetheBlochIonGasModel::G4BetheBlochIonGasModel(), G4BetheBlochModel(), G4BetheBlochNoDeltaModel::G4BetheBlochNoDeltaModel(), and operator=().

◆ ~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 219 of file G4BetheBlochModel.cc.

225{
226 return Z*ComputeCrossSectionPerElectron(p,kinEnergy,cutEnergy,maxEnergy);
227}
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 187 of file G4BetheBlochModel.cc.

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

Referenced by ComputeCrossSectionPerAtom(), and CrossSectionPerVolume().

◆ ComputeDEDXPerVolume()

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

Reimplemented from G4VEmModel.

Reimplemented in G4BetheBlochNoDeltaModel.

Definition at line 248 of file G4BetheBlochModel.cc.

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

Referenced by G4BetheBlochNoDeltaModel::ComputeDEDXPerVolume().

◆ CorrectionsAlongStep()

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

Reimplemented from G4VEmModel.

Definition at line 342 of file G4BetheBlochModel.cc.

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

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Reimplemented in G4BetheBlochNoDeltaModel.

Definition at line 231 of file G4BetheBlochModel.cc.

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

◆ 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 129 of file G4BetheBlochModel.cc.

132{
133 // this method is called only for ions, so no check if it is an ion
134 if(isAlpha) { return 1.0; }
135 chargeSquare = corr->EffectiveChargeSquareRatio(p, mat, kinEnergy);
136 return chargeSquare;
137}

◆ GetParticleCharge()

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

Reimplemented from G4VEmModel.

Definition at line 141 of file G4BetheBlochModel.cc.

144{
145 // this method is called only for ions, so no check if it is an ion
146 return corr->GetParticleCharge(p, mat, kineticEnergy);
147}

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 93 of file G4BetheBlochModel.cc.

95{
96 if(p != particle) { SetupParameters(p); }
97
98 // always false before the run
100
101 // initialisation once
102 if(nullptr == fParticleChange) {
103 const G4String& pname = particle->GetParticleName();
104 if(G4EmParameters::Instance()->UseICRU90Data() &&
105 (pname == "proton" || pname == "GenericIon" || pname == "alpha")) {
106 fICRU90 = nist->GetICRU90StoppingData();
107 }
108 if (pname == "GenericIon") {
109 isIon = true;
110 } else if (pname == "alpha") {
111 isAlpha = true;
112 } else if (particle->GetPDGCharge() > 1.1*CLHEP::eplus) {
113 isIon = true;
114 }
115
116 fParticleChange = GetParticleChangeForLoss();
117 if(UseAngularGeneratorFlag() && nullptr == GetAngularDistribution()) {
118 SetAngularDistribution(new G4DeltaAngle());
119 }
120 }
121 // initialisation for each new run
122 if(IsMaster() && nullptr != fICRU90) {
123 fICRU90->Initialise();
124 }
125}
static G4EmParameters * Instance()
G4VEmAngularDistribution * GetAngularDistribution()
G4bool IsMaster() const
void SetDeexcitationFlag(G4bool val)
void SetAngularDistribution(G4VEmAngularDistribution *)
G4bool UseAngularGeneratorFlag() const
G4ParticleChangeForLoss * GetParticleChangeForLoss()

◆ MaxSecondaryEnergy()

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

Reimplemented from G4VEmModel.

Definition at line 490 of file G4BetheBlochModel.cc.

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

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

◆ MinEnergyCut()

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

Reimplemented from G4VEmModel.

Definition at line 178 of file G4BetheBlochModel.cc.

180{
182}

◆ 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 378 of file G4BetheBlochModel.cc.

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

◆ 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: