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

#include <G4eBremsstrahlungRelModel.hh>

+ Inheritance diagram for G4eBremsstrahlungRelModel:

Public Member Functions

 G4eBremsstrahlungRelModel (const G4ParticleDefinition *p=0, const G4String &nam="eBremLPM")
 
virtual ~G4eBremsstrahlungRelModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double tkin, G4double Z, G4double, G4double cutEnergy, G4double maxEnergy=DBL_MAX)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double)
 
void SetLPMconstant (G4double val)
 
G4double LPMconstant () const
 
- 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 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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., 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 *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
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)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
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
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void ForceBuildTable (G4bool val)
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Member Functions

virtual G4double ComputeDXSectionPerAtom (G4double gammaEnergy)
 
void SetCurrentElement (const G4double)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Protected Attributes

G4NistManagernist
 
const G4ParticleDefinitionparticle
 
G4ParticleDefinitiontheGamma
 
G4ParticleChangeForLossfParticleChange
 
G4double bremFactor
 
G4double particleMass
 
G4double kinEnergy
 
G4double totalEnergy
 
G4double currentZ
 
G4double densityFactor
 
G4double densityCorr
 
G4bool isElectron
 
- Protected Attributes inherited from G4VEmModel
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

Detailed Description

Definition at line 60 of file G4eBremsstrahlungRelModel.hh.

Constructor & Destructor Documentation

◆ G4eBremsstrahlungRelModel()

G4eBremsstrahlungRelModel::G4eBremsstrahlungRelModel ( const G4ParticleDefinition p = 0,
const G4String nam = "eBremLPM" 
)

Definition at line 83 of file G4eBremsstrahlungRelModel.cc.

85 : G4VEmModel(nam),
86 particle(0),
87 bremFactor(fine_structure_const*classic_electr_radius*classic_electr_radius*16./3.),
88 isElectron(true),
89 fMigdalConstant(classic_electr_radius*electron_Compton_length*electron_Compton_length*4.0*pi),
90 fLPMconstant(fine_structure_const*electron_mass_c2*electron_mass_c2/(4.*pi*hbarc)*0.5),
91 fXiLPM(0), fPhiLPM(0), fGLPM(0),
92 use_completescreening(false),isInitialised(false)
93{
96
97 lowKinEnergy = 0.1*GeV;
98 SetLowEnergyLimit(lowKinEnergy);
99
101
102 SetLPMFlag(true);
103 //SetAngularDistribution(new G4ModifiedTsai());
105
106 particleMass = kinEnergy = totalEnergy = currentZ = z13 = z23 = lnZ = Fel
107 = Finel = fCoulomb = fMax = densityFactor = densityCorr = lpmEnergy
108 = xiLPM = phiLPM = gLPM = klpm = kp = 0.0;
109
110 energyThresholdLPM = 1.e39;
111
112 InitialiseConstants();
113 if(p) { SetParticle(p); }
114}
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4NistManager * Instance()
void SetLPMFlag(G4bool val)
Definition: G4VEmModel.hh:634
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:515
const G4ParticleDefinition * particle
G4ParticleChangeForLoss * fParticleChange

◆ ~G4eBremsstrahlungRelModel()

G4eBremsstrahlungRelModel::~G4eBremsstrahlungRelModel ( )
virtual

Definition at line 129 of file G4eBremsstrahlungRelModel.cc.

130{
131}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4eBremsstrahlungRelModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition p,
G4double  tkin,
G4double  Z,
G4double  ,
G4double  cutEnergy,
G4double  maxEnergy = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 260 of file G4eBremsstrahlungRelModel.cc.

266{
267 if(!particle) { SetParticle(p); }
268 if(kineticEnergy < lowKinEnergy) { return 0.0; }
269 G4double cut = std::min(cutEnergy, kineticEnergy);
270 G4double tmax = std::min(maxEnergy, kineticEnergy);
271
272 if(cut >= tmax) { return 0.0; }
273
275
276 G4double cross = ComputeXSectionPerAtom(cut);
277
278 // allow partial integration
279 if(tmax < kinEnergy) { cross -= ComputeXSectionPerAtom(tmax); }
280
281 cross *= Z*Z*bremFactor;
282
283 return cross;
284}
double G4double
Definition: G4Types.hh:64

◆ ComputeDEDXPerVolume()

G4double G4eBremsstrahlungRelModel::ComputeDEDXPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 190 of file G4eBremsstrahlungRelModel.cc.

195{
196 if(!particle) { SetParticle(p); }
197 if(kineticEnergy < lowKinEnergy) { return 0.0; }
198 G4double cut = std::min(cutEnergy, kineticEnergy);
199 if(cut == 0.0) { return 0.0; }
200
201 SetupForMaterial(particle, material,kineticEnergy);
202
203 const G4ElementVector* theElementVector = material->GetElementVector();
204 const G4double* theAtomicNumDensityVector = material->GetAtomicNumDensityVector();
205
206 G4double dedx = 0.0;
207
208 // loop for elements in the material
209 for (size_t i=0; i<material->GetNumberOfElements(); i++) {
210
211 G4VEmModel::SetCurrentElement((*theElementVector)[i]);
212 SetCurrentElement((*theElementVector)[i]->GetZ());
213
214 dedx += theAtomicNumDensityVector[i]*currentZ*currentZ*ComputeBremLoss(cut);
215 }
216 dedx *= bremFactor;
217
218
219 return dedx;
220}
std::vector< G4Element * > G4ElementVector
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:189
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:215
void SetCurrentElement(const G4Element *)
Definition: G4VEmModel.hh:384
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double)

◆ ComputeDXSectionPerAtom()

G4double G4eBremsstrahlungRelModel::ComputeDXSectionPerAtom ( G4double  gammaEnergy)
protectedvirtual

Reimplemented in G4SeltzerBergerModel.

Definition at line 428 of file G4eBremsstrahlungRelModel.cc.

433{
434
435 if(gammaEnergy < 0.0) { return 0.0; }
436
437 G4double y = gammaEnergy/totalEnergy;
438
439 G4double main=0.,secondTerm=0.;
440
441 if (use_completescreening|| currentZ<5) {
442 // ** form factors complete screening case **
443 main = (3./4.*y*y - y + 1.) * ( (Fel-fCoulomb) + Finel/currentZ );
444 secondTerm = (1.-y)/12.*(1.+1./currentZ);
445 }
446 else {
447 // ** intermediate screening using Thomas-Fermi FF from Tsai only valid for Z>=5**
448 G4double dd=100.*electron_mass_c2*y/(totalEnergy-gammaEnergy);
449 G4double gg=dd/z13;
450 G4double eps=dd/z23;
451 G4double phi1=Phi1(gg,currentZ), phi1m2=Phi1M2(gg,currentZ);
452 G4double psi1=Psi1(eps,currentZ), psi1m2=Psi1M2(eps,currentZ);
453
454 main = (3./4.*y*y - y + 1.) * ( (0.25*phi1-1./3.*lnZ-fCoulomb) + (0.25*psi1-2./3.*lnZ)/currentZ );
455 secondTerm = (1.-y)/8.*(phi1m2+psi1m2/currentZ);
456 }
457 G4double cross = main+secondTerm;
458 return cross;
459}

Referenced by SampleSecondaries().

◆ Initialise()

void G4eBremsstrahlungRelModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector cuts 
)
virtual

Implements G4VEmModel.

Reimplemented in G4SeltzerBergerModel.

Definition at line 172 of file G4eBremsstrahlungRelModel.cc.

174{
175 if(p) { SetParticle(p); }
176
177 lowKinEnergy = LowEnergyLimit();
178
179 currentZ = 0.;
180
182
183 if(isInitialised) { return; }
185 isInitialised = true;
186}
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:123
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:95

Referenced by G4SeltzerBergerModel::Initialise().

◆ LPMconstant()

G4double G4eBremsstrahlungRelModel::LPMconstant ( ) const
inline

Definition at line 250 of file G4eBremsstrahlungRelModel.hh.

251{
252 return fLPMconstant;
253}

◆ SampleSecondaries()

void G4eBremsstrahlungRelModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  vdp,
const G4MaterialCutsCouple couple,
const G4DynamicParticle dp,
G4double  cutEnergy,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Reimplemented in G4SeltzerBergerModel.

Definition at line 463 of file G4eBremsstrahlungRelModel.cc.

469{
470 G4double kineticEnergy = dp->GetKineticEnergy();
471 if(kineticEnergy < lowKinEnergy) { return; }
472 G4double cut = std::min(cutEnergy, kineticEnergy);
473 G4double emax = std::min(maxEnergy, kineticEnergy);
474 if(cut >= emax) { return; }
475
476 SetupForMaterial(particle, couple->GetMaterial(), kineticEnergy);
477
478 const G4Element* elm =
479 SelectRandomAtom(couple,particle,kineticEnergy,cut,emax);
480 SetCurrentElement(elm->GetZ());
481
482 kinEnergy = kineticEnergy;
483 totalEnergy = kineticEnergy + particleMass;
485
486 //G4double fmax= fMax;
487 G4bool highe = true;
488 if(totalEnergy < energyThresholdLPM) { highe = false; }
489
490 G4double xmin = log(cut*cut + densityCorr);
491 G4double xmax = log(emax*emax + densityCorr);
492 G4double gammaEnergy, f, x;
493
494 do {
495 x = exp(xmin + G4UniformRand()*(xmax - xmin)) - densityCorr;
496 if(x < 0.0) { x = 0.0; }
497 gammaEnergy = sqrt(x);
498 if(highe) { f = ComputeRelDXSectionPerAtom(gammaEnergy); }
499 else { f = ComputeDXSectionPerAtom(gammaEnergy); }
500
501 if ( f > fMax ) {
502 G4cout << "### G4eBremsstrahlungRelModel Warning: Majoranta exceeded! "
503 << f << " > " << fMax
504 << " Egamma(MeV)= " << gammaEnergy
505 << " Ee(MeV)= " << kineticEnergy
506 << " " << GetName()
507 << G4endl;
508 }
509
510 } while (f < fMax*G4UniformRand());
511
512 //
513 // angles of the emitted gamma. ( Z - axis along the parent particle)
514 // use general interface
515 //
516
517 G4ThreeVector gammaDirection =
520 couple->GetMaterial());
521
522 // create G4DynamicParticle object for the Gamma
523 G4DynamicParticle* gamma = new G4DynamicParticle(theGamma,gammaDirection,
524 gammaEnergy);
525 vdp->push_back(gamma);
526
527 G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2));
528 G4ThreeVector direction = (totMomentum*dp->GetMomentumDirection()
529 - gammaEnergy*gammaDirection).unit();
530
531 // energy of primary
532 G4double finalE = kineticEnergy - gammaEnergy;
533
534 // stop tracking and create new secondary instead of primary
535 if(gammaEnergy > SecondaryThreshold()) {
538 G4DynamicParticle* el =
540 direction, finalE);
541 vdp->push_back(el);
542
543 // continue tracking
544 } else {
547 }
548}
@ fStopAndKill
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
const G4Material * GetMaterial() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:508
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:459
const G4String & GetName() const
Definition: G4VEmModel.hh:655
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:557
void ProposeTrackStatus(G4TrackStatus status)
virtual G4double ComputeDXSectionPerAtom(G4double gammaEnergy)
int G4lrint(double ad)
Definition: templates.hh:163

◆ SetCurrentElement()

void G4eBremsstrahlungRelModel::SetCurrentElement ( const G4double  Z)
inlineprotected

Definition at line 185 of file G4eBremsstrahlungRelModel.hh.

186{
187 if(Z != currentZ) {
188 currentZ = Z;
189
190 G4int iz = G4int(Z);
191 z13 = nist->GetZ13(iz);
192 z23 = z13*z13;
193 lnZ = nist->GetLOGZ(iz);
194
195 if (iz <= 4) {
196 Fel = Fel_light[iz];
197 Finel = Finel_light[iz] ;
198 }
199 else {
200 Fel = facFel - lnZ/3. ;
201 Finel = facFinel - 2.*lnZ/3. ;
202 }
203
204 fCoulomb = GetCurrentElement()->GetfCoulomb();
205 fMax = Fel-fCoulomb + Finel/currentZ + (1.+1./currentZ)/12.;
206 }
207}
int G4int
Definition: G4Types.hh:66
G4double GetfCoulomb() const
Definition: G4Element.hh:201
G4double GetZ13(G4double Z)
G4double GetLOGZ(G4int Z)
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:391

Referenced by ComputeCrossSectionPerAtom(), ComputeDEDXPerVolume(), SampleSecondaries(), and G4SeltzerBergerModel::SampleSecondaries().

◆ SetLPMconstant()

void G4eBremsstrahlungRelModel::SetLPMconstant ( G4double  val)
inline

Definition at line 242 of file G4eBremsstrahlungRelModel.hh.

243{
244 fLPMconstant = val;
245}

◆ SetupForMaterial()

void G4eBremsstrahlungRelModel::SetupForMaterial ( const G4ParticleDefinition ,
const G4Material mat,
G4double  kineticEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 145 of file G4eBremsstrahlungRelModel.cc.

148{
149 densityFactor = mat->GetElectronDensity()*fMigdalConstant;
150 lpmEnergy = mat->GetRadlen()*fLPMconstant;
151
152 // Threshold for LPM effect (i.e. below which LPM hidden by density effect)
153 if (LPMFlag()) {
154 energyThresholdLPM=sqrt(densityFactor)*lpmEnergy;
155 } else {
156 energyThresholdLPM=1.e39; // i.e. do not use LPM effect
157 }
158 // calculate threshold for density effect
159 kinEnergy = kineticEnergy;
160 totalEnergy = kineticEnergy + particleMass;
162
163 // define critical gamma energies (important for integration/dicing)
164 klpm=totalEnergy*totalEnergy/lpmEnergy;
165 kp=sqrt(densityCorr);
166
167}
G4double GetElectronDensity() const
Definition: G4Material.hh:216
G4double GetRadlen() const
Definition: G4Material.hh:219
G4bool LPMFlag() const
Definition: G4VEmModel.hh:564

Referenced by ComputeDEDXPerVolume(), SampleSecondaries(), and G4SeltzerBergerModel::SampleSecondaries().

Member Data Documentation

◆ bremFactor

G4double G4eBremsstrahlungRelModel::bremFactor
protected

◆ currentZ

◆ densityCorr

G4double G4eBremsstrahlungRelModel::densityCorr
protected

◆ densityFactor

G4double G4eBremsstrahlungRelModel::densityFactor
protected

◆ fParticleChange

G4ParticleChangeForLoss* G4eBremsstrahlungRelModel::fParticleChange
protected

◆ isElectron

G4bool G4eBremsstrahlungRelModel::isElectron
protected

◆ kinEnergy

◆ nist

G4NistManager* G4eBremsstrahlungRelModel::nist
protected

Definition at line 127 of file G4eBremsstrahlungRelModel.hh.

Referenced by G4eBremsstrahlungRelModel(), and SetCurrentElement().

◆ particle

◆ particleMass

◆ theGamma

G4ParticleDefinition* G4eBremsstrahlungRelModel::theGamma
protected

◆ totalEnergy


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