Geant4 11.1.1
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=nullptr, const G4String &nam="eBremParam")
 
 ~G4eBremParametrizedModel () override
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *) override
 
G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double tkin, G4double Z, G4double, G4double cutEnergy, G4double maxEnergy=DBL_MAX) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy) override
 
void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double) override
 
G4eBremParametrizedModeloperator= (const G4eBremParametrizedModel &right)=delete
 
 G4eBremParametrizedModel (const G4eBremParametrizedModel &)=delete
 
- 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 void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
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 GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
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 MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
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 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) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
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)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

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
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
 
size_t currentCoupleIndex = 0
 
size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

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

Constructor & Destructor Documentation

◆ G4eBremParametrizedModel() [1/2]

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

Definition at line 98 of file G4eBremParametrizedModel.cc.

100 : G4VEmModel(nam),
101 particle(nullptr),
102 fMigdalConstant(classic_electr_radius*electron_Compton_length*electron_Compton_length*4.0*pi),
103 bremFactor(fine_structure_const*classic_electr_radius*classic_electr_radius*16./3.),
104 isInitialised(false),
105 isElectron(true)
106{
108
109 minThreshold = 0.1*keV;
110 lowKinEnergy = 10.*MeV;
111 SetLowEnergyLimit(lowKinEnergy);
112
114
116
119
120 InitialiseConstants();
121 if(nullptr != p) { SetParticle(p); }
122}
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4NistManager * Instance()
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:607
const G4ParticleDefinition * particle

◆ ~G4eBremParametrizedModel()

G4eBremParametrizedModel::~G4eBremParametrizedModel ( )
overridedefault

◆ G4eBremParametrizedModel() [2/2]

G4eBremParametrizedModel::G4eBremParametrizedModel ( const G4eBremParametrizedModel )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 263 of file G4eBremParametrizedModel.cc.

269{
270 if(!particle) { SetParticle(p); }
271 if(kineticEnergy < lowKinEnergy) { return 0.0; }
272 G4double cut = std::min(cutEnergy, kineticEnergy);
273 G4double tmax = std::min(maxEnergy, kineticEnergy);
274
275 if(cut >= tmax) { return 0.0; }
276
277 SetCurrentElement(Z);
278
279 G4double cross = ComputeXSectionPerAtom(cut);
280
281 // allow partial integration
282 if(tmax < kinEnergy) { cross -= ComputeXSectionPerAtom(tmax); }
283
284 cross *= Z*Z*bremFactor;
285
286 return cross;
287}
double G4double
Definition: G4Types.hh:83
const G4int Z[17]

◆ ComputeDEDXPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 197 of file G4eBremParametrizedModel.cc.

202{
203 if(!particle) { SetParticle(p); }
204 if(kineticEnergy < lowKinEnergy) { return 0.0; }
205 G4double cut = std::min(cutEnergy, kineticEnergy);
206 if(cut == 0.0) { return 0.0; }
207
208 SetupForMaterial(particle, material,kineticEnergy);
209
210 const G4ElementVector* theElementVector = material->GetElementVector();
211 const G4double* theAtomicNumDensityVector = material->GetAtomicNumDensityVector();
212
213 G4double dedx = 0.0;
214
215 // loop for elements in the material
216 for (size_t i=0; i<material->GetNumberOfElements(); i++) {
217
218 G4VEmModel::SetCurrentElement((*theElementVector)[i]);
219 SetCurrentElement((*theElementVector)[i]->GetZ());
220
221 dedx += theAtomicNumDensityVector[i]*currentZ*currentZ*ComputeBremLoss(cut);
222 }
223 dedx *= bremFactor;
224
225 return dedx;
226}
std::vector< const G4Element * > G4ElementVector
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:185
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:211
void SetCurrentElement(const G4Element *)
Definition: G4VEmModel.hh:493
void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double) override

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 171 of file G4eBremParametrizedModel.cc.

173{
174 if(p) { SetParticle(p); }
175
176 lowKinEnergy = LowEnergyLimit();
177
178 currentZ = 0.;
179
180 if(IsMaster()) { InitialiseElementSelectors(p, cuts); }
181
182 if(isInitialised) { return; }
184 isInitialised = true;
185}
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4bool IsMaster() const
Definition: G4VEmModel.hh:725
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:139
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:109
G4ParticleChangeForLoss * fParticleChange

◆ InitialiseLocal()

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

Reimplemented from G4VEmModel.

Definition at line 189 of file G4eBremParametrizedModel.cc.

191{
193}
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:831
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:823

◆ MinEnergyCut()

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

Reimplemented from G4VEmModel.

Definition at line 148 of file G4eBremParametrizedModel.cc.

150{
151 return minThreshold;
152}

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 470 of file G4eBremParametrizedModel.cc.

476{
477 G4double kineticEnergy = dp->GetKineticEnergy();
478 if(kineticEnergy < lowKinEnergy) { return; }
479 G4double cut = std::min(cutEnergy, kineticEnergy);
480 G4double emax = std::min(maxEnergy, kineticEnergy);
481 if(cut >= emax) { return; }
482
483 SetupForMaterial(particle, couple->GetMaterial(),kineticEnergy);
484
485 const G4Element* elm = SelectTargetAtom(couple,particle,kineticEnergy,
486 dp->GetLogKineticEnergy(),cut,emax);
487 SetCurrentElement(elm->GetZ());
488
489 kinEnergy = kineticEnergy;
490 totalEnergy = kineticEnergy + particleMass;
492
493 G4double xmin = G4Log(cut*cut + densityCorr);
494 G4double xmax = G4Log(emax*emax + densityCorr);
495 G4double gammaEnergy, f, x;
496
497 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
498
499 do {
500 x = G4Exp(xmin + rndmEngine->flat()*(xmax - xmin)) - densityCorr;
501 if(x < 0.0) x = 0.0;
502 gammaEnergy = sqrt(x);
503 f = ComputeDXSectionPerAtom(gammaEnergy);
504
505 if ( f > fMax ) {
506 G4cout << "### G4eBremParametrizedModel Warning: Majoranta exceeded! "
507 << f << " > " << fMax
508 << " Egamma(MeV)= " << gammaEnergy
509 << " E(mEV)= " << kineticEnergy
510 << G4endl;
511 }
512
513 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
514 } while (f < fMax*rndmEngine->flat());
515
516 //
517 // angles of the emitted gamma. ( Z - axis along the parent particle)
518 // use general interface
519 //
520 G4ThreeVector gammaDirection =
523 couple->GetMaterial());
524
525 // create G4DynamicParticle object for the Gamma
526 auto gamma = new G4DynamicParticle(theGamma,gammaDirection, gammaEnergy);
527 vdp->push_back(gamma);
528
529 G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2));
530 G4ThreeVector direction = (totMomentum*dp->GetMomentumDirection()
531 - gammaEnergy*gammaDirection).unit();
532
533 // energy of primary
534 G4double finalE = kineticEnergy - gammaEnergy;
535
536 // stop tracking and create new secondary instead of primary
537 if(gammaEnergy > SecondaryThreshold()) {
540 auto el =
542 direction, finalE);
543 vdp->push_back(el);
544
545 // continue tracking
546 } else {
549 }
550}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
G4double G4Log(G4double x)
Definition: G4Log.hh:227
@ fStopAndKill
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
virtual double flat()=0
const G4ThreeVector & GetMomentumDirection() const
G4double GetLogKineticEnergy() 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:600
const G4Element * SelectTargetAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:577
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:669
void ProposeTrackStatus(G4TrackStatus status)
double flat()
Definition: G4AblaRandom.cc:48
int G4lrint(double ad)
Definition: templates.hh:134

◆ SetupForMaterial()

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

Reimplemented from G4VEmModel.

Definition at line 156 of file G4eBremParametrizedModel.cc.

159{
160 densityFactor = mat->GetElectronDensity()*fMigdalConstant;
161
162 // calculate threshold for density effect
163 kinEnergy = kineticEnergy;
164 totalEnergy = kineticEnergy + particleMass;
166}
G4double GetElectronDensity() const
Definition: G4Material.hh:212

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

◆ facFinel

G4double G4eBremParametrizedModel::facFinel
protected

Definition at line 143 of file G4eBremParametrizedModel.hh.

◆ fCoulomb

G4double G4eBremParametrizedModel::fCoulomb
protected

Definition at line 144 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

◆ Fel

G4double G4eBremParametrizedModel::Fel
protected

Definition at line 142 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

◆ Finel

G4double G4eBremParametrizedModel::Finel
protected

Definition at line 142 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

◆ fMax

G4double G4eBremParametrizedModel::fMax
protected

Definition at line 144 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel(), and SampleSecondaries().

◆ fParticleChange

G4ParticleChangeForLoss* G4eBremParametrizedModel::fParticleChange
protected

Definition at line 128 of file G4eBremParametrizedModel.hh.

Referenced by Initialise(), and SampleSecondaries().

◆ isElectron

G4bool G4eBremParametrizedModel::isElectron
protected

Definition at line 156 of file G4eBremParametrizedModel.hh.

◆ kinEnergy

G4double G4eBremParametrizedModel::kinEnergy
protected

◆ lnZ

G4double G4eBremParametrizedModel::lnZ
protected

Definition at line 139 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

◆ minThreshold

G4double G4eBremParametrizedModel::minThreshold
protected

Definition at line 132 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel(), and MinEnergyCut().

◆ nist

G4NistManager* G4eBremParametrizedModel::nist
protected

Definition at line 125 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 127 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 130 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 130 of file G4eBremParametrizedModel.hh.

◆ z13

G4double G4eBremParametrizedModel::z13
protected

Definition at line 139 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

◆ z23

G4double G4eBremParametrizedModel::z23
protected

Definition at line 139 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().


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