Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LindhardSorensenIonModel Class Reference

#include <G4LindhardSorensenIonModel.hh>

+ Inheritance diagram for G4LindhardSorensenIonModel:

Public Member Functions

 G4LindhardSorensenIonModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="LindhardSorensen")
 
 ~G4LindhardSorensenIonModel () override
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
 
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
 
G4LindhardSorensenIonModeloperator= (const G4LindhardSorensenIonModel &right)=delete
 
 G4LindhardSorensenIonModel (const G4LindhardSorensenIonModel &)=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)
 
void SetMasterThread (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)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Member Functions

G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kinEnergy) override
 
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 62 of file G4LindhardSorensenIonModel.hh.

Constructor & Destructor Documentation

◆ G4LindhardSorensenIonModel() [1/2]

G4LindhardSorensenIonModel::G4LindhardSorensenIonModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "LindhardSorensen" )
explicit

Definition at line 72 of file G4LindhardSorensenIonModel.cc.

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}
G4double G4Log(G4double x)
Definition G4Log.hh:227
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
static G4NistManager * Instance()
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67

◆ ~G4LindhardSorensenIonModel()

G4LindhardSorensenIonModel::~G4LindhardSorensenIonModel ( )
override

Definition at line 87 of file G4LindhardSorensenIonModel.cc.

87 {
88 if(isFirst) {
89 delete lsdata;
90 delete fIonData;
91 lsdata = nullptr;
92 fIonData = nullptr;
93 }
94}

◆ G4LindhardSorensenIonModel() [2/2]

G4LindhardSorensenIonModel::G4LindhardSorensenIonModel ( const G4LindhardSorensenIonModel & )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4LindhardSorensenIonModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition * p,
G4double kineticEnergy,
G4double Z,
G4double A,
G4double cutEnergy,
G4double maxEnergy )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 203 of file G4LindhardSorensenIonModel.cc.

209{
210 return Z*ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
211}
G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)

◆ ComputeCrossSectionPerElectron()

G4double G4LindhardSorensenIonModel::ComputeCrossSectionPerElectron ( const G4ParticleDefinition * p,
G4double kineticEnergy,
G4double cutEnergy,
G4double maxEnergy )

Definition at line 183 of file G4LindhardSorensenIonModel.cc.

188{
189 // take into account formfactor
190 G4double tmax = MaxSecondaryEnergy(p, kinEnergy);
191 G4double emax = std::min(tmax, maxKinEnergy);
192 G4double escaled = kinEnergy*pRatio;
193 G4double cross = (escaled <= fElimit)
194 ? fBraggModel->ComputeCrossSectionPerElectron(p,kinEnergy,cutEnergy,emax)
195 : fBBModel->ComputeCrossSectionPerElectron(p,kinEnergy,cutEnergy,emax);
196 // G4cout << "LS: e= " << kinEnergy << " tmin= " << cutEnergy
197 // << " tmax= " << maxEnergy << " cross= " << cross << G4endl;
198 return cross;
199}
double G4double
Definition G4Types.hh:83
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override

Referenced by ComputeCrossSectionPerAtom(), and CrossSectionPerVolume().

◆ ComputeDEDXPerVolume()

G4double G4LindhardSorensenIonModel::ComputeDEDXPerVolume ( const G4Material * mat,
const G4ParticleDefinition * p,
G4double kineticEnergy,
G4double cutEnergy )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 229 of file G4LindhardSorensenIonModel.cc.

233{
234 // formfactor is taken into account in CorrectionsAlongStep(..)
235 G4double tmax = MaxSecondaryEnergy(p, kinEnergy);
236 G4double cutEnergy = std::min(std::min(cut,tmax), tlimit);
237
238 G4double escaled = kinEnergy*pRatio;
239 G4double dedx = (escaled <= fElimit)
240 ? fBraggModel->ComputeDEDXPerVolume(mat, p, kinEnergy, cutEnergy)
241 : fBBModel->ComputeDEDXPerVolume(mat, p, kinEnergy, cutEnergy);
242
243 //G4cout << "E(MeV)=" << kinEnergy/MeV << " dedx=" << dedx
244 // << " " << mat->GetName() << " Ecut(MeV)=" << cutEnergy << G4endl;
245 return dedx;
246}
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override

◆ CorrectionsAlongStep()

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

Reimplemented from G4VEmModel.

Definition at line 250 of file G4LindhardSorensenIonModel.cc.

255{
256 // no correction at the last step
257 const G4double preKinEnergy = dp->GetKineticEnergy();
258 if(eloss >= preKinEnergy) { return; }
259
260 const G4ParticleDefinition* p = dp->GetDefinition();
261 SetParticle(p);
262 const G4Material* mat = couple->GetMaterial();
263 const G4double eDensity = mat->GetElectronDensity();
264 const G4double e = std::max(preKinEnergy - eloss*0.5, preKinEnergy*0.5);
265 const G4double tmax = MaxSecondaryEnergy(p, e);
266 const G4double escaled = e*pRatio;
267 const G4double tau = e/mass;
268 const G4double q2 = corr->EffectiveChargeSquareRatio(p, mat, e);
269 const G4int Z = p->GetAtomicNumber();
270
271 G4double res = 0.0;
272 if(escaled <= fElimit) {
273 // data from ICRU73 or ICRU90
274 if(Z > 2 && Z <= 80) {
275 res = fIonData->GetDEDX(mat, Z, escaled, G4Log(escaled));
276 /*
277 G4cout << " GetDEDX for Z=" << Z << " in " << mat->GetName()
278 << " Escaled=" << escaled << " E="
279 << e << " dEdx=" << res << G4endl;
280 */
281 }
282 if(res > 0.0) {
283 auto pcuts = couple->GetProductionCuts();
284 G4double cut = (nullptr == pcuts) ? tmax : pcuts->GetProductionCut(1);
285 if(cut < tmax) {
286 const G4double x = cut/tmax;
287 res += (G4Log(x)*(tau + 1.)*(tau + 1.)/(tau * (tau + 2.0)) + 1.0 - x)
288 *q2*CLHEP::twopi_mc2_rcl2*eDensity;
289 }
290 res *= length;
291 } else {
292 // simplified correction
293 res = eloss*q2/chargeSquare;
294 }
295 } else {
296 // Lindhard-Sorensen model
297 const G4double gam = tau + 1.0;
298 const G4double beta2 = tau * (tau+2.0)/(gam*gam);
299 G4double deltaL0 = 2.0*corr->BarkasCorrection(p, mat, e)*(charge-1.)/charge;
300 G4double deltaL = lsdata->GetDeltaL(Zin, gam);
301
302 res = eloss +
303 CLHEP::twopi_mc2_rcl2*q2*eDensity*(deltaL+deltaL0)*length/beta2;
304 /*
305 G4cout << " E(GeV)=" << preKinEnergy/GeV << " eloss(MeV)=" << eloss
306 << " L= " << eloss*beta2/(twopi_mc2_rcl2*q2*eDensity*length)
307 << " dL0= " << deltaL0
308 << " dL= " << deltaL << " dE(MeV)=" << res - eloss << G4endl;
309 */
310 }
311 if(res > preKinEnergy || 2*res < eloss) { res = eloss; }
312 /*
313 G4cout << " G4LindhardSorensenIonModel::CorrectionsAlongStep: E(GeV)="
314 << preKinEnergy/GeV << " eloss(MeV)=" << eloss
315 << " res(MeV)=" << res << G4endl;
316 */
317 eloss = res;
318}
int G4int
Definition G4Types.hh:85
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, const G4double kineticEnergy)
G4double BarkasCorrection(const G4ParticleDefinition *, const G4Material *, const G4double kineticEnergy, const G4bool isInitialized=false)
G4double GetDEDX(const G4Material *, const G4int Z, const G4double e, const G4double loge) const
G4double GetDeltaL(G4int Z, G4double gamma) const
const G4Material * GetMaterial() const
G4ProductionCuts * GetProductionCuts() const
G4double GetElectronDensity() const
G4int GetAtomicNumber() const

◆ CrossSectionPerVolume()

G4double G4LindhardSorensenIonModel::CrossSectionPerVolume ( const G4Material * material,
const G4ParticleDefinition * p,
G4double kineticEnergy,
G4double cutEnergy,
G4double maxEnergy )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 215 of file G4LindhardSorensenIonModel.cc.

221{
222 return material->GetElectronDensity()
223 *ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
224}

◆ GetChargeSquareRatio()

G4double G4LindhardSorensenIonModel::GetChargeSquareRatio ( const G4ParticleDefinition * p,
const G4Material * mat,
G4double kineticEnergy )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 131 of file G4LindhardSorensenIonModel.cc.

134{
135 chargeSquare = corr->EffectiveChargeSquareRatio(p, mat, kinEnergy);
136 return chargeSquare;
137}

◆ GetParticleCharge()

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

Reimplemented from G4VEmModel.

Definition at line 142 of file G4LindhardSorensenIonModel.cc.

145{
146 return corr->GetParticleCharge(p,mat,kinEnergy);
147}
G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, const G4double kineticEnergy)

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 98 of file G4LindhardSorensenIonModel.cc.

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}
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4VEmAngularDistribution * GetAngularDistribution()
void SetDeexcitationFlag(G4bool val)
void SetAngularDistribution(G4VEmAngularDistribution *)
G4bool UseAngularGeneratorFlag() const
G4ParticleChangeForLoss * GetParticleChangeForLoss()

◆ MaxSecondaryEnergy()

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

Reimplemented from G4VEmModel.

Definition at line 440 of file G4LindhardSorensenIonModel.cc.

442{
443 // here particle type is checked for any method
444 SetParticle(pd);
445 G4double tau = kinEnergy/mass;
446 return 2.0*CLHEP::electron_mass_c2*tau*(tau + 2.) /
447 (1. + 2.0*(tau + 1.)*eRatio + eRatio*eRatio);
448}

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

◆ MinEnergyCut()

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

Reimplemented from G4VEmModel.

Definition at line 174 of file G4LindhardSorensenIonModel.cc.

176{
178}
G4double GetMeanExcitationEnergy() const
G4IonisParamMat * GetIonisation() const

◆ operator=()

G4LindhardSorensenIonModel & G4LindhardSorensenIonModel::operator= ( const G4LindhardSorensenIonModel & right)
delete

◆ SampleSecondaries()

void G4LindhardSorensenIonModel::SampleSecondaries ( std::vector< G4DynamicParticle * > * vdp,
const G4MaterialCutsCouple * couple,
const G4DynamicParticle * dp,
G4double tmin,
G4double maxEnergy )
overridevirtual

Implements G4VEmModel.

Definition at line 322 of file G4LindhardSorensenIonModel.cc.

328{
329 G4double kineticEnergy = dp->GetKineticEnergy();
330 // take into account formfactor
331 const G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(), kineticEnergy);
332 const G4double minKinEnergy = std::min(cut, tmax);
333 const G4double maxKinEnergy = std::min(maxEnergy, tmax);
334 if(minKinEnergy >= maxKinEnergy) { return; }
335
336 //G4cout << "G4LindhardSorensenIonModel::SampleSecondaries Emin= "
337 // << minKinEnergy << " Emax= " << maxKinEnergy << G4endl;
338
339 G4double totEnergy = kineticEnergy + mass;
340 G4double etot2 = totEnergy*totEnergy;
341 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
342
343 G4double deltaKinEnergy, f;
344 G4double f1 = 0.0;
345 G4double fmax = 1.0;
346 if( 0.0 < spin ) { fmax += 0.5*maxKinEnergy*maxKinEnergy/etot2; }
347
348 CLHEP::HepRandomEngine* rndmEngineMod = G4Random::getTheEngine();
349 G4double rndm[2];
350
351 // sampling without nuclear size effect
352 do {
353 rndmEngineMod->flatArray(2, rndm);
354 deltaKinEnergy = minKinEnergy*maxKinEnergy
355 /(minKinEnergy*(1.0 - rndm[0]) + maxKinEnergy*rndm[0]);
356
357 f = 1.0 - beta2*deltaKinEnergy/tmax;
358 if( 0.0 < spin ) {
359 f1 = 0.5*deltaKinEnergy*deltaKinEnergy/etot2;
360 f += f1;
361 }
362
363 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
364 } while( fmax*rndm[1] > f);
365
366 // projectile formfactor - suppresion of high energy
367 // delta-electron production at high energy
368
369 G4double x = formfact*deltaKinEnergy;
370 if(x > 1.e-6) {
371
372 G4double x1 = 1.0 + x;
373 G4double grej = 1.0/(x1*x1);
374 if( 0.0 < spin ) {
375 G4double x2 = 0.5*electron_mass_c2*deltaKinEnergy/(mass*mass);
376 grej *= (1.0 + magMoment2*(x2 - f1/f)/(1.0 + x2));
377 }
378 if(grej > 1.1) {
379 G4cout << "### G4LindhardSorensenIonModel WARNING: grej= " << grej
380 << " " << dp->GetDefinition()->GetParticleName()
381 << " Ekin(MeV)= " << kineticEnergy
382 << " delEkin(MeV)= " << deltaKinEnergy
383 << G4endl;
384 }
385 if(rndmEngineMod->flat() > grej) { return; }
386 }
387
388 G4ThreeVector deltaDirection;
389
391
392 const G4Material* mat = couple->GetMaterial();
394
395 deltaDirection =
396 GetAngularDistribution()->SampleDirection(dp, deltaKinEnergy, Z, mat);
397
398 } else {
399
400 G4double deltaMomentum =
401 std::sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
402 G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
403 (deltaMomentum * dp->GetTotalMomentum());
404 cost = std::min(cost, 1.0);
405 G4double sint = std::sqrt((1.0 - cost)*(1.0 + cost));
406
407 G4double phi = CLHEP::twopi*rndmEngineMod->flat();
408
409 deltaDirection.set(sint*std::cos(phi),sint*std::sin(phi), cost) ;
410 deltaDirection.rotateUz(dp->GetMomentumDirection());
411 }
412 /*
413 G4cout << "### G4LindhardSorensenIonModel "
414 << dp->GetDefinition()->GetParticleName()
415 << " Ekin(MeV)= " << kineticEnergy
416 << " delEkin(MeV)= " << deltaKinEnergy
417 << " tmin(MeV)= " << minKinEnergy
418 << " tmax(MeV)= " << maxKinEnergy
419 << " dir= " << dp->GetMomentumDirection()
420 << " dirDelta= " << deltaDirection
421 << G4endl;
422 */
423 // create G4DynamicParticle object for delta ray
424 auto delta = new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
425
426 vdp->push_back(delta);
427
428 // Change kinematics of primary particle
429 kineticEnergy -= deltaKinEnergy;
430 G4ThreeVector finalP = dp->GetMomentum() - delta->GetMomentum();
431 finalP = finalP.unit();
432
433 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
434 fParticleChange->SetProposedMomentumDirection(finalP);
435}
#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
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
const G4String & GetParticleName() const
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4int SelectRandomAtomNumber(const G4Material *) const

◆ SetChargeSquareRatio()

void G4LindhardSorensenIonModel::SetChargeSquareRatio ( G4double val)
inlineprotected

Definition at line 184 of file G4LindhardSorensenIonModel.hh.

185{
186 chargeSquare = val;
187}

The documentation for this class was generated from the following files: