76G4double G4LivermoreBremsstrahlungModel::ylimit[] = {0.0};
77G4double G4LivermoreBremsstrahlungModel::expnumlim = -12.;
80static const G4double alpha = CLHEP::twopi*CLHEP::fine_structure_const;
81static const G4double epeaklimit= 300*CLHEP::MeV;
82static const G4double elowlimit = 20*CLHEP::keV;
98 for(
size_t i=0; i<101; ++i) {
121 for(
size_t i=0; i<numOfElm; ++i) {
122 G4int Z = (*theElmTable)[i]->GetZasInt();
124 else if(
Z > 100) {
Z = 100; }
127 if(!dataSB[
Z]) { ReadData(
Z, path); }
138 return "/livermore/brem/br";
143void G4LivermoreBremsstrahlungModel::ReadData(
G4int Z,
const char* path)
145 if(dataSB[
Z]) {
return; }
146 const char* datadir = path;
151 G4Exception(
"G4LivermoreBremsstrahlungModel::ReadData()",
"em0006",
156 std::ostringstream ost;
158 std::ifstream fin(ost.str().c_str());
159 if( !fin.is_open()) {
161 ed <<
"Bremsstrahlung data file <" << ost.str().c_str()
162 <<
"> is not opened!";
163 G4Exception(
"G4LivermoreBremsstrahlungModel::ReadData()",
"em0003",
165 "G4LEDATA version should be G4EMLOW8.0 or later.");
174 ylimit[
Z] = v->
Value(0.97, emaxlog, idx, idy);
177 ed <<
"Bremsstrahlung data file <" << ost.str().c_str()
178 <<
"> is not retrieved!";
179 G4Exception(
"G4LivermoreBremsstrahlungModel::ReadData()",
"em0005",
181 "G4LEDATA version should be G4EMLOW8.0 or later.");
211 if(xxx < expnumlim) { cross = 0.0; }
212 else { cross *=
G4Exp(xxx); }
225 std::vector<G4DynamicParticle*>* vdp,
232 G4double cut = std::min(cutEnergy, kineticEnergy);
233 G4double emax = std::min(maxEnergy, kineticEnergy);
234 if(cut >= emax) {
return; }
262 ((kineticEnergy > epeaklimit) || (kineticEnergy < elowlimit))) {
263 G4double ylim = std::min(ylimit[
Z],1.1*dataSB[
Z]->
Value(0.97,y,idx,idy));
264 if(ylim > vmax) { vmax = ylim; }
266 if(x0 < 0.05) { vmax *= 1.2; }
271 if(x < 0.0) { x = 0.0; }
272 gammaEnergy = sqrt(x);
273 G4double x1 = gammaEnergy/kineticEnergy;
274 v = dataSB[
Z]->
Value(x1, y, idx, idy);
281 G4double e2 = kineticEnergy - gammaEnergy;
286 if(xxx < expnumlim) { v = 0.0; }
287 else { v *=
G4Exp(xxx); }
290 if (v > 1.05*vmax && nwarn < 5) {
293 ed <<
"### G4LivermoreBremsstrahlungModel Warning: Majoranta exceeded! "
294 << v <<
" > " << vmax <<
" by " << v/vmax
295 <<
" Egamma(MeV)= " << gammaEnergy
296 <<
" Ee(MeV)= " << kineticEnergy
300 ed <<
"\n ### G4LivermoreBremsstrahlungModel Warnings stopped";
302 G4Exception(
"G4LivermoreBremsstrahlungModel::SampleScattering",
"em0044",
320 vdp->push_back(gamma);
323 - gammaEnergy*gammaDirection).unit();
333 G4double finalE = kineticEnergy - gammaEnergy;
357 G4AutoLock l(&LivermoreBremsstrahlungModelMutex);
358 if(!dataSB[
Z]) { ReadData(
Z); }
std::vector< G4Element * > G4ElementTable
const char * G4FindDataDir(const char *)
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
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
static G4ElementTable * GetElementTable()
static size_t GetNumberOfElements()
G4LivermoreBremsstrahlungModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="LowEnBrem")
virtual ~G4LivermoreBremsstrahlungModel()
G4String DirectoryPath() const
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void InitialiseForElement(const G4ParticleDefinition *, G4int Z) override
G4double ComputeDXSectionPerAtom(G4double gammaEnergy) override
const G4Material * GetMaterial() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
const G4String & GetParticleName() const
G4bool Retrieve(std::ifstream &fIn)
G4double Value(G4double x, G4double y, std::size_t &lastidx, std::size_t &lastidy) const
void SetBicubicInterpolation(G4bool)
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4VEmAngularDistribution * GetAngularDistribution()
void SetLPMFlag(G4bool val)
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
void SetLowEnergyLimit(G4double)
void SetAngularDistribution(G4VEmAngularDistribution *)
G4double SecondaryThreshold() const
void ProposeTrackStatus(G4TrackStatus status)
G4double fPrimaryParticleMass
G4double fPrimaryKinEnergy
const G4ParticleDefinition * fPrimaryParticle
static const G4double gBremFactor
G4double fPrimaryTotalEnergy
G4ParticleDefinition * fGammaParticle
void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double) override
G4ParticleChangeForLoss * fParticleChange
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override