Geant4 10.7.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//
29// GEANT4 Class header file
30//
31//
32// File name: G4LindhardSorensenIonModel
33//
34// Author: Alexander Bagulya & Vladimir Ivanchenko
35//
36// Creation date: 16.04.2018
37//
38//
39// -------------------------------------------------------------------
40//
41
42
43//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
44//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
45
47#include "Randomize.hh"
49#include "G4SystemOfUnits.hh"
50#include "G4Electron.hh"
51#include "G4LossTableManager.hh"
52#include "G4EmCorrections.hh"
54#include "G4Log.hh"
55#include "G4DeltaAngle.hh"
57#include "G4BraggIonModel.hh"
58
59//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
60
61using namespace std;
62
63G4LindhardSorensenData* G4LindhardSorensenIonModel::lsdata = nullptr;
64std::vector<G4float>* G4LindhardSorensenIonModel::fact[] = {nullptr};
65
67 const G4String& nam)
68 : G4VEmModel(nam),
69 particle(nullptr),
70 tlimit(DBL_MAX),
71 twoln10(2.0*G4Log(10.0))
72{
73 fParticleChange = nullptr;
74 theElectron = G4Electron::Electron();
75 SetParticle(theElectron);
78 fBraggIonModel = new G4BraggIonModel();
79 SetLowEnergyLimit(2.0*MeV);
80}
81
82//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
83
85{}
86
87//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
88
90 const G4DataVector& ptr)
91{
92 fBraggIonModel->Initialise(p, ptr);
93 SetParticle(p);
94 //G4cout << "G4LindhardSorensenIonModel::Initialise for "
95 // << p->GetParticleName() << G4endl;
96
97 // always false before the run
99
100 if(nullptr == fParticleChange) {
101 fParticleChange = GetParticleChangeForLoss();
102 if(UseAngularGeneratorFlag() && nullptr == GetAngularDistribution()) {
104 }
105 }
106 if(IsMaster() && nullptr == lsdata) {
107 lsdata = new G4LindhardSorensenData();
108 }
109}
110
111//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
112
115 const G4Material*,
116 G4double)
117{
118 return chargeSquare;
119}
120
121//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
122
125 const G4Material*,
126 G4double)
127{
128 return charge;
129}
130
131//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
132
133void G4LindhardSorensenIonModel::SetupParameters()
134{
135 mass = particle->GetPDGMass();
136 spin = particle->GetPDGSpin();
137 charge = particle->GetPDGCharge()*inveplus;
138 Zin = G4lrint(charge);
139 chargeSquare = charge*charge;
140 ratio = electron_mass_c2/mass;
141 static const G4double aMag = 1./(0.5*eplus*hbar_Planck*c_squared);
142 G4double magmom = particle->GetPDGMagneticMoment()*mass*aMag;
143 magMoment2 = magmom*magmom - 1.0;
144 if(Zin <= 1) {
145 formfact = (spin == 0.0 && mass < GeV) ? 1.181e-6 : 1.548e-6;
146 } else {
147 G4double x = nist->GetA27(Zin);
148 formfact = 3.969e-6*x*x;
149 }
150 tlimit = std::sqrt(0.414/formfact +
151 electron_mass_c2*electron_mass_c2) - electron_mass_c2;
152}
153
154//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
155
157 const G4MaterialCutsCouple* couple)
158{
160}
161
162//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
163
166 const G4ParticleDefinition* p,
167 G4double kineticEnergy,
168 G4double cutEnergy,
169 G4double maxKinEnergy)
170{
171 G4double cross = 0.0;
172 // take into account formfactor
173 G4double tmax = std::min(MaxSecondaryEnergy(p, kineticEnergy),tlimit);
174 G4double maxEnergy = std::min(tmax,maxKinEnergy);
175 if(cutEnergy < maxEnergy) {
176
177 G4double totEnergy = kineticEnergy + mass;
178 G4double energy2 = totEnergy*totEnergy;
179 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
180
181 cross = (maxEnergy - cutEnergy)/(cutEnergy*maxEnergy)
182 - beta2*G4Log(maxEnergy/cutEnergy)/tmax;
183
184 // +term for spin=1/2 particle
185 if( 0.0 < spin ) { cross += 0.5*(maxEnergy - cutEnergy)/energy2; }
186
187 cross *= twopi_mc2_rcl2*chargeSquare/beta2;
188 }
189
190 // G4cout << "BB: e= " << kineticEnergy << " tmin= " << cutEnergy
191 // << " tmax= " << tmax << " cross= " << cross << G4endl;
192
193 return cross;
194}
195
196//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
197
199 const G4ParticleDefinition* p,
200 G4double kineticEnergy,
202 G4double cutEnergy,
203 G4double maxEnergy)
204{
205 return Z*ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
206}
207
208//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
209
211 const G4Material* material,
212 const G4ParticleDefinition* p,
213 G4double kineticEnergy,
214 G4double cutEnergy,
215 G4double maxEnergy)
216{
217 return material->GetElectronDensity()
218 *ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
219}
220
221//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
222
225 const G4ParticleDefinition* p,
226 G4double kineticEnergy,
227 G4double cut)
228{
229 // formfactor is taken into account in CorrectionsAlongStep(..)
230 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
231 G4double cutEnergy = std::min(cut,tmax);
232
233 G4double tau = kineticEnergy/mass;
234 G4double gam = tau + 1.0;
235 G4double bg2 = tau * (tau+2.0);
236 G4double beta2 = bg2/(gam*gam);
237
238 G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
239 G4double eexc2 = eexc*eexc;
240
241 G4double eDensity = material->GetElectronDensity();
242
243 G4double dedx = G4Log(2.0*electron_mass_c2*bg2*cutEnergy/eexc2)
244 - (1.0 + cutEnergy/tmax)*beta2;
245
246 if(0.0 < spin) {
247 G4double del = 0.5*cutEnergy/(kineticEnergy + mass);
248 dedx += del*del;
249 }
250
251 // density correction
252 G4double x = G4Log(bg2)/twoln10;
253 dedx -= material->GetIonisation()->DensityCorrection(x);
254
255 // shell correction
256 dedx -= 2.0*corr->ShellCorrection(p,material,kineticEnergy);
257
258 //High order correction different for hadrons and ions
259 dedx += 2.0*corr->BarkasCorrection(p,material,kineticEnergy);
260 dedx = std::max(dedx, 0.0);
261
262 // now compute the total ionization loss
263 dedx *= twopi_mc2_rcl2*chargeSquare*eDensity/beta2;
264
265 //G4cout << "E(MeV)= " << kineticEnergy/MeV << " dedx= " << dedx
266 // << " " << material->GetName() << G4endl;
267
268 return dedx;
269}
270
271//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
272
273void
275 const G4DynamicParticle* dp,
276 G4double& eloss,
277 G4double&,
278 G4double length)
279{
280 const G4ParticleDefinition* p = dp->GetDefinition();
281 SetParticle(p);
282 G4double eDensity = couple->GetMaterial()->GetElectronDensity();
283 G4double preKinEnergy = dp->GetKineticEnergy();
284 G4double e = preKinEnergy - eloss*0.5;
286
287 G4double tau = e/mass;
288 G4double gam = tau + 1.0;
289 G4double beta2 = tau * (tau+2.0)/(gam*gam);
290 G4double deltaL0 =
291 2.0*corr->BarkasCorrection (p, couple->GetMaterial(), e)*(charge-1.)/charge;
292 G4double deltaL = lsdata->GetDeltaL(Zin, gam);
293
294 G4double elossnew =
295 eloss + twopi_mc2_rcl2*chargeSquare*eDensity*(deltaL+deltaL0)*length/beta2;
296 /*
297 G4cout << "G4LindhardSorensenIonModel::CorrectionsAlongStep: E(GeV)= "
298 << preKinEnergy/GeV << " eloss(MeV)= " << eloss
299 << " L= " << eloss*beta2/(twopi_mc2_rcl2*chargeSquare*eDensity*length)
300 << " dL0= " << deltaL0
301 << " dL= " << deltaL << G4endl;
302 */
303 if(elossnew > preKinEnergy) { elossnew = preKinEnergy; }
304 else if(elossnew < 0.0) { elossnew = eloss*0.5; }
305
306 eloss = elossnew;
307}
308
309//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
310
312 vector<G4DynamicParticle*>* vdp,
313 const G4MaterialCutsCouple* couple,
314 const G4DynamicParticle* dp,
315 G4double minKinEnergy,
316 G4double maxEnergy)
317{
318 G4double kineticEnergy = dp->GetKineticEnergy();
319 // take into account formfactor
320 G4double tmax =
321 std::min(MaxSecondaryEnergy(dp->GetDefinition(),kineticEnergy),tlimit);
322
323 G4double maxKinEnergy = std::min(maxEnergy,tmax);
324 if(minKinEnergy >= maxKinEnergy) { return; }
325
326 //G4cout << "G4LindhardSorensenIonModel::SampleSecondaries Emin= "
327 // << minKinEnergy << " Emax= " << maxKinEnergy << G4endl;
328
329 G4double totEnergy = kineticEnergy + mass;
330 G4double etot2 = totEnergy*totEnergy;
331 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
332
333 G4double deltaKinEnergy, f;
334 G4double f1 = 0.0;
335 G4double fmax = 1.0;
336 if( 0.0 < spin ) { fmax += 0.5*maxKinEnergy*maxKinEnergy/etot2; }
337
338 CLHEP::HepRandomEngine* rndmEngineMod = G4Random::getTheEngine();
339 G4double rndm[2];
340
341 // sampling without nuclear size effect
342 do {
343 rndmEngineMod->flatArray(2, rndm);
344 deltaKinEnergy = minKinEnergy*maxKinEnergy
345 /(minKinEnergy*(1.0 - rndm[0]) + maxKinEnergy*rndm[0]);
346
347 f = 1.0 - beta2*deltaKinEnergy/tmax;
348 if( 0.0 < spin ) {
349 f1 = 0.5*deltaKinEnergy*deltaKinEnergy/etot2;
350 f += f1;
351 }
352
353 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
354 } while( fmax*rndm[1] > f);
355
356 // projectile formfactor - suppresion of high energy
357 // delta-electron production at high energy
358
359 G4double x = formfact*deltaKinEnergy*(deltaKinEnergy + 2*electron_mass_c2);
360 if(x > 1.e-6) {
361
362 G4double x1 = 1.0 + x;
363 G4double grej = 1.0/(x1*x1);
364 if( 0.0 < spin ) {
365 G4double x2 = 0.5*electron_mass_c2*deltaKinEnergy/(mass*mass);
366 grej *= (1.0 + magMoment2*(x2 - f1/f)/(1.0 + x2));
367 }
368 if(grej > 1.1) {
369 G4cout << "### G4LindhardSorensenIonModel WARNING: grej= " << grej
370 << " " << dp->GetDefinition()->GetParticleName()
371 << " Ekin(MeV)= " << kineticEnergy
372 << " delEkin(MeV)= " << deltaKinEnergy
373 << G4endl;
374 }
375 if(rndmEngineMod->flat() > grej) { return; }
376 }
377
378 G4ThreeVector deltaDirection;
379
381
382 const G4Material* mat = couple->GetMaterial();
384
385 deltaDirection =
386 GetAngularDistribution()->SampleDirection(dp, deltaKinEnergy, Z, mat);
387
388 } else {
389
390 G4double deltaMomentum =
391 sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
392 G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
393 (deltaMomentum * dp->GetTotalMomentum());
394 if(cost > 1.0) { cost = 1.0; }
395 G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
396
397 G4double phi = twopi*rndmEngineMod->flat();
398
399 deltaDirection.set(sint*cos(phi),sint*sin(phi), cost) ;
400 deltaDirection.rotateUz(dp->GetMomentumDirection());
401 }
402 /*
403 G4cout << "### G4LindhardSorensenIonModel "
404 << dp->GetDefinition()->GetParticleName()
405 << " Ekin(MeV)= " << kineticEnergy
406 << " delEkin(MeV)= " << deltaKinEnergy
407 << " tmin(MeV)= " << minKinEnergy
408 << " tmax(MeV)= " << maxKinEnergy
409 << " dir= " << dp->GetMomentumDirection()
410 << " dirDelta= " << deltaDirection
411 << G4endl;
412 */
413 // create G4DynamicParticle object for delta ray
414 G4DynamicParticle* delta =
415 new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
416
417 vdp->push_back(delta);
418
419 // Change kinematics of primary particle
420 kineticEnergy -= deltaKinEnergy;
421 G4ThreeVector finalP = dp->GetMomentum() - delta->GetMomentum();
422 finalP = finalP.unit();
423
424 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
425 fParticleChange->SetProposedMomentumDirection(finalP);
426}
427
428//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
429
432 G4double kinEnergy)
433{
434 // here particle type is checked for any method
435 SetParticle(pd);
436 G4double tau = kinEnergy/mass;
437 G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
438 (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
439 // formfactor is not taken into account
440 return tmax;
441}
442
443//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4double G4Log(G4double x)
Definition: G4Log.hh:226
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#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
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
G4double GetTotalMomentum() const
static G4Electron * Electron()
Definition: G4Electron.cc:93
G4double BarkasCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double ShellCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double GetMeanExcitationEnergy() const
G4double DensityCorrection(G4double x)
G4double GetDeltaL(G4int Z, G4double gamma) const
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override
void CorrectionsAlongStep(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double &eloss, G4double &, G4double length) 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 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
Definition: G4Material.hh:224
G4double GetElectronDensity() const
Definition: G4Material.hh:215
G4double GetA27(G4int Z) const
static G4NistManager * Instance()
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
G4double GetPDGMagneticMoment() const
G4double GetPDGCharge() const
const G4String & GetParticleName() const
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
virtual void SetParticleAndCharge(const G4ParticleDefinition *, G4double q2)
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:604
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:611
G4double inveplus
Definition: G4VEmModel.hh:446
G4bool IsMaster() const
Definition: G4VEmModel.hh:736
G4int SelectRandomAtomNumber(const G4Material *)
Definition: G4VEmModel.cc:315
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:813
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:618
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:708
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:118
int G4lrint(double ad)
Definition: templates.hh:134
#define DBL_MAX
Definition: templates.hh:62