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

#include <G4MuBremsstrahlungModel.hh>

+ Inheritance diagram for G4MuBremsstrahlungModel:

Public Member Functions

 G4MuBremsstrahlungModel (const G4ParticleDefinition *p=0, const G4String &nam="MuBrem")
 
virtual ~G4MuBremsstrahlungModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
void SetLowestKineticEnergy (G4double e)
 
- 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 Member Functions

G4double ComputMuBremLoss (G4double Z, G4double tkin, G4double cut)
 
G4double ComputeMicroscopicCrossSection (G4double tkin, G4double Z, G4double cut)
 
virtual G4double ComputeDMicroscopicCrossSection (G4double tkin, G4double Z, G4double gammaEnergy)
 
void SetParticle (const G4ParticleDefinition *)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Protected Attributes

const G4ParticleDefinitionparticle
 
G4NistManagernist
 
G4double mass
 
G4double rmass
 
G4double cc
 
G4double coeff
 
G4double sqrte
 
G4double bh
 
G4double bh1
 
G4double btf
 
G4double btf1
 
- Protected Attributes inherited from G4VEmModel
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

Detailed Description

Definition at line 71 of file G4MuBremsstrahlungModel.hh.

Constructor & Destructor Documentation

◆ G4MuBremsstrahlungModel()

G4MuBremsstrahlungModel::G4MuBremsstrahlungModel ( const G4ParticleDefinition p = 0,
const G4String nam = "MuBrem" 
)

Definition at line 82 of file G4MuBremsstrahlungModel.cc.

84 : G4VEmModel(nam),
85 particle(0),
86 sqrte(sqrt(exp(1.))),
87 bh(202.4),
88 bh1(446.),
89 btf(183.),
90 btf1(1429.),
91 fParticleChange(0),
92 lowestKinEnergy(1.0*GeV),
93 minThreshold(1.0*keV)
94{
95 theGamma = G4Gamma::Gamma();
97
98 mass = rmass = cc = coeff = 1.0;
99
100 fDN[0] = 0.0;
101 for(G4int i=1; i<93; ++i) {
102 G4double dn = 1.54*nist->GetA27(i);
103 fDN[i] = dn;
104 if(1 < i) {
105 fDN[i] /= std::pow(dn, 1./G4double(i));
106 }
107 }
108
109 if(p) { SetParticle(p); }
110}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void SetParticle(const G4ParticleDefinition *)
const G4ParticleDefinition * particle
G4double GetA27(G4int Z)
static G4NistManager * Instance()

◆ ~G4MuBremsstrahlungModel()

G4MuBremsstrahlungModel::~G4MuBremsstrahlungModel ( )
virtual

Definition at line 114 of file G4MuBremsstrahlungModel.cc.

115{
116 size_t n = partialSumSigma.size();
117 if(n > 0) {
118 for(size_t i=0; i<n; i++) {
119 delete partialSumSigma[i];
120 }
121 }
122}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4MuBremsstrahlungModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kineticEnergy,
G4double  Z,
G4double  A,
G4double  cutEnergy,
G4double  maxEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 351 of file G4MuBremsstrahlungModel.cc.

357{
358 G4double cross = 0.0;
359 if (kineticEnergy <= lowestKinEnergy) return cross;
360 G4double tmax = std::min(maxEnergy, kineticEnergy);
361 G4double cut = std::min(cutEnergy, kineticEnergy);
362 if(cut < minThreshold) cut = minThreshold;
363 if (cut >= tmax) return cross;
364
365 cross = ComputeMicroscopicCrossSection (kineticEnergy, Z, cut);
366 if(tmax < kineticEnergy) {
367 cross -= ComputeMicroscopicCrossSection(kineticEnergy, Z, tmax);
368 }
369 return cross;
370}
G4double ComputeMicroscopicCrossSection(G4double tkin, G4double Z, G4double cut)

◆ ComputeDEDXPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 184 of file G4MuBremsstrahlungModel.cc.

189{
190 G4double dedx = 0.0;
191 if (kineticEnergy <= lowestKinEnergy) return dedx;
192
193 G4double tmax = kineticEnergy;
194 G4double cut = std::min(cutEnergy,tmax);
195 if(cut < minThreshold) cut = minThreshold;
196
197 const G4ElementVector* theElementVector = material->GetElementVector();
198 const G4double* theAtomicNumDensityVector =
199 material->GetAtomicNumDensityVector();
200
201 // loop for elements in the material
202 for (size_t i=0; i<material->GetNumberOfElements(); i++) {
203
204 G4double loss =
205 ComputMuBremLoss((*theElementVector)[i]->GetZ(), kineticEnergy, cut);
206
207 dedx += loss*theAtomicNumDensityVector[i];
208 }
209 // G4cout << "BR e= " << kineticEnergy << " dedx= " << dedx << G4endl;
210 if(dedx < 0.) dedx = 0.;
211 return dedx;
212}
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
G4double ComputMuBremLoss(G4double Z, G4double tkin, G4double cut)

◆ ComputeDMicroscopicCrossSection()

G4double G4MuBremsstrahlungModel::ComputeDMicroscopicCrossSection ( G4double  tkin,
G4double  Z,
G4double  gammaEnergy 
)
protectedvirtual

Reimplemented in G4hBremsstrahlungModel.

Definition at line 296 of file G4MuBremsstrahlungModel.cc.

301{
302 G4double dxsection = 0.;
303
304 if( gammaEnergy > tkin) return dxsection ;
305
306 G4double E = tkin + mass ;
307 G4double v = gammaEnergy/E ;
308 G4double delta = 0.5*mass*mass*v/(E-gammaEnergy) ;
309 G4double rab0=delta*sqrte ;
310
311 G4int iz = G4int(Z);
312 if(iz < 1) iz = 1;
313 else if(iz > 92) iz = 92;
314
315 G4double z13 = 1.0/nist->GetZ13(iz);
316 G4double dnstar = fDN[iz];
317
318 G4double b,b1;
319
320 if(1 == iz) {
321 b = bh;
322 b1 = bh1;
323 } else {
324 b = btf;
325 b1 = btf1;
326 }
327
328 // nucleus contribution logarithm
329 G4double rab1=b*z13;
330 G4double fn=log(rab1/(dnstar*(electron_mass_c2+rab0*rab1))*
331 (mass+delta*(dnstar*sqrte-2.))) ;
332 if(fn <0.) fn = 0. ;
333 // electron contribution logarithm
334 G4double epmax1=E/(1.+0.5*mass*rmass/E) ;
335 G4double fe=0.;
336 if(gammaEnergy<epmax1)
337 {
338 G4double rab2=b1*z13*z13 ;
339 fe=log(rab2*mass/((1.+delta*rmass/(electron_mass_c2*sqrte))*
340 (electron_mass_c2+rab0*rab2))) ;
341 if(fe<0.) fe=0. ;
342 }
343
344 dxsection = coeff*(1.-v*(1. - 0.75*v))*Z*(fn*Z + fe)/gammaEnergy;
345
346 return dxsection;
347}
G4fissionEvent * fe
G4double GetZ13(G4double Z)

Referenced by ComputeMicroscopicCrossSection(), ComputMuBremLoss(), and SampleSecondaries().

◆ ComputeMicroscopicCrossSection()

G4double G4MuBremsstrahlungModel::ComputeMicroscopicCrossSection ( G4double  tkin,
G4double  Z,
G4double  cut 
)
protected

Definition at line 253 of file G4MuBremsstrahlungModel.cc.

257{
258 G4double totalEnergy = tkin + mass;
259 G4double ak1 = 2.3;
260 G4int k2 = 4;
261 G4double xgi[]={0.03377,0.16940,0.38069,0.61931,0.83060,0.96623};
262 G4double wgi[]={0.08566,0.18038,0.23396,0.23396,0.18038,0.08566};
263 G4double cross = 0.;
264
265 if(cut >= tkin) return cross;
266
267 G4double vcut = cut/totalEnergy;
268 G4double vmax = tkin/totalEnergy;
269
270 G4double aaa = log(vcut);
271 G4double bbb = log(vmax);
272 G4int kkk = (G4int)((bbb-aaa)/ak1)+k2 ;
273 G4double hhh = (bbb-aaa)/G4double(kkk);
274
275 G4double aa = aaa;
276
277 for(G4int l=0; l<kkk; l++)
278 {
279 for(G4int i=0; i<6; i++)
280 {
281 G4double ep = exp(aa + xgi[i]*hhh)*totalEnergy;
282 cross += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
283 }
284 aa += hhh;
285 }
286
287 cross *=hhh;
288
289 //G4cout << "BR e= " << tkin<< " cross= " << cross/barn << G4endl;
290
291 return cross;
292}
virtual G4double ComputeDMicroscopicCrossSection(G4double tkin, G4double Z, G4double gammaEnergy)

Referenced by ComputeCrossSectionPerAtom().

◆ ComputMuBremLoss()

G4double G4MuBremsstrahlungModel::ComputMuBremLoss ( G4double  Z,
G4double  tkin,
G4double  cut 
)
protected

Definition at line 216 of file G4MuBremsstrahlungModel.cc.

218{
219 G4double totalEnergy = mass + tkin;
220 G4double ak1 = 0.05;
221 G4int k2=5;
222 G4double xgi[]={0.03377,0.16940,0.38069,0.61931,0.83060,0.96623};
223 G4double wgi[]={0.08566,0.18038,0.23396,0.23396,0.18038,0.08566};
224 G4double loss = 0.;
225
226 G4double vcut = cut/totalEnergy;
227 G4double vmax = tkin/totalEnergy;
228
229 G4double aaa = 0.;
230 G4double bbb = vcut;
231 if(vcut>vmax) bbb=vmax ;
232 G4int kkk = (G4int)((bbb-aaa)/ak1)+k2 ;
233 G4double hhh=(bbb-aaa)/float(kkk) ;
234
235 G4double aa = aaa;
236 for(G4int l=0; l<kkk; l++)
237 {
238 for(G4int i=0; i<6; i++)
239 {
240 G4double ep = (aa + xgi[i]*hhh)*totalEnergy;
241 loss += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
242 }
243 aa += hhh;
244 }
245
246 loss *=hhh*totalEnergy ;
247
248 return loss;
249}

Referenced by ComputeDEDXPerVolume().

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 134 of file G4MuBremsstrahlungModel.cc.

136{
137 if(p) { SetParticle(p); }
138
139 // partial cross section is computed for fixed energy
140 G4double fixedEnergy = 0.5*HighEnergyLimit();
141
142 const G4ProductionCutsTable* theCoupleTable=
144 if(theCoupleTable) {
145 G4int numOfCouples = theCoupleTable->GetTableSize();
146
147 G4int nn = partialSumSigma.size();
148 G4int nc = cuts.size();
149
150 // do we need to perform initialisation?
151 if(nn == numOfCouples) { return; }
152
153 // clear old data
154 if(nn > 0) {
155 for (G4int ii=0; ii<nn; ii++){
156 G4DataVector* a = partialSumSigma[ii];
157 if ( a ) { delete a; }
158 }
159 partialSumSigma.clear();
160 }
161 // fill new data
162 if (numOfCouples>0) {
163 for (G4int i=0; i<numOfCouples; i++) {
164 G4double cute = DBL_MAX;
165
166 // protection for usage with extrapolator
167 if(i < nc) { cute = cuts[i]; }
168
169 const G4MaterialCutsCouple* couple =
170 theCoupleTable->GetMaterialCutsCouple(i);
171 const G4Material* material = couple->GetMaterial();
172 G4DataVector* dv = ComputePartialSumSigma(material,fixedEnergy,cute);
173 partialSumSigma.push_back(dv);
174 }
175 }
176 }
177
178 // define pointer to G4ParticleChange
179 if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
180}
const G4Material * GetMaterial() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:95
#define DBL_MAX
Definition: templates.hh:83

◆ MinEnergyCut()

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

Definition at line 126 of file G4MuBremsstrahlungModel.cc.

128{
129 return minThreshold;
130}

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 403 of file G4MuBremsstrahlungModel.cc.

409{
410 G4double kineticEnergy = dp->GetKineticEnergy();
411 // check against insufficient energy
412 G4double tmax = std::min(kineticEnergy, maxEnergy);
413 G4double tmin = std::min(kineticEnergy, minEnergy);
414 if(tmin < minThreshold) tmin = minThreshold;
415 if(tmin >= tmax) return;
416
417 // ===== sampling of energy transfer ======
418
419 G4ParticleMomentum partDirection = dp->GetMomentumDirection();
420
421 // select randomly one element constituing the material
422 const G4Element* anElement = SelectRandomAtom(couple);
423 G4double Z = anElement->GetZ();
424
425 G4double totalEnergy = kineticEnergy + mass;
426 G4double totalMomentum = sqrt(kineticEnergy*(kineticEnergy + 2.0*mass));
427
428 G4double func1 = tmin*
429 ComputeDMicroscopicCrossSection(kineticEnergy,Z,tmin);
430
431 G4double lnepksi, epksi;
432 G4double func2;
433
434 G4double xmin = log(tmin/MeV);
435 G4double xmax = log(kineticEnergy/tmin);
436
437 do {
438 lnepksi = xmin + G4UniformRand()*xmax;
439 epksi = MeV*exp(lnepksi);
440 func2 = epksi*ComputeDMicroscopicCrossSection(kineticEnergy,Z,epksi);
441
442 } while(func2 < func1*G4UniformRand());
443
444 G4double gEnergy = epksi;
445
446 // ===== sample angle =====
447
448 G4double gam = totalEnergy/mass;
449 G4double rmax = gam*std::min(1.0, totalEnergy/gEnergy - 1.0);
450 G4double rmax2= rmax*rmax;
451 G4double x = G4UniformRand()*rmax2/(1.0 + rmax2);
452
453 G4double theta = sqrt(x/(1.0 - x))/gam;
454 G4double sint = sin(theta);
455 G4double phi = twopi * G4UniformRand() ;
456 G4double dirx = sint*cos(phi), diry = sint*sin(phi), dirz = cos(theta) ;
457
458 G4ThreeVector gDirection(dirx, diry, dirz);
459 gDirection.rotateUz(partDirection);
460
461 partDirection *= totalMomentum;
462 partDirection -= gEnergy*gDirection;
463 partDirection = partDirection.unit();
464
465 // primary change
466 kineticEnergy -= gEnergy;
467 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
468 fParticleChange->SetProposedMomentumDirection(partDirection);
469
470 // save secondary
471 G4DynamicParticle* aGamma =
472 new G4DynamicParticle(theGamma,gDirection,gEnergy);
473 vdp->push_back(aGamma);
474}
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector unit() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)

◆ SetLowestKineticEnergy()

void G4MuBremsstrahlungModel::SetLowestKineticEnergy ( G4double  e)
inline

Definition at line 160 of file G4MuBremsstrahlungModel.hh.

161{
162 lowestKinEnergy = e;
163}

◆ SetParticle()

void G4MuBremsstrahlungModel::SetParticle ( const G4ParticleDefinition p)
inlineprotected

Definition at line 168 of file G4MuBremsstrahlungModel.hh.

169{
170 if(!particle) {
171 particle = p;
173 rmass=mass/CLHEP::electron_mass_c2 ;
174 cc=CLHEP::classic_electr_radius/rmass ;
175 coeff= 16.*CLHEP::fine_structure_const*cc*cc/3. ;
176 }
177}

Referenced by G4MuBremsstrahlungModel(), and Initialise().

Member Data Documentation

◆ bh

G4double G4MuBremsstrahlungModel::bh
protected

◆ bh1

G4double G4MuBremsstrahlungModel::bh1
protected

Definition at line 141 of file G4MuBremsstrahlungModel.hh.

Referenced by ComputeDMicroscopicCrossSection().

◆ btf

G4double G4MuBremsstrahlungModel::btf
protected

◆ btf1

G4double G4MuBremsstrahlungModel::btf1
protected

Definition at line 143 of file G4MuBremsstrahlungModel.hh.

Referenced by ComputeDMicroscopicCrossSection().

◆ cc

G4double G4MuBremsstrahlungModel::cc
protected

Definition at line 137 of file G4MuBremsstrahlungModel.hh.

Referenced by G4MuBremsstrahlungModel(), and SetParticle().

◆ coeff

◆ mass

◆ nist

◆ particle

const G4ParticleDefinition* G4MuBremsstrahlungModel::particle
protected

◆ rmass

G4double G4MuBremsstrahlungModel::rmass
protected

◆ sqrte

G4double G4MuBremsstrahlungModel::sqrte
protected

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