Geant4 10.7.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=nullptr, const G4String &nam="eBremParam")
 
virtual ~G4eBremParametrizedModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *) 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 SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy) override
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double) override
 
- 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 *, 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 G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
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)
 
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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
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 *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
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
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

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()

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

Definition at line 98 of file G4eBremParametrizedModel.cc.

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

◆ ~G4eBremParametrizedModel()

G4eBremParametrizedModel::~G4eBremParametrizedModel ( )
virtual

Definition at line 134 of file G4eBremParametrizedModel.cc.

135{
136}

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

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

◆ ComputeDEDXPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 199 of file G4eBremParametrizedModel.cc.

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

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 173 of file G4eBremParametrizedModel.cc.

175{
176 if(p) { SetParticle(p); }
177
178 lowKinEnergy = LowEnergyLimit();
179
180 currentZ = 0.;
181
182 if(IsMaster()) { InitialiseElementSelectors(p, cuts); }
183
184 if(isInitialised) { return; }
186 isInitialised = true;
187}
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
G4bool IsMaster() const
Definition: G4VEmModel.hh:736
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:118
G4ParticleChangeForLoss * fParticleChange

◆ InitialiseLocal()

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

Reimplemented from G4VEmModel.

Definition at line 191 of file G4eBremParametrizedModel.cc.

193{
195}
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:842
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:834

◆ MinEnergyCut()

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

Reimplemented from G4VEmModel.

Definition at line 150 of file G4eBremParametrizedModel.cc.

152{
153 return minThreshold;
154}

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 472 of file G4eBremParametrizedModel.cc.

478{
479 G4double kineticEnergy = dp->GetKineticEnergy();
480 if(kineticEnergy < lowKinEnergy) { return; }
481 G4double cut = std::min(cutEnergy, kineticEnergy);
482 G4double emax = std::min(maxEnergy, kineticEnergy);
483 if(cut >= emax) { return; }
484
485 SetupForMaterial(particle, couple->GetMaterial(),kineticEnergy);
486
487 const G4Element* elm = SelectTargetAtom(couple,particle,kineticEnergy,
488 dp->GetLogKineticEnergy(),cut,emax);
489 SetCurrentElement(elm->GetZ());
490
491 kinEnergy = kineticEnergy;
492 totalEnergy = kineticEnergy + particleMass;
494
495 G4double xmin = G4Log(cut*cut + densityCorr);
496 G4double xmax = G4Log(emax*emax + densityCorr);
497 G4double gammaEnergy, f, x;
498
499 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
500
501 do {
502 x = G4Exp(xmin + rndmEngine->flat()*(xmax - xmin)) - densityCorr;
503 if(x < 0.0) x = 0.0;
504 gammaEnergy = sqrt(x);
505 f = ComputeDXSectionPerAtom(gammaEnergy);
506
507 if ( f > fMax ) {
508 G4cout << "### G4eBremParametrizedModel Warning: Majoranta exceeded! "
509 << f << " > " << fMax
510 << " Egamma(MeV)= " << gammaEnergy
511 << " E(mEV)= " << kineticEnergy
512 << G4endl;
513 }
514
515 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
516 } while (f < fMax*rndmEngine->flat());
517
518 //
519 // angles of the emitted gamma. ( Z - axis along the parent particle)
520 // use general interface
521 //
522 G4ThreeVector gammaDirection =
525 couple->GetMaterial());
526
527 // create G4DynamicParticle object for the Gamma
528 G4DynamicParticle* gamma = new G4DynamicParticle(theGamma,gammaDirection,
529 gammaEnergy);
530 vdp->push_back(gamma);
531
532 G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2));
533 G4ThreeVector direction = (totMomentum*dp->GetMomentumDirection()
534 - gammaEnergy*gammaDirection).unit();
535
536 // energy of primary
537 G4double finalE = kineticEnergy - gammaEnergy;
538
539 // stop tracking and create new secondary instead of primary
540 if(gammaEnergy > SecondaryThreshold()) {
543 G4DynamicParticle* el =
545 direction, finalE);
546 vdp->push_back(el);
547
548 // continue tracking
549 } else {
552 }
553}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double G4Log(G4double x)
Definition: G4Log.hh:226
@ 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:130
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:611
const G4Element * SelectTargetAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:587
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:680
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 158 of file G4eBremParametrizedModel.cc.

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

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

◆ facFinel

G4double G4eBremParametrizedModel::facFinel
protected

Definition at line 145 of file G4eBremParametrizedModel.hh.

◆ fCoulomb

G4double G4eBremParametrizedModel::fCoulomb
protected

Definition at line 146 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

◆ Fel

G4double G4eBremParametrizedModel::Fel
protected

Definition at line 144 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

◆ Finel

G4double G4eBremParametrizedModel::Finel
protected

Definition at line 144 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

◆ fMax

G4double G4eBremParametrizedModel::fMax
protected

Definition at line 146 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel(), and SampleSecondaries().

◆ fParticleChange

G4ParticleChangeForLoss* G4eBremParametrizedModel::fParticleChange
protected

Definition at line 130 of file G4eBremParametrizedModel.hh.

Referenced by Initialise(), and SampleSecondaries().

◆ isElectron

G4bool G4eBremParametrizedModel::isElectron
protected

Definition at line 148 of file G4eBremParametrizedModel.hh.

◆ kinEnergy

G4double G4eBremParametrizedModel::kinEnergy
protected

◆ lnZ

G4double G4eBremParametrizedModel::lnZ
protected

Definition at line 141 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

◆ minThreshold

G4double G4eBremParametrizedModel::minThreshold
protected

Definition at line 134 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel(), and MinEnergyCut().

◆ nist

G4NistManager* G4eBremParametrizedModel::nist
protected

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

◆ z13

G4double G4eBremParametrizedModel::z13
protected

Definition at line 141 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

◆ z23

G4double G4eBremParametrizedModel::z23
protected

Definition at line 141 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().


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