80G4double G4SeltzerBergerModel::gYLimitData[] = { 0.0 };
82G4SBBremTable* G4SeltzerBergerModel::gSBSamplingTable =
nullptr;
85const G4double G4SeltzerBergerModel::gBremFactor
86 = 16. * CLHEP::fine_structure_const * CLHEP::classic_electr_radius
87 * CLHEP::classic_electr_radius/3.;
90const G4double G4SeltzerBergerModel::gMigdalConstant
91 = 4. * CLHEP::pi * CLHEP::classic_electr_radius
92 * CLHEP::electron_Compton_length * CLHEP::electron_Compton_length;
94static constexpr G4double twoMass = 2* CLHEP::electron_mass_c2;
95static constexpr G4double kAlpha = CLHEP::twopi*CLHEP::fine_structure_const;
96static std::once_flag applyOnce;
104 1.98550718e-02, 1.01666761e-01, 2.37233795e-01, 4.08282679e-01,
105 5.91717321e-01, 7.62766205e-01, 8.98333239e-01, 9.80144928e-01
108 5.06142681e-02, 1.11190517e-01, 1.56853323e-01, 1.81341892e-01,
109 1.81341892e-01, 1.56853323e-01, 1.11190517e-01, 5.06142681e-02
116 fGammaParticle(
G4Gamma::Gamma()),
117 fLowestKinEnergy(1.0*
CLHEP::keV)
121 if (fPrimaryParticle != p) { SetParticle(p); }
128 for (std::size_t iz = 0; iz < gMaxZet; ++iz) {
129 if (gSBDCSData[iz]) {
130 delete gSBDCSData[iz];
131 gSBDCSData[iz] =
nullptr;
134 if (gSBSamplingTable) {
135 delete gSBSamplingTable;
136 gSBSamplingTable =
nullptr;
145 if (fPrimaryParticle != p) {
152 std::call_once(applyOnce, [
this]() { isInitializer =
true; });
159 for (
auto const & elm : *elemTable) {
160 G4int Z = std::max(1,std::min(elm->GetZasInt(), gMaxZet-1));
162 if (gSBDCSData[Z] ==
nullptr) ReadData(Z);
166 if (fIsUseSamplingTables) {
167 if (
nullptr == gSBSamplingTable) {
170 gSBSamplingTable->Initialize(std::max(fLowestKinEnergy,
LowEnergyLimit()),
184 if (
nullptr != trmodel) {
185 trmodel->Initialise(p, cuts);
186 fIsScatOffElectron =
true;
198 fPrimaryParticle = p;
202void G4SeltzerBergerModel::ReadData(
G4int Z) {
204 if (gSBDCSData[Z] !=
nullptr)
return;
206 if (gSBDCSData[Z] ==
nullptr) {
207 std::ostringstream ost;
209 std::ifstream fin(ost.str().c_str());
210 if (!fin.is_open()) {
212 ed <<
"Bremsstrahlung data file <" << ost.str().c_str()
213 <<
"> is not opened!";
215 ed,
"G4LEDATA version should be G4EMLOW6.23 or later.");
220 auto v =
new G4Physics2DVector();
221 if (v->Retrieve(fin)) {
222 v->SetBicubicInterpolation(fIsUseBicubicInterpolation);
224 gYLimitData[Z] = v->Value(0.97, emaxlog, fIndx, fIndy);
228 ed <<
"Bremsstrahlung data file <" << ost.str().c_str()
229 <<
"> is not retrieved!";
231 ed,
"G4LEDATA version should be G4EMLOW6.23 or later.");
242 return std::max(fLowestKinEnergy, cut);
254 fPrimaryKinEnergy = kinEnergy;
255 fPrimaryTotalEnergy = kinEnergy + CLHEP::electron_mass_c2;
256 fDensityCorr = fDensityFactor*fPrimaryTotalEnergy*fPrimaryTotalEnergy;
268 if (
nullptr == fPrimaryParticle) {
271 if (kineticEnergy <= fLowestKinEnergy) {
275 G4double tmax = std::min(cutEnergy, kineticEnergy);
284 const std::size_t numberOfElements = theElemVector->size();
287 for (std::size_t ie = 0; ie < numberOfElements; ++ie) {
289 G4int Z = (*theElemVector)[ie]->GetZasInt();
290 fCurrentIZ = std::min(Z, gMaxZet);
291 dedx += (Z*Z)*theAtomNumDensVector[ie]*ComputeBremLoss(tmax);
295 return std::max(dedx, 0.);
317 const G4double alphaMax = tmax/fPrimaryTotalEnergy;
323 for (
G4int l = 0; l < nSub; ++l) {
324 for (
G4int igl = 0; igl < 8; ++igl) {
326 const G4double k = (alpha_i+gXGL[igl]*delta)*fPrimaryTotalEnergy;
328 const G4double dcs = ComputeDXSectionPerAtom(k);
330 dedxInteg += gWGL[igl]*dcs/(1.0+fDensityCorr/(k*k));
336 dedxInteg *= delta*fPrimaryTotalEnergy;
337 return std::max(dedxInteg,0.);
351 if (
nullptr == fPrimaryParticle) {
354 if (kineticEnergy <= fLowestKinEnergy) {
358 const G4double tmin = std::min(cut, kineticEnergy);
359 const G4double tmax = std::min(maxEnergy, kineticEnergy);
364 fCurrentIZ = std::min(
G4lrint(Z), gMaxZet);
367 crossSection = ComputeXSectionPerAtom(tmin);
372 if (tmax < kineticEnergy) {
373 crossSection -= ComputeXSectionPerAtom(tmax);
376 crossSection *= Z*Z*gBremFactor;
377 return std::max(crossSection, 0.);
399 const G4double alphaMax =
G4Log(fPrimaryKinEnergy/fPrimaryTotalEnergy);
400 const G4int nSub = (
G4int)(0.45*(alphaMax-alphaMin))+4;
404 for (
G4int l = 0; l < nSub; ++l) {
405 for (
G4int igl = 0; igl < 8; ++igl) {
407 const G4double k =
G4Exp(alpha_i+gXGL[igl]*delta)*fPrimaryTotalEnergy;
409 const G4double dcs = ComputeDXSectionPerAtom(k);
411 xSection += gWGL[igl]*dcs/(1.0+fDensityCorr/(k*k));
419 return std::max(xSection, 0.);
425 if (gammaEnergy < 0.0 || fPrimaryKinEnergy <= 0.0) {
429 const G4double x = gammaEnergy/fPrimaryKinEnergy;
434 fCurrentIZ = std::max(std::min(fCurrentIZ, gMaxZet-1), 1);
435 if (
nullptr == gSBDCSData[fCurrentIZ]) {
437 ReadData(fCurrentIZ);
441 const G4double pt2 = fPrimaryKinEnergy*(fPrimaryKinEnergy + twoMass);
442 const G4double invb2 = fPrimaryTotalEnergy*fPrimaryTotalEnergy/pt2;
443 G4double val = gSBDCSData[fCurrentIZ]->Value(x,y,fIndx,fIndy);
444 dxsec = val*invb2*CLHEP::millibarn/gBremFactor;
447 const G4double invbeta1 = std::sqrt(invb2);
448 const G4double e2 = fPrimaryKinEnergy - gammaEnergy;
451 (e2 + CLHEP::electron_mass_c2)/std::sqrt(e2*(e2 + twoMass));
452 const G4double dum0 = kAlpha*fCurrentIZ*(invbeta1-invbeta2);
453 if (dum0 < gExpNumLimit) {
456 dxsec *=
G4Exp(dum0);
474 const G4double tmin = std::min(cutEnergy, kinEnergy);
475 const G4double tmax = std::min(maxEnergy, kinEnergy);
482 logKinEnergy, tmin, tmax);
483 fCurrentIZ = std::max(std::min(elm->
GetZasInt(), gMaxZet-1), 1);
485 const G4double totMomentum = std::sqrt(kinEnergy*(kinEnergy + twoMass));
493 const G4double gammaEnergy = !fIsUseSamplingTables
494 ? SampleEnergyTransfer(kinEnergy, logKinEnergy, tmin, tmax)
495 : gSBSamplingTable->SampleEnergyTransfer(kinEnergy, logKinEnergy, tmin,
496 fDensityCorr, fCurrentIZ, couple->
GetIndex(), fIsElectron);
498 if (gammaEnergy <= 0.) {
505 fPrimaryTotalEnergy-gammaEnergy, fCurrentIZ, couple->
GetMaterial());
508 vdp->push_back(gamma);
513 const G4double finalE = kinEnergy - gammaEnergy;
538G4SeltzerBergerModel::SampleEnergyTransfer(
const G4double kinEnergy,
551 if (
nullptr == gSBDCSData[fCurrentIZ]) {
552 ReadData(fCurrentIZ);
554 vmax = gSBDCSData[fCurrentIZ]->Value(x0, y, fIndx, fIndy)*1.02;
556 static const G4double kEPeakLim = 300.*CLHEP::MeV;
557 static const G4double kELowLim = 20.*CLHEP::keV;
559 if (fIsElectron && x0 < 0.97 &&
560 ((kinEnergy>kEPeakLim) || (kinEnergy<kELowLim))) {
561 G4double ylim = std::min(gYLimitData[fCurrentIZ],
562 1.1*gSBDCSData[fCurrentIZ]->
Value(0.97,y,fIndx,fIndy));
563 vmax = std::max(vmax, ylim);
570 static const G4int kNCountMax = 100;
571 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
574 for (
G4int nn = 0;
nn < kNCountMax; ++
nn) {
577 std::sqrt(std::max(
G4Exp(xmin + rndm[0]*xrange)-fDensityCorr,0.));
578 v = gSBDCSData[fCurrentIZ]->Value(gammaEnergy/kinEnergy, y, fIndx, fIndy);
581 const G4double e1 = kinEnergy - tmin;
583 (e1 + CLHEP::electron_mass_c2)/std::sqrt(e1*(e1 + twoMass));
584 const G4double e2 = kinEnergy-gammaEnergy;
586 (e2 + CLHEP::electron_mass_c2)/std::sqrt(e2*(e2 + twoMass));
587 const G4double dum0 = kAlpha*fCurrentIZ*(invbeta1-invbeta2);
588 if (dum0 < gExpNumLimit) {
594 if (v > 1.05*vmax && fNumWarnings < 11) {
597 ed <<
"### G4SeltzerBergerModel Warning: Majoranta exceeded! "
598 << v <<
" > " << vmax <<
" by " << v/vmax
600 <<
" Egamma(MeV)= " << gammaEnergy
601 <<
" Ee(MeV)= " << kinEnergy
602 <<
" Z= " << fCurrentIZ <<
" " << fPrimaryParticle->GetParticleName();
604 if (10 == fNumWarnings) {
605 ed <<
"\n ### G4SeltzerBergerModel Warnings stopped";
607 G4Exception(
"G4SeltzerBergerModel::SampleScattering",
"em0044",
610 if (v >= vmax*rndm[1]) {
G4TemplateAutoLock< G4Mutex > G4AutoLock
std::vector< const G4Element * > G4ElementVector
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
G4double G4Log(G4double x)
#define G4MUTEX_INITIALIZER
CLHEP::Hep3Vector G4ThreeVector
virtual void flatArray(const int size, double *vect)=0
const G4ThreeVector & GetMomentumDirection() const
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
static const G4ElementTable * GetElementTable()
static G4EmParameters * Instance()
G4bool EnableSamplingTable() const
const G4String & GetDirLEDATA() const
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
const G4double * GetAtomicNumDensityVector() const
G4double GetElectronDensity() const
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double ekin, G4double zet, G4double, G4double cutEnergy, G4double maxEnergy=DBL_MAX) override
G4ParticleChangeForLoss * fParticleChange
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double) override
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy) override
G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double cutEnergy) override
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double ekin, G4double cutEnergy) override
~G4SeltzerBergerModel() override
G4SeltzerBergerModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="eBremSB")
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
G4VEmAngularDistribution * GetAngularDistribution()
G4double LowEnergyLimit() const
std::vector< G4EmElementSelector * > * GetElementSelectors()
G4VEmModel * GetTripletModel()
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
G4double HighEnergyLimit() const
void SetCurrentElement(const G4Element *)
void SetLowEnergyLimit(G4double)
G4VEmModel(const G4String &nam)
void SetAngularDistribution(G4VEmAngularDistribution *)
const G4Element * SelectTargetAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
G4double SecondaryThreshold() const
G4ParticleChangeForLoss * GetParticleChangeForLoss()