Geant4 10.7.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=nullptr, const G4String &nam="MuBrem")
 
 ~G4MuBremsstrahlungModel ()=default
 
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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
void SetLowestKineticEnergy (G4double e)
 
virtual 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 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 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
 
G4ParticleDefinitiontheGamma
 
G4ParticleChangeForLossfParticleChange
 
G4double lowestKinEnergy
 
G4double minThreshold
 
- 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 [6]
 
static const G4double wgi [6]
 
static G4double fDN [93] = {0.0}
 

Detailed Description

Definition at line 70 of file G4MuBremsstrahlungModel.hh.

Constructor & Destructor Documentation

◆ G4MuBremsstrahlungModel() [1/2]

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

Definition at line 91 of file G4MuBremsstrahlungModel.cc.

93 : G4VEmModel(nam),
94 particle(nullptr),
95 sqrte(sqrt(G4Exp(1.))),
96 bh(202.4),
97 bh1(446.),
98 btf(183.),
99 btf1(1429.),
100 fParticleChange(nullptr),
101 lowestKinEnergy(1.0*GeV),
102 minThreshold(0.9*keV)
103{
106
107 lowestKinEnergy = 1.*GeV;
108
109 mass = rmass = cc = coeff = 1.0;
110
111 if(0.0 == fDN[1]) {
112 for(G4int i=1; i<93; ++i) {
113 G4double dn = 1.54*nist->GetA27(i);
114 fDN[i] = dn;
115 if(1 < i) {
116 fDN[i] /= std::pow(dn, 1./G4double(i));
117 }
118 }
119 }
121
122 if(p) { SetParticle(p); }
123}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
void SetParticle(const G4ParticleDefinition *)
G4ParticleDefinition * theGamma
G4ParticleChangeForLoss * fParticleChange
const G4ParticleDefinition * particle
G4double GetA27(G4int Z) const
static G4NistManager * Instance()
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:618

◆ ~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 336 of file G4MuBremsstrahlungModel.cc.

342{
343 G4double cross = 0.0;
344 if (kineticEnergy <= lowestKinEnergy) return cross;
345 G4double tmax = std::min(maxEnergy, kineticEnergy);
346 G4double cut = std::min(cutEnergy, kineticEnergy);
347 if(cut < minThreshold) cut = minThreshold;
348 if (cut >= tmax) return cross;
349
350 cross = ComputeMicroscopicCrossSection (kineticEnergy, Z, cut);
351 if(tmax < kineticEnergy) {
352 cross -= ComputeMicroscopicCrossSection(kineticEnergy, Z, tmax);
353 }
354 return cross;
355}
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 169 of file G4MuBremsstrahlungModel.cc.

174{
175 G4double dedx = 0.0;
176 if (kineticEnergy <= lowestKinEnergy) { return dedx; }
177
178 G4double tmax = kineticEnergy;
179 G4double cut = std::min(cutEnergy,tmax);
180 if(cut < minThreshold) { cut = minThreshold; }
181
182 const G4ElementVector* theElementVector = material->GetElementVector();
183 const G4double* theAtomicNumDensityVector =
184 material->GetAtomicNumDensityVector();
185
186 // loop for elements in the material
187 for (size_t i=0; i<material->GetNumberOfElements(); i++) {
188
189 G4double loss =
190 ComputMuBremLoss((*theElementVector)[i]->GetZ(), kineticEnergy, cut);
191
192 dedx += loss*theAtomicNumDensityVector[i];
193 }
194 // G4cout << "BR e= " << kineticEnergy << " dedx= " << dedx << G4endl;
195 if(dedx < 0.) dedx = 0.;
196 return dedx;
197}
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
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 281 of file G4MuBremsstrahlungModel.cc.

286{
287 G4double dxsection = 0.;
288
289 if(gammaEnergy > tkin) { return dxsection; }
290
291 G4double E = tkin + mass ;
292 G4double v = gammaEnergy/E ;
293 G4double delta = 0.5*mass*mass*v/(E-gammaEnergy) ;
294 G4double rab0 = delta*sqrte ;
295
296 G4int iz = G4lrint(Z);
297 if(iz < 1) { iz = 1; }
298 else if(iz > 92) { iz = 92; }
299
300 G4double z13 = 1.0/nist->GetZ13(iz);
301 G4double dnstar = fDN[iz];
302
303 G4double b,b1;
304
305 if(1 == iz) {
306 b = bh;
307 b1 = bh1;
308 } else {
309 b = btf;
310 b1 = btf1;
311 }
312
313 // nucleus contribution logarithm
314 G4double rab1=b*z13;
315 G4double fn=G4Log(rab1/(dnstar*(electron_mass_c2+rab0*rab1))*
316 (mass+delta*(dnstar*sqrte-2.))) ;
317 if(fn <0.) { fn = 0.; }
318 // electron contribution logarithm
319 G4double epmax1=E/(1.+0.5*mass*rmass/E) ;
320 G4double fe=0.;
321 if(gammaEnergy<epmax1)
322 {
323 G4double rab2=b1*z13*z13 ;
324 fe=G4Log(rab2*mass/((1.+delta*rmass/(electron_mass_c2*sqrte))*
325 (electron_mass_c2+rab0*rab2))) ;
326 if(fe<0.) { fe=0.; }
327 }
328
329 dxsection = coeff*(1.-v*(1. - 0.75*v))*Z*(fn*Z + fe)/gammaEnergy;
330
331 return dxsection;
332}
G4fissionEvent * fe
G4double G4Log(G4double x)
Definition: G4Log.hh:226
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 238 of file G4MuBremsstrahlungModel.cc.

242{
243 G4double totalEnergy = tkin + mass;
244 static const G4double ak1 = 2.3;
245 static const G4int k2 = 4;
246 G4double cross = 0.;
247
248 if(cut >= tkin) return cross;
249
250 G4double vcut = cut/totalEnergy;
251 G4double vmax = tkin/totalEnergy;
252
253 G4double aaa = G4Log(vcut);
254 G4double bbb = G4Log(vmax);
255 G4int kkk = (G4int)((bbb-aaa)/ak1)+k2 ;
256 if(kkk < 1) { kkk = 1; }
257
258 G4double hhh = (bbb-aaa)/G4double(kkk);
259
260 G4double aa = aaa;
261
262 for(G4int l=0; l<kkk; l++)
263 {
264 for(G4int i=0; i<6; i++)
265 {
266 G4double ep = G4Exp(aa + xgi[i]*hhh)*totalEnergy;
267 cross += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
268 }
269 aa += hhh;
270 }
271
272 cross *=hhh;
273
274 //G4cout << "BR e= " << tkin<< " cross= " << cross/barn << G4endl;
275
276 return cross;
277}
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 201 of file G4MuBremsstrahlungModel.cc.

203{
204 G4double totalEnergy = mass + tkin;
205 static const G4double ak1 = 0.05;
206 static const G4int k2=5;
207 G4double loss = 0.;
208
209 G4double vcut = cut/totalEnergy;
210 G4double vmax = tkin/totalEnergy;
211
212 G4double aaa = 0.;
213 G4double bbb = vcut;
214 if(vcut>vmax) { bbb = vmax; }
215 G4int kkk = (G4int)((bbb-aaa)/ak1)+k2;
216 if(kkk < 1) { kkk = 1; }
217
218 G4double hhh=(bbb-aaa)/G4double(kkk);
219
220 G4double aa = aaa;
221 for(G4int l=0; l<kkk; l++)
222 {
223 for(G4int i=0; i<6; i++)
224 {
225 G4double ep = (aa + xgi[i]*hhh)*totalEnergy;
226 loss += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
227 }
228 aa += hhh;
229 }
230
231 loss *=hhh*totalEnergy ;
232
233 return loss;
234}

Referenced by ComputeDEDXPerVolume().

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 144 of file G4MuBremsstrahlungModel.cc.

146{
147 if(p) { SetParticle(p); }
148
149 // define pointer to G4ParticleChange
151
152 if(IsMaster() && p == particle && lowestKinEnergy < HighEnergyLimit()) {
154 }
155}
G4bool IsMaster() const
Definition: G4VEmModel.hh:736
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:118

◆ InitialiseLocal()

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

Reimplemented from G4VEmModel.

Definition at line 159 of file G4MuBremsstrahlungModel.cc.

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

◆ MinEnergyCut()

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

Reimplemented from G4VEmModel.

Definition at line 127 of file G4MuBremsstrahlungModel.cc.

129{
130 return minThreshold;
131}

◆ MinPrimaryEnergy()

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

Reimplemented from G4VEmModel.

Definition at line 135 of file G4MuBremsstrahlungModel.cc.

138{
139 return std::max(lowestKinEnergy,cut);
140}

◆ 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 359 of file G4MuBremsstrahlungModel.cc.

365{
366 G4double kineticEnergy = dp->GetKineticEnergy();
367 // check against insufficient energy
368 G4double tmax = std::min(kineticEnergy, maxEnergy);
369 G4double tmin = std::min(kineticEnergy, minEnergy);
370 if(tmin < minThreshold) tmin = minThreshold;
371 if(tmin >= tmax) return;
372
373 // ===== sampling of energy transfer ======
374
375 G4ParticleMomentum partDirection = dp->GetMomentumDirection();
376
377 // select randomly one element constituing the material
378 const G4Element* anElement = SelectRandomAtom(couple,particle,kineticEnergy);
379 G4double Z = anElement->GetZ();
380
381 G4double func1 = tmin*
382 ComputeDMicroscopicCrossSection(kineticEnergy,Z,tmin);
383
384 G4double lnepksi, epksi;
385 G4double func2;
386
387 G4double xmin = G4Log(tmin/CLHEP::MeV);
388 G4double xmax = G4Log(kineticEnergy/tmin);
389
390 do {
391 lnepksi = xmin + G4UniformRand()*xmax;
392 epksi = CLHEP::MeV*G4Exp(lnepksi);
393 func2 = epksi*ComputeDMicroscopicCrossSection(kineticEnergy,Z,epksi);
394
395 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
396 } while(func2 < func1*G4UniformRand());
397
398 G4double gEnergy = epksi;
399
400 // ===== sample angle =====
401
402 //
403 // angles of the emitted gamma using general interface
404
405 G4ThreeVector gamDir =
406 GetAngularDistribution()->SampleDirection(dp, gEnergy, Z,
407 couple->GetMaterial());
408 // create G4DynamicParticle object for the Gamma
409 G4DynamicParticle* gamma = new G4DynamicParticle(theGamma, gamDir, gEnergy);
410 vdp->push_back(gamma);
411 // compute post-interaction kinematics of primary e-/e+ based on
412 // energy-momentum conservation
413 const G4double totMomentum =
414 std::sqrt(kineticEnergy*(kineticEnergy + 2.0*mass));
415 G4ThreeVector dir =
416 (totMomentum*dp->GetMomentumDirection()-gEnergy*gamDir).unit();
417 const G4double finalE = kineticEnergy - gEnergy;
418 // if secondary gamma energy is higher than threshold(very high by default)
419 // then stop tracking the primary particle and create new secondary e-/e+
420 // instead of the primary one
421 if (gEnergy > SecondaryThreshold()) {
424 G4DynamicParticle* newdp = new G4DynamicParticle(particle, dir, finalE);
425 vdp->push_back(newdp);
426 } else { // continue tracking the primary e-/e+ otherwise
429 }
430}
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:52
const G4ThreeVector & GetMomentumDirection() 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 * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:570
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:680
void ProposeTrackStatus(G4TrackStatus status)

◆ 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 142 of file G4MuBremsstrahlungModel.hh.

Referenced by ComputeDMicroscopicCrossSection().

◆ btf

G4double G4MuBremsstrahlungModel::btf
protected

◆ btf1

G4double G4MuBremsstrahlungModel::btf1
protected

Definition at line 144 of file G4MuBremsstrahlungModel.hh.

Referenced by ComputeDMicroscopicCrossSection().

◆ cc

G4double G4MuBremsstrahlungModel::cc
protected

Definition at line 138 of file G4MuBremsstrahlungModel.hh.

Referenced by G4MuBremsstrahlungModel(), and SetParticle().

◆ coeff

◆ fDN

G4double G4MuBremsstrahlungModel::fDN = {0.0}
staticprotected

◆ fParticleChange

G4ParticleChangeForLoss* G4MuBremsstrahlungModel::fParticleChange
protected

Definition at line 147 of file G4MuBremsstrahlungModel.hh.

Referenced by Initialise(), and SampleSecondaries().

◆ lowestKinEnergy

◆ mass

◆ minThreshold

G4double G4MuBremsstrahlungModel::minThreshold
protected

◆ nist

◆ particle

◆ rmass

G4double G4MuBremsstrahlungModel::rmass
protected

◆ sqrte

G4double G4MuBremsstrahlungModel::sqrte
protected

◆ theGamma

G4ParticleDefinition* G4MuBremsstrahlungModel::theGamma
protected

Definition at line 146 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 153 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 152 of file G4MuBremsstrahlungModel.hh.

Referenced by ComputeMicroscopicCrossSection(), and ComputMuBremLoss().


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