Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LindhardSorensenIonModel.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: G4LindhardSorensenIonModel
32//
33// Author: Alexander Bagulya & Vladimir Ivanchenko
34//
35// Creation date: 16.04.2018
36//
37//
38// -------------------------------------------------------------------
39//
40
41//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
42//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
43
45#include "Randomize.hh"
47#include "G4SystemOfUnits.hh"
48#include "G4Electron.hh"
49#include "G4LossTableManager.hh"
50#include "G4EmCorrections.hh"
52#include "G4Log.hh"
53#include "G4DeltaAngle.hh"
55#include "G4BraggModel.hh"
56#include "G4BetheBlochModel.hh"
57#include "G4IonICRU73Data.hh"
58#include "G4AutoLock.hh"
59
60//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
61
62using namespace std;
63
64G4LindhardSorensenData* G4LindhardSorensenIonModel::lsdata = nullptr;
65G4IonICRU73Data* G4LindhardSorensenIonModel::fIonData = nullptr;
66
67namespace
68{
69 G4Mutex ionXSMutex = G4MUTEX_INITIALIZER;
70}
71
73 const G4String& nam)
74 : G4VEmModel(nam),
75 twoln10(2.0*G4Log(10.0))
76{
77 theElectron = G4Electron::Electron();
80 fBraggModel = new G4BraggModel();
81 fBBModel = new G4BetheBlochModel();
82 fElimit = 2.0*CLHEP::MeV;
83}
84
85//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
86
88 if(isFirst) {
89 delete lsdata;
90 delete fIonData;
91 lsdata = nullptr;
92 fIonData = nullptr;
93 }
94}
95
96//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
97
99 const G4DataVector& ptr)
100{
101 fBraggModel->Initialise(p, ptr);
102 fBBModel->Initialise(p, ptr);
103 SetParticle(p);
104
105 // always false before the run
106 SetDeexcitationFlag(false);
107
108 if(nullptr == fParticleChange) {
109 fParticleChange = GetParticleChangeForLoss();
110 if(UseAngularGeneratorFlag() && nullptr == GetAngularDistribution()) {
112 }
113 }
114 if(nullptr == lsdata) {
115 G4AutoLock l(&ionXSMutex);
116 if(nullptr == lsdata) {
117 isFirst = true;
118 lsdata = new G4LindhardSorensenData();
119 fIonData = new G4IonICRU73Data();
120 }
121 l.unlock();
122 }
123 if(isFirst) {
124 fIonData->Initialise();
125 }
126}
127
128//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
129
132 const G4Material* mat,
133 G4double kinEnergy)
134{
135 chargeSquareRatio = corr->EffectiveChargeSquareRatio(p, mat, kinEnergy);
136 fBraggModel->SetChargeSquareRatio(chargeSquareRatio);
137 fBBModel->SetChargeSquareRatio(chargeSquareRatio);
138 return chargeSquareRatio;
139}
140
141//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
142
145 const G4Material* mat,
146 G4double kinEnergy)
147{
148 return corr->GetParticleCharge(p,mat,kinEnergy);
149}
150
151//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
152
153void G4LindhardSorensenIonModel::SetupParameters()
154{
155 mass = particle->GetPDGMass();
156 spin = particle->GetPDGSpin();
157 charge = particle->GetPDGCharge()*inveplus;
158 Zin = G4lrint(std::abs(charge));
159 eRatio = CLHEP::electron_mass_c2/mass;
160 pRatio = CLHEP::proton_mass_c2/mass;
161 const G4double aMag =
162 1./(0.5*CLHEP::eplus*CLHEP::hbar_Planck*CLHEP::c_squared);
163 G4double magmom = particle->GetPDGMagneticMoment()*mass*aMag;
164 magMoment2 = magmom*magmom - 1.0;
165 G4double x = 0.8426*CLHEP::GeV;
166 if(spin == 0.0 && mass < GeV) { x = 0.736*CLHEP::GeV; }
167 else if (Zin > 1) { x /= nist->GetA27(Zin); }
168
169 formfact = 2.0*CLHEP::electron_mass_c2/(x*x);
170 tlimit = 2.0/formfact;
171}
172
173//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
174
180
181//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
182
185 const G4ParticleDefinition* p,
186 G4double kinEnergy,
187 G4double cutEnergy,
188 G4double maxKinEnergy)
189{
190 // take into account formfactor
191 G4double tmax = MaxSecondaryEnergy(p, kinEnergy);
192 G4double emax = std::min(tmax, maxKinEnergy);
193 G4double escaled = kinEnergy*pRatio;
194 G4double cross = (escaled <= fElimit)
195 ? fBraggModel->ComputeCrossSectionPerElectron(p,kinEnergy,cutEnergy,emax)
196 : fBBModel->ComputeCrossSectionPerElectron(p,kinEnergy,cutEnergy,emax);
197 // G4cout << "LS: e= " << kinEnergy << " tmin= " << cutEnergy
198 // << " tmax= " << maxEnergy << " cross= " << cross << G4endl;
199 return cross;
200}
201
202//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
203
205 const G4ParticleDefinition* p,
206 G4double kineticEnergy,
208 G4double cutEnergy,
209 G4double maxEnergy)
210{
211 return Z*ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
212}
213
214//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
215
217 const G4Material* material,
218 const G4ParticleDefinition* p,
219 G4double kineticEnergy,
220 G4double cutEnergy,
221 G4double maxEnergy)
222{
223 return material->GetElectronDensity()
224 *ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
225}
226
227//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
228
231 const G4ParticleDefinition* p,
232 G4double kinEnergy,
233 G4double cut)
234{
235 // formfactor is taken into account in CorrectionsAlongStep(..)
236 G4double tmax = MaxSecondaryEnergy(p, kinEnergy);
237 G4double cutEnergy = std::min(std::min(cut,tmax), tlimit);
238
239 G4double escaled = kinEnergy*pRatio;
240 G4double dedx = (escaled <= fElimit)
241 ? fBraggModel->ComputeDEDXPerVolume(mat, p, kinEnergy, cutEnergy)
242 : fBBModel->ComputeDEDXPerVolume(mat, p, kinEnergy, cutEnergy);
243
244 //G4cout << "E(MeV)=" << kinEnergy/MeV << " dedx=" << dedx
245 // << " " << mat->GetName() << " Ecut(MeV)=" << cutEnergy << G4endl;
246 return dedx;
247}
248
249//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
250
252 const G4Material* mat,
253 const G4ParticleDefinition* p,
254 const G4double preKinEnergy,
255 const G4double cut,
256 const G4double& length,
257 G4double& eloss)
258{
259 if(eloss >= preKinEnergy) { return; }
260
261 SetParticle(p);
262 const G4double eDensity = mat->GetElectronDensity();
263 const G4double e = std::max(preKinEnergy - eloss*0.5, preKinEnergy*0.5);
264 const G4double tmax = MaxSecondaryEnergy(p, e);
265 const G4double escaled = e*pRatio;
266 const G4double tau = e/mass;
267 const G4double q2 = corr->EffectiveChargeSquareRatio(p, mat, e);
268 const G4int Z = p->GetAtomicNumber();
269
270 G4double res = 0.0;
271 if (escaled <= fElimit) {
272 // data from ICRU73 or ICRU90
273 if(Z > 2 && Z <= 80) {
274 res = fIonData->GetDEDX(mat, Z, escaled, G4Log(escaled));
275 /*
276 G4cout << " GetDEDX for Z=" << Z << " in " << mat->GetName()
277 << " Escaled=" << escaled << " E="
278 << e << " dEdx=" << res << G4endl;
279 */
280 }
281 if (res > 0.0) {
282 if (cut < tmax) {
283 const G4double x = cut/tmax;
284 res += (G4Log(x)*(tau + 1.)*(tau + 1.)/(tau * (tau + 2.0)) + 1.0 - x)
285 *q2*CLHEP::twopi_mc2_rcl2*eDensity;
286 }
287 res *= length;
288 } else {
289 // simplified correction
290 res = eloss*q2/chargeSquareRatio;
291 }
292 } else {
293 // Lindhard-Sorensen model
294 const G4double gam = tau + 1.0;
295 const G4double beta2 = tau * (tau+2.0)/(gam*gam);
296 G4double deltaL0 = 2.0*corr->BarkasCorrection(p, mat, e)*(charge-1.)/charge;
297 G4double deltaL = lsdata->GetDeltaL(Zin, gam);
298
299 res = eloss +
300 CLHEP::twopi_mc2_rcl2*q2*eDensity*(deltaL+deltaL0)*length/beta2;
301 /*
302 G4cout << " E(GeV)=" << preKinEnergy/GeV << " eloss(MeV)=" << eloss
303 << " L= " << eloss*beta2/(twopi_mc2_rcl2*q2*eDensity*length)
304 << " dL0= " << deltaL0
305 << " dL= " << deltaL << " dE(MeV)=" << res - eloss << G4endl;
306 */
307 }
308 if(res > preKinEnergy || 2*res < eloss) { res = eloss; }
309 /*
310 G4cout << " G4LindhardSorensenIonModel::CorrectionsAlongStep: E(GeV)="
311 << preKinEnergy/GeV << " eloss(MeV)=" << eloss
312 << " res(MeV)=" << res << G4endl;
313 */
314 eloss = res;
315}
316
317//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
318
320 std::vector<G4DynamicParticle*>* vdp,
321 const G4MaterialCutsCouple* couple,
322 const G4DynamicParticle* dp,
323 G4double cut,
324 G4double maxEnergy)
325{
326 G4double kineticEnergy = dp->GetKineticEnergy();
327 // take into account formfactor
328 const G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(), kineticEnergy);
329 const G4double minKinEnergy = std::min(cut, tmax);
330 const G4double maxKinEnergy = std::min(maxEnergy, tmax);
331 if(minKinEnergy >= maxKinEnergy) { return; }
332
333 //G4cout << "G4LindhardSorensenIonModel::SampleSecondaries Emin= "
334 // << minKinEnergy << " Emax= " << maxKinEnergy << G4endl;
335
336 G4double totEnergy = kineticEnergy + mass;
337 G4double etot2 = totEnergy*totEnergy;
338 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
339
340 G4double deltaKinEnergy, f;
341 G4double f1 = 0.0;
342 G4double fmax = 1.0;
343 if( 0.0 < spin ) { fmax += 0.5*maxKinEnergy*maxKinEnergy/etot2; }
344
345 CLHEP::HepRandomEngine* rndmEngineMod = G4Random::getTheEngine();
346 G4double rndm[2];
347
348 // sampling without nuclear size effect
349 do {
350 rndmEngineMod->flatArray(2, rndm);
351 deltaKinEnergy = minKinEnergy*maxKinEnergy
352 /(minKinEnergy*(1.0 - rndm[0]) + maxKinEnergy*rndm[0]);
353
354 f = 1.0 - beta2*deltaKinEnergy/tmax;
355 if( 0.0 < spin ) {
356 f1 = 0.5*deltaKinEnergy*deltaKinEnergy/etot2;
357 f += f1;
358 }
359
360 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
361 } while( fmax*rndm[1] > f);
362
363 // projectile formfactor - suppresion of high energy
364 // delta-electron production at high energy
365
366 G4double x = formfact*deltaKinEnergy;
367 if(x > 1.e-6) {
368
369 G4double x1 = 1.0 + x;
370 G4double grej = 1.0/(x1*x1);
371 if( 0.0 < spin ) {
372 G4double x2 = 0.5*electron_mass_c2*deltaKinEnergy/(mass*mass);
373 grej *= (1.0 + magMoment2*(x2 - f1/f)/(1.0 + x2));
374 }
375 if(grej > 1.1) {
376 G4cout << "### G4LindhardSorensenIonModel WARNING: grej= " << grej
377 << " " << dp->GetDefinition()->GetParticleName()
378 << " Ekin(MeV)= " << kineticEnergy
379 << " delEkin(MeV)= " << deltaKinEnergy
380 << G4endl;
381 }
382 if(rndmEngineMod->flat() > grej) { return; }
383 }
384
385 G4ThreeVector deltaDirection;
386
388
389 const G4Material* mat = couple->GetMaterial();
391
392 deltaDirection =
393 GetAngularDistribution()->SampleDirection(dp, deltaKinEnergy, Z, mat);
394
395 } else {
396
397 G4double deltaMomentum =
398 std::sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
399 G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
400 (deltaMomentum * dp->GetTotalMomentum());
401 cost = std::min(cost, 1.0);
402 G4double sint = std::sqrt((1.0 - cost)*(1.0 + cost));
403
404 G4double phi = CLHEP::twopi*rndmEngineMod->flat();
405
406 deltaDirection.set(sint*std::cos(phi),sint*std::sin(phi), cost) ;
407 deltaDirection.rotateUz(dp->GetMomentumDirection());
408 }
409 /*
410 G4cout << "### G4LindhardSorensenIonModel "
411 << dp->GetDefinition()->GetParticleName()
412 << " Ekin(MeV)= " << kineticEnergy
413 << " delEkin(MeV)= " << deltaKinEnergy
414 << " tmin(MeV)= " << minKinEnergy
415 << " tmax(MeV)= " << maxKinEnergy
416 << " dir= " << dp->GetMomentumDirection()
417 << " dirDelta= " << deltaDirection
418 << G4endl;
419 */
420 // create G4DynamicParticle object for delta ray
421 auto delta = new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
422
423 vdp->push_back(delta);
424
425 // Change kinematics of primary particle
426 kineticEnergy -= deltaKinEnergy;
427 G4ThreeVector finalP = dp->GetMomentum() - delta->GetMomentum();
428 finalP = finalP.unit();
429
430 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
431 fParticleChange->SetProposedMomentumDirection(finalP);
432}
433
434//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
435
438 G4double kinEnergy)
439{
440 // here particle type is checked for any method
441 SetParticle(pd);
442 G4double tau = kinEnergy/mass;
443 return 2.0*CLHEP::electron_mass_c2*tau*(tau + 2.) /
444 (1. + 2.0*(tau + 1.)*eRatio + eRatio*eRatio);
445}
446
447//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4TemplateAutoLock< G4Mutex > G4AutoLock
G4double G4Log(G4double x)
Definition G4Log.hh:169
#define G4MUTEX_INITIALIZER
std::mutex G4Mutex
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
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
G4double GetTotalMomentum() const
static G4Electron * Electron()
Definition G4Electron.cc:91
G4double GetMeanExcitationEnergy() const
void CorrectionsAlongStep(const G4Material *, const G4ParticleDefinition *, const G4double kinEnergy, const G4double cutEnergy, const G4double &length, G4double &eloss) override
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
G4double GetChargeSquareRatio(const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4LindhardSorensenIonModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="LindhardSorensen")
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 GetParticleCharge(const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
const G4Material * GetMaterial() const
G4IonisParamMat * GetIonisation() const
G4double GetElectronDensity() const
static G4NistManager * Instance()
G4double GetPDGMagneticMoment() const
G4int GetAtomicNumber() 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
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