Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4SeltzerBergerModel Class Reference

#include <G4SeltzerBergerModel.hh>

Inheritance diagram for G4SeltzerBergerModel:
Collaboration diagram for G4SeltzerBergerModel:

Public Member Functions

 G4SeltzerBergerModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="eBremSB")
 
virtual ~G4SeltzerBergerModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z) override
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy) override
 
void SetBicubicInterpolationFlag (G4bool)
 
- Public Member Functions inherited from G4eBremsstrahlungRelModel
 G4eBremsstrahlungRelModel (const G4ParticleDefinition *p=0, const G4String &nam="eBremLPM")
 
virtual ~G4eBremsstrahlungRelModel ()
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double tkin, G4double Z, G4double, G4double cutEnergy, G4double maxEnergy=DBL_MAX) override
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double) override
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut) override
 
void SetLPMconstant (G4double val)
 
G4double LPMconstant () const
 
void SetLowestKinEnergy (G4double)
 
G4double LowestKinEnergy () const
 
- 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 *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector
< G4EmElementSelector * > * 
GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual 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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
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
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
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 SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 

Protected Member Functions

virtual G4double ComputeDXSectionPerAtom (G4double gammaEnergy) override
 
virtual G4String DirectoryPath () const
 
- Protected Member Functions inherited from G4eBremsstrahlungRelModel
void SetCurrentElement (G4int)
 
- 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 densityFactor
 
G4double densityCorr
 
G4int currentZ
 
G4bool isElectron
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 61 of file G4SeltzerBergerModel.hh.

Constructor & Destructor Documentation

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

Definition at line 78 of file G4SeltzerBergerModel.cc.

80  : G4eBremsstrahlungRelModel(p,nam),useBicubicInterpolation(false)
81 {
84  SetLPMFlag(false);
85  nwarn = 0;
86  idx = idy = 0;
87 }
G4eBremsstrahlungRelModel(const G4ParticleDefinition *p=0, const G4String &nam="eBremLPM")
void SetLPMFlag(G4bool val)
Definition: G4VEmModel.hh:773
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:731
static constexpr double keV
Definition: G4SIunits.hh:216

Here is the call graph for this function:

G4SeltzerBergerModel::~G4SeltzerBergerModel ( )
virtual

Definition at line 91 of file G4SeltzerBergerModel.cc.

92 {
93  if(IsMaster()) {
94  for(size_t i=0; i<101; ++i) {
95  if(dataSB[i]) {
96  delete dataSB[i];
97  dataSB[i] = nullptr;
98  }
99  }
100  }
101 }
G4bool IsMaster() const
Definition: G4VEmModel.hh:717

Here is the call graph for this function:

Member Function Documentation

G4double G4SeltzerBergerModel::ComputeDXSectionPerAtom ( G4double  gammaEnergy)
overrideprotectedvirtual

Reimplemented from G4eBremsstrahlungRelModel.

Definition at line 189 of file G4SeltzerBergerModel.cc.

190 {
191 
192  if(gammaEnergy < 0.0 || kinEnergy <= 0.0) { return 0.0; }
193  G4double x = gammaEnergy/kinEnergy;
195  G4int Z = G4lrint(currentZ);
196 
197  //G4cout << "G4SeltzerBergerModel::ComputeDXSectionPerAtom Z= " << Z
198  // << " x= " << x << " y= " << y << " " << dataSB[Z] << G4endl;
199  if(nullptr == dataSB[Z]) { InitialiseForElement(0, Z); }
200  /*
201  G4ExceptionDescription ed;
202  ed << "Bremsstrahlung data for Z= " << Z
203  << " are not initialized!";
204  G4Exception("G4SeltzerBergerModel::ComputeDXSectionPerAtom()","em0005",
205  FatalException, ed,
206  "G4LEDATA version should be G4EMLOW6.23 or later.");
207  }
208  */
209  G4double invb2 =
211  G4double cross = dataSB[Z]->Value(x,y,idx,idy)*invb2*millibarn/bremFactor;
212 
213  if(!isElectron) {
214  G4double invbeta1 = sqrt(invb2);
215  G4double e2 = kinEnergy - gammaEnergy;
216  if(e2 > 0.0) {
217  G4double invbeta2 = (e2 + particleMass)/sqrt(e2*(e2 + 2*particleMass));
219  G4double xxx = alpha*currentZ*(invbeta1 - invbeta2);
220  if(xxx < expnumlim) { cross = 0.0; }
221  else { cross *= G4Exp(xxx); }
222  } else {
223  cross = 0.0;
224  }
225  }
226 
227  return cross;
228 }
int G4int
Definition: G4Types.hh:78
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z) override
G4double Value(G4double x, G4double y, size_t &lastidx, size_t &lastidy) const
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
int G4lrint(double ad)
Definition: templates.hh:163
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76
static constexpr double fine_structure_const
static const G4double alpha
static constexpr double twopi
Definition: SystemOfUnits.h:55
static constexpr double millibarn
Definition: G4SIunits.hh:106

Here is the call graph for this function:

G4String G4SeltzerBergerModel::DirectoryPath ( ) const
protectedvirtual

Definition at line 134 of file G4SeltzerBergerModel.cc.

135 {
136  return "/brem_SB/br";
137 }
void G4SeltzerBergerModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector cuts 
)
overridevirtual

Reimplemented from G4eBremsstrahlungRelModel.

Reimplemented in G4ePolarizedBremsstrahlungModel.

Definition at line 105 of file G4SeltzerBergerModel.cc.

107 {
108  // Access to elements
109  if(IsMaster()) {
110 
111  // check environment variable
112  // Build the complete string identifying the file with the data set
113  char* path = getenv("G4LEDATA");
114 
115  const G4ElementTable* theElmTable = G4Element::GetElementTable();
116  size_t numOfElm = G4Element::GetNumberOfElements();
117  if(numOfElm > 0) {
118  for(size_t i=0; i<numOfElm; ++i) {
119  G4int Z = G4lrint(((*theElmTable)[i])->GetZ());
120  if(Z < 1) { Z = 1; }
121  else if(Z > 100) { Z = 100; }
122  //G4cout << "Z= " << Z << G4endl;
123  // Initialisation
124  if(nullptr == dataSB[Z]) { ReadData(Z, path); }
125  }
126  }
127  }
128 
130 }
G4bool IsMaster() const
Definition: G4VEmModel.hh:717
int G4int
Definition: G4Types.hh:78
static size_t GetNumberOfElements()
Definition: G4Element.cc:405
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
int G4lrint(double ad)
Definition: templates.hh:163
std::vector< G4Element * > G4ElementTable
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:398

Here is the call graph for this function:

Here is the caller graph for this function:

void G4SeltzerBergerModel::InitialiseForElement ( const G4ParticleDefinition ,
G4int  Z 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 379 of file G4SeltzerBergerModel.cc.

381 {
382  G4AutoLock l(&SeltzerBergerModelMutex);
383  // G4cout << "G4SeltzerBergerModel::InitialiseForElement Z= " << Z << G4endl;
384  if(nullptr == dataSB[Z]) { ReadData(Z); }
385 }

Here is the caller graph for this function:

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

Reimplemented from G4eBremsstrahlungRelModel.

Reimplemented in G4ePolarizedBremsstrahlungModel.

Definition at line 233 of file G4SeltzerBergerModel.cc.

238 {
239  G4double kineticEnergy = dp->GetKineticEnergy();
240  G4double cut = std::min(cutEnergy, kineticEnergy);
241  G4double emax = std::min(maxEnergy, kineticEnergy);
242  if(cut >= emax) { return; }
243 
244  SetupForMaterial(particle, couple->GetMaterial(), kineticEnergy);
245 
246  const G4Element* elm =
247  SelectRandomAtom(couple,particle,kineticEnergy,cut,emax);
248  SetCurrentElement(elm->GetZasInt());
249 
250  totalEnergy = kineticEnergy + particleMass;
251  densityCorr = densityFactor*totalEnergy*totalEnergy;
252  G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2));
253  /*
254  G4cout << "G4SeltzerBergerModel::SampleSecondaries E(MeV)= "
255  << kineticEnergy/MeV
256  << " Z= " << Z << " cut(MeV)= " << cut/MeV
257  << " emax(MeV)= " << emax/MeV << " corr= " << densityCorr << G4endl;
258  */
259  G4double xmin = G4Log(cut*cut + densityCorr);
260  G4double xmax = G4Log(emax*emax + densityCorr);
261  G4double y = G4Log(kineticEnergy/MeV);
262 
263  G4double gammaEnergy, v;
264 
265  // majoranta
266  G4double x0 = cut/kineticEnergy;
267  G4double vmax;
268  if(currentZ <= 92) {
269  vmax = dataSB[currentZ]->Value(x0, y, idx, idy)*1.02;
270  } else {
271  idx = idy = 0;
272  vmax = dataSB[currentZ]->Value(x0, y, idx, idy)*1.2;
273  }
274 
275  static const G4double epeaklimit= 300*CLHEP::MeV;
276  static const G4double elowlimit = 20*CLHEP::keV;
277 
278  // majoranta corrected for e-
279  if(isElectron && x0 < 0.97 &&
280  ((kineticEnergy > epeaklimit) || (kineticEnergy < elowlimit))) {
281  G4double ylim = std::min(ylimit[currentZ],1.1*dataSB[currentZ]->Value(0.97,y,idx,idy));
282  if(ylim > vmax) { vmax = ylim; }
283  }
284  if(x0 < 0.05) { vmax *= 1.2; }
285 
286  //G4cout<<"y= "<<y<<" xmin= "<<xmin<<" xmax= "<<xmax
287  //<<" vmax= "<<vmax<<G4endl;
288  static const G4int ncountmax = 100;
289  CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
290  G4double rndm[2];
291 
292  for(G4int nn=0; nn<ncountmax; ++nn) {
293  rndmEngine->flatArray(2, rndm);
294  G4double x = G4Exp(xmin + rndm[0]*(xmax - xmin)) - densityCorr;
295  if(x < 0.0) { x = 0.0; }
296  gammaEnergy = sqrt(x);
297  G4double x1 = gammaEnergy/kineticEnergy;
298  v = dataSB[currentZ]->Value(x1, y, idx, idy);
299 
300  // correction for positrons
301  if(!isElectron) {
302  G4double e1 = kineticEnergy - cut;
303  G4double invbeta1 = (e1 + particleMass)/sqrt(e1*(e1 + 2*particleMass));
304  G4double e2 = kineticEnergy - gammaEnergy;
305  G4double invbeta2 = (e2 + particleMass)/sqrt(e2*(e2 + 2*particleMass));
306  G4double xxx = twopi*fine_structure_const*currentZ*(invbeta1 - invbeta2);
307 
308  if(xxx < expnumlim) { v = 0.0; }
309  else { v *= G4Exp(xxx); }
310  }
311 
312  if (v > 1.05*vmax && nwarn < 5) {
313  ++nwarn;
315  ed << "### G4SeltzerBergerModel Warning: Majoranta exceeded! "
316  << v << " > " << vmax << " by " << v/vmax
317  << " Niter= " << nn
318  << " Egamma(MeV)= " << gammaEnergy
319  << " Ee(MeV)= " << kineticEnergy
320  << " Z= " << currentZ << " " << particle->GetParticleName();
321 
322  if ( 20 == nwarn ) {
323  ed << "\n ### G4SeltzerBergerModel Warnings stopped";
324  }
325  G4Exception("G4SeltzerBergerModel::SampleScattering","em0044",
326  JustWarning, ed,"");
327 
328  }
329  if(v >= vmax*rndm[1]) { break; }
330  }
331 
332  //
333  // angles of the emitted gamma. ( Z - axis along the parent particle)
334  // use general interface
335  //
336 
337  G4ThreeVector gammaDirection =
338  GetAngularDistribution()->SampleDirection(dp, totalEnergy-gammaEnergy,
339  currentZ, couple->GetMaterial());
340 
341  // create G4DynamicParticle object for the Gamma
342  G4DynamicParticle* gamma =
343  new G4DynamicParticle(theGamma,gammaDirection,gammaEnergy);
344  vdp->push_back(gamma);
345 
346  G4ThreeVector direction = (totMomentum*dp->GetMomentumDirection()
347  - gammaEnergy*gammaDirection).unit();
348 
349  /*
350  G4cout << "### G4SBModel: v= "
351  << " Eg(MeV)= " << gammaEnergy
352  << " Ee(MeV)= " << kineticEnergy
353  << " DirE " << direction << " DirG " << gammaDirection
354  << G4endl;
355  */
356  // energy of primary
357  G4double finalE = kineticEnergy - gammaEnergy;
358 
359  // stop tracking and create new secondary instead of primary
360  if(gammaEnergy > SecondaryThreshold()) {
363  G4DynamicParticle* el =
364  new G4DynamicParticle(const_cast<G4ParticleDefinition*>(particle),
365  direction, finalE);
366  vdp->push_back(el);
367 
368  // continue tracking
369  } else {
372  }
373 }
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:668
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double GetKineticEnergy() const
static constexpr double keV
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:616
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
const G4ParticleDefinition * particle
static constexpr double twopi
Definition: G4SIunits.hh:76
static constexpr double electron_mass_c2
G4double Value(G4double x, G4double y, size_t &lastidx, size_t &lastidy) const
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
const G4ThreeVector & GetMomentumDirection() const
static constexpr double MeV
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static const G4double emax
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double) override
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static constexpr double MeV
Definition: G4SIunits.hh:214
static const G4double epeaklimit
G4ParticleChangeForLoss * fParticleChange
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
static const G4double elowlimit
static constexpr double fine_structure_const
virtual void flatArray(const int size, double *vect)=0
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:541
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:377
const G4Material * GetMaterial() const

Here is the call graph for this function:

Here is the caller graph for this function:

void G4SeltzerBergerModel::SetBicubicInterpolationFlag ( G4bool  val)
inline

Definition at line 106 of file G4SeltzerBergerModel.hh.

107 {
108  useBicubicInterpolation = val;
109 }

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