Geant4  10.02.p03
G4MuBremsstrahlungModel Class Reference

#include <G4MuBremsstrahlungModel.hh>

Inheritance diagram for G4MuBremsstrahlungModel:
Collaboration 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 void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
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)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double)
 
- 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, const G4ParticleDefinition *, G4double)
 
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 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 *> *)
 
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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
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
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
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 SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () 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)
 

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
G4ParticleChangeForLoss * GetParticleChangeForLoss ()
 
G4ParticleChangeForGamma * GetParticleChangeForGamma ()
 
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
G4ElementDatafElementData
 
G4VParticleChange * pParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Private Member Functions

G4MuBremsstrahlungModeloperator= (const G4MuBremsstrahlungModel &right)
 
 G4MuBremsstrahlungModel (const G4MuBremsstrahlungModel &)
 

Private Attributes

G4ParticleDefinitiontheGamma
 
G4ParticleChangeForLoss * fParticleChange
 
G4double lowestKinEnergy
 
G4double minThreshold
 

Static Private Attributes

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

Additional Inherited Members

- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 72 of file G4MuBremsstrahlungModel.hh.

Constructor & Destructor Documentation

◆ G4MuBremsstrahlungModel() [1/2]

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

Definition at line 91 of file G4MuBremsstrahlungModel.cc.

93  : G4VEmModel(nam),
94  particle(0),
95  sqrte(sqrt(G4Exp(1.))),
96  bh(202.4),
97  bh1(446.),
98  btf(183.),
99  btf1(1429.),
100  fParticleChange(0),
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  }
120 
121  if(p) { SetParticle(p); }
122 }
G4ParticleDefinition * theGamma
G4double GetA27(G4int Z)
void SetParticle(const G4ParticleDefinition *)
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:69
const G4ParticleDefinition * particle
static const double GeV
Definition: G4SIunits.hh:214
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4ParticleChangeForLoss * fParticleChange
static const double keV
Definition: G4SIunits.hh:213
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ ~G4MuBremsstrahlungModel()

G4MuBremsstrahlungModel::~G4MuBremsstrahlungModel ( )
virtual

Definition at line 126 of file G4MuBremsstrahlungModel.cc.

127 {}

◆ G4MuBremsstrahlungModel() [2/2]

G4MuBremsstrahlungModel::G4MuBremsstrahlungModel ( const G4MuBremsstrahlungModel )
private

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

346 {
347  G4double cross = 0.0;
348  if (kineticEnergy <= lowestKinEnergy) return cross;
349  G4double tmax = std::min(maxEnergy, kineticEnergy);
350  G4double cut = std::min(cutEnergy, kineticEnergy);
351  if(cut < minThreshold) cut = minThreshold;
352  if (cut >= tmax) return cross;
353 
354  cross = ComputeMicroscopicCrossSection (kineticEnergy, Z, cut);
355  if(tmax < kineticEnergy) {
356  cross -= ComputeMicroscopicCrossSection(kineticEnergy, Z, tmax);
357  }
358  return cross;
359 }
G4double ComputeMicroscopicCrossSection(G4double tkin, G4double Z, G4double cut)
Float_t Z
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ ComputeDEDXPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 173 of file G4MuBremsstrahlungModel.cc.

178 {
179  G4double dedx = 0.0;
180  if (kineticEnergy <= lowestKinEnergy) { return dedx; }
181 
182  G4double tmax = kineticEnergy;
183  G4double cut = std::min(cutEnergy,tmax);
184  if(cut < minThreshold) { cut = minThreshold; }
185 
186  const G4ElementVector* theElementVector = material->GetElementVector();
187  const G4double* theAtomicNumDensityVector =
188  material->GetAtomicNumDensityVector();
189 
190  // loop for elements in the material
191  for (size_t i=0; i<material->GetNumberOfElements(); i++) {
192 
193  G4double loss =
194  ComputMuBremLoss((*theElementVector)[i]->GetZ(), kineticEnergy, cut);
195 
196  dedx += loss*theAtomicNumDensityVector[i];
197  }
198  // G4cout << "BR e= " << kineticEnergy << " dedx= " << dedx << G4endl;
199  if(dedx < 0.) dedx = 0.;
200  return dedx;
201 }
std::vector< G4Element * > G4ElementVector
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:216
G4double ComputMuBremLoss(G4double Z, G4double tkin, G4double cut)
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ ComputeDMicroscopicCrossSection()

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

Reimplemented in G4hBremsstrahlungModel.

Definition at line 285 of file G4MuBremsstrahlungModel.cc.

290 {
291  G4double dxsection = 0.;
292 
293  if(gammaEnergy > tkin) { return dxsection; }
294 
295  G4double E = tkin + mass ;
296  G4double v = gammaEnergy/E ;
297  G4double delta = 0.5*mass*mass*v/(E-gammaEnergy) ;
298  G4double rab0 = delta*sqrte ;
299 
300  G4int iz = G4lrint(Z);
301  if(iz < 1) { iz = 1; }
302  else if(iz > 92) { iz = 92; }
303 
304  G4double z13 = 1.0/nist->GetZ13(iz);
305  G4double dnstar = fDN[iz];
306 
307  G4double b,b1;
308 
309  if(1 == iz) {
310  b = bh;
311  b1 = bh1;
312  } else {
313  b = btf;
314  b1 = btf1;
315  }
316 
317  // nucleus contribution logarithm
318  G4double rab1=b*z13;
319  G4double fn=G4Log(rab1/(dnstar*(electron_mass_c2+rab0*rab1))*
320  (mass+delta*(dnstar*sqrte-2.))) ;
321  if(fn <0.) { fn = 0.; }
322  // electron contribution logarithm
323  G4double epmax1=E/(1.+0.5*mass*rmass/E) ;
324  G4double fe=0.;
325  if(gammaEnergy<epmax1)
326  {
327  G4double rab2=b1*z13*z13 ;
328  fe=G4Log(rab2*mass/((1.+delta*rmass/(electron_mass_c2*sqrte))*
329  (electron_mass_c2+rab0*rab2))) ;
330  if(fe<0.) { fe=0.; }
331  }
332 
333  dxsection = coeff*(1.-v*(1. - 0.75*v))*Z*(fn*Z + fe)/gammaEnergy;
334 
335  return dxsection;
336 }
G4double GetZ13(G4double Z)
int G4int
Definition: G4Types.hh:78
Float_t Z
G4double iz
Definition: TRTMaterials.hh:39
float electron_mass_c2
Definition: hepunit.py:274
G4double G4Log(G4double x)
Definition: G4Log.hh:230
int G4lrint(double ad)
Definition: templates.hh:163
static const G4double b1
double G4double
Definition: G4Types.hh:76
G4fissionEvent * fe
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeMicroscopicCrossSection()

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

Definition at line 242 of file G4MuBremsstrahlungModel.cc.

246 {
247  G4double totalEnergy = tkin + mass;
248  static const G4double ak1 = 2.3;
249  static const G4int k2 = 4;
250  G4double cross = 0.;
251 
252  if(cut >= tkin) return cross;
253 
254  G4double vcut = cut/totalEnergy;
255  G4double vmax = tkin/totalEnergy;
256 
257  G4double aaa = G4Log(vcut);
258  G4double bbb = G4Log(vmax);
259  G4int kkk = (G4int)((bbb-aaa)/ak1)+k2 ;
260  if(kkk < 1) { kkk = 1; }
261 
262  G4double hhh = (bbb-aaa)/G4double(kkk);
263 
264  G4double aa = aaa;
265 
266  for(G4int l=0; l<kkk; l++)
267  {
268  for(G4int i=0; i<6; i++)
269  {
270  G4double ep = G4Exp(aa + xgi[i]*hhh)*totalEnergy;
271  cross += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
272  }
273  aa += hhh;
274  }
275 
276  cross *=hhh;
277 
278  //G4cout << "BR e= " << tkin<< " cross= " << cross/barn << G4endl;
279 
280  return cross;
281 }
int G4int
Definition: G4Types.hh:78
virtual G4double ComputeDMicroscopicCrossSection(G4double tkin, G4double Z, G4double gammaEnergy)
static const G4double xgi[6]
static const G4double wgi[6]
Float_t Z
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
double G4double
Definition: G4Types.hh:76
TH1F * hhh
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputMuBremLoss()

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

Definition at line 205 of file G4MuBremsstrahlungModel.cc.

207 {
208  G4double totalEnergy = mass + tkin;
209  static const G4double ak1 = 0.05;
210  static const G4int k2=5;
211  G4double loss = 0.;
212 
213  G4double vcut = cut/totalEnergy;
214  G4double vmax = tkin/totalEnergy;
215 
216  G4double aaa = 0.;
217  G4double bbb = vcut;
218  if(vcut>vmax) { bbb = vmax; }
219  G4int kkk = (G4int)((bbb-aaa)/ak1)+k2;
220  if(kkk < 1) { kkk = 1; }
221 
222  G4double hhh=(bbb-aaa)/G4double(kkk);
223 
224  G4double aa = aaa;
225  for(G4int l=0; l<kkk; l++)
226  {
227  for(G4int i=0; i<6; i++)
228  {
229  G4double ep = (aa + xgi[i]*hhh)*totalEnergy;
230  loss += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
231  }
232  aa += hhh;
233  }
234 
235  loss *=hhh*totalEnergy ;
236 
237  return loss;
238 }
int G4int
Definition: G4Types.hh:78
virtual G4double ComputeDMicroscopicCrossSection(G4double tkin, G4double Z, G4double gammaEnergy)
static const G4double xgi[6]
static const G4double wgi[6]
Float_t Z
double G4double
Definition: G4Types.hh:76
TH1F * hhh
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 148 of file G4MuBremsstrahlungModel.cc.

150 {
151  if(p) { SetParticle(p); }
152 
153  // define pointer to G4ParticleChange
155 
156  if(IsMaster() && p == particle && lowestKinEnergy < HighEnergyLimit()) {
157  InitialiseElementSelectors(p, cuts);
158  }
159 }
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:120
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148
G4bool IsMaster() const
Definition: G4VEmModel.hh:718
void SetParticle(const G4ParticleDefinition *)
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
const G4ParticleDefinition * particle
G4ParticleChangeForLoss * fParticleChange
Here is the call graph for this function:
Here is the caller graph for this function:

◆ InitialiseLocal()

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

Reimplemented from G4VEmModel.

Definition at line 163 of file G4MuBremsstrahlungModel.cc.

165 {
166  if(p == particle && lowestKinEnergy < HighEnergyLimit()) {
168  }
169 }
void SetElementSelectors(std::vector< G4EmElementSelector *> *)
Definition: G4VEmModel.hh:810
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
const G4ParticleDefinition * particle
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:802
Here is the call graph for this function:

◆ MinEnergyCut()

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

Reimplemented from G4VEmModel.

Definition at line 131 of file G4MuBremsstrahlungModel.cc.

133 {
134  return minThreshold;
135 }

◆ MinPrimaryEnergy()

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

Reimplemented from G4VEmModel.

Definition at line 139 of file G4MuBremsstrahlungModel.cc.

142 {
143  return std::max(lowestKinEnergy,cut);
144 }

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 363 of file G4MuBremsstrahlungModel.cc.

369 {
370  G4double kineticEnergy = dp->GetKineticEnergy();
371  // check against insufficient energy
372  G4double tmax = std::min(kineticEnergy, maxEnergy);
373  G4double tmin = std::min(kineticEnergy, minEnergy);
374  if(tmin < minThreshold) tmin = minThreshold;
375  if(tmin >= tmax) return;
376 
377  // ===== sampling of energy transfer ======
378 
379  G4ParticleMomentum partDirection = dp->GetMomentumDirection();
380 
381  // select randomly one element constituing the material
382  const G4Element* anElement = SelectRandomAtom(couple,particle,kineticEnergy);
383  G4double Z = anElement->GetZ();
384 
385  G4double totalEnergy = kineticEnergy + mass;
386  G4double totalMomentum = sqrt(kineticEnergy*(kineticEnergy + 2.0*mass));
387 
388  G4double func1 = tmin*
389  ComputeDMicroscopicCrossSection(kineticEnergy,Z,tmin);
390 
391  G4double lnepksi, epksi;
392  G4double func2;
393 
394  G4double xmin = G4Log(tmin/MeV);
395  G4double xmax = G4Log(kineticEnergy/tmin);
396 
397  do {
398  lnepksi = xmin + G4UniformRand()*xmax;
399  epksi = MeV*G4Exp(lnepksi);
400  func2 = epksi*ComputeDMicroscopicCrossSection(kineticEnergy,Z,epksi);
401 
402  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
403  } while(func2 < func1*G4UniformRand());
404 
405  G4double gEnergy = epksi;
406 
407  // ===== sample angle =====
408 
409  G4double gam = totalEnergy/mass;
410  G4double rmax = gam*std::min(1.0, totalEnergy/gEnergy - 1.0);
411  G4double rmax2= rmax*rmax;
412  G4double x = G4UniformRand()*rmax2/(1.0 + rmax2);
413 
414  G4double theta = sqrt(x/(1.0 - x))/gam;
415  G4double sint = sin(theta);
416  G4double phi = twopi * G4UniformRand() ;
417  G4double dirx = sint*cos(phi), diry = sint*sin(phi), dirz = cos(theta) ;
418 
419  G4ThreeVector gDirection(dirx, diry, dirz);
420  gDirection.rotateUz(partDirection);
421 
422  partDirection *= totalMomentum;
423  partDirection -= gEnergy*gDirection;
424  partDirection = partDirection.unit();
425 
426  // primary change
427  kineticEnergy -= gEnergy;
428  fParticleChange->SetProposedKineticEnergy(kineticEnergy);
429  fParticleChange->SetProposedMomentumDirection(partDirection);
430 
431  // save secondary
432  G4DynamicParticle* aGamma =
433  new G4DynamicParticle(theGamma,gDirection,gEnergy);
434  vdp->push_back(aGamma);
435 }
static const double MeV
Definition: G4SIunits.hh:211
int func1(int i)
Definition: XFunc.cc:40
G4ParticleDefinition * theGamma
int func2(int i)
Definition: XFunc.cc:51
virtual G4double ComputeDMicroscopicCrossSection(G4double tkin, G4double Z, G4double gammaEnergy)
G4double GetKineticEnergy() const
#define G4UniformRand()
Definition: Randomize.hh:97
Float_t Z
const G4ParticleDefinition * particle
Hep3Vector unit() const
static const double twopi
Definition: G4SIunits.hh:75
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4ParticleChangeForLoss * fParticleChange
const G4ThreeVector & GetMomentumDirection() const
double G4double
Definition: G4Types.hh:76
G4double GetZ() const
Definition: G4Element.hh:131
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:544
Here is the call graph for this function:

◆ SetLowestKineticEnergy()

void G4MuBremsstrahlungModel::SetLowestKineticEnergy ( G4double  e)
inline

Definition at line 163 of file G4MuBremsstrahlungModel.hh.

164 {
165  lowestKinEnergy = e;
166 }

◆ SetParticle()

void G4MuBremsstrahlungModel::SetParticle ( const G4ParticleDefinition p)
inlineprotected

Definition at line 171 of file G4MuBremsstrahlungModel.hh.

172 {
173  if(!particle) {
174  particle = p;
175  mass = particle->GetPDGMass();
179  }
180 }
static const double fine_structure_const
static const double classic_electr_radius
const G4ParticleDefinition * particle
static const double electron_mass_c2
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ bh

G4double G4MuBremsstrahlungModel::bh
protected

Definition at line 142 of file G4MuBremsstrahlungModel.hh.

◆ bh1

G4double G4MuBremsstrahlungModel::bh1
protected

Definition at line 143 of file G4MuBremsstrahlungModel.hh.

◆ btf

G4double G4MuBremsstrahlungModel::btf
protected

Definition at line 144 of file G4MuBremsstrahlungModel.hh.

◆ btf1

G4double G4MuBremsstrahlungModel::btf1
protected

Definition at line 145 of file G4MuBremsstrahlungModel.hh.

◆ cc

G4double G4MuBremsstrahlungModel::cc
protected

Definition at line 139 of file G4MuBremsstrahlungModel.hh.

◆ coeff

G4double G4MuBremsstrahlungModel::coeff
protected

Definition at line 140 of file G4MuBremsstrahlungModel.hh.

◆ fDN

G4double G4MuBremsstrahlungModel::fDN = {0.0}
staticprivate

Definition at line 158 of file G4MuBremsstrahlungModel.hh.

◆ fParticleChange

G4ParticleChangeForLoss* G4MuBremsstrahlungModel::fParticleChange
private

Definition at line 150 of file G4MuBremsstrahlungModel.hh.

◆ lowestKinEnergy

G4double G4MuBremsstrahlungModel::lowestKinEnergy
private

Definition at line 152 of file G4MuBremsstrahlungModel.hh.

◆ mass

G4double G4MuBremsstrahlungModel::mass
protected

Definition at line 137 of file G4MuBremsstrahlungModel.hh.

◆ minThreshold

G4double G4MuBremsstrahlungModel::minThreshold
private

Definition at line 153 of file G4MuBremsstrahlungModel.hh.

◆ nist

G4NistManager* G4MuBremsstrahlungModel::nist
protected

Definition at line 136 of file G4MuBremsstrahlungModel.hh.

◆ particle

const G4ParticleDefinition* G4MuBremsstrahlungModel::particle
protected

Definition at line 135 of file G4MuBremsstrahlungModel.hh.

◆ rmass

G4double G4MuBremsstrahlungModel::rmass
protected

Definition at line 138 of file G4MuBremsstrahlungModel.hh.

◆ sqrte

G4double G4MuBremsstrahlungModel::sqrte
protected

Definition at line 141 of file G4MuBremsstrahlungModel.hh.

◆ theGamma

G4ParticleDefinition* G4MuBremsstrahlungModel::theGamma
private

Definition at line 149 of file G4MuBremsstrahlungModel.hh.

◆ wgi

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

Definition at line 156 of file G4MuBremsstrahlungModel.hh.

◆ xgi

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

Definition at line 155 of file G4MuBremsstrahlungModel.hh.


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