Geant4 11.4.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=nullptr, const G4String &nam="eBremLPM")
 ~G4eBremsstrahlungRelModel () override
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double ekin, G4double cutEnergy) override
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double ekin, G4double zet, G4double, G4double cutEnergy, G4double maxEnergy=DBL_MAX) override
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy) override
void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double) override
G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cutEnergy) override
Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
virtual ~G4VEmModel ()
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
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 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 G4Material *, const G4ParticleDefinition *, const G4double kinEnergy, const G4double cutEnergy, const G4double &length, G4double &eloss)
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
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)
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)
void SetMasterThread (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 = true
G4bool fIsScatOffElectron = false
G4bool fIsLPMActive = false
G4int fCurrentIZ = 0
const G4ParticleDefinitionfPrimaryParticle = nullptr
G4ParticleDefinitionfGammaParticle = nullptr
G4ParticleChangeForLossfParticleChange = nullptr
G4double fPrimaryParticleMass = 0.
G4double fPrimaryKinEnergy = 0.
G4double fPrimaryTotalEnergy = 0.
G4double fDensityFactor = 0.
G4double fDensityCorr = 0.
G4double fLowestKinEnergy
G4double fNucTerm = 0.
G4double fSumTerm = 0.
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

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 = nullptr,
const G4String & nam = "eBremLPM" )
explicit

Definition at line 126 of file G4eBremsstrahlungRelModel.cc.

128: G4VEmModel(nam)
129{
131 //
132 fLowestKinEnergy = 1.0*CLHEP::MeV;
134 //
135 fLPMEnergyThreshold = 1.e+39;
136 fLPMEnergy = 0.;
137 SetAngularDistribution(new G4ModifiedTsai());
138 //
139 if (nullptr != p) {
140 SetParticle(p);
141 }
142}
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
void SetLowEnergyLimit(G4double)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67
void SetAngularDistribution(G4VEmAngularDistribution *)
void SetParticle(const G4ParticleDefinition *p)

Referenced by G4LivermoreBremsstrahlungModel::G4LivermoreBremsstrahlungModel().

◆ ~G4eBremsstrahlungRelModel()

G4eBremsstrahlungRelModel::~G4eBremsstrahlungRelModel ( )
override

Definition at line 144 of file G4eBremsstrahlungRelModel.cc.

145{
146 if (fIsInitializer) {
147 // clear ElementData container
148 for (auto const & ptr : gElementData) { delete ptr; }
149 gElementData.clear();
150 }
151}

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 320 of file G4eBremsstrahlungRelModel.cc.

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

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

◆ ComputeDXSectionPerAtom()

G4double G4eBremsstrahlungRelModel::ComputeDXSectionPerAtom ( G4double gammaEnergy)
protectedvirtual

Reimplemented in G4LivermoreBremsstrahlungModel.

Definition at line 476 of file G4eBremsstrahlungRelModel.cc.

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

Referenced by SampleSecondaries().

◆ Initialise()

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

Implements G4VEmModel.

Reimplemented in G4LivermoreBremsstrahlungModel.

Definition at line 153 of file G4eBremsstrahlungRelModel.cc.

155{
156 // parameters in each thread
157 if (fPrimaryParticle != p) {
158 SetParticle(p);
159 }
160 fUseLPM = G4EmParameters::Instance()->LPM();
161 fCurrentIZ = 0;
162
163 // init static element data and precompute LPM functions only once
164 std::call_once(applyOnce, [this]() { fIsInitializer = true; });
165
166 // for all treads and derived classes
167 if (fIsInitializer || gElementData.empty()) {
168 G4AutoLock l(&theBremRelMutex);
169 if (gElementData.empty()) {
170 gElementData.resize(gMaxZet+1, nullptr);
171 }
172 InitialiseElementData();
173 l.unlock();
174 }
175
176 // element selectors are initialized in the master thread
177 if (IsMaster()) {
179 }
180 // initialisation in all threads
181 if (nullptr == fParticleChange) {
183 }
184 if (GetTripletModel()) {
185 GetTripletModel()->Initialise(p, cuts);
186 fIsScatOffElectron = true;
187 }
188}
G4TemplateAutoLock< G4Mutex > G4AutoLock
static G4EmParameters * Instance()
G4bool LPM() const
G4bool IsMaster() const
G4VEmModel * GetTripletModel()
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)=0
G4ParticleChangeForLoss * GetParticleChangeForLoss()
G4ParticleChangeForLoss * fParticleChange

Referenced by G4LivermoreBremsstrahlungModel::Initialise().

◆ InitialiseLocal()

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

Reimplemented from G4VEmModel.

Definition at line 190 of file G4eBremsstrahlungRelModel.cc.

192{
194}
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
std::vector< G4EmElementSelector * > * GetElementSelectors()

◆ MinPrimaryEnergy()

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

Reimplemented from G4VEmModel.

Definition at line 227 of file G4eBremsstrahlungRelModel.cc.

230{
231 return std::max(fLowestKinEnergy, cut);
232}

◆ SampleSecondaries()

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

Implements G4VEmModel.

Reimplemented in G4LivermoreBremsstrahlungModel.

Definition at line 540 of file G4eBremsstrahlungRelModel.cc.

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

Definition at line 206 of file G4eBremsstrahlungRelModel.cc.

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

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

Member Data Documentation

◆ fCurrentIZ

◆ fDensityCorr

G4double G4eBremsstrahlungRelModel::fDensityCorr = 0.
protected

◆ fDensityFactor

G4double G4eBremsstrahlungRelModel::fDensityFactor = 0.
protected

Definition at line 152 of file G4eBremsstrahlungRelModel.hh.

Referenced by SetupForMaterial().

◆ fGammaParticle

G4ParticleDefinition* G4eBremsstrahlungRelModel::fGammaParticle = nullptr
protected

◆ fIsElectron

G4bool G4eBremsstrahlungRelModel::fIsElectron = true
protected

◆ fIsLPMActive

G4bool G4eBremsstrahlungRelModel::fIsLPMActive = false
protected

Definition at line 142 of file G4eBremsstrahlungRelModel.hh.

Referenced by SampleSecondaries(), and SetupForMaterial().

◆ fIsScatOffElectron

G4bool G4eBremsstrahlungRelModel::fIsScatOffElectron = false
protected

◆ fLowestKinEnergy

G4double G4eBremsstrahlungRelModel::fLowestKinEnergy
protected

Definition at line 154 of file G4eBremsstrahlungRelModel.hh.

Referenced by G4eBremsstrahlungRelModel(), and MinPrimaryEnergy().

◆ fNucTerm

G4double G4eBremsstrahlungRelModel::fNucTerm = 0.
protected

Definition at line 156 of file G4eBremsstrahlungRelModel.hh.

Referenced by ComputeDXSectionPerAtom(), and SampleSecondaries().

◆ fParticleChange

G4ParticleChangeForLoss* G4eBremsstrahlungRelModel::fParticleChange = nullptr
protected

◆ fPrimaryKinEnergy

G4double G4eBremsstrahlungRelModel::fPrimaryKinEnergy = 0.
protected

◆ fPrimaryParticle

const G4ParticleDefinition* G4eBremsstrahlungRelModel::fPrimaryParticle = nullptr
protected

◆ fPrimaryParticleMass

G4double G4eBremsstrahlungRelModel::fPrimaryParticleMass = 0.
protected

◆ fPrimaryTotalEnergy

◆ fSumTerm

G4double G4eBremsstrahlungRelModel::fSumTerm = 0.
protected

Definition at line 157 of file G4eBremsstrahlungRelModel.hh.

Referenced by ComputeDXSectionPerAtom(), and SampleSecondaries().

◆ 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 167 of file G4eBremsstrahlungRelModel.hh.

Referenced by ComputeCrossSectionPerAtom(), ComputeDEDXPerVolume(), and G4LivermoreBremsstrahlungModel::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 168 of file G4eBremsstrahlungRelModel.hh.

Referenced by SetupForMaterial().


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