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

#include <G4SeltzerBergerModel.hh>

+ Inheritance diagram for G4SeltzerBergerModel:

Public Member Functions

 G4SeltzerBergerModel (const G4ParticleDefinition *p=0, const G4String &nam="eBremSB")
 
virtual ~G4SeltzerBergerModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy)
 
void SetBicubicInterpolationFlag (G4bool)
 
- Public Member Functions inherited from G4eBremsstrahlungRelModel
 G4eBremsstrahlungRelModel (const G4ParticleDefinition *p=0, const G4String &nam="eBremLPM")
 
virtual ~G4eBremsstrahlungRelModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double tkin, G4double Z, G4double, G4double cutEnergy, G4double maxEnergy=DBL_MAX)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double)
 
void SetLPMconstant (G4double val)
 
G4double LPMconstant () const
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., 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 *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
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)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void ForceBuildTable (G4bool val)
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Member Functions

virtual G4double ComputeDXSectionPerAtom (G4double gammaEnergy)
 
- Protected Member Functions inherited from G4eBremsstrahlungRelModel
virtual G4double ComputeDXSectionPerAtom (G4double gammaEnergy)
 
void SetCurrentElement (const G4double)
 
- 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
G4NistManagernist
 
const G4ParticleDefinitionparticle
 
G4ParticleDefinitiontheGamma
 
G4ParticleChangeForLossfParticleChange
 
G4double bremFactor
 
G4double particleMass
 
G4double kinEnergy
 
G4double totalEnergy
 
G4double currentZ
 
G4double densityFactor
 
G4double densityCorr
 
G4bool isElectron
 
- Protected Attributes inherited from G4VEmModel
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

Detailed Description

Definition at line 61 of file G4SeltzerBergerModel.hh.

Constructor & Destructor Documentation

◆ G4SeltzerBergerModel()

G4SeltzerBergerModel::G4SeltzerBergerModel ( const G4ParticleDefinition p = 0,
const G4String nam = "eBremSB" 
)

Definition at line 78 of file G4SeltzerBergerModel.cc.

80 : G4eBremsstrahlungRelModel(p,nam),useBicubicInterpolation(false)
81{
83 SetLPMFlag(false);
84 nwarn = 0;
85}
void SetLPMFlag(G4bool val)
Definition: G4VEmModel.hh:634
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592

◆ ~G4SeltzerBergerModel()

G4SeltzerBergerModel::~G4SeltzerBergerModel ( )
virtual

Definition at line 89 of file G4SeltzerBergerModel.cc.

90{
91 for(size_t i=0; i<101; ++i) {
92 if(dataSB[i]) {
93 delete dataSB[i];
94 dataSB[i] = 0;
95 }
96 }
97}

Member Function Documentation

◆ ComputeDXSectionPerAtom()

G4double G4SeltzerBergerModel::ComputeDXSectionPerAtom ( G4double  gammaEnergy)
protectedvirtual

Reimplemented from G4eBremsstrahlungRelModel.

Definition at line 175 of file G4SeltzerBergerModel.cc.

176{
177
178 if(gammaEnergy < 0.0 || kinEnergy <= 0.0) { return 0.0; }
179 G4double x = gammaEnergy/kinEnergy;
180 G4double y = log(kinEnergy/MeV);
181 G4int Z = G4int(currentZ);
182 //G4cout << "G4SeltzerBergerModel::ComputeDXSectionPerAtom Z= " << Z
183 // << " x= " << x << " y= " << y << " " << dataSB[Z] << G4endl;
184 if(!dataSB[Z]) { ReadData(Z); }
186 G4double cross = dataSB[Z]->Value(x,y)*invb2*millibarn/bremFactor;
187
188 if(!isElectron) {
189 G4double invbeta1 = sqrt(invb2);
190 G4double e2 = kinEnergy - gammaEnergy;
191 if(e2 > 0.0) {
192 G4double invbeta2 = (e2 + particleMass)/sqrt(e2*(e2 + 2*particleMass));
193 G4double xxx = twopi*fine_structure_const*currentZ*(invbeta1 - invbeta2);
194 if(xxx < expnumlim) { cross = 0.0; }
195 else { cross *= exp(xxx); }
196 } else {
197 cross = 0.0;
198 }
199 }
200
201 return cross;
202}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
G4double Value(G4double x, G4double y)

◆ Initialise()

void G4SeltzerBergerModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector cuts 
)
virtual

Reimplemented from G4eBremsstrahlungRelModel.

Definition at line 101 of file G4SeltzerBergerModel.cc.

103{
104 // check environment variable
105 // Build the complete string identifying the file with the data set
106 char* path = getenv("G4LEDATA");
107
108 // Access to elements
109 const G4ElementTable* theElmTable = G4Element::GetElementTable();
110 size_t numOfElm = G4Element::GetNumberOfElements();
111 if(numOfElm > 0) {
112 for(size_t i=0; i<numOfElm; ++i) {
113 G4int Z = G4int(((*theElmTable)[i])->GetZ());
114 if(Z < 1) { Z = 1; }
115 else if(Z > 100) { Z = 100; }
116 //G4cout << "Z= " << Z << G4endl;
117 // Initialisation
118 if(!dataSB[Z]) { ReadData(Z, path); }
119 }
120 }
121
123}
std::vector< G4Element * > G4ElementTable
static size_t GetNumberOfElements()
Definition: G4Element.cc:406
static const G4ElementTable * GetElementTable()
Definition: G4Element.cc:399
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)

◆ SampleSecondaries()

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

Reimplemented from G4eBremsstrahlungRelModel.

Definition at line 207 of file G4SeltzerBergerModel.cc.

212{
213 G4double kineticEnergy = dp->GetKineticEnergy();
214 G4double cut = std::min(cutEnergy, kineticEnergy);
215 G4double emax = std::min(maxEnergy, kineticEnergy);
216 if(cut >= emax) { return; }
217
218 SetupForMaterial(particle, couple->GetMaterial(), kineticEnergy);
219
220 const G4Element* elm =
221 SelectRandomAtom(couple,particle,kineticEnergy,cut,emax);
222 SetCurrentElement(elm->GetZ());
223 G4int Z = G4int(currentZ);
224
225 totalEnergy = kineticEnergy + particleMass;
227 G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2));
228 /*
229 G4cout << "G4SeltzerBergerModel::SampleSecondaries E(MeV)= "
230 << kineticEnergy/MeV
231 << " Z= " << Z << " cut(MeV)= " << cut/MeV
232 << " emax(MeV)= " << emax/MeV << " corr= " << densityCorr << G4endl;
233 */
234 G4double xmin = log(cut*cut + densityCorr);
235 G4double xmax = log(emax*emax + densityCorr);
236 G4double y = log(kineticEnergy/MeV);
237
238 G4double gammaEnergy, v;
239
240 // majoranta
241 G4double x0 = cut/kineticEnergy;
242 G4double vmax = dataSB[Z]->Value(x0, y)*1.02;
243 // G4double invbeta1 = 0;
244
245 const G4double epeaklimit= 300*MeV;
246 const G4double elowlimit = 10*keV;
247
248 // majoranta corrected for e-
249 if(isElectron && x0 < 0.97 &&
250 ((kineticEnergy > epeaklimit) || (kineticEnergy < elowlimit))) {
251 G4double ylim = std::min(ylimit[Z],1.1*dataSB[Z]->Value(0.97, y));
252 if(ylim > vmax) { vmax = ylim; }
253 }
254 if(x0 < 0.05) { vmax *= 1.2; }
255
256 //G4cout<<"y= "<<y<<" xmin= "<<xmin<<" xmax= "<<xmax<<" vmax= "<<vmax<<G4endl;
257 // G4int ncount = 0;
258 do {
259 //++ncount;
260 G4double x = exp(xmin + G4UniformRand()*(xmax - xmin)) - densityCorr;
261 if(x < 0.0) { x = 0.0; }
262 gammaEnergy = sqrt(x);
263 G4double x1 = gammaEnergy/kineticEnergy;
264 v = dataSB[Z]->Value(x1, y);
265
266 // correction for positrons
267 if(!isElectron) {
268 G4double e1 = kineticEnergy - cut;
269 G4double invbeta1 = (e1 + particleMass)/sqrt(e1*(e1 + 2*particleMass));
270 G4double e2 = kineticEnergy - gammaEnergy;
271 G4double invbeta2 = (e2 + particleMass)/sqrt(e2*(e2 + 2*particleMass));
272 G4double xxx = twopi*fine_structure_const*currentZ*(invbeta1 - invbeta2);
273
274 if(xxx < expnumlim) { v = 0.0; }
275 else { v *= exp(xxx); }
276 }
277
278 if (v > 1.05*vmax && nwarn < 20) {
279 ++nwarn;
280 G4cout << "### G4SeltzerBergerModel Warning: Majoranta exceeded! "
281 << v << " > " << vmax << " by " << v/vmax
282 << " Egamma(MeV)= " << gammaEnergy
283 << " Ee(MeV)= " << kineticEnergy
284 << " Z= " << Z << " " << particle->GetParticleName()
285 //<< " ncount= " << ncount
286 << G4endl;
287
288 if ( 20 == nwarn ) {
289 G4cout << "### G4SeltzerBergerModel Warnings will not be printed anymore"
290 << G4endl;
291 }
292 }
293 } while (v < vmax*G4UniformRand());
294
295 //
296 // angles of the emitted gamma. ( Z - axis along the parent particle)
297 // use general interface
298 //
299
300 G4ThreeVector gammaDirection =
302 Z, couple->GetMaterial());
303
304 // create G4DynamicParticle object for the Gamma
305 G4DynamicParticle* gamma =
306 new G4DynamicParticle(theGamma,gammaDirection,gammaEnergy);
307 vdp->push_back(gamma);
308
309 G4ThreeVector direction = (totMomentum*dp->GetMomentumDirection()
310 - gammaEnergy*gammaDirection).unit();
311
312 /*
313 G4cout << "### G4SBModel: v= "
314 << " Eg(MeV)= " << gammaEnergy
315 << " Ee(MeV)= " << kineticEnergy
316 << " DirE " << direction << " DirG " << gammaDirection
317 << G4endl;
318 */
319 // energy of primary
320 G4double finalE = kineticEnergy - gammaEnergy;
321
322 // stop tracking and create new secondary instead of primary
323 if(gammaEnergy > SecondaryThreshold()) {
326 G4DynamicParticle* el =
328 direction, finalE);
329 vdp->push_back(el);
330
331 // continue tracking
332 } else {
335 }
336}
@ fStopAndKill
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
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()
Definition: G4VEmModel.hh:508
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:286
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:459
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:557
void ProposeTrackStatus(G4TrackStatus status)
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double)
const G4ParticleDefinition * particle
G4ParticleChangeForLoss * fParticleChange

◆ SetBicubicInterpolationFlag()

void G4SeltzerBergerModel::SetBicubicInterpolationFlag ( G4bool  val)
inline

Definition at line 100 of file G4SeltzerBergerModel.hh.

101{
102 useBicubicInterpolation = val;
103}

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