Geant4 11.2.2
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=nullptr, const G4String &nam="MuBrem")
 
 ~G4MuBremsstrahlungModel ()=default
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *) override
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
 
G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
void SetLowestKineticEnergy (G4double e)
 
G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double) override
 
G4MuBremsstrahlungModeloperator= (const G4MuBremsstrahlungModel &right)=delete
 
 G4MuBremsstrahlungModel (const G4MuBremsstrahlungModel &)=delete
 
- 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 level, const G4ParticleDefinition *, G4double kineticEnergy)
 
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 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 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 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)
 
void SetLPMFlag (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

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 = nullptr
 
G4ParticleDefinitiontheGamma = nullptr
 
G4ParticleChangeForLossfParticleChange = nullptr
 
G4NistManagernist = nullptr
 
G4double mass = 1.0
 
G4double rmass = 1.0
 
G4double cc = 1.0
 
G4double coeff = 1.0
 
G4double sqrte
 
G4double bh = 202.4
 
G4double bh1 = 446.
 
G4double btf = 183.
 
G4double btf1 = 1429.
 
G4double lowestKinEnergy
 
G4double minThreshold
 
- 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
 
std::size_t currentCoupleIndex = 0
 
std::size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Static Protected Attributes

static const G4double xgi [6]
 
static const G4double wgi [6]
 
static G4double fDN [93] = {0.0}
 

Detailed Description

Definition at line 69 of file G4MuBremsstrahlungModel.hh.

Constructor & Destructor Documentation

◆ G4MuBremsstrahlungModel() [1/2]

G4MuBremsstrahlungModel::G4MuBremsstrahlungModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "MuBrem" )
explicit

Definition at line 83 of file G4MuBremsstrahlungModel.cc.

85 : G4VEmModel(nam),
86 sqrte(std::sqrt(G4Exp(1.))),
87 lowestKinEnergy(0.1*CLHEP::GeV),
88 minThreshold(0.9*CLHEP::keV)
89{
92
94
95 if (nullptr != p) { SetParticle(p); }
96 if (0.0 == fDN[1]) {
97 for (G4int i=1; i<93; ++i) {
98 G4double dn = 1.54*nist->GetA27(i);
99 fDN[i] = dn;
100 if(1 < i) {
101 fDN[i] /= std::pow(dn, 1./G4double(i));
102 }
103 }
104 }
105}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:180
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
void SetParticle(const G4ParticleDefinition *)
G4ParticleDefinition * theGamma
G4double GetA27(G4int Z) const
static G4NistManager * Instance()
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67
void SetAngularDistribution(G4VEmAngularDistribution *)

◆ ~G4MuBremsstrahlungModel()

G4MuBremsstrahlungModel::~G4MuBremsstrahlungModel ( )
default

◆ G4MuBremsstrahlungModel() [2/2]

G4MuBremsstrahlungModel::G4MuBremsstrahlungModel ( const G4MuBremsstrahlungModel & )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 311 of file G4MuBremsstrahlungModel.cc.

317{
318 G4double cross = 0.0;
319 if (kineticEnergy <= lowestKinEnergy) return cross;
320 G4double tmax = std::min(maxEnergy, kineticEnergy);
321 G4double cut = std::min(cutEnergy, kineticEnergy);
322 if (cut < minThreshold) cut = minThreshold;
323 if (cut >= tmax) return cross;
324
325 cross = ComputeMicroscopicCrossSection (kineticEnergy, Z, cut);
326 if(tmax < kineticEnergy) {
327 cross -= ComputeMicroscopicCrossSection(kineticEnergy, Z, tmax);
328 }
329 return cross;
330}
G4double ComputeMicroscopicCrossSection(G4double tkin, G4double Z, G4double cut)

◆ ComputeDEDXPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 163 of file G4MuBremsstrahlungModel.cc.

168{
169 G4double dedx = 0.0;
170 if (kineticEnergy <= lowestKinEnergy) { return dedx; }
171
172 G4double cut = std::max(cutEnergy, minThreshold);
173 cut = std::min(cut, kineticEnergy);
174
175 const G4ElementVector* theElementVector = material->GetElementVector();
176 const G4double* theAtomicNumDensityVector =
177 material->GetAtomicNumDensityVector();
178
179 // loop for elements in the material
180 for (size_t i=0; i<material->GetNumberOfElements(); ++i) {
181 G4double loss =
182 ComputMuBremLoss((*theElementVector)[i]->GetZ(), kineticEnergy, cut);
183 dedx += loss*theAtomicNumDensityVector[i];
184 }
185 // G4cout << "BR e= " << kineticEnergy << " dedx= " << dedx << G4endl;
186 dedx = std::max(dedx, 0.);
187 return dedx;
188}
std::vector< const G4Element * > G4ElementVector
const G4ElementVector * GetElementVector() const
const G4double * GetAtomicNumDensityVector() const
std::size_t GetNumberOfElements() const
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 259 of file G4MuBremsstrahlungModel.cc.

264{
265 G4double dxsection = 0.;
266 if(gammaEnergy > tkin) { return dxsection; }
267
268 G4double E = tkin + mass ;
269 G4double v = gammaEnergy/E ;
270 G4double delta = 0.5*mass*mass*v/(E-gammaEnergy) ;
271 G4double rab0 = delta*sqrte ;
272
273 G4int iz = G4lrint(Z);
274 if(iz < 1) { iz = 1; }
275 else if(iz > 92) { iz = 92; }
276
277 G4double z13 = 1.0/nist->GetZ13(iz);
278 G4double dnstar = fDN[iz];
279
280 G4double b,b1;
281 if(1 == iz) {
282 b = bh;
283 b1 = bh1;
284 } else {
285 b = btf;
286 b1 = btf1;
287 }
288
289 // nucleus contribution logarithm
290 G4double rab1 = b*z13;
291 G4double fn = G4Log(rab1/(dnstar*(CLHEP::electron_mass_c2+rab0*rab1))*
292 (mass + delta*(dnstar*sqrte-2.)));
293 fn = std::max(fn, 0.);
294 // electron contribution logarithm
295 G4double epmax1 = E/(1.+0.5*mass*rmass/E);
296 G4double fe = 0.;
297 if(gammaEnergy < epmax1) {
298 G4double rab2 = b1*z13*z13;
299 fe = G4Log(rab2*mass/((1.+delta*rmass/(CLHEP::electron_mass_c2*sqrte))*
300 (CLHEP::electron_mass_c2+rab0*rab2)));
301 fe = std::max(fe, 0.);
302 }
303
304 dxsection = coeff*(1.-v*(1. - 0.75*v))*Z*(fn*Z + fe)/gammaEnergy;
305 dxsection = std::max(dxsection, 0.0);
306 return dxsection;
307}
G4fissionEvent * fe
G4double G4Log(G4double x)
Definition G4Log.hh:227
G4double GetZ13(G4double Z) const
int G4lrint(double ad)
Definition templates.hh:134

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

◆ ComputeMicroscopicCrossSection()

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

Definition at line 221 of file G4MuBremsstrahlungModel.cc.

225{
226 G4double totalEnergy = tkin + mass;
227 static const G4double ak1 = 2.3;
228 static const G4int k2 = 4;
229 G4double cross = 0.;
230
231 if(cut >= tkin) return cross;
232
233 G4double vcut = cut/totalEnergy;
234 G4double vmax = tkin/totalEnergy;
235
236 G4double aaa = G4Log(vcut);
237 G4double bbb = G4Log(vmax);
238 G4int kkk = (G4int)((bbb-aaa)/ak1) + k2 ;
239 if(kkk > 8) { kkk = 8; }
240 else if (kkk < 1) { kkk = 1; }
241 G4double hhh = (bbb-aaa)/(G4double)(kkk);
242 G4double aa = aaa;
243
244 for(G4int l=0; l<kkk; ++l) {
245 for(G4int i=0; i<6; ++i) {
246 G4double ep = G4Exp(aa + xgi[i]*hhh)*totalEnergy;
247 cross += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
248 }
249 aa += hhh;
250 }
251
252 cross *= hhh;
253 //G4cout << "BR e= " << tkin<< " cross= " << cross/barn << G4endl;
254 return cross;
255}
static const G4double xgi[6]
virtual G4double ComputeDMicroscopicCrossSection(G4double tkin, G4double Z, G4double gammaEnergy)
static const G4double wgi[6]

Referenced by ComputeCrossSectionPerAtom().

◆ ComputMuBremLoss()

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

Definition at line 192 of file G4MuBremsstrahlungModel.cc.

194{
195 G4double totalEnergy = mass + tkin;
196 static const G4double ak1 = 0.05;
197 static const G4int k2 = 5;
198 G4double loss = 0.;
199
200 G4double vcut = cut/totalEnergy;
201 G4int kkk = (G4int)(vcut/ak1) + k2;
202 if (kkk > 8) { kkk = 8; }
203 else if (kkk < 1) { kkk = 1; }
204 G4double hhh = vcut/(G4double)(kkk);
205
206 G4double aa = 0.;
207 for(G4int l=0; l<kkk; ++l) {
208 for(G4int i=0; i<6; ++i) {
209 G4double ep = (aa + xgi[i]*hhh)*totalEnergy;
210 loss += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
211 }
212 aa += hhh;
213 }
214
215 loss *= hhh*totalEnergy;
216 return loss;
217}

Referenced by ComputeDEDXPerVolume().

◆ Initialise()

void G4MuBremsstrahlungModel::Initialise ( const G4ParticleDefinition * p,
const G4DataVector & cuts )
overridevirtual

Implements G4VEmModel.

Definition at line 139 of file G4MuBremsstrahlungModel.cc.

141{
142 SetParticle(p);
143 if(nullptr == fParticleChange) {
145 }
146 if(IsMaster() && p == particle && lowestKinEnergy < HighEnergyLimit()) {
148 }
149}
G4ParticleChangeForLoss * fParticleChange
const G4ParticleDefinition * particle
G4bool IsMaster() const
G4double HighEnergyLimit() const
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
G4ParticleChangeForLoss * GetParticleChangeForLoss()

◆ InitialiseLocal()

void G4MuBremsstrahlungModel::InitialiseLocal ( const G4ParticleDefinition * p,
G4VEmModel * masterModel )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 153 of file G4MuBremsstrahlungModel.cc.

155{
158 }
159}
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
std::vector< G4EmElementSelector * > * GetElementSelectors()

◆ MinEnergyCut()

G4double G4MuBremsstrahlungModel::MinEnergyCut ( const G4ParticleDefinition * ,
const G4MaterialCutsCouple *  )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 109 of file G4MuBremsstrahlungModel.cc.

111{
112 return minThreshold;
113}

◆ MinPrimaryEnergy()

G4double G4MuBremsstrahlungModel::MinPrimaryEnergy ( const G4Material * ,
const G4ParticleDefinition * ,
G4double cut )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 117 of file G4MuBremsstrahlungModel.cc.

120{
121 return std::max(lowestKinEnergy, cut);
122}

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 334 of file G4MuBremsstrahlungModel.cc.

340{
341 G4double kineticEnergy = dp->GetKineticEnergy();
342 // check against insufficient energy
343 G4double tmax = std::min(kineticEnergy, maxEnergy);
344 G4double tmin = std::min(kineticEnergy, minEnergy);
345 tmin = std::max(tmin, minThreshold);
346 if(tmin >= tmax) return;
347
348 // ===== sampling of energy transfer ======
349
350 G4ParticleMomentum partDirection = dp->GetMomentumDirection();
351
352 // select randomly one element constituing the material
353 const G4Element* anElement = SelectRandomAtom(couple,particle,kineticEnergy);
354 G4double Z = anElement->GetZ();
355 G4double func1 = tmin*
356 ComputeDMicroscopicCrossSection(kineticEnergy, Z, tmin);
357
358 G4double gEnergy;
359 G4double func2;
360
361 G4double xmin = G4Log(tmin/minThreshold);
362 G4double xmax = G4Log(tmax/tmin);
363
364 do {
365 gEnergy = minThreshold*G4Exp(xmin + G4UniformRand()*xmax);
366 func2 = gEnergy*ComputeDMicroscopicCrossSection(kineticEnergy, Z, gEnergy);
367
368 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
369 } while(func2 < func1*G4UniformRand());
370
371 // angles of the emitted gamma using general interface
372 G4ThreeVector gamDir =
373 GetAngularDistribution()->SampleDirection(dp, gEnergy, Z,
374 couple->GetMaterial());
375 // create G4DynamicParticle object for the Gamma
376 G4DynamicParticle* gamma = new G4DynamicParticle(theGamma, gamDir, gEnergy);
377 vdp->push_back(gamma);
378
379 // compute post-interaction kinematics of primary e-/e+ based on
380 // energy-momentum conservation
381 const G4double totMomentum =
382 std::sqrt(kineticEnergy*(kineticEnergy + 2.0*mass));
383 G4ThreeVector dir =
384 (totMomentum*dp->GetMomentumDirection() - gEnergy*gamDir).unit();
385 const G4double finalE = kineticEnergy - gEnergy;
386
387 // if secondary gamma energy is higher than threshold(very high by default)
388 // then stop tracking the primary particle and create new secondary e-/e+
389 // instead of the primary one
390 if (gEnergy > SecondaryThreshold()) {
393 G4DynamicParticle* newdp = new G4DynamicParticle(particle, dir, finalE);
394 vdp->push_back(newdp);
395 } else {
396 // continue tracking the primary e-/e+ otherwise
399 }
400}
@ fStopAndKill
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector unit() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition G4Element.hh:119
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()
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
G4double SecondaryThreshold() const
void ProposeTrackStatus(G4TrackStatus status)

◆ SetLowestKineticEnergy()

void G4MuBremsstrahlungModel::SetLowestKineticEnergy ( G4double e)
inline

Definition at line 156 of file G4MuBremsstrahlungModel.hh.

157{
158 lowestKinEnergy = e;
159}

◆ SetParticle()

void G4MuBremsstrahlungModel::SetParticle ( const G4ParticleDefinition * p)
protected

Definition at line 126 of file G4MuBremsstrahlungModel.cc.

127{
128 if(nullptr == particle) {
129 particle = p;
131 rmass = mass/CLHEP::electron_mass_c2 ;
132 cc = CLHEP::classic_electr_radius/rmass ;
133 coeff = 16.*CLHEP::fine_structure_const*cc*cc/3. ;
134 }
135}

Referenced by G4MuBremsstrahlungModel(), and Initialise().

Member Data Documentation

◆ bh

G4double G4MuBremsstrahlungModel::bh = 202.4
protected

◆ bh1

G4double G4MuBremsstrahlungModel::bh1 = 446.
protected

Definition at line 142 of file G4MuBremsstrahlungModel.hh.

Referenced by ComputeDMicroscopicCrossSection().

◆ btf

G4double G4MuBremsstrahlungModel::btf = 183.
protected

◆ btf1

G4double G4MuBremsstrahlungModel::btf1 = 1429.
protected

Definition at line 144 of file G4MuBremsstrahlungModel.hh.

Referenced by ComputeDMicroscopicCrossSection().

◆ cc

G4double G4MuBremsstrahlungModel::cc = 1.0
protected

Definition at line 138 of file G4MuBremsstrahlungModel.hh.

Referenced by SetParticle().

◆ coeff

G4double G4MuBremsstrahlungModel::coeff = 1.0
protected

◆ fDN

G4double G4MuBremsstrahlungModel::fDN = {0.0}
staticprotected

◆ fParticleChange

G4ParticleChangeForLoss* G4MuBremsstrahlungModel::fParticleChange = nullptr
protected

Definition at line 133 of file G4MuBremsstrahlungModel.hh.

Referenced by Initialise(), and SampleSecondaries().

◆ lowestKinEnergy

G4double G4MuBremsstrahlungModel::lowestKinEnergy
protected

◆ mass

◆ minThreshold

G4double G4MuBremsstrahlungModel::minThreshold
protected

◆ nist

G4NistManager* G4MuBremsstrahlungModel::nist = nullptr
protected

◆ particle

const G4ParticleDefinition* G4MuBremsstrahlungModel::particle = nullptr
protected

◆ rmass

G4double G4MuBremsstrahlungModel::rmass = 1.0
protected

Definition at line 137 of file G4MuBremsstrahlungModel.hh.

Referenced by ComputeDMicroscopicCrossSection(), and SetParticle().

◆ sqrte

G4double G4MuBremsstrahlungModel::sqrte
protected

◆ theGamma

G4ParticleDefinition* G4MuBremsstrahlungModel::theGamma = nullptr
protected

Definition at line 132 of file G4MuBremsstrahlungModel.hh.

Referenced by G4MuBremsstrahlungModel(), and SampleSecondaries().

◆ wgi

const G4double G4MuBremsstrahlungModel::wgi
staticprotected
Initial value:
=
{0.08566,0.18038,0.23396,0.23396,0.18038,0.08566}

Definition at line 79 of file G4MuBremsstrahlungModel.hh.

Referenced by ComputeMicroscopicCrossSection(), and ComputMuBremLoss().

◆ xgi

const G4double G4MuBremsstrahlungModel::xgi
staticprotected
Initial value:
=
{0.03377,0.16940,0.38069,0.61931,0.83060,0.96623}

Definition at line 77 of file G4MuBremsstrahlungModel.hh.

Referenced by ComputeMicroscopicCrossSection(), and ComputMuBremLoss().


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