Geant4  10.02.p03
G4eBremParametrizedModel Class Reference

#include <G4eBremParametrizedModel.hh>

Inheritance diagram for G4eBremParametrizedModel:
Collaboration diagram for G4eBremParametrizedModel:

Public Member Functions

 G4eBremParametrizedModel (const G4ParticleDefinition *p=0, const G4String &nam="eBremParam")
 
virtual ~G4eBremParametrizedModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
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)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int, const G4ParticleDefinition *, G4double)
 
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 MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
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 *> *)
 
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=0)
 
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 Attributes

G4NistManagernist
 
const G4ParticleDefinitionparticle
 
G4ParticleDefinitiontheGamma
 
G4ParticleChangeForLoss * fParticleChange
 
G4double minThreshold
 
G4double particleMass
 
G4double kinEnergy
 
G4double totalEnergy
 
G4double currentZ
 
G4double z13
 
G4double z23
 
G4double lnZ
 
G4double densityFactor
 
G4double densityCorr
 
G4double Fel
 
G4double Finel
 
G4double facFel
 
G4double facFinel
 
G4double fMax
 
G4double fCoulomb
 
G4bool isElectron
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChange * pParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Static Protected Attributes

static const G4double xgi [8]
 
static const G4double wgi [8]
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Private Member Functions

void InitialiseConstants ()
 
G4double ComputeBremLoss (G4double cutEnergy)
 
G4double ComputeXSectionPerAtom (G4double cutEnergy)
 
G4double ComputeDXSectionPerAtom (G4double gammaEnergy)
 
G4double ComputeParametrizedDXSectionPerAtom (G4double kineticEnergy, G4double gammaEnergy, G4double Z)
 
void SetParticle (const G4ParticleDefinition *p)
 
G4double ScreenFunction1 (G4double ScreenVariable)
 
G4double ScreenFunction2 (G4double ScreenVariable)
 
void SetCurrentElement (const G4double)
 
G4eBremParametrizedModeloperator= (const G4eBremParametrizedModel &right)
 
 G4eBremParametrizedModel (const G4eBremParametrizedModel &)
 

Private Attributes

G4double lowKinEnergy
 
G4double fMigdalConstant
 
G4double bremFactor
 
G4bool isInitialised
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLoss * GetParticleChangeForLoss ()
 
G4ParticleChangeForGamma * GetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Detailed Description

Definition at line 61 of file G4eBremParametrizedModel.hh.

Constructor & Destructor Documentation

◆ G4eBremParametrizedModel() [1/2]

G4eBremParametrizedModel::G4eBremParametrizedModel ( const G4ParticleDefinition p = 0,
const G4String nam = "eBremParam" 
)

Definition at line 100 of file G4eBremParametrizedModel.cc.

102  : G4VEmModel(nam),
103  particle(0),
104  isElectron(true),
107  isInitialised(false)
108 {
110 
111  minThreshold = 0.1*keV;
112  lowKinEnergy = 10.*MeV;
114 
116 
118 
120  = densityFactor = densityCorr = fMax = fCoulomb = 0.;
121 
123  if(p) { SetParticle(p); }
124 }
static const double MeV
Definition: G4SIunits.hh:211
static G4NistManager * Instance()
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:69
int fine_structure_const
Definition: hepunit.py:287
float electron_Compton_length
Definition: hepunit.py:289
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
int classic_electr_radius
Definition: hepunit.py:288
static const double pi
Definition: G4SIunits.hh:74
void SetParticle(const G4ParticleDefinition *p)
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:624
const G4ParticleDefinition * particle
static const double keV
Definition: G4SIunits.hh:213
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
Here is the call graph for this function:

◆ ~G4eBremParametrizedModel()

G4eBremParametrizedModel::~G4eBremParametrizedModel ( )
virtual

Definition at line 136 of file G4eBremParametrizedModel.cc.

137 {
138 }

◆ G4eBremParametrizedModel() [2/2]

G4eBremParametrizedModel::G4eBremParametrizedModel ( const G4eBremParametrizedModel )
private

Member Function Documentation

◆ ComputeBremLoss()

G4double G4eBremParametrizedModel::ComputeBremLoss ( G4double  cutEnergy)
private

Definition at line 234 of file G4eBremParametrizedModel.cc.

235 {
236  G4double loss = 0.0;
237 
238  // number of intervals and integration step
239  G4double vcut = cut/totalEnergy;
240  G4int n = (G4int)(20*vcut) + 3;
241  G4double delta = vcut/G4double(n);
242 
243  G4double e0 = 0.0;
244  G4double xs;
245 
246  // integration
247  for(G4int l=0; l<n; l++) {
248 
249  for(G4int i=0; i<8; i++) {
250 
251  G4double eg = (e0 + xgi[i]*delta)*totalEnergy;
252 
253  xs = ComputeDXSectionPerAtom(eg);
254 
255  loss += wgi[i]*xs/(1.0 + densityCorr/(eg*eg));
256  }
257  e0 += delta;
258  }
259 
260  loss *= delta*totalEnergy;
261 
262  return loss;
263 }
static const G4double xgi[8]
int G4int
Definition: G4Types.hh:78
Char_t n[5]
G4double ComputeDXSectionPerAtom(G4double gammaEnergy)
double G4double
Definition: G4Types.hh:76
static const G4double wgi[8]
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeCrossSectionPerAtom()

G4double G4eBremParametrizedModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition p,
G4double  tkin,
G4double  Z,
G4double  ,
G4double  cutEnergy,
G4double  maxEnergy = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 267 of file G4eBremParametrizedModel.cc.

273 {
274  if(!particle) { SetParticle(p); }
275  if(kineticEnergy < lowKinEnergy) { return 0.0; }
276  G4double cut = std::min(cutEnergy, kineticEnergy);
277  G4double tmax = std::min(maxEnergy, kineticEnergy);
278 
279  if(cut >= tmax) { return 0.0; }
280 
282 
283  G4double cross = ComputeXSectionPerAtom(cut);
284 
285  // allow partial integration
286  if(tmax < kinEnergy) { cross -= ComputeXSectionPerAtom(tmax); }
287 
288  cross *= Z*Z*bremFactor;
289 
290  return cross;
291 }
G4double ComputeXSectionPerAtom(G4double cutEnergy)
Float_t Z
void SetParticle(const G4ParticleDefinition *p)
const G4ParticleDefinition * particle
void SetCurrentElement(const G4double)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ ComputeDEDXPerVolume()

G4double G4eBremParametrizedModel::ComputeDEDXPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 201 of file G4eBremParametrizedModel.cc.

206 {
207  if(!particle) { SetParticle(p); }
208  if(kineticEnergy < lowKinEnergy) { return 0.0; }
209  G4double cut = std::min(cutEnergy, kineticEnergy);
210  if(cut == 0.0) { return 0.0; }
211 
212  SetupForMaterial(particle, material,kineticEnergy);
213 
214  const G4ElementVector* theElementVector = material->GetElementVector();
215  const G4double* theAtomicNumDensityVector = material->GetAtomicNumDensityVector();
216 
217  G4double dedx = 0.0;
218 
219  // loop for elements in the material
220  for (size_t i=0; i<material->GetNumberOfElements(); i++) {
221 
222  G4VEmModel::SetCurrentElement((*theElementVector)[i]);
223  SetCurrentElement((*theElementVector)[i]->GetZ());
224 
225  dedx += theAtomicNumDensityVector[i]*currentZ*currentZ*ComputeBremLoss(cut);
226  }
227  dedx *= bremFactor;
228 
229  return dedx;
230 }
std::vector< G4Element * > G4ElementVector
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:216
void SetParticle(const G4ParticleDefinition *p)
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
const G4ParticleDefinition * particle
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double)
void SetCurrentElement(const G4double)
double G4double
Definition: G4Types.hh:76
void SetCurrentElement(const G4Element *)
Definition: G4VEmModel.hh:459
G4double ComputeBremLoss(G4double cutEnergy)
Here is the call graph for this function:

◆ ComputeDXSectionPerAtom()

G4double G4eBremParametrizedModel::ComputeDXSectionPerAtom ( G4double  gammaEnergy)
private

Definition at line 441 of file G4eBremParametrizedModel.cc.

442 {
443 
444  if(gammaEnergy < 0.0) { return 0.0; }
445 
446  G4double y = gammaEnergy/totalEnergy;
447 
448  G4double main=0.;
449  //secondTerm=0.;
450 
451  // ** form factors complete screening case **
452  // only valid for high energies (and if LPM suppression does not play a role)
453  main = (3./4.*y*y - y + 1.) * ( (Fel-fCoulomb) + Finel/currentZ );
454  // secondTerm = (1.-y)/12.*(1.+1./currentZ);
455 
456  std::cout<<" F1(0) "<<ScreenFunction1(0.) <<std::endl;
457  std::cout<<" F1(0) "<<ScreenFunction2(0.) <<std::endl;
458  std::cout<<"Ekin = "<<kinEnergy<<std::endl;
459  std::cout<<"Z = "<<currentZ<<std::endl;
460  std::cout<<"main = "<<main<<std::endl;
461  std::cout<<" y = "<<y<<std::endl;
462  std::cout<<" Fel-fCoulomb "<< (Fel-fCoulomb) <<std::endl;
463 
465  std::cout<<"main2 = "<<main2<<std::endl;
466  std::cout<<"main2tot = "<<main2 * ( (Fel-fCoulomb) + Finel/currentZ )/(Fel-fCoulomb);
467 
468  G4double cross = main2; //main+secondTerm;
469  return cross;
470 }
G4double ScreenFunction2(G4double ScreenVariable)
int main(int argc, char **argv)
Definition: genwindef.cpp:30
G4double ScreenFunction1(G4double ScreenVariable)
Double_t y
G4double ComputeParametrizedDXSectionPerAtom(G4double kineticEnergy, G4double gammaEnergy, G4double Z)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeParametrizedDXSectionPerAtom()

G4double G4eBremParametrizedModel::ComputeParametrizedDXSectionPerAtom ( G4double  kineticEnergy,
G4double  gammaEnergy,
G4double  Z 
)
private

Definition at line 363 of file G4eBremParametrizedModel.cc.

366 {
368  G4double FZ = lnZ* (4.- 0.55*lnZ);
369  G4double Z3 = z13;
370  G4double ZZ = z13*nist->GetZ13(G4lrint(Z)+1);
371 
372  totalEnergy = kineticEnergy + electron_mass_c2;
373 
374  // G4double x, epsil, greject, migdal, grejmax, q;
375  G4double epsil, greject;
376  G4double U = G4Log(kineticEnergy/electron_mass_c2);
377  G4double U2 = U*U;
378 
379  // precalculated parameters
380  G4double ah, bh;
381 
382  if (kineticEnergy > tlow) {
383 
384  G4double ah1 = ah10 + ZZ* (ah11 + ZZ* ah12);
385  G4double ah2 = ah20 + ZZ* (ah21 + ZZ* ah22);
386  G4double ah3 = ah30 + ZZ* (ah31 + ZZ* ah32);
387 
388  G4double bh1 = bh10 + ZZ* (bh11 + ZZ* bh12);
389  G4double bh2 = bh20 + ZZ* (bh21 + ZZ* bh22);
390  G4double bh3 = bh30 + ZZ* (bh31 + ZZ* bh32);
391 
392  ah = 1. + (ah1*U2 + ah2*U + ah3) / (U2*U);
393  bh = 0.75 + (bh1*U2 + bh2*U + bh3) / (U2*U);
394 
395  // limit of the screening variable
396  G4double screenfac =
397  136.*electron_mass_c2/(Z3*totalEnergy);
398 
399  epsil = gammaEnergy/totalEnergy; // epsil = x*kineticEnergy/totalEnergy;
400  G4double screenvar = screenfac*epsil/(1.0-epsil);
401  G4double F1 = max(ScreenFunction1(screenvar) - FZ ,0.);
402  G4double F2 = max(ScreenFunction2(screenvar) - FZ ,0.);
403 
404 
405  greject = (F1 - epsil* (ah*F1 - bh*epsil*F2))/8.; // 1./(42.392 - FZ);
406 
407  std::cout << " yy = "<<epsil<<std::endl;
408  std::cout << " F1/(...) "<<F1/(42.392 - FZ)<<std::endl;
409  std::cout << " F2/(...) "<<F2/(42.392 - FZ)<<std::endl;
410  std::cout << " (42.392 - FZ) " << (42.392 - FZ) <<std::endl;
411 
412  } else {
413 
414  G4double al0 = al00 + ZZ* (al01 + ZZ* al02);
415  G4double al1 = al10 + ZZ* (al11 + ZZ* al12);
416  G4double al2 = al20 + ZZ* (al21 + ZZ* al22);
417 
418  G4double bl0 = bl00 + ZZ* (bl01 + ZZ* bl02);
419  G4double bl1 = bl10 + ZZ* (bl11 + ZZ* bl12);
420  G4double bl2 = bl20 + ZZ* (bl21 + ZZ* bl22);
421 
422  ah = al0 + al1*U + al2*U2;
423  bh = bl0 + bl1*U + bl2*U2;
424 
425  G4double x=gammaEnergy/kineticEnergy;
426  greject=(1. + x* (ah + bh*x));
427 
428  /*
429  // Compute the maximum of the rejection function
430  grejmax = max(1. + xmin* (ah + bh*xmin), 1.+ah+bh);
431  G4double xm = -ah/(2.*bh);
432  if ( xmin < xm && xm < xmax) grejmax = max(grejmax, 1.+ xm* (ah + bh*xm));
433  */
434  }
435 
436  return greject;
437 }
static const G4double bh12
G4double ScreenFunction2(G4double ScreenVariable)
static const G4double al00
static const G4double bh10
static const G4double bh32
static const G4double al12
static const G4double bl02
static const G4double ah11
static const G4double bl11
G4double GetZ13(G4double Z)
static const G4double bh22
G4double ScreenFunction1(G4double ScreenVariable)
static const G4double ah31
static const G4double al01
static const G4double bl21
static const G4double ah22
static const G4double al11
static const G4double bh20
static const G4double bl00
static const G4double al22
static const G4double bl01
static const G4double al20
static const G4double al02
static const G4double bh11
Float_t Z
static const G4double tlow
float electron_mass_c2
Definition: hepunit.py:274
static const G4double bh30
G4double G4Log(G4double x)
Definition: G4Log.hh:230
static const G4double ah12
static const G4double ah32
static const G4double al10
static const G4double ah20
int G4lrint(double ad)
Definition: templates.hh:163
static const G4double bh21
static const G4double bl22
static const G4double bh31
static const G4double bl10
static const G4double bl12
void SetCurrentElement(const G4double)
double G4double
Definition: G4Types.hh:76
static const G4double al21
static const G4double ah21
static const G4double ah30
static const G4double bl20
static const G4double ah10
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeXSectionPerAtom()

G4double G4eBremParametrizedModel::ComputeXSectionPerAtom ( G4double  cutEnergy)
private

Definition at line 295 of file G4eBremParametrizedModel.cc.

296 {
297  G4double cross = 0.0;
298 
299  // number of intervals and integration step
300  G4double vcut = G4Log(cut/totalEnergy);
302  G4int n = (G4int)(0.45*(vmax - vcut)) + 4;
303  // n=1; // integration test
304  G4double delta = (vmax - vcut)/G4double(n);
305 
306  G4double e0 = vcut;
307  G4double xs;
308 
309  // integration
310  for(G4int l=0; l<n; l++) {
311 
312  for(G4int i=0; i<8; i++) {
313 
314  G4double eg = G4Exp(e0 + xgi[i]*delta)*totalEnergy;
315 
316  xs = ComputeDXSectionPerAtom(eg);
317 
318  cross += wgi[i]*xs/(1.0 + densityCorr/(eg*eg));
319  }
320  e0 += delta;
321  }
322 
323  cross *= delta;
324 
325  return cross;
326 }
static const G4double xgi[8]
int G4int
Definition: G4Types.hh:78
Char_t n[5]
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4double ComputeDXSectionPerAtom(G4double gammaEnergy)
double G4double
Definition: G4Types.hh:76
static const G4double wgi[8]
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 175 of file G4eBremParametrizedModel.cc.

177 {
178  if(p) { SetParticle(p); }
179 
181 
182  currentZ = 0.;
183 
184  if(IsMaster()) { InitialiseElementSelectors(p, cuts); }
185 
186  if(isInitialised) { return; }
188  isInitialised = true;
189 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:120
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148
G4ParticleChangeForLoss * fParticleChange
G4bool IsMaster() const
Definition: G4VEmModel.hh:718
void SetParticle(const G4ParticleDefinition *p)
Here is the call graph for this function:

◆ InitialiseConstants()

void G4eBremParametrizedModel::InitialiseConstants ( )
private

Definition at line 128 of file G4eBremParametrizedModel.cc.

129 {
130  facFel = G4Log(184.15);
131  facFinel = G4Log(1194.);
132 }
G4double G4Log(G4double x)
Definition: G4Log.hh:230
Here is the call graph for this function:
Here is the caller graph for this function:

◆ InitialiseLocal()

void G4eBremParametrizedModel::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 193 of file G4eBremParametrizedModel.cc.

195 {
197 }
void SetElementSelectors(std::vector< G4EmElementSelector *> *)
Definition: G4VEmModel.hh:810
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:802
Here is the call graph for this function:

◆ MinEnergyCut()

G4double G4eBremParametrizedModel::MinEnergyCut ( const G4ParticleDefinition ,
const G4MaterialCutsCouple  
)
virtual

Reimplemented from G4VEmModel.

Definition at line 152 of file G4eBremParametrizedModel.cc.

154 {
155  return minThreshold;
156 }

◆ operator=()

G4eBremParametrizedModel& G4eBremParametrizedModel::operator= ( const G4eBremParametrizedModel right)
private

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 474 of file G4eBremParametrizedModel.cc.

480 {
481  G4double kineticEnergy = dp->GetKineticEnergy();
482  if(kineticEnergy < lowKinEnergy) { return; }
483  G4double cut = std::min(cutEnergy, kineticEnergy);
484  G4double emax = std::min(maxEnergy, kineticEnergy);
485  if(cut >= emax) { return; }
486 
487  SetupForMaterial(particle, couple->GetMaterial(),kineticEnergy);
488 
489  const G4Element* elm =
490  SelectRandomAtom(couple,particle,kineticEnergy,cut,emax);
491  SetCurrentElement(elm->GetZ());
492 
493  kinEnergy = kineticEnergy;
494  totalEnergy = kineticEnergy + particleMass;
496 
497  G4double xmin = G4Log(cut*cut + densityCorr);
498  G4double xmax = G4Log(emax*emax + densityCorr);
499  G4double gammaEnergy, f, x;
500 
501  CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
502 
503  do {
504  x = G4Exp(xmin + rndmEngine->flat()*(xmax - xmin)) - densityCorr;
505  if(x < 0.0) x = 0.0;
506  gammaEnergy = sqrt(x);
507  f = ComputeDXSectionPerAtom(gammaEnergy);
508 
509  if ( f > fMax ) {
510  G4cout << "### G4eBremParametrizedModel Warning: Majoranta exceeded! "
511  << f << " > " << fMax
512  << " Egamma(MeV)= " << gammaEnergy
513  << " E(mEV)= " << kineticEnergy
514  << G4endl;
515  }
516 
517  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
518  } while (f < fMax*rndmEngine->flat());
519 
520  //
521  // angles of the emitted gamma. ( Z - axis along the parent particle)
522  // use general interface
523  //
524  G4ThreeVector gammaDirection =
525  GetAngularDistribution()->SampleDirection(dp, totalEnergy-gammaEnergy,
526  G4lrint(currentZ),
527  couple->GetMaterial());
528 
529  // create G4DynamicParticle object for the Gamma
530  G4DynamicParticle* gamma = new G4DynamicParticle(theGamma,gammaDirection,
531  gammaEnergy);
532  vdp->push_back(gamma);
533 
534  G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2));
535  G4ThreeVector direction = (totMomentum*dp->GetMomentumDirection()
536  - gammaEnergy*gammaDirection).unit();
537 
538  // energy of primary
539  G4double finalE = kineticEnergy - gammaEnergy;
540 
541  // stop tracking and create new secondary instead of primary
542  if(gammaEnergy > SecondaryThreshold()) {
543  fParticleChange->ProposeTrackStatus(fStopAndKill);
544  fParticleChange->SetProposedKineticEnergy(0.0);
545  G4DynamicParticle* el =
546  new G4DynamicParticle(const_cast<G4ParticleDefinition*>(particle),
547  direction, finalE);
548  vdp->push_back(el);
549 
550  // continue tracking
551  } else {
552  fParticleChange->SetProposedMomentumDirection(direction);
553  fParticleChange->SetProposedKineticEnergy(finalE);
554  }
555 }
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:669
const G4Material * GetMaterial() const
G4ParticleChangeForLoss * fParticleChange
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:617
virtual double flat()=0
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
float electron_mass_c2
Definition: hepunit.py:274
double flat()
Definition: G4AblaRandom.cc:47
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
int G4lrint(double ad)
Definition: templates.hh:163
const G4ThreeVector & GetMomentumDirection() const
const G4ParticleDefinition * particle
#define G4endl
Definition: G4ios.hh:61
G4double ComputeDXSectionPerAtom(G4double gammaEnergy)
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double)
void SetCurrentElement(const G4double)
double G4double
Definition: G4Types.hh:76
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:544
Here is the call graph for this function:

◆ ScreenFunction1()

G4double G4eBremParametrizedModel::ScreenFunction1 ( G4double  ScreenVariable)
private

Definition at line 332 of file G4eBremParametrizedModel.cc.

333 {
334  G4double screenVal;
335 
336  if (ScreenVariable > 1.)
337  screenVal = 42.24 - 8.368*G4Log(ScreenVariable+0.952);
338  else
339  screenVal = 42.392 - ScreenVariable* (7.796 - 1.961*ScreenVariable);
340 
341  return screenVal;
342 }
G4double G4Log(G4double x)
Definition: G4Log.hh:230
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ScreenFunction2()

G4double G4eBremParametrizedModel::ScreenFunction2 ( G4double  ScreenVariable)
private

Definition at line 348 of file G4eBremParametrizedModel.cc.

349 {
350  G4double screenVal;
351 
352  if (ScreenVariable > 1.)
353  screenVal = 42.24 - 8.368*G4Log(ScreenVariable+0.952);
354  else
355  screenVal = 41.734 - ScreenVariable* (6.484 - 1.250*ScreenVariable);
356 
357  return screenVal;
358 }
G4double G4Log(G4double x)
Definition: G4Log.hh:230
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetCurrentElement()

void G4eBremParametrizedModel::SetCurrentElement ( const G4double  Z)
inlineprivate

Definition at line 164 of file G4eBremParametrizedModel.hh.

165 {
166  std::cout<<"SetCurrentElement Z="<<Z<<std::endl;
167  if(Z != currentZ) {
168  currentZ = Z;
169 
170  G4int iz = G4int(Z);
171  z13 = nist->GetZ13(iz);
172  z23 = z13*z13;
173  lnZ = nist->GetLOGZ(iz);
174 
175  Fel = facFel - lnZ/3. ;
176  Finel = facFinel - 2.*lnZ/3. ;
177 
179  fMax = Fel-fCoulomb + Finel/currentZ + (1.+1./currentZ)/12.;
180 
181  }
182 }
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:467
G4double GetZ13(G4double Z)
int G4int
Definition: G4Types.hh:78
Float_t Z
G4double iz
Definition: TRTMaterials.hh:39
G4double GetfCoulomb() const
Definition: G4Element.hh:190
G4double GetLOGZ(G4int Z)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetParticle()

void G4eBremParametrizedModel::SetParticle ( const G4ParticleDefinition p)
private

Definition at line 142 of file G4eBremParametrizedModel.cc.

143 {
144  particle = p;
145  particleMass = p->GetPDGMass();
146  if(p == G4Electron::Electron()) { isElectron = true; }
147  else { isElectron = false;}
148 }
const G4ParticleDefinition * particle
static G4Electron * Electron()
Definition: G4Electron.cc:94
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetupForMaterial()

void G4eBremParametrizedModel::SetupForMaterial ( const G4ParticleDefinition ,
const G4Material mat,
G4double  kineticEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 160 of file G4eBremParametrizedModel.cc.

163 {
165 
166  // calculate threshold for density effect
167  kinEnergy = kineticEnergy;
168  totalEnergy = kineticEnergy + particleMass;
170 }
G4double GetElectronDensity() const
Definition: G4Material.hh:217
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ bremFactor

G4double G4eBremParametrizedModel::bremFactor
private

Definition at line 157 of file G4eBremParametrizedModel.hh.

◆ currentZ

G4double G4eBremParametrizedModel::currentZ
protected

Definition at line 142 of file G4eBremParametrizedModel.hh.

◆ densityCorr

G4double G4eBremParametrizedModel::densityCorr
protected

Definition at line 145 of file G4eBremParametrizedModel.hh.

◆ densityFactor

G4double G4eBremParametrizedModel::densityFactor
protected

Definition at line 144 of file G4eBremParametrizedModel.hh.

◆ facFel

G4double G4eBremParametrizedModel::facFel
protected

Definition at line 147 of file G4eBremParametrizedModel.hh.

◆ facFinel

G4double G4eBremParametrizedModel::facFinel
protected

Definition at line 147 of file G4eBremParametrizedModel.hh.

◆ fCoulomb

G4double G4eBremParametrizedModel::fCoulomb
protected

Definition at line 148 of file G4eBremParametrizedModel.hh.

◆ Fel

G4double G4eBremParametrizedModel::Fel
protected

Definition at line 146 of file G4eBremParametrizedModel.hh.

◆ Finel

G4double G4eBremParametrizedModel::Finel
protected

Definition at line 146 of file G4eBremParametrizedModel.hh.

◆ fMax

G4double G4eBremParametrizedModel::fMax
protected

Definition at line 148 of file G4eBremParametrizedModel.hh.

◆ fMigdalConstant

G4double G4eBremParametrizedModel::fMigdalConstant
private

Definition at line 156 of file G4eBremParametrizedModel.hh.

◆ fParticleChange

G4ParticleChangeForLoss* G4eBremParametrizedModel::fParticleChange
protected

Definition at line 132 of file G4eBremParametrizedModel.hh.

◆ isElectron

G4bool G4eBremParametrizedModel::isElectron
protected

Definition at line 150 of file G4eBremParametrizedModel.hh.

◆ isInitialised

G4bool G4eBremParametrizedModel::isInitialised
private

Definition at line 159 of file G4eBremParametrizedModel.hh.

◆ kinEnergy

G4double G4eBremParametrizedModel::kinEnergy
protected

Definition at line 140 of file G4eBremParametrizedModel.hh.

◆ lnZ

G4double G4eBremParametrizedModel::lnZ
protected

Definition at line 143 of file G4eBremParametrizedModel.hh.

◆ lowKinEnergy

G4double G4eBremParametrizedModel::lowKinEnergy
private

Definition at line 155 of file G4eBremParametrizedModel.hh.

◆ minThreshold

G4double G4eBremParametrizedModel::minThreshold
protected

Definition at line 136 of file G4eBremParametrizedModel.hh.

◆ nist

G4NistManager* G4eBremParametrizedModel::nist
protected

Definition at line 129 of file G4eBremParametrizedModel.hh.

◆ particle

const G4ParticleDefinition* G4eBremParametrizedModel::particle
protected

Definition at line 130 of file G4eBremParametrizedModel.hh.

◆ particleMass

G4double G4eBremParametrizedModel::particleMass
protected

Definition at line 139 of file G4eBremParametrizedModel.hh.

◆ theGamma

G4ParticleDefinition* G4eBremParametrizedModel::theGamma
protected

Definition at line 131 of file G4eBremParametrizedModel.hh.

◆ totalEnergy

G4double G4eBremParametrizedModel::totalEnergy
protected

Definition at line 141 of file G4eBremParametrizedModel.hh.

◆ wgi

const G4double G4eBremParametrizedModel::wgi
staticprotected
Initial value:
={ 0.0506, 0.1112, 0.1569, 0.1813,
0.1813, 0.1569, 0.1112, 0.0506 }

Definition at line 134 of file G4eBremParametrizedModel.hh.

◆ xgi

const G4double G4eBremParametrizedModel::xgi
staticprotected
Initial value:
={ 0.0199, 0.1017, 0.2372, 0.4083,
0.5917, 0.7628, 0.8983, 0.9801 }

Definition at line 134 of file G4eBremParametrizedModel.hh.

◆ z13

G4double G4eBremParametrizedModel::z13
protected

Definition at line 143 of file G4eBremParametrizedModel.hh.

◆ z23

G4double G4eBremParametrizedModel::z23
protected

Definition at line 143 of file G4eBremParametrizedModel.hh.


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