Geant4 10.7.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 &) override
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double ekin, G4double cutEnergy) override
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double ekin, G4double zet, G4double, G4double cutEnergy, G4double maxEnergy=DBL_MAX) override
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy) override
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double) override
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cutEnergy) 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 ComputeDXSectionPerAtom (G4double gammaEnergy)
 
void SetParticle (const G4ParticleDefinition *p)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Protected Attributes

G4bool fIsElectron
 
G4bool fIsScatOffElectron
 
G4bool fIsLPMActive
 
G4int fCurrentIZ
 
G4double fPrimaryParticleMass
 
G4double fPrimaryKinEnergy
 
G4double fPrimaryTotalEnergy
 
G4double fDensityFactor
 
G4double fDensityCorr
 
G4double fLowestKinEnergy
 
G4double fNucTerm
 
G4double fSumTerm
 
const G4ParticleDefinitionfPrimaryParticle
 
G4ParticleDefinitionfGammaParticle
 
G4ParticleChangeForLossfParticleChange
 
- 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
 

Static Protected Attributes

static const G4double gBremFactor
 
static const G4double gMigdalConstant
 

Detailed Description

Definition at line 59 of file G4eBremsstrahlungRelModel.hh.

Constructor & Destructor Documentation

◆ G4eBremsstrahlungRelModel()

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

Definition at line 116 of file G4eBremsstrahlungRelModel.cc.

118: G4VEmModel(nam), fIsElectron(true), fIsScatOffElectron(false),
119 fIsLPMActive(false), fPrimaryParticle(nullptr), fIsUseCompleteScreening(false)
120{
121 fCurrentIZ = 0;
122 //
126 fDensityFactor = 0.;
127 fDensityCorr = 0.;
128 fNucTerm = 0.;
129 fSumTerm = 0.;
130 //
131 fPrimaryParticle = nullptr;
133 fParticleChange = nullptr;
134 //
135 fLowestKinEnergy = 1.0*MeV;
137 //
138 fLPMEnergyThreshold = 1.e+39;
139 fLPMEnergy = 0.;
140
141 SetLPMFlag(true);
142 //
144 //SetAngularDistribution(new G4DipBustGenerator());
145 //
146 if (p) {
147 SetParticle(p);
148 }
149}
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
void SetLPMFlag(G4bool val)
Definition: G4VEmModel.hh:806
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:618
void SetParticle(const G4ParticleDefinition *p)
const G4ParticleDefinition * fPrimaryParticle
G4ParticleDefinition * fGammaParticle
G4ParticleChangeForLoss * fParticleChange

◆ ~G4eBremsstrahlungRelModel()

G4eBremsstrahlungRelModel::~G4eBremsstrahlungRelModel ( )
virtual

Definition at line 151 of file G4eBremsstrahlungRelModel.cc.

152{
153 if (IsMaster()) {
154 // clear ElementData container
155 for (size_t iz = 0; iz < gElementData.size(); ++iz) {
156 if (gElementData[iz]) {
157 delete gElementData[iz];
158 }
159 }
160 gElementData.clear();
161 // clear LPMFunctions (if any)
162 if (LPMFlag()) {
163 gLPMFuncs.fLPMFuncG.clear();
164 gLPMFuncs.fLPMFuncPhi.clear();
165 gLPMFuncs.fIsInitialized = false;
166 }
167 }
168}
G4bool IsMaster() const
Definition: G4VEmModel.hh:736
G4bool LPMFlag() const
Definition: G4VEmModel.hh:687

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4eBremsstrahlungRelModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition p,
G4double  ekin,
G4double  zet,
G4double  ,
G4double  cutEnergy,
G4double  maxEnergy = DBL_MAX 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 325 of file G4eBremsstrahlungRelModel.cc.

332{
333 G4double crossSection = 0.0;
334 if (!fPrimaryParticle) {
335 SetParticle(p);
336 }
337 if (kineticEnergy < LowEnergyLimit()) {
338 return crossSection;
339 }
340 // min/max kinetic energy limits of the DCS integration:
341 const G4double tmin = std::min(cut, kineticEnergy);
342 const G4double tmax = std::min(maxEnergy, kineticEnergy);
343 // zero restricted x-section if e- kinetic energy is below gamma cut
344 if (tmin >= tmax) {
345 return crossSection;
346 }
347 fCurrentIZ = std::min(G4lrint(Z), gMaxZet);
348 // integrate numerically (dependent part of) the DCS between the kin. limits:
349 // a. integrate between tmin and kineticEnergy of the e-
350 crossSection = ComputeXSectionPerAtom(tmin);
351 // allow partial integration: only if maxEnergy < kineticEnergy
352 // b. integrate between tmax and kineticEnergy (tmax=maxEnergy in this case)
353 // (so the result in this case is the integral of DCS between tmin and
354 // maxEnergy)
355 if (tmax < kineticEnergy) {
356 crossSection -= ComputeXSectionPerAtom(tmax);
357 }
358 // multiply with the constant factors: 16\alpha r_0^2/3 Z^2
359 crossSection *= Z*Z*gBremFactor;
360 return std::max(crossSection, 0.);
361}
double G4double
Definition: G4Types.hh:83
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
int G4lrint(double ad)
Definition: templates.hh:134

◆ ComputeDEDXPerVolume()

G4double G4eBremsstrahlungRelModel::ComputeDEDXPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  ekin,
G4double  cutEnergy 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 241 of file G4eBremsstrahlungRelModel.cc.

245{
246 G4double dedx = 0.0;
247 if (!fPrimaryParticle) {
248 SetParticle(p);
249 }
250 if (kineticEnergy < LowEnergyLimit()) {
251 return dedx;
252 }
253 // maximum value of the dE/dx integral (the minimum is 0 of course)
254 G4double tmax = std::min(cutEnergy, kineticEnergy);
255 if (tmax == 0.0) {
256 return dedx;
257 }
258 // sets kinematical and material related variables
259 SetupForMaterial(fPrimaryParticle, material,kineticEnergy);
260 // get element compositions of the material
261 const G4ElementVector* theElemVector = material->GetElementVector();
262 const G4double* theAtomNumDensVector = material->GetAtomicNumDensityVector();
263 const size_t numberOfElements = theElemVector->size();
264 // loop over the elements of the material and compute their contributions to
265 // the restricted dE/dx by numerical integration of the dependent part of DCS
266 for (size_t ie = 0; ie < numberOfElements; ++ie) {
267 G4VEmModel::SetCurrentElement((*theElemVector)[ie]);
268 //SetCurrentElement((*theElementVector)[i]->GetZasInt());
269 const G4double zet = (*theElemVector)[ie]->GetZ();
270 fCurrentIZ = std::min(G4lrint(zet), gMaxZet);
271 dedx += theAtomNumDensVector[ie]*zet*zet*ComputeBremLoss(tmax);
272 }
273 // apply the constant factor C/Z = 16\alpha r_0^2/3
274 dedx *= gBremFactor;
275 return std::max(dedx,0.);
276}
std::vector< G4Element * > G4ElementVector
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:214
void SetCurrentElement(const G4Element *)
Definition: G4VEmModel.hh:487
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double) override

◆ ComputeDXSectionPerAtom()

G4double G4eBremsstrahlungRelModel::ComputeDXSectionPerAtom ( G4double  gammaEnergy)
protectedvirtual

Reimplemented in G4LivermoreBremsstrahlungModel, and G4SeltzerBergerModel.

Definition at line 481 of file G4eBremsstrahlungRelModel.cc.

482{
483 G4double dxsec = 0.0;
484 if (gammaEnergy < 0.) {
485 return dxsec;
486 }
487 const G4double y = gammaEnergy/fPrimaryTotalEnergy;
488 const G4double onemy = 1.-y;
489 const G4double dum0 = onemy+0.75*y*y;
490 const ElementData* elDat = gElementData[fCurrentIZ];
491 // use complete screening and L_el, L_inel from Dirac-Fock model instead of TF
492 if (fCurrentIZ < 5 || fIsUseCompleteScreening) {
493 dxsec = dum0*elDat->fZFactor1;
494 dxsec += onemy*elDat->fZFactor2;
495 if (fIsScatOffElectron) {
496 fSumTerm = dxsec;
497 fNucTerm = dum0*elDat->fZFactor11+onemy/12.;
498 }
499 } else {
500 // use Tsai's analytical approx. (Tsai Eqs. [3.38-3.41]) to the 'universal'
501 // numerical screening functions computed by using the TF model of the atom
502 const G4double invZ = 1./(G4double)fCurrentIZ;
503 const G4double Fz = elDat->fFz;
504 const G4double logZ = elDat->fLogZ;
505 const G4double dum1 = y/(fPrimaryTotalEnergy-gammaEnergy);
506 const G4double gamma = dum1*elDat->fGammaFactor;
507 const G4double epsilon = dum1*elDat->fEpsilonFactor;
508 // evaluate the screening functions
509 G4double phi1, phi1m2, psi1, psi1m2;
510 ComputeScreeningFunctions(phi1, phi1m2, psi1, psi1m2, gamma, epsilon);
511 dxsec = dum0*((0.25*phi1-Fz) + (0.25*psi1-2.*logZ/3.)*invZ);
512 dxsec += 0.125*onemy*(phi1m2 + psi1m2*invZ);
513 if (fIsScatOffElectron) {
514 fSumTerm = dxsec;
515 fNucTerm = dum0*(0.25*phi1-Fz) + 0.125*onemy*phi1m2;
516 }
517 }
518 return std::max(dxsec,0.0);
519}
double epsilon(double density, double temperature)

Referenced by SampleSecondaries().

◆ Initialise()

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

Implements G4VEmModel.

Reimplemented in G4LivermoreBremsstrahlungModel, G4ePolarizedBremsstrahlungModel, and G4SeltzerBergerModel.

Definition at line 170 of file G4eBremsstrahlungRelModel.cc.

172{
173 if (p) {
174 SetParticle(p);
175 }
176 fCurrentIZ = 0;
177 // init element data and precompute LPM functions (only if lpmflag is true)
178 if (IsMaster()) {
179 InitialiseElementData();
180 if (LPMFlag()) { InitLPMFunctions(); }
183 }
184 }
186 if (GetTripletModel()) {
187 GetTripletModel()->Initialise(p, cuts);
188 fIsScatOffElectron = true;
189 }
190}
G4VEmModel * GetTripletModel()
Definition: G4VEmModel.hh:628
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)=0
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:118

Referenced by G4LivermoreBremsstrahlungModel::Initialise().

◆ InitialiseLocal()

void G4eBremsstrahlungRelModel::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 192 of file G4eBremsstrahlungRelModel.cc.

194{
197 }
198}
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:842
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:834

◆ MinPrimaryEnergy()

G4double G4eBremsstrahlungRelModel::MinPrimaryEnergy ( const G4Material ,
const G4ParticleDefinition ,
G4double  cutEnergy 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 231 of file G4eBremsstrahlungRelModel.cc.

234{
235 return std::max(fLowestKinEnergy, cut);
236}

◆ SampleSecondaries()

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

Implements G4VEmModel.

Reimplemented in G4LivermoreBremsstrahlungModel, G4SeltzerBergerModel, and G4ePolarizedBremsstrahlungModel.

Definition at line 545 of file G4eBremsstrahlungRelModel.cc.

550{
551 const G4double kineticEnergy = dp->GetKineticEnergy();
552// const G4double logKineticEnergy = dp->GetLogKineticEnergy();
553 if (kineticEnergy < LowEnergyLimit()) {
554 return;
555 }
556 // min, max kinetic energy limits
557 const G4double tmin = std::min(cutEnergy, kineticEnergy);
558 const G4double tmax = std::min(maxEnergy, kineticEnergy);
559 if (tmin >= tmax) {
560 return;
561 }
562 //
563 SetupForMaterial(fPrimaryParticle, couple->GetMaterial(), kineticEnergy);
564 const G4Element* elm = SelectTargetAtom(couple,fPrimaryParticle,kineticEnergy,
565 dp->GetLogKineticEnergy(),tmin,tmax);
566 //
567 fCurrentIZ = elm->GetZasInt();
568 const ElementData* elDat = gElementData[fCurrentIZ];
569 const G4double funcMax = elDat->fZFactor1+elDat->fZFactor2;
570 // get the random engine
571 G4double rndm[2];
572 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
573 // min max of the transformed variable: x(k) = ln(k^2+k_p^2) that is in [ln(k_c^2+k_p^2), ln(E_k^2+k_p^2)]
574 const G4double xmin = G4Log(tmin*tmin+fDensityCorr);
575 const G4double xrange = G4Log(tmax*tmax+fDensityCorr)-xmin;
576 G4double gammaEnergy, funcVal;
577 do {
578 rndmEngine->flatArray(2, rndm);
579 gammaEnergy = std::sqrt(std::max(G4Exp(xmin+rndm[0]*xrange)-fDensityCorr, 0.0));
580 funcVal = fIsLPMActive
581 ? ComputeRelDXSectionPerAtom(gammaEnergy)
582 : ComputeDXSectionPerAtom(gammaEnergy);
583 // cross-check of proper function maximum in the rejection
584// if (funcVal > funcMax) {
585// G4cout << "### G4eBremsstrahlungRelModel Warning: Majoranta exceeded! "
586// << funcVal << " > " << funcMax
587// << " Egamma(MeV)= " << gammaEnergy
588// << " Ee(MeV)= " << kineticEnergy
589// << " " << GetName()
590// << G4endl;
591// }
592 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
593 } while (funcVal < funcMax*rndm[1]);
594 //
595 // scattering off nucleus or off e- by triplet model
596 if (fIsScatOffElectron && rndmEngine->flat()*fSumTerm>fNucTerm) {
597 GetTripletModel()->SampleSecondaries(vdp, couple, dp, cutEnergy, maxEnergy);
598 return;
599 }
600 //
601 // angles of the emitted gamma. ( Z - axis along the parent particle)
602 // use general interface
603 G4ThreeVector gamDir =
605 fCurrentIZ, couple->GetMaterial());
606 // create G4DynamicParticle object for the Gamma
608 gammaEnergy);
609 vdp->push_back(gamma);
610 // compute post-interaction kinematics of primary e-/e+ based on
611 // energy-momentum conservation
612 const G4double totMomentum = std::sqrt(kineticEnergy*(
613 fPrimaryTotalEnergy + CLHEP::electron_mass_c2));
614 G4ThreeVector dir =
615 (totMomentum*dp->GetMomentumDirection()-gammaEnergy*gamDir).unit();
616 const G4double finalE = kineticEnergy-gammaEnergy;
617 // if secondary gamma energy is higher than threshold(very high by default)
618 // then stop tracking the primary particle and create new secondary e-/e+
619 // instead of the primary one
620 if (gammaEnergy > SecondaryThreshold()) {
624 const_cast<G4ParticleDefinition*>(fPrimaryParticle), dir, finalE);
625 vdp->push_back(el);
626 } else { // continue tracking the primary e-/e+ otherwise
629 }
630}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double G4Log(G4double x)
Definition: G4Log.hh:226
@ fStopAndKill
virtual double flat()=0
virtual void flatArray(const int size, double *vect)=0
const G4ThreeVector & GetMomentumDirection() const
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const
G4int GetZasInt() 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:611
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
const G4Element * SelectTargetAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:587
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:680
void ProposeTrackStatus(G4TrackStatus status)
virtual G4double ComputeDXSectionPerAtom(G4double gammaEnergy)

◆ SetParticle()

void G4eBremsstrahlungRelModel::SetParticle ( const G4ParticleDefinition p)
protected

◆ SetupForMaterial()

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

Reimplemented from G4VEmModel.

Reimplemented in G4SeltzerBergerModel.

Definition at line 210 of file G4eBremsstrahlungRelModel.cc.

213{
215 fLPMEnergy = gLPMconstant*mat->GetRadlen();
216 // threshold for LPM effect (i.e. below which LPM hidden by density effect)
217 if (LPMFlag()) {
218 fLPMEnergyThreshold = std::sqrt(fDensityFactor)*fLPMEnergy;
219 } else {
220 fLPMEnergyThreshold = 1.e+39; // i.e. do not use LPM effect
221 }
222 // calculate threshold for density effect: k_p = sqrt(fDensityCorr)
223 fPrimaryKinEnergy = kineticEnergy;
226 // set activation flag for LPM effects in the DCS
227 fIsLPMActive = (fPrimaryTotalEnergy>fLPMEnergyThreshold);
228}
G4double GetElectronDensity() const
Definition: G4Material.hh:215
G4double GetRadlen() const
Definition: G4Material.hh:218

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

Member Data Documentation

◆ fCurrentIZ

◆ fDensityCorr

◆ fDensityFactor

G4double G4eBremsstrahlungRelModel::fDensityFactor
protected

◆ fGammaParticle

◆ fIsElectron

◆ fIsLPMActive

G4bool G4eBremsstrahlungRelModel::fIsLPMActive
protected

◆ fIsScatOffElectron

G4bool G4eBremsstrahlungRelModel::fIsScatOffElectron
protected

◆ fLowestKinEnergy

G4double G4eBremsstrahlungRelModel::fLowestKinEnergy
protected

◆ fNucTerm

G4double G4eBremsstrahlungRelModel::fNucTerm
protected

◆ fParticleChange

◆ fPrimaryKinEnergy

◆ fPrimaryParticle

◆ fPrimaryParticleMass

◆ fPrimaryTotalEnergy

◆ fSumTerm

G4double G4eBremsstrahlungRelModel::fSumTerm
protected

◆ gBremFactor

const G4double G4eBremsstrahlungRelModel::gBremFactor
staticprotected
Initial value:
= 16. * CLHEP::fine_structure_const * CLHEP::classic_electr_radius
* CLHEP::classic_electr_radius/3.

Definition at line 162 of file G4eBremsstrahlungRelModel.hh.

Referenced by ComputeCrossSectionPerAtom(), ComputeDEDXPerVolume(), G4LivermoreBremsstrahlungModel::ComputeDXSectionPerAtom(), and G4SeltzerBergerModel::ComputeDXSectionPerAtom().

◆ gMigdalConstant

const G4double G4eBremsstrahlungRelModel::gMigdalConstant
staticprotected
Initial value:
= 4. * CLHEP::pi * CLHEP::classic_electr_radius
* CLHEP::electron_Compton_length * CLHEP::electron_Compton_length

Definition at line 163 of file G4eBremsstrahlungRelModel.hh.

Referenced by SetupForMaterial(), and G4SeltzerBergerModel::SetupForMaterial().


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