Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LivermoreBremsstrahlungModel Class Reference

#include <G4LivermoreBremsstrahlungModel.hh>

+ Inheritance diagram for G4LivermoreBremsstrahlungModel:

Public Member Functions

 G4LivermoreBremsstrahlungModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="LowEnBrem")
 
virtual ~G4LivermoreBremsstrahlungModel ()
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
void InitialiseForElement (const G4ParticleDefinition *, G4int Z) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy) override
 
void SetBicubicInterpolationFlag (G4bool)
 
G4LivermoreBremsstrahlungModeloperator= (const G4LivermoreBremsstrahlungModel &right)=delete
 
 G4LivermoreBremsstrahlungModel (const G4LivermoreBremsstrahlungModel &)=delete
 
- Public Member Functions inherited from G4eBremsstrahlungRelModel
 G4eBremsstrahlungRelModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="eBremLPM")
 
 ~G4eBremsstrahlungRelModel () 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 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 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)
 
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 *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
void SetLPMFlag (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Member Functions

G4double ComputeDXSectionPerAtom (G4double gammaEnergy) override
 
G4String DirectoryPath () const
 
- Protected Member Functions inherited from G4eBremsstrahlungRelModel
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 *)
 

Additional Inherited Members

- Protected Attributes inherited from G4eBremsstrahlungRelModel
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 inherited from G4eBremsstrahlungRelModel
static const G4double gBremFactor
 
static const G4double gMigdalConstant
 

Detailed Description

Definition at line 60 of file G4LivermoreBremsstrahlungModel.hh.

Constructor & Destructor Documentation

◆ G4LivermoreBremsstrahlungModel() [1/2]

G4LivermoreBremsstrahlungModel::G4LivermoreBremsstrahlungModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "LowEnBrem" )
explicit

Definition at line 85 of file G4LivermoreBremsstrahlungModel.cc.

87 : G4eBremsstrahlungRelModel(p,nam),useBicubicInterpolation(false)
88{
89 SetLowEnergyLimit(10.0*eV);
91}
void SetLowEnergyLimit(G4double)
void SetAngularDistribution(G4VEmAngularDistribution *)
G4eBremsstrahlungRelModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="eBremLPM")

◆ ~G4LivermoreBremsstrahlungModel()

G4LivermoreBremsstrahlungModel::~G4LivermoreBremsstrahlungModel ( )
virtual

Definition at line 95 of file G4LivermoreBremsstrahlungModel.cc.

96{
97 if(IsMaster()) {
98 for(size_t i=0; i<101; ++i) {
99 if(dataSB[i]) {
100 delete dataSB[i];
101 dataSB[i] = nullptr;
102 }
103 }
104 }
105}
G4bool IsMaster() const

◆ G4LivermoreBremsstrahlungModel() [2/2]

G4LivermoreBremsstrahlungModel::G4LivermoreBremsstrahlungModel ( const G4LivermoreBremsstrahlungModel & )
delete

Member Function Documentation

◆ ComputeDXSectionPerAtom()

G4double G4LivermoreBremsstrahlungModel::ComputeDXSectionPerAtom ( G4double gammaEnergy)
overrideprotectedvirtual

Reimplemented from G4eBremsstrahlungRelModel.

Definition at line 189 of file G4LivermoreBremsstrahlungModel.cc.

190{
191 if(gammaEnergy < 0.0 || fPrimaryKinEnergy <= 0.0) { return 0.0; }
192 G4double x = gammaEnergy/fPrimaryKinEnergy;
194 G4int Z = fCurrentIZ;
195
196 //G4cout << "G4LivermoreBremsstrahlungModel::ComputeDXSectionPerAtom Z= " << Z
197 // << " x= " << x << " y= " << y << " " << dataSB[Z] << G4endl;
198 if(!dataSB[Z]) { InitialiseForElement(0, Z); }
199
202 G4double cross = dataSB[Z]->Value(x,y,idx,idy)*invb2*millibarn/gBremFactor;
203
204 if(!fIsElectron) {
205 G4double invbeta1 = sqrt(invb2);
206 G4double e2 = fPrimaryKinEnergy - gammaEnergy;
207 if(e2 > 0.0) {
208 G4double invbeta2 = (e2 + fPrimaryParticleMass)
209 /sqrt(e2*(e2 + 2.*fPrimaryParticleMass));
210 G4double xxx = alpha*fCurrentIZ*(invbeta1 - invbeta2);
211 if(xxx < expnumlim) { cross = 0.0; }
212 else { cross *= G4Exp(xxx); }
213 } else {
214 cross = 0.0;
215 }
216 }
217
218 return cross;
219}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:180
G4double G4Log(G4double x)
Definition G4Log.hh:227
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
void InitialiseForElement(const G4ParticleDefinition *, G4int Z) override
G4double Value(G4double x, G4double y, std::size_t &lastidx, std::size_t &lastidy) const

◆ DirectoryPath()

G4String G4LivermoreBremsstrahlungModel::DirectoryPath ( ) const
protected

Definition at line 136 of file G4LivermoreBremsstrahlungModel.cc.

137{
138 return "/livermore/brem/br";
139}

◆ Initialise()

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

Reimplemented from G4eBremsstrahlungRelModel.

Definition at line 109 of file G4LivermoreBremsstrahlungModel.cc.

111{
112 // Access to elements
113 if(IsMaster()) {
114 // check environment variable
115 // Build the complete string identifying the file with the data set
116 const char* path = G4FindDataDir("G4LEDATA");
117
118 const G4ElementTable* theElmTable = G4Element::GetElementTable();
119 size_t numOfElm = G4Element::GetNumberOfElements();
120 if(numOfElm > 0) {
121 for(size_t i=0; i<numOfElm; ++i) {
122 G4int Z = (*theElmTable)[i]->GetZasInt();
123 if(Z < 1) { Z = 1; }
124 else if(Z > 100) { Z = 100; }
125 //G4cout << "Z= " << Z << G4endl;
126 // Initialisation
127 if(!dataSB[Z]) { ReadData(Z, path); }
128 }
129 }
130 }
132}
std::vector< G4Element * > G4ElementTable
const char * G4FindDataDir(const char *)
static G4ElementTable * GetElementTable()
Definition G4Element.cc:389
static size_t GetNumberOfElements()
Definition G4Element.cc:393
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override

◆ InitialiseForElement()

void G4LivermoreBremsstrahlungModel::InitialiseForElement ( const G4ParticleDefinition * ,
G4int Z )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 353 of file G4LivermoreBremsstrahlungModel.cc.

356{
357 G4AutoLock l(&LivermoreBremsstrahlungModelMutex);
358 if(!dataSB[Z]) { ReadData(Z); }
359 l.unlock();
360}

Referenced by ComputeDXSectionPerAtom().

◆ operator=()

G4LivermoreBremsstrahlungModel & G4LivermoreBremsstrahlungModel::operator= ( const G4LivermoreBremsstrahlungModel & right)
delete

◆ SampleSecondaries()

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

Reimplemented from G4eBremsstrahlungRelModel.

Definition at line 224 of file G4LivermoreBremsstrahlungModel.cc.

230{
231 G4double kineticEnergy = dp->GetKineticEnergy();
232 G4double cut = std::min(cutEnergy, kineticEnergy);
233 G4double emax = std::min(maxEnergy, kineticEnergy);
234 if(cut >= emax) { return; }
235 // sets total energy, kinetic energy and density correction
236 SetupForMaterial(fPrimaryParticle, couple->GetMaterial(), kineticEnergy);
237
238 const G4Element* elm =
239 SelectRandomAtom(couple,fPrimaryParticle,kineticEnergy,cut,emax);
240 fCurrentIZ = elm->GetZasInt();
241 G4int Z = fCurrentIZ;
242
243 G4double totMomentum = sqrt(kineticEnergy*(fPrimaryTotalEnergy+electron_mass_c2));
244 /*
245 G4cout << "G4LivermoreBremsstrahlungModel::SampleSecondaries E(MeV)= "
246 << kineticEnergy/MeV
247 << " Z= " << Z << " cut(MeV)= " << cut/MeV
248 << " emax(MeV)= " << emax/MeV << " corr= " << fDensityCorr << G4endl;
249 */
250 G4double xmin = G4Log(cut*cut + fDensityCorr);
251 G4double xmax = G4Log(emax*emax + fDensityCorr);
252 G4double y = G4Log(kineticEnergy/MeV);
253
254 G4double gammaEnergy, v;
255
256 // majoranta
257 G4double x0 = cut/kineticEnergy;
258 G4double vmax = dataSB[Z]->Value(x0, y, idx, idy)*1.02;
259
260 // majoranta corrected for e-
261 if(fIsElectron && x0 < 0.97 &&
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; }
265 }
266 if(x0 < 0.05) { vmax *= 1.2; }
267
268 do {
269 //++ncount;
270 G4double x = G4Exp(xmin + G4UniformRand()*(xmax - xmin)) - fDensityCorr;
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);
275
276 // correction for positrons
277 if(!fIsElectron) {
278 G4double e1 = kineticEnergy - cut;
279 G4double invbeta1 = (e1 + fPrimaryParticleMass)
280 /sqrt(e1*(e1 + 2*fPrimaryParticleMass));
281 G4double e2 = kineticEnergy - gammaEnergy;
282 G4double invbeta2 = (e2 + fPrimaryParticleMass)
283 /sqrt(e2*(e2 + 2*fPrimaryParticleMass));
284 G4double xxx = twopi*fine_structure_const*fCurrentIZ*(invbeta1 - invbeta2);
285
286 if(xxx < expnumlim) { v = 0.0; }
287 else { v *= G4Exp(xxx); }
288 }
289
290 if (v > 1.05*vmax && nwarn < 5) {
291 ++nwarn;
293 ed << "### G4LivermoreBremsstrahlungModel Warning: Majoranta exceeded! "
294 << v << " > " << vmax << " by " << v/vmax
295 << " Egamma(MeV)= " << gammaEnergy
296 << " Ee(MeV)= " << kineticEnergy
297 << " Z= " << Z << " " << fPrimaryParticle->GetParticleName();
298
299 if ( 20 == nwarn ) {
300 ed << "\n ### G4LivermoreBremsstrahlungModel Warnings stopped";
301 }
302 G4Exception("G4LivermoreBremsstrahlungModel::SampleScattering","em0044",
303 JustWarning, ed,"");
304
305 }
306 } while (v < vmax*G4UniformRand());
307
308 //
309 // angles of the emitted gamma. ( Z - axis along the parent particle)
310 // use general interface
311 //
312
313 G4ThreeVector gammaDirection =
315 Z, couple->GetMaterial());
316
317 // create G4DynamicParticle object for the Gamma
318 G4DynamicParticle* gamma =
319 new G4DynamicParticle(fGammaParticle,gammaDirection,gammaEnergy);
320 vdp->push_back(gamma);
321
322 G4ThreeVector direction = (totMomentum*dp->GetMomentumDirection()
323 - gammaEnergy*gammaDirection).unit();
324
325 /*
326 G4cout << "### G4SBModel: v= "
327 << " Eg(MeV)= " << gammaEnergy
328 << " Ee(MeV)= " << kineticEnergy
329 << " DirE " << direction << " DirG " << gammaDirection
330 << G4endl;
331 */
332 // energy of primary
333 G4double finalE = kineticEnergy - gammaEnergy;
334
335 // stop tracking and create new secondary instead of primary
336 if(gammaEnergy > SecondaryThreshold()) {
341 direction, finalE);
342 vdp->push_back(el);
343
344 // continue tracking
345 } else {
348 }
349}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
@ fStopAndKill
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector unit() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4int GetZasInt() const
Definition G4Element.hh:120
const G4Material * GetMaterial() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
const G4String & GetParticleName() const
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4VEmAngularDistribution * GetAngularDistribution()
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)
G4double SecondaryThreshold() const
void ProposeTrackStatus(G4TrackStatus status)
const G4ParticleDefinition * fPrimaryParticle
void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double) override
G4ParticleChangeForLoss * fParticleChange

◆ SetBicubicInterpolationFlag()

void G4LivermoreBremsstrahlungModel::SetBicubicInterpolationFlag ( G4bool val)
inline

Definition at line 98 of file G4LivermoreBremsstrahlungModel.hh.

99{
100 useBicubicInterpolation = val;
101}

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