Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BetheBlochModel.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// -------------------------------------------------------------------
27//
28// GEANT4 Class header file
29//
30//
31// File name: G4BetheBlochModel
32//
33// Author: Vladimir Ivanchenko on base of Laszlo Urban code
34//
35// Creation date: 03.01.2002
36//
37// Modifications:
38//
39// 04-12-02 Fix problem of G4DynamicParticle constructor (V.Ivanchenko)
40// 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
41// 27-01-03 Make models region aware (V.Ivanchenko)
42// 13-02-03 Add name (V.Ivanchenko)
43// 24-03-05 Add G4EmCorrections (V.Ivanchenko)
44// 11-04-05 Major optimisation of internal interfaces (V.Ivanchenko)
45// 11-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
46// 12-02-06 move G4LossTableManager::Instance()->EmCorrections()
47// in constructor (mma)
48// 12-08-08 Added methods GetParticleCharge, GetChargeSquareRatio,
49// CorrectionsAlongStep needed for ions(V.Ivanchenko)
50//
51// -------------------------------------------------------------------
52//
53
54//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
55//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
56
57#include "G4BetheBlochModel.hh"
58#include "Randomize.hh"
60#include "G4SystemOfUnits.hh"
61#include "G4NistManager.hh"
62#include "G4Electron.hh"
63#include "G4LossTableManager.hh"
64#include "G4EmCorrections.hh"
65#include "G4EmParameters.hh"
68#include "G4Log.hh"
69#include "G4DeltaAngle.hh"
70#include <vector>
71
72//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
73
75 const G4String& nam)
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}
86
87//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
88
90
91//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
92
94 const G4DataVector&)
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()) {
119 }
120 }
121 // initialisation for each new run
122 if(IsMaster() && nullptr != fICRU90) {
123 fICRU90->Initialise();
124 }
125}
126
127//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
128
130 const G4Material* mat,
131 G4double kinEnergy)
132{
133 // this method is called only for ions, so no check if it is an ion
134 chargeSquareRatio = corr->EffectiveChargeSquareRatio(p, mat, kinEnergy);
135 return chargeSquareRatio;
136}
137
138//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
139
141 const G4Material* mat,
142 G4double kineticEnergy)
143{
144 // this method is called only for ions, so no check if it is an ion
145 return corr->GetParticleCharge(p, mat, kineticEnergy);
146}
147
148//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
149
150void G4BetheBlochModel::SetupParameters(const G4ParticleDefinition* p)
151{
152 particle = p;
153 mass = particle->GetPDGMass();
154 spin = particle->GetPDGSpin();
155 G4double q = particle->GetPDGCharge()*inveplus;
156 ratio = electron_mass_c2/mass;
157 constexpr G4double aMag = 1./(0.5*eplus*CLHEP::hbar_Planck*CLHEP::c_squared);
158 G4double magmom = particle->GetPDGMagneticMoment()*mass*aMag;
159 magMoment2 = magmom*magmom - 1.0;
160 formfact = 0.0;
161 tlimit = DBL_MAX;
162 if(particle->GetLeptonNumber() == 0) {
163 G4double x = 0.8426*CLHEP::GeV;
164 if(spin == 0.0 && mass < CLHEP::GeV) { x = 0.736*CLHEP::GeV; }
165 else if (mass > CLHEP::GeV) {
166 G4int iz = G4lrint(std::abs(q));
167 if(iz > 1) { x /= nist->GetA27(iz); }
168 }
169 formfact = 2.0*CLHEP::electron_mass_c2/(x*x);
170 tlimit = 2.0/formfact;
171 }
172}
173
174//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
175
181
182//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
183
186 G4double kineticEnergy,
187 G4double cut,
188 G4double maxKinEnergy)
189{
190 G4double cross = 0.0;
191 const G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
192 const G4double cutEnergy = std::min(std::min(cut,tmax), tlimit);
193 const 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*chargeSquareRatio/beta2;
207 }
208
209 // G4cout << "BB: e= " << kineticEnergy << " tmin= " << cutEnergy
210 // << " tmax= " << tmax << " cross= " << cross << G4endl;
211
212 return cross;
213}
214
215//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
216
218 const G4ParticleDefinition* p,
219 G4double kinEnergy,
221 G4double cutEnergy,
222 G4double maxEnergy)
223{
224 return Z*ComputeCrossSectionPerElectron(p,kinEnergy,cutEnergy,maxEnergy);
225}
226
227//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
228
230 const G4Material* mat,
231 const G4ParticleDefinition* p,
232 G4double kinEnergy,
233 G4double cutEnergy,
234 G4double maxEnergy)
235{
236 chargeSquareRatio = GetChargeSquareRatio(p, mat, kinEnergy);
237 G4double sigma = mat->GetElectronDensity()
238 *ComputeCrossSectionPerElectron(p,kinEnergy,cutEnergy,maxEnergy);
239 return sigma;
240}
241
242//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
243
245 const G4ParticleDefinition* p,
246 G4double kineticEnergy,
247 G4double cut)
248{
249 const G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
250 // projectile formfactor limit energy loss
251 const G4double cutEnergy = std::min(std::min(cut,tmax), tlimit);
252 chargeSquareRatio = GetChargeSquareRatio(p, material, kineticEnergy);
253
254 G4double tau = kineticEnergy/mass;
255 G4double gam = tau + 1.0;
256 G4double bg2 = tau * (tau+2.0);
257 G4double beta2 = bg2/(gam*gam);
258 G4double xc = cutEnergy/tmax;
259
260 G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
261 G4double eexc2 = eexc*eexc;
262
263 G4double eDensity = material->GetElectronDensity();
264
265 // added ICRU90 stopping data for limited list of materials
266 /*
267 G4cout << "### DEDX ICRI90:" << (nullptr != fICRU90)
268 << " Ekin=" << kineticEnergy
269 << " " << p->GetParticleName()
270 << " q2=" << chargeSquare
271 << " inside " << material->GetName() << G4endl;
272 */
273 if(nullptr != fICRU90 && kineticEnergy < fProtonTlimit) {
274 if(material != currentMaterial) {
275 currentMaterial = material;
276 baseMaterial = material->GetBaseMaterial()
277 ? material->GetBaseMaterial() : material;
278 iICRU90 = fICRU90->GetIndex(baseMaterial);
279 }
280 if(iICRU90 >= 0) {
281 G4double dedx = 0.0;
282 // only for alpha
283 if(isAlpha) {
284 if(kineticEnergy <= fAlphaTlimit) {
285 dedx = fICRU90->GetElectronicDEDXforAlpha(iICRU90, kineticEnergy);
286 } else {
287 const G4double e = kineticEnergy*CLHEP::proton_mass_c2/mass;
288 dedx = fICRU90->GetElectronicDEDXforProton(iICRU90, e)*chargeSquareRatio;
289 }
290 } else {
291 const G4double e = kineticEnergy*CLHEP::proton_mass_c2/mass;
292 dedx = fICRU90->GetElectronicDEDXforProton(iICRU90, e)*chargeSquareRatio;
293 }
294 dedx *= material->GetDensity();
295 if(cutEnergy < tmax) {
296 dedx += (G4Log(xc) + (1.0 - xc)*beta2)*CLHEP::twopi_mc2_rcl2
297 *(eDensity*chargeSquareRatio/beta2);
298 }
299 //G4cout << " iICRU90=" << iICRU90 << " dedx=" << dedx << G4endl;
300 if(dedx > 0.0) { return dedx; }
301 }
302 }
303 // general Bethe-Bloch formula
304 G4double dedx = G4Log(2.0*CLHEP::electron_mass_c2*bg2*cutEnergy/eexc2)
305 - (1.0 + xc)*beta2;
306
307 if(0.0 < spin) {
308 G4double del = 0.5*cutEnergy/(kineticEnergy + mass);
309 dedx += del*del;
310 }
311
312 // density correction
313 G4double x = G4Log(bg2)/twoln10;
314 dedx -= material->GetIonisation()->DensityCorrection(x);
315
316 // shell correction
317 dedx -= 2.0*corr->ShellCorrection(p,material,kineticEnergy);
318
319 // now compute the total ionization loss
320 dedx *= CLHEP::twopi_mc2_rcl2*chargeSquareRatio*eDensity/beta2;
321
322 //High order correction different for hadrons and ions
323 dedx += corr->HighOrderCorrections(p,material,kineticEnergy,cutEnergy);
324 dedx = std::max(dedx, 0.0);
325
326 /*
327 G4cout << "E(MeV)= " << kineticEnergy/CLHEP::MeV << " dedx= " << dedx
328 << " " << material->GetName() << G4endl;
329 */
330 return dedx;
331}
332
333//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
334
336 const G4ParticleDefinition* p,
337 const G4double preKinEnergy,
338 const G4double,
339 const G4double& /*length*/,
340 G4double& eloss)
341{
342 // no correction for alpha
343 if(isAlpha) { return; }
344
345 // no correction at the last step or at small step
346 if(eloss >= preKinEnergy || eloss < preKinEnergy*0.05) { return; }
347
348 // corrections for all charged particles with Q > 1
349 if(p != particle) { SetupParameters(p); }
350 if(!isIon) { return; }
351
352 // effective energy and charge at a step
353 const G4double e = std::max(preKinEnergy - eloss*0.5, preKinEnergy*0.5);
354 const G4double q2 = GetChargeSquareRatio(p, mat, e);
355 const G4double qfactor = q2/chargeSquareRatio;
356
357 /*
358 G4cout << "G4BetheBlochModel::CorrectionsAlongStep: Epre(MeV)="
359 << preKinEnergy << " Eeff(MeV)=" << e
360 << " eloss=" << eloss << " elossnew=" << eloss*qfactor
361 << " qfactor=" << qfactor << " Qpre=" << q20
362 << p->GetParticleName() <<G4endl;
363 */
364 eloss *= qfactor;
365}
366
367//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
368
369void G4BetheBlochModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
370 const G4MaterialCutsCouple* couple,
371 const G4DynamicParticle* dp,
372 G4double cut,
373 G4double maxEnergy)
374{
375 G4double kinEnergy = dp->GetKineticEnergy();
376 const G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(), kinEnergy);
377 const G4double minKinEnergy = std::min(cut, tmax);
378 const G4double maxKinEnergy = std::min(maxEnergy, tmax);
379 if(minKinEnergy >= maxKinEnergy) { return; }
380
381 //G4cout << "G4BetheBlochModel::SampleSecondaries Emin= " << minKinEnergy
382 // << " Emax= " << maxKinEnergy << G4endl;
383
384 const G4double totEnergy = kinEnergy + mass;
385 const G4double etot2 = totEnergy*totEnergy;
386 const G4double beta2 = kinEnergy*(kinEnergy + 2.0*mass)/etot2;
387
388 G4double deltaKinEnergy, f;
389 G4double f1 = 0.0;
390 G4double fmax = 1.0;
391 if( 0.0 < spin ) { fmax += 0.5*maxKinEnergy*maxKinEnergy/etot2; }
392
393 CLHEP::HepRandomEngine* rndmEngineMod = G4Random::getTheEngine();
394 G4double rndm[2];
395
396 // sampling without nuclear size effect
397 do {
398 rndmEngineMod->flatArray(2, rndm);
399 deltaKinEnergy = minKinEnergy*maxKinEnergy
400 /(minKinEnergy*(1.0 - rndm[0]) + maxKinEnergy*rndm[0]);
401
402 f = 1.0 - beta2*deltaKinEnergy/tmax;
403 if( 0.0 < spin ) {
404 f1 = 0.5*deltaKinEnergy*deltaKinEnergy/etot2;
405 f += f1;
406 }
407
408 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
409 } while( fmax*rndm[1] > f);
410
411 // projectile formfactor - suppresion of high energy
412 // delta-electron production at high energy
413
414 G4double x = formfact*deltaKinEnergy;
415 if(x > 1.e-6) {
416
417 G4double x1 = 1.0 + x;
418 G4double grej = 1.0/(x1*x1);
419 if( 0.0 < spin ) {
420 G4double x2 = 0.5*electron_mass_c2*deltaKinEnergy/(mass*mass);
421 grej *= (1.0 + magMoment2*(x2 - f1/f)/(1.0 + x2));
422 }
423 if(grej > 1.1) {
424 G4cout << "### G4BetheBlochModel WARNING: grej= " << grej
425 << " " << dp->GetDefinition()->GetParticleName()
426 << " Ekin(MeV)= " << kinEnergy
427 << " delEkin(MeV)= " << deltaKinEnergy
428 << G4endl;
429 }
430 if(rndmEngineMod->flat() > grej) { return; }
431 }
432
433 G4ThreeVector deltaDirection;
434
436 const G4Material* mat = couple->GetMaterial();
437 deltaDirection =
438 GetAngularDistribution()->SampleDirection(dp, deltaKinEnergy,
440 mat);
441 } else {
442
443 G4double deltaMomentum =
444 std::sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
445 G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
446 (deltaMomentum * dp->GetTotalMomentum());
447 cost = std::min(cost, 1.0);
448 const G4double sint = std::sqrt((1.0 - cost)*(1.0 + cost));
449 const G4double phi = twopi*rndmEngineMod->flat();
450
451 deltaDirection.set(sint*std::cos(phi),sint*std::sin(phi), cost) ;
452 deltaDirection.rotateUz(dp->GetMomentumDirection());
453 }
454 /*
455 G4cout << "### G4BetheBlochModel "
456 << dp->GetDefinition()->GetParticleName()
457 << " Ekin(MeV)= " << kinEnergy
458 << " delEkin(MeV)= " << deltaKinEnergy
459 << " tmin(MeV)= " << minKinEnergy
460 << " tmax(MeV)= " << maxKinEnergy
461 << " dir= " << dp->GetMomentumDirection()
462 << " dirDelta= " << deltaDirection
463 << G4endl;
464 */
465 // create G4DynamicParticle object for delta ray
466 auto delta = new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
467
468 vdp->push_back(delta);
469
470 // Change kinematics of primary particle
471 kinEnergy -= deltaKinEnergy;
472 G4ThreeVector finalP = dp->GetMomentum() - delta->GetMomentum();
473 finalP = finalP.unit();
474
475 fParticleChange->SetProposedKineticEnergy(kinEnergy);
476 fParticleChange->SetProposedMomentumDirection(finalP);
477}
478
479//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
480
482 G4double kinEnergy)
483{
484 // here particle type is checked for the case,
485 // when this model is shared between particles
486 if(pd != particle) { SetupParameters(pd); }
487 G4double tau = kinEnergy/mass;
488 return 2.0*CLHEP::electron_mass_c2*tau*(tau + 2.) /
489 (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
490}
491
492//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4double G4Log(G4double x)
Definition G4Log.hh:169
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#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
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double GetParticleCharge(const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
G4double GetChargeSquareRatio(const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
~G4BetheBlochModel() override
G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override
G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
void CorrectionsAlongStep(const G4Material *, const G4ParticleDefinition *, const G4double kinEnergy, const G4double cutEnergy, const G4double &length, G4double &eloss) override
G4BetheBlochModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="BetheBloch")
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
G4double GetTotalMomentum() const
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4EmParameters * Instance()
G4double DensityCorrection(G4double x) const
G4double GetMeanExcitationEnergy() const
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
const G4Material * GetMaterial() const
G4double GetDensity() const
const G4Material * GetBaseMaterial() const
G4IonisParamMat * GetIonisation() const
G4double GetElectronDensity() const
std::size_t GetIndex() const
static G4NistManager * Instance()
G4double GetPDGMagneticMoment() const
const G4String & GetParticleName() const
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4VEmAngularDistribution * GetAngularDistribution()
G4double inveplus
G4int SelectRandomAtomNumber(const G4Material *) const
G4bool IsMaster() const
void SetLowEnergyLimit(G4double)
void SetDeexcitationFlag(G4bool val)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67
void SetAngularDistribution(G4VEmAngularDistribution *)
G4bool UseAngularGeneratorFlag() const
G4ParticleChangeForLoss * GetParticleChangeForLoss()
int G4lrint(double ad)
Definition templates.hh:134
#define DBL_MAX
Definition templates.hh:62