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

#include <G4eBremParametrizedModel.hh>

+ Inheritance 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 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 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 Attributes

G4NistManagernist
 
const G4ParticleDefinitionparticle
 
G4ParticleDefinitiontheGamma
 
G4ParticleChangeForLossfParticleChange
 
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
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

Static Protected Attributes

static const G4double xgi [8]
 
static const G4double wgi [8]
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
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()

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

Definition at line 74 of file G4eBremParametrizedModel.cc.

76 : G4VEmModel(nam),
77 particle(0),
78 isElectron(true),
79 fMigdalConstant(classic_electr_radius*electron_Compton_length*electron_Compton_length*4.0*pi),
80 bremFactor(fine_structure_const*classic_electr_radius*classic_electr_radius*16./3.),
81 isInitialised(false)
82{
84
85 minThreshold = 0.1*keV;
86 lowKinEnergy = 10.*MeV;
87 SetLowEnergyLimit(lowKinEnergy);
88
90
92
95
96 InitialiseConstants();
97 if(p) { SetParticle(p); }
98}
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4NistManager * Instance()
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:515
const G4ParticleDefinition * particle

◆ ~G4eBremParametrizedModel()

G4eBremParametrizedModel::~G4eBremParametrizedModel ( )
virtual

Definition at line 110 of file G4eBremParametrizedModel.cc.

111{
112}

Member Function Documentation

◆ 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 234 of file G4eBremParametrizedModel.cc.

240{
241 if(!particle) { SetParticle(p); }
242 if(kineticEnergy < lowKinEnergy) { return 0.0; }
243 G4double cut = std::min(cutEnergy, kineticEnergy);
244 G4double tmax = std::min(maxEnergy, kineticEnergy);
245
246 if(cut >= tmax) { return 0.0; }
247
248 SetCurrentElement(Z);
249
250 G4double cross = ComputeXSectionPerAtom(cut);
251
252 // allow partial integration
253 if(tmax < kinEnergy) { cross -= ComputeXSectionPerAtom(tmax); }
254
255 cross *= Z*Z*bremFactor;
256
257 return cross;
258}
double G4double
Definition: G4Types.hh:64

◆ ComputeDEDXPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 167 of file G4eBremParametrizedModel.cc.

172{
173 if(!particle) { SetParticle(p); }
174 if(kineticEnergy < lowKinEnergy) { return 0.0; }
175 G4double cut = std::min(cutEnergy, kineticEnergy);
176 if(cut == 0.0) { return 0.0; }
177
178 SetupForMaterial(particle, material,kineticEnergy);
179
180 const G4ElementVector* theElementVector = material->GetElementVector();
181 const G4double* theAtomicNumDensityVector = material->GetAtomicNumDensityVector();
182
183 G4double dedx = 0.0;
184
185 // loop for elements in the material
186 for (size_t i=0; i<material->GetNumberOfElements(); i++) {
187
188 G4VEmModel::SetCurrentElement((*theElementVector)[i]);
189 SetCurrentElement((*theElementVector)[i]->GetZ());
190
191 dedx += theAtomicNumDensityVector[i]*currentZ*currentZ*ComputeBremLoss(cut);
192 }
193 dedx *= bremFactor;
194
195
196 return dedx;
197}
std::vector< G4Element * > G4ElementVector
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:189
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:215
void SetCurrentElement(const G4Element *)
Definition: G4VEmModel.hh:384
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double)

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 149 of file G4eBremParametrizedModel.cc.

151{
152 if(p) { SetParticle(p); }
153
154 lowKinEnergy = LowEnergyLimit();
155
156 currentZ = 0.;
157
159
160 if(isInitialised) { return; }
162 isInitialised = true;
163}
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:123
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:95
G4ParticleChangeForLoss * fParticleChange

◆ MinEnergyCut()

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

Definition at line 126 of file G4eBremParametrizedModel.cc.

128{
129 return minThreshold;
130}

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 466 of file G4eBremParametrizedModel.cc.

472{
473 G4double kineticEnergy = dp->GetKineticEnergy();
474 if(kineticEnergy < lowKinEnergy) { return; }
475 G4double cut = std::min(cutEnergy, kineticEnergy);
476 G4double emax = std::min(maxEnergy, kineticEnergy);
477 if(cut >= emax) { return; }
478
479 SetupForMaterial(particle, couple->GetMaterial(),kineticEnergy);
480
481 const G4Element* elm =
482 SelectRandomAtom(couple,particle,kineticEnergy,cut,emax);
483 SetCurrentElement(elm->GetZ());
484
485 kinEnergy = kineticEnergy;
486 totalEnergy = kineticEnergy + particleMass;
488
489 G4double xmin = log(cut*cut + densityCorr);
490 G4double xmax = log(emax*emax + densityCorr);
491 G4double gammaEnergy, f, x;
492
493 do {
494 x = exp(xmin + G4UniformRand()*(xmax - xmin)) - densityCorr;
495 if(x < 0.0) x = 0.0;
496 gammaEnergy = sqrt(x);
497 f = ComputeDXSectionPerAtom(gammaEnergy);
498
499 if ( f > fMax ) {
500 G4cout << "### G4eBremParametrizedModel Warning: Majoranta exceeded! "
501 << f << " > " << fMax
502 << " Egamma(MeV)= " << gammaEnergy
503 << " E(mEV)= " << kineticEnergy
504 << G4endl;
505 }
506
507 } while (f < fMax*G4UniformRand());
508
509 //
510 // angles of the emitted gamma. ( Z - axis along the parent particle)
511 // use general interface
512 //
513 G4ThreeVector gammaDirection =
516 couple->GetMaterial());
517
518 // create G4DynamicParticle object for the Gamma
519 G4DynamicParticle* gamma = new G4DynamicParticle(theGamma,gammaDirection,
520 gammaEnergy);
521 vdp->push_back(gamma);
522
523 G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2));
524 G4ThreeVector direction = (totMomentum*dp->GetMomentumDirection()
525 - gammaEnergy*gammaDirection).unit();
526
527 // energy of primary
528 G4double finalE = kineticEnergy - gammaEnergy;
529
530 // stop tracking and create new secondary instead of primary
531 if(gammaEnergy > SecondaryThreshold()) {
534 G4DynamicParticle* el =
536 direction, finalE);
537 vdp->push_back(el);
538
539 // continue tracking
540 } else {
543 }
544}
@ 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)
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:508
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)
int G4lrint(double ad)
Definition: templates.hh:163

◆ SetupForMaterial()

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

Reimplemented from G4VEmModel.

Definition at line 134 of file G4eBremParametrizedModel.cc.

137{
138 densityFactor = mat->GetElectronDensity()*fMigdalConstant;
139
140 // calculate threshold for density effect
141 kinEnergy = kineticEnergy;
142 totalEnergy = kineticEnergy + particleMass;
144}
G4double GetElectronDensity() const
Definition: G4Material.hh:216

Referenced by ComputeDEDXPerVolume(), and SampleSecondaries().

Member Data Documentation

◆ currentZ

G4double G4eBremParametrizedModel::currentZ
protected

◆ densityCorr

G4double G4eBremParametrizedModel::densityCorr
protected

◆ densityFactor

G4double G4eBremParametrizedModel::densityFactor
protected

◆ facFel

G4double G4eBremParametrizedModel::facFel
protected

Definition at line 136 of file G4eBremParametrizedModel.hh.

◆ facFinel

G4double G4eBremParametrizedModel::facFinel
protected

Definition at line 136 of file G4eBremParametrizedModel.hh.

◆ fCoulomb

G4double G4eBremParametrizedModel::fCoulomb
protected

Definition at line 137 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

◆ Fel

G4double G4eBremParametrizedModel::Fel
protected

Definition at line 135 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

◆ Finel

G4double G4eBremParametrizedModel::Finel
protected

Definition at line 135 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

◆ fMax

G4double G4eBremParametrizedModel::fMax
protected

Definition at line 137 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel(), and SampleSecondaries().

◆ fParticleChange

G4ParticleChangeForLoss* G4eBremParametrizedModel::fParticleChange
protected

Definition at line 121 of file G4eBremParametrizedModel.hh.

Referenced by Initialise(), and SampleSecondaries().

◆ isElectron

G4bool G4eBremParametrizedModel::isElectron
protected

Definition at line 139 of file G4eBremParametrizedModel.hh.

◆ kinEnergy

G4double G4eBremParametrizedModel::kinEnergy
protected

◆ lnZ

G4double G4eBremParametrizedModel::lnZ
protected

Definition at line 132 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

◆ minThreshold

G4double G4eBremParametrizedModel::minThreshold
protected

Definition at line 125 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel(), and MinEnergyCut().

◆ nist

G4NistManager* G4eBremParametrizedModel::nist
protected

Definition at line 118 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

◆ particle

const G4ParticleDefinition* G4eBremParametrizedModel::particle
protected

◆ particleMass

G4double G4eBremParametrizedModel::particleMass
protected

◆ theGamma

G4ParticleDefinition* G4eBremParametrizedModel::theGamma
protected

Definition at line 120 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel(), and SampleSecondaries().

◆ totalEnergy

G4double G4eBremParametrizedModel::totalEnergy
protected

◆ 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 123 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 123 of file G4eBremParametrizedModel.hh.

◆ z13

G4double G4eBremParametrizedModel::z13
protected

Definition at line 132 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

◆ z23

G4double G4eBremParametrizedModel::z23
protected

Definition at line 132 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().


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