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

#include <G4AtimaEnergyLossModel.hh>

+ Inheritance diagram for G4AtimaEnergyLossModel:

Public Member Functions

 G4AtimaEnergyLossModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="Atima")
 
virtual ~G4AtimaEnergyLossModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
 
virtual G4double ComputeCrossSectionPerElectron (const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double) override
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &, G4double &, G4double) override
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
- 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

virtual 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 57 of file G4AtimaEnergyLossModel.hh.

Constructor & Destructor Documentation

◆ G4AtimaEnergyLossModel()

G4AtimaEnergyLossModel::G4AtimaEnergyLossModel ( const G4ParticleDefinition p = nullptr,
const G4String nam = "Atima" 
)
explicit

Definition at line 74 of file G4AtimaEnergyLossModel.cc.

76 : G4VEmModel(nam),
77 particle(nullptr),
78 tlimit(DBL_MAX),
79 isIon(false)
80{
81 g4calc = G4Pow::GetInstance();
82 fParticleChange = nullptr;
83 theElectron = G4Electron::Electron();
84 SetParticle(theElectron);
87 SetLowEnergyLimit(2.0*MeV);
88 MLN10 = 2.30258509299;
89 atomic_mass_unit = 931.4940954; // MeV/c^2
90 dedx_constant = 0.3070749187; //4*pi*Na*me*c^2*r_e^2 //MeV cm^2
91 electron_mass = 0.510998928; // MeV/c^2
92 fine_structure = 1.0/137.035999139;
93 domega2dx_constant = dedx_constant*electron_mass; //4*pi*Na*me*c^2*r_e^2 //MeV^2 cm^2
94 if(tableE[0] == 0.0) {
95 G4double logmin = 0.;
96 G4double logmax = 5.;
97 stepE = (logmax-logmin)/199.;
98 for(G4int i=0; i<200; ++i){
99 tableE[i] = G4Exp(MLN10*(logmin + ((G4double)i)*stepE));
100 }
101 }
102}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
static G4NistManager * Instance()
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764
#define DBL_MAX
Definition: templates.hh:62

◆ ~G4AtimaEnergyLossModel()

G4AtimaEnergyLossModel::~G4AtimaEnergyLossModel ( )
virtual

Definition at line 106 of file G4AtimaEnergyLossModel.cc.

107{}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 226 of file G4AtimaEnergyLossModel.cc.

232{
234 (p,kineticEnergy,cutEnergy,maxEnergy);
235 return cross;
236}
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)

◆ ComputeCrossSectionPerElectron()

G4double G4AtimaEnergyLossModel::ComputeCrossSectionPerElectron ( const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
)
virtual

Definition at line 195 of file G4AtimaEnergyLossModel.cc.

199{
200 G4double cross = 0.0;
201 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
202 G4double maxEnergy = std::min(tmax,maxKinEnergy);
203 if(cutEnergy < maxEnergy) {
204
205 G4double totEnergy = kineticEnergy + mass;
206 G4double energy2 = totEnergy*totEnergy;
207 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
208
209 cross = (maxEnergy - cutEnergy)/(cutEnergy*maxEnergy)
210 - beta2*G4Log(maxEnergy/cutEnergy)/tmax;
211
212 // +term for spin=1/2 particle
213 if( 0.0 < spin ) { cross += 0.5*(maxEnergy - cutEnergy)/energy2; }
214
215 cross *= twopi_mc2_rcl2*chargeSquare/beta2;
216 }
217
218 // G4cout << "BB: e= " << kineticEnergy << " tmin= " << cutEnergy
219 // << " tmax= " << tmax << " cross= " << cross << G4endl;
220
221 return cross;
222}
G4double G4Log(G4double x)
Definition: G4Log.hh:226
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override

Referenced by ComputeCrossSectionPerAtom(), and CrossSectionPerVolume().

◆ ComputeDEDXPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 255 of file G4AtimaEnergyLossModel.cc.

259{
260 //Call to ATIMA Stopping Power function
261 G4double zt = material->GetIonisation()->GetZeffective();
262 zt = std::min(zt,93.);
263 G4double at = nist->GetAtomicMassAmu(G4lrint(zt));
264 G4double dedx = StoppingPower(p->GetPDGMass(), p->GetPDGCharge(),
265 kineticEnergy/(MeV), at, zt) *material->GetDensity()/(g/cm3);
266 dedx = std::max(dedx, 0.0);
267
268 // G4cout << "E(MeV)= " << kineticEnergy/MeV << " dedx= " << dedx
269 // << " " << material->GetName() << G4endl;
270 return dedx;
271}
G4double GetZeffective() const
G4double GetDensity() const
Definition: G4Material.hh:178
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:224
G4double GetAtomicMassAmu(const G4String &symb) const
G4double GetPDGCharge() const
int G4lrint(double ad)
Definition: templates.hh:134

Referenced by CorrectionsAlongStep().

◆ CorrectionsAlongStep()

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

Reimplemented from G4VEmModel.

Definition at line 275 of file G4AtimaEnergyLossModel.cc.

280{
281 if(isIon) {
282 const G4ParticleDefinition* p = dp->GetDefinition();
283 const G4Material* mat = couple->GetMaterial();
284 G4double cutEnergy = DBL_MAX;
286 G4double kineticEnergy = dp->GetKineticEnergy();
287 eloss = ComputeDEDXPerVolume(mat, p, kineticEnergy, cutEnergy)*length/(cm);
288 }
289
290}
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double) override
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4Material * GetMaterial() const
virtual void SetParticleAndCharge(const G4ParticleDefinition *, G4double q2)
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:604

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 240 of file G4AtimaEnergyLossModel.cc.

246{
247 G4double eDensity = material->GetElectronDensity();
249 (p,kineticEnergy,cutEnergy,maxEnergy);
250 return cross;
251}
G4double GetElectronDensity() const
Definition: G4Material.hh:215

◆ GetChargeSquareRatio() [1/2]

G4double G4AtimaEnergyLossModel::GetChargeSquareRatio ( ) const
inlineprotected

Definition at line 208 of file G4AtimaEnergyLossModel.hh.

209{
210 return chargeSquare;
211}

◆ GetChargeSquareRatio() [2/2]

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

Reimplemented from G4VEmModel.

Definition at line 134 of file G4AtimaEnergyLossModel.cc.

137{
138 // this method is called only for ions
139 G4double q2 = corr->EffectiveChargeSquareRatio(p,mat,kineticEnergy);
140 corrFactor = q2*corr->EffectiveChargeCorrection(p,mat,kineticEnergy);
141 return corrFactor;
142}
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double EffectiveChargeCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)

◆ GetParticleCharge()

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

Reimplemented from G4VEmModel.

Definition at line 146 of file G4AtimaEnergyLossModel.cc.

149{
150 //G4cout<<"G4AtimaEnergyLossModel::GetParticleCharge e= "<<kineticEnergy <<
151 // " q= " << corr->GetParticleCharge(p,mat,kineticEnergy) <<G4endl;
152 // this method is called only for ions, so no check if it is an ion
153 return corr->GetParticleCharge(p,mat,kineticEnergy);
154}
G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)

◆ Initialise()

void G4AtimaEnergyLossModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector  
)
overridevirtual

Implements G4VEmModel.

Definition at line 111 of file G4AtimaEnergyLossModel.cc.

113{
114 SetGenericIon(p);
115 SetParticle(p);
116
117 //G4cout << "G4AtimaEnergyLossModel::Initialise for " << p->GetParticleName()
118 // << " isIon= " << isIon
119 // << G4endl;
120
121 // always false before the run
122 SetDeexcitationFlag(false);
123
124 if(nullptr == fParticleChange) {
125 fParticleChange = GetParticleChangeForLoss();
128 }
129 }
130}
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:611
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 G4AtimaEnergyLossModel::MaxSecondaryEnergy ( const G4ParticleDefinition pd,
G4double  kinEnergy 
)
overrideprotectedvirtual

Reimplemented from G4VEmModel.

Definition at line 410 of file G4AtimaEnergyLossModel.cc.

412{
413 // here particle type is checked for any method
414 SetParticle(pd);
415 G4double tau = kinEnergy/mass;
416 G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
417 (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
418 return std::min(tmax,tlimit);
419}

Referenced by ComputeCrossSectionPerElectron(), and SampleSecondaries().

◆ MinEnergyCut()

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

Reimplemented from G4VEmModel.

Definition at line 186 of file G4AtimaEnergyLossModel.cc.

188{
190}
G4double GetMeanExcitationEnergy() const

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 294 of file G4AtimaEnergyLossModel.cc.

299{
300 G4double kineticEnergy = dp->GetKineticEnergy();
301 G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(),kineticEnergy);
302
303 G4double maxKinEnergy = std::min(maxEnergy,tmax);
304 if(minKinEnergy >= maxKinEnergy) { return; }
305
306 //G4cout << "G4AtimaEnergyLossModel::SampleSecondaries Emin= " << minKinEnergy
307 // << " Emax= " << maxKinEnergy << G4endl;
308
309 G4double totEnergy = kineticEnergy + mass;
310 G4double etot2 = totEnergy*totEnergy;
311 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
312
313 G4double deltaKinEnergy, f;
314 G4double f1 = 0.0;
315 G4double fmax = 1.0;
316 if( 0.0 < spin ) { fmax += 0.5*maxKinEnergy*maxKinEnergy/etot2; }
317
318 CLHEP::HepRandomEngine* rndmEngineMod = G4Random::getTheEngine();
319 G4double rndm[2];
320
321 // sampling without nuclear size effect
322 do {
323 rndmEngineMod->flatArray(2, rndm);
324 deltaKinEnergy = minKinEnergy*maxKinEnergy
325 /(minKinEnergy*(1.0 - rndm[0]) + maxKinEnergy*rndm[0]);
326
327 f = 1.0 - beta2*deltaKinEnergy/tmax;
328 if( 0.0 < spin ) {
329 f1 = 0.5*deltaKinEnergy*deltaKinEnergy/etot2;
330 f += f1;
331 }
332
333 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
334 } while( fmax*rndm[1] > f);
335
336 // projectile formfactor - suppresion of high energy
337 // delta-electron production at high energy
338
339 G4double x = formfact*deltaKinEnergy*(deltaKinEnergy + 2*electron_mass_c2);
340 if(x > 1.e-6) {
341
342 G4double x10 = 1.0 + x;
343 G4double grej = 1.0/(x10*x10);
344 if( 0.0 < spin ) {
345 G4double x2 = 0.5*electron_mass_c2*deltaKinEnergy/(mass*mass);
346 grej *= (1.0 + magMoment2*(x2 - f1/f)/(1.0 + x2));
347 }
348 if(grej > 1.1) {
349 G4cout << "### G4AtimaEnergyLossModel WARNING: grej= " << grej
350 << " " << dp->GetDefinition()->GetParticleName()
351 << " Ekin(MeV)= " << kineticEnergy
352 << " delEkin(MeV)= " << deltaKinEnergy
353 << G4endl;
354 }
355 if(rndmEngineMod->flat() > grej) { return; }
356 }
357
358 G4ThreeVector deltaDirection;
359
361
362 const G4Material* mat = couple->GetMaterial();
364
365 deltaDirection =
366 GetAngularDistribution()->SampleDirection(dp, deltaKinEnergy, Z, mat);
367
368 } else {
369
370 G4double deltaMomentum =
371 std::sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
372 G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
373 (deltaMomentum * dp->GetTotalMomentum());
374 if(cost > 1.0) { cost = 1.0; }
375 G4double sint = std::sqrt((1.0 - cost)*(1.0 + cost));
376
377 G4double phi = twopi*rndmEngineMod->flat();
378
379 deltaDirection.set(sint*cos(phi),sint*sin(phi), cost) ;
380 deltaDirection.rotateUz(dp->GetMomentumDirection());
381 }
382 /*
383 G4cout << "### G4AtimaEnergyLossModel "
384 << dp->GetDefinition()->GetParticleName()
385 << " Ekin(MeV)= " << kineticEnergy
386 << " delEkin(MeV)= " << deltaKinEnergy
387 << " tmin(MeV)= " << minKinEnergy
388 << " tmax(MeV)= " << maxKinEnergy
389 << " dir= " << dp->GetMomentumDirection()
390 << " dirDelta= " << deltaDirection
391 << G4endl;
392 */
393 // create G4DynamicParticle object for delta ray
394 G4DynamicParticle* delta =
395 new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
396
397 vdp->push_back(delta);
398
399 // Change kinematics of primary particle
400 kineticEnergy -= deltaKinEnergy;
401 G4ThreeVector finalP = dp->GetMomentum() - delta->GetMomentum();
402 finalP = finalP.unit();
403
404 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
405 fParticleChange->SetProposedMomentumDirection(finalP);
406}
#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 G4AtimaEnergyLossModel::SetChargeSquareRatio ( G4double  val)
inlineprotected

Definition at line 215 of file G4AtimaEnergyLossModel.hh.

216{
217 chargeSquare = val;
218}

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