Geant4 11.1.1
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 "G4BraggIonModel.hh"
56#include "G4BetheBlochModel.hh"
57#include "G4IonICRU73Data.hh"
58
59//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
60
61using namespace std;
62
63G4LindhardSorensenData* G4LindhardSorensenIonModel::lsdata = nullptr;
64G4IonICRU73Data* G4LindhardSorensenIonModel::fIonData = nullptr;
65std::vector<G4float>* G4LindhardSorensenIonModel::fact[] = {nullptr};
66
68 const G4String& nam)
69 : G4VEmModel(nam),
70 particle(nullptr),
71 twoln10(2.0*G4Log(10.0))
72{
73 fParticleChange = nullptr;
74 theElectron = G4Electron::Electron();
77 fBraggModel = new G4BraggIonModel();
78 fBBModel = new G4BetheBlochModel();
79 fElimit = 2.0*CLHEP::MeV;
80}
81
82//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
83
85
86//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
87
89 const G4DataVector& ptr)
90{
91 fBraggModel->Initialise(p, ptr);
92 fBBModel->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()) {
107 if(nullptr == lsdata) {
108 lsdata = new G4LindhardSorensenData();
109 }
110 if(nullptr == fIonData) {
111 fIonData = new G4IonICRU73Data();
112 }
113 fIonData->Initialise();
114 }
115}
116
117//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
118
121 const G4Material* mat,
122 G4double kinEnergy)
123{
124 return corr->EffectiveChargeSquareRatio(p,mat,kinEnergy);
125}
126
127//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
128
131 const G4Material* mat,
132 G4double kinEnergy)
133{
134 return corr->GetParticleCharge(p,mat,kinEnergy);
135}
136
137//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
138
139void G4LindhardSorensenIonModel::SetupParameters()
140{
141 mass = particle->GetPDGMass();
142 spin = particle->GetPDGSpin();
143 charge = particle->GetPDGCharge()*inveplus;
144 Zin = G4lrint(std::abs(charge));
145 chargeSquare = charge*charge;
146 eRatio = CLHEP::electron_mass_c2/mass;
147 pRatio = CLHEP::proton_mass_c2/mass;
148 const G4double aMag =
149 1./(0.5*CLHEP::eplus*CLHEP::hbar_Planck*CLHEP::c_squared);
150 G4double magmom = particle->GetPDGMagneticMoment()*mass*aMag;
151 magMoment2 = magmom*magmom - 1.0;
152 G4double x = 0.8426*CLHEP::GeV;
153 if(spin == 0.0 && mass < GeV) { x = 0.736*CLHEP::GeV; }
154 else if (Zin > 1) { x /= nist->GetA27(Zin); }
155
156 formfact = 2.0*CLHEP::electron_mass_c2/(x*x);
157 tlimit = 2.0/formfact;
158}
159
160//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
161
163 const G4MaterialCutsCouple* couple)
164{
166}
167
168//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
169
172 const G4ParticleDefinition* p,
173 G4double kinEnergy,
174 G4double cutEnergy,
175 G4double maxKinEnergy)
176{
177 // take into account formfactor
178 G4double tmax = MaxSecondaryEnergy(p, kinEnergy);
179 G4double emax = std::min(tmax, maxKinEnergy);
180 G4double escaled = kinEnergy*pRatio;
181 G4double cross = (escaled <= fElimit)
182 ? fBraggModel->ComputeCrossSectionPerElectron(p,kinEnergy,cutEnergy,emax)
183 : fBBModel->ComputeCrossSectionPerElectron(p,kinEnergy,cutEnergy,emax);
184 // G4cout << "LS: e= " << kinEnergy << " tmin= " << cutEnergy
185 // << " tmax= " << maxEnergy << " cross= " << cross << G4endl;
186 return cross;
187}
188
189//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
190
192 const G4ParticleDefinition* p,
193 G4double kineticEnergy,
195 G4double cutEnergy,
196 G4double maxEnergy)
197{
198 return Z*ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
199}
200
201//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
202
204 const G4Material* material,
205 const G4ParticleDefinition* p,
206 G4double kineticEnergy,
207 G4double cutEnergy,
208 G4double maxEnergy)
209{
210 return material->GetElectronDensity()
211 *ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
212}
213
214//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
215
218 const G4ParticleDefinition* p,
219 G4double kinEnergy,
220 G4double cut)
221{
222 // formfactor is taken into account in CorrectionsAlongStep(..)
223 G4double tmax = MaxSecondaryEnergy(p, kinEnergy);
224 G4double cutEnergy = std::min(std::min(cut,tmax), tlimit);
225
226 G4double escaled = kinEnergy*pRatio;
227 G4double dedx = (escaled <= fElimit)
228 ? fBraggModel->ComputeDEDXPerVolume(mat, p, kinEnergy, cutEnergy)
229 : fBBModel->ComputeDEDXPerVolume(mat, p, kinEnergy, cutEnergy);
230
231 //G4cout << "E(MeV)=" << kinEnergy/MeV << " dedx=" << dedx
232 // << " " << material->GetName() << " Ecut(MeV)=" << cutEnergy << G4endl;
233 return dedx;
234}
235
236//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
237
239 const G4MaterialCutsCouple* couple,
240 const G4DynamicParticle* dp,
241 const G4double& length,
242 G4double& eloss)
243{
244 // no correction at the last step
245 const G4double preKinEnergy = dp->GetKineticEnergy();
246 if(eloss >= preKinEnergy) { return; }
247
248 const G4ParticleDefinition* p = dp->GetDefinition();
249 SetParticle(p);
250 const G4Material* mat = couple->GetMaterial();
251 const G4double eDensity = mat->GetElectronDensity();
252 const G4double e = std::max(preKinEnergy - eloss*0.5, preKinEnergy*0.5);
253 const G4double tmax = MaxSecondaryEnergy(p, e);
254 const G4double escaled = e*pRatio;
255 const G4double tau = e/mass;
256
257 const G4double q20 = corr->EffectiveChargeSquareRatio(p, mat, preKinEnergy);
258 const G4double q2 = corr->EffectiveChargeSquareRatio(p, mat, e);
259 const G4int Z = p->GetAtomicNumber();
260
261 G4double res;
262 if(escaled <= fElimit) {
263 // data from ICRU73 or ICRU90
264 res = fIonData->GetDEDX(mat, Z, escaled, G4Log(escaled));
265 /*
266 G4cout << "GetDEDX for Z=" << Z << " in " << mat->GetName()
267 << " Escaled=" << escaled << " E="
268 << e << " dEdx=" << res << G4endl;
269 */
270 if(res > 0.0) {
271 auto pcuts = couple->GetProductionCuts();
272 G4double cut = (nullptr == pcuts) ? tmax : pcuts->GetProductionCut(1);
273 if(cut < tmax) {
274 const G4double x = cut/tmax;
275 res += (G4Log(x)*(tau + 1.)*(tau + 1.)/(tau * (tau + 2.0)) + 1.0 - x)
276 *q2*CLHEP::twopi_mc2_rcl2*eDensity;
277 }
278 res *= length;
279 } else {
280 // simplified correction
281 res = eloss*q2/q20;
282 }
283 } else {
284 // Lindhard-Sorensen model
285 const G4double gam = tau + 1.0;
286 const G4double beta2 = tau * (tau+2.0)/(gam*gam);
287 G4double deltaL0 = 2.0*corr->BarkasCorrection(p, mat, e)*(charge-1.)/charge;
288 G4double deltaL = lsdata->GetDeltaL(Zin, gam);
289
290 res = eloss +
291 CLHEP::twopi_mc2_rcl2*q2*eDensity*(deltaL+deltaL0)*length/beta2;
292 /*
293 G4cout << "G4LindhardSorensenIonModel::CorrectionsAlongStep: E(GeV)= "
294 << preKinEnergy/GeV << " eloss(MeV)= " << eloss
295 << " L= " << eloss*beta2/(twopi_mc2_rcl2*chargeSquare*eDensity*length)
296 << " dL0= " << deltaL0
297 << " dL= " << deltaL << G4endl;
298 */
299 }
300 if(res > preKinEnergy) { res = preKinEnergy; }
301 else if(res < 0.0) { res = eloss; }
302 /*
303 G4cout << "G4LindhardSorensenIonModel::CorrectionsAlongStep: E(GeV)="
304 << preKinEnergy/GeV << " eloss(MeV)=" << eloss
305 << " res(MeV)=" << res << G4endl;
306 */
307 eloss = res;
308}
309
310//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
311
313 vector<G4DynamicParticle*>* vdp,
314 const G4MaterialCutsCouple* couple,
315 const G4DynamicParticle* dp,
316 G4double minKinEnergy,
317 G4double maxEnergy)
318{
319 G4double kineticEnergy = dp->GetKineticEnergy();
320 // take into account formfactor
321 G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(),kineticEnergy);
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;
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 std::sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
392 G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
393 (deltaMomentum * dp->GetTotalMomentum());
394 cost = std::min(cost, 1.0);
395 G4double sint = std::sqrt((1.0 - cost)*(1.0 + cost));
396
397 G4double phi = CLHEP::twopi*rndmEngineMod->flat();
398
399 deltaDirection.set(sint*std::cos(phi),sint*std::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 auto delta = new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
415
416 vdp->push_back(delta);
417
418 // Change kinematics of primary particle
419 kineticEnergy -= deltaKinEnergy;
420 G4ThreeVector finalP = dp->GetMomentum() - delta->GetMomentum();
421 finalP = finalP.unit();
422
423 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
424 fParticleChange->SetProposedMomentumDirection(finalP);
425}
426
427//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
428
431 G4double kinEnergy)
432{
433 // here particle type is checked for any method
434 SetParticle(pd);
435 G4double tau = kinEnergy/mass;
436 return 2.0*CLHEP::electron_mass_c2*tau*(tau + 2.) /
437 (1. + 2.0*(tau + 1.)*eRatio + eRatio*eRatio);
438}
439
440//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4double G4Log(G4double x)
Definition: G4Log.hh:227
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
#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
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 Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) 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 EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double BarkasCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double GetDEDX(const G4Material *, const G4int Z, const G4double e, const G4double loge) const
G4double GetMeanExcitationEnergy() const
G4double GetDeltaL(G4int Z, G4double gamma) const
~G4LindhardSorensenIonModel() 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 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")
void CorrectionsAlongStep(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, const G4double &length, G4double &eloss) 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 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
G4ProductionCuts * GetProductionCuts() const
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:221
G4double GetElectronDensity() const
Definition: G4Material.hh:212
G4double GetA27(G4int Z) const
static G4NistManager * Instance()
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
G4double GetPDGMagneticMoment() const
G4int GetAtomicNumber() const
G4double GetPDGCharge() const
const G4String & GetParticleName() const
G4double GetProductionCut(G4int index) const
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:600
G4double inveplus
Definition: G4VEmModel.hh:427
G4int SelectRandomAtomNumber(const G4Material *) const
Definition: G4VEmModel.cc:253
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
int G4lrint(double ad)
Definition: templates.hh:134