Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 G4MaterialCutsCouple *, const G4DynamicParticle *, 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 60 of file G4eBremsstrahlungRelModel.hh.

Constructor & Destructor Documentation

◆ G4eBremsstrahlungRelModel()

G4eBremsstrahlungRelModel::G4eBremsstrahlungRelModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "eBremLPM" )
explicit

Definition at line 144 of file G4eBremsstrahlungRelModel.cc.

146: G4VEmModel(nam), fLPMFuncs(gLPMFuncs()), fElementData(gElementData())
147{
149 //
150 fLowestKinEnergy = 1.0*CLHEP::MeV;
152 //
153 fLPMEnergyThreshold = 1.e+39;
154 fLPMEnergy = 0.;
155 SetAngularDistribution(new G4ModifiedTsai());
156 //
157 if (nullptr != p) {
158 SetParticle(p);
159 }
160}
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 162 of file G4eBremsstrahlungRelModel.cc.

163{
164 if (fIsInitializer) {
165 // clear ElementData container
166 for (auto const & ptr : *fElementData) { delete ptr; }
167 fElementData->clear();
168 // clear LPMFunctions (if any)
169 if (fLPMFuncs->fIsInitialized) {
170 fLPMFuncs->fLPMFuncG.clear();
171 fLPMFuncs->fLPMFuncPhi.clear();
172 fLPMFuncs->fIsInitialized = false;
173 }
174 }
175}

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

352{
353 G4double crossSection = 0.0;
354 if (nullptr == fPrimaryParticle) {
355 SetParticle(p);
356 }
357 if (kineticEnergy < LowEnergyLimit()) {
358 return crossSection;
359 }
360 // min/max kinetic energy limits of the DCS integration:
361 const G4double tmin = std::min(cut, kineticEnergy);
362 const G4double tmax = std::min(maxEnergy, kineticEnergy);
363 // zero restricted x-section if e- kinetic energy is below gamma cut
364 if (tmin >= tmax) {
365 return crossSection;
366 }
367 fCurrentIZ = std::min(G4lrint(Z), gMaxZet);
368 // integrate numerically (dependent part of) the DCS between the kin. limits:
369 // a. integrate between tmin and kineticEnergy of the e-
370 crossSection = ComputeXSectionPerAtom(tmin);
371 // allow partial integration: only if maxEnergy < kineticEnergy
372 // b. integrate between tmax and kineticEnergy (tmax=maxEnergy in this case)
373 // (so the result in this case is the integral of DCS between tmin and
374 // maxEnergy)
375 if (tmax < kineticEnergy) {
376 crossSection -= ComputeXSectionPerAtom(tmax);
377 }
378 // multiply with the constant factors: 16\alpha r_0^2/3 Z^2
379 crossSection *= Z*Z*gBremFactor;
380 return std::max(crossSection, 0.);
381}
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 262 of file G4eBremsstrahlungRelModel.cc.

266{
267 G4double dedx = 0.0;
268 if (nullptr == fPrimaryParticle) {
269 SetParticle(p);
270 }
271 if (kineticEnergy < LowEnergyLimit()) {
272 return dedx;
273 }
274 // maximum value of the dE/dx integral (the minimum is 0 of course)
275 G4double tmax = std::min(cutEnergy, kineticEnergy);
276 if (tmax == 0.0) {
277 return dedx;
278 }
279 // sets kinematical and material related variables
280 SetupForMaterial(fPrimaryParticle, material,kineticEnergy);
281 // get element compositions of the material
282 const G4ElementVector* theElemVector = material->GetElementVector();
283 const G4double* theAtomNumDensVector = material->GetAtomicNumDensityVector();
284 const std::size_t numberOfElements = theElemVector->size();
285 // loop over the elements of the material and compute their contributions to
286 // the restricted dE/dx by numerical integration of the dependent part of DCS
287 for (std::size_t ie = 0; ie < numberOfElements; ++ie) {
288 G4VEmModel::SetCurrentElement((*theElemVector)[ie]);
289 G4int zet = (*theElemVector)[ie]->GetZasInt();
290 fCurrentIZ = std::min(zet, gMaxZet);
291 dedx += (zet*zet)*theAtomNumDensVector[ie]*ComputeBremLoss(tmax);
292 }
293 // apply the constant factor C/Z = 16\alpha r_0^2/3
294 dedx *= gBremFactor;
295 return std::max(dedx,0.);
296}
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 501 of file G4eBremsstrahlungRelModel.cc.

502{
503 G4double dxsec = 0.0;
504 if (gammaEnergy < 0.) {
505 return dxsec;
506 }
507 const G4double y = gammaEnergy/fPrimaryTotalEnergy;
508 const G4double onemy = 1.-y;
509 const G4double dum0 = onemy+0.75*y*y;
510 const ElementData* elDat = (*fElementData)[fCurrentIZ];
511 // use complete screening and L_el, L_inel from Dirac-Fock model instead of TF
512 if (fCurrentIZ < 5 || fIsUseCompleteScreening) {
513 dxsec = dum0*elDat->fZFactor1;
514 dxsec += onemy*elDat->fZFactor2;
515 if (fIsScatOffElectron) {
516 fSumTerm = dxsec;
517 fNucTerm = dum0*elDat->fZFactor11+onemy/12.;
518 }
519 } else {
520 // use Tsai's analytical approx. (Tsai Eqs. [3.38-3.41]) to the 'universal'
521 // numerical screening functions computed by using the TF model of the atom
522 const G4double invZ = 1./(G4double)fCurrentIZ;
523 const G4double Fz = elDat->fFz;
524 const G4double logZ = elDat->fLogZ;
525 const G4double dum1 = y/(fPrimaryTotalEnergy-gammaEnergy);
526 const G4double gamma = dum1*elDat->fGammaFactor;
527 const G4double epsilon = dum1*elDat->fEpsilonFactor;
528 // evaluate the screening functions
529 G4double phi1, phi1m2, psi1, psi1m2;
530 ComputeScreeningFunctions(phi1, phi1m2, psi1, psi1m2, gamma, epsilon);
531 dxsec = dum0*((0.25*phi1-Fz) + (0.25*psi1-2.*logZ/3.)*invZ);
532 dxsec += 0.125*onemy*(phi1m2 + psi1m2*invZ);
533 if (fIsScatOffElectron) {
534 fSumTerm = dxsec;
535 fNucTerm = dum0*(0.25*phi1-Fz) + 0.125*onemy*phi1m2;
536 }
537 }
538 return std::max(dxsec,0.0);
539}
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 177 of file G4eBremsstrahlungRelModel.cc.

179{
180 // parameters in each thread
181 if (fPrimaryParticle != p) {
182 SetParticle(p);
183 }
184 fUseLPM = G4EmParameters::Instance()->LPM();
185 fCurrentIZ = 0;
186
187 // init static element data and precompute LPM functions only once
188 std::call_once(applyOnce, [this]() { fIsInitializer = true; });
189
190 // for all treads and derived classes
191 if (fIsInitializer || fElementData->empty()) {
192 G4AutoLock l(&theBremRelMutex);
193 if (fElementData->empty()) {
194 fElementData->resize(gMaxZet+1, nullptr);
195 }
196 InitialiseElementData();
197 InitLPMFunctions();
198 l.unlock();
199 }
200
201 // element selectors are initialized in the master thread
202 if (IsMaster()) {
204 }
205 // initialisation in all threads
206 if (nullptr == fParticleChange) {
208 }
209 if (GetTripletModel()) {
210 GetTripletModel()->Initialise(p, cuts);
211 fIsScatOffElectron = true;
212 }
213}
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 215 of file G4eBremsstrahlungRelModel.cc.

217{
219}
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 252 of file G4eBremsstrahlungRelModel.cc.

255{
256 return std::max(fLowestKinEnergy, cut);
257}

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

570{
571 const G4double kineticEnergy = dp->GetKineticEnergy();
572 if (kineticEnergy < LowEnergyLimit()) {
573 return;
574 }
575 // min, max kinetic energy limits
576 const G4double tmin = std::min(cutEnergy, kineticEnergy);
577 const G4double tmax = std::min(maxEnergy, kineticEnergy);
578 if (tmin >= tmax) {
579 return;
580 }
581 //
582 SetupForMaterial(fPrimaryParticle, couple->GetMaterial(), kineticEnergy);
583 const G4Element* elm = SelectTargetAtom(couple,fPrimaryParticle,kineticEnergy,
584 dp->GetLogKineticEnergy(),tmin,tmax);
585 //
586 fCurrentIZ = elm->GetZasInt();
587 const ElementData* elDat = (*fElementData)[fCurrentIZ];
588 const G4double funcMax = elDat->fZFactor1+elDat->fZFactor2;
589 // get the random engine
590 G4double rndm[2];
591 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
592 // 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)]
593 const G4double xmin = G4Log(tmin*tmin+fDensityCorr);
594 const G4double xrange = G4Log(tmax*tmax+fDensityCorr)-xmin;
595 G4double gammaEnergy, funcVal;
596 do {
597 rndmEngine->flatArray(2, rndm);
598 gammaEnergy = std::sqrt(std::max(G4Exp(xmin+rndm[0]*xrange)-fDensityCorr, 0.0));
599 funcVal = fIsLPMActive
600 ? ComputeRelDXSectionPerAtom(gammaEnergy)
601 : ComputeDXSectionPerAtom(gammaEnergy);
602 // cross-check of proper function maximum in the rejection
603// if (funcVal > funcMax) {
604// G4cout << "### G4eBremsstrahlungRelModel Warning: Majoranta exceeded! "
605// << funcVal << " > " << funcMax
606// << " Egamma(MeV)= " << gammaEnergy
607// << " Ee(MeV)= " << kineticEnergy
608// << " " << GetName()
609// << G4endl;
610// }
611 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
612 } while (funcVal < funcMax*rndm[1]);
613 //
614 // scattering off nucleus or off e- by triplet model
615 if (fIsScatOffElectron && rndmEngine->flat()*fSumTerm>fNucTerm) {
616 GetTripletModel()->SampleSecondaries(vdp, couple, dp, cutEnergy, maxEnergy);
617 return;
618 }
619 //
620 // angles of the emitted gamma. ( Z - axis along the parent particle)
621 // use general interface
622 G4ThreeVector gamDir =
624 fCurrentIZ, couple->GetMaterial());
625 // create G4DynamicParticle object for the Gamma
626 auto gamma = new G4DynamicParticle(fGammaParticle, gamDir, gammaEnergy);
627 vdp->push_back(gamma);
628 // compute post-interaction kinematics of primary e-/e+ based on
629 // energy-momentum conservation
630 const G4double totMomentum = std::sqrt(kineticEnergy*(
631 fPrimaryTotalEnergy + CLHEP::electron_mass_c2));
632 G4ThreeVector dir =
633 (totMomentum*dp->GetMomentumDirection()-gammaEnergy*gamDir).unit();
634 const G4double finalE = kineticEnergy-gammaEnergy;
635 // if secondary gamma energy is higher than threshold(very high by default)
636 // then stop tracking the primary particle and create new secondary e-/e+
637 // instead of the primary one
638 if (gammaEnergy > SecondaryThreshold()) {
639 fParticleChange->ProposeTrackStatus(fStopAndKill);
640 fParticleChange->SetProposedKineticEnergy(0.0);
641 auto el = new G4DynamicParticle(
642 const_cast<G4ParticleDefinition*>(fPrimaryParticle), dir, finalE);
643 vdp->push_back(el);
644 } else { // continue tracking the primary e-/e+ otherwise
645 fParticleChange->SetProposedMomentumDirection(dir);
646 fParticleChange->SetProposedKineticEnergy(finalE);
647 }
648}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:180
G4double G4Log(G4double x)
Definition G4Log.hh:227
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 231 of file G4eBremsstrahlungRelModel.cc.

234{
236 fLPMEnergy = gLPMconstant*mat->GetRadlen();
237 // threshold for LPM effect (i.e. below which LPM hidden by density effect)
238 if (fUseLPM) {
239 fLPMEnergyThreshold = std::sqrt(fDensityFactor)*fLPMEnergy;
240 } else {
241 fLPMEnergyThreshold = 1.e+39; // i.e. do not use LPM effect
242 }
243 // calculate threshold for density effect: k_p = sqrt(fDensityCorr)
244 fPrimaryKinEnergy = kineticEnergy;
247 // set activation flag for LPM effects in the DCS
248 fIsLPMActive = (fPrimaryTotalEnergy>fLPMEnergyThreshold);
249}
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 164 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 154 of file G4eBremsstrahlungRelModel.hh.

Referenced by SampleSecondaries(), and SetupForMaterial().

◆ fIsScatOffElectron

G4bool G4eBremsstrahlungRelModel::fIsScatOffElectron = false
protected

◆ fLowestKinEnergy

G4double G4eBremsstrahlungRelModel::fLowestKinEnergy
protected

Definition at line 166 of file G4eBremsstrahlungRelModel.hh.

Referenced by G4eBremsstrahlungRelModel(), and MinPrimaryEnergy().

◆ fNucTerm

G4double G4eBremsstrahlungRelModel::fNucTerm = 0.
protected

Definition at line 168 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 169 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 179 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 180 of file G4eBremsstrahlungRelModel.hh.

Referenced by SetupForMaterial().


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