Geant4 10.7.0
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, G4double &eloss, G4double &, G4double length) 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 Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
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 G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual 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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
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 LPMFlag () 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 SetLPMFlag (G4bool val)
 
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 *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Member Functions

G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kinEnergy) override
 
G4double GetChargeSquareRatio () const
 
void SetChargeSquareRatio (G4double val)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Additional Inherited Members

- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

Detailed Description

Definition at line 63 of file G4LindhardSorensenIonModel.hh.

Constructor & Destructor Documentation

◆ G4LindhardSorensenIonModel() [1/2]

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

Definition at line 66 of file G4LindhardSorensenIonModel.cc.

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}
G4double G4Log(G4double x)
Definition: G4Log.hh:226
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
static G4NistManager * Instance()
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764
#define DBL_MAX
Definition: templates.hh:62

◆ ~G4LindhardSorensenIonModel()

G4LindhardSorensenIonModel::~G4LindhardSorensenIonModel ( )
override

Definition at line 84 of file G4LindhardSorensenIonModel.cc.

85{}

◆ 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 198 of file G4LindhardSorensenIonModel.cc.

204{
205 return Z*ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
206}
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 165 of file G4LindhardSorensenIonModel.cc.

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}
double G4double
Definition: G4Types.hh:83
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override

Referenced by ComputeCrossSectionPerAtom(), and CrossSectionPerVolume().

◆ ComputeDEDXPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 224 of file G4LindhardSorensenIonModel.cc.

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}
G4double BarkasCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double ShellCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double GetMeanExcitationEnergy() const
G4double DensityCorrection(G4double x)
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:224
G4double GetElectronDensity() const
Definition: G4Material.hh:215

◆ CorrectionsAlongStep()

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

Reimplemented from G4VEmModel.

Definition at line 274 of file G4LindhardSorensenIonModel.cc.

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}
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetDeltaL(G4int Z, G4double gamma) const
const G4Material * GetMaterial() const
virtual void SetParticleAndCharge(const G4ParticleDefinition *, G4double q2)
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:604

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 210 of file G4LindhardSorensenIonModel.cc.

216{
217 return material->GetElectronDensity()
218 *ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
219}

◆ GetChargeSquareRatio() [1/2]

G4double G4LindhardSorensenIonModel::GetChargeSquareRatio ( ) const
inlineprotected

Definition at line 188 of file G4LindhardSorensenIonModel.hh.

189{
190 return chargeSquare;
191}

◆ GetChargeSquareRatio() [2/2]

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

Reimplemented from G4VEmModel.

Definition at line 114 of file G4LindhardSorensenIonModel.cc.

117{
118 return chargeSquare;
119}

◆ GetParticleCharge()

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

Reimplemented from G4VEmModel.

Definition at line 124 of file G4LindhardSorensenIonModel.cc.

127{
128 return charge;
129}

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 89 of file G4LindhardSorensenIonModel.cc.

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}
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:611
G4bool IsMaster() const
Definition: G4VEmModel.hh:736
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

◆ MaxSecondaryEnergy()

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

Reimplemented from G4VEmModel.

Definition at line 431 of file G4LindhardSorensenIonModel.cc.

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}

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

◆ MinEnergyCut()

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

Reimplemented from G4VEmModel.

Definition at line 156 of file G4LindhardSorensenIonModel.cc.

158{
160}

◆ 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 311 of file G4LindhardSorensenIonModel.cc.

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}
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
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 *)
Definition: G4VEmModel.cc:315

◆ SetChargeSquareRatio()

void G4LindhardSorensenIonModel::SetChargeSquareRatio ( G4double  val)
inlineprotected

Definition at line 195 of file G4LindhardSorensenIonModel.hh.

196{
197 chargeSquare = val;
198}

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