Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/muons/src/G4MuBremsstrahlungModel.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  1 //
  2 // ********************************************************************
  3 // * License and Disclaimer                                           *
  4 // *                                                                  *
  5 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
  6 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
  7 // * conditions of the Geant4 Software License,  included in the file *
  8 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
  9 // * include a list of copyright holders.                             *
 10 // *                                                                  *
 11 // * Neither the authors of this software system, nor their employing *
 12 // * institutes,nor the agencies providing financial support for this *
 13 // * work  make  any representation or  warranty, express or implied, *
 14 // * regarding  this  software system or assume any liability for its *
 15 // * use.  Please see the license in the file  LICENSE  and URL above *
 16 // * for the full disclaimer and the limitation of liability.         *
 17 // *                                                                  *
 18 // * This  code  implementation is the result of  the  scientific and *
 19 // * technical work of the GEANT4 collaboration.                      *
 20 // * By using,  copying,  modifying or  distributing the software (or *
 21 // * any work based  on the software)  you  agree  to acknowledge its *
 22 // * use  in  resulting  scientific  publications,  and indicate your *
 23 // * acceptance of all terms of the Geant4 Software license.          *
 24 // ********************************************************************
 25 //
 26 //
 27 // -------------------------------------------------------------------
 28 //
 29 // GEANT4 Class file
 30 //
 31 //
 32 // File name:     G4MuBremsstrahlungModel
 33 //
 34 // Author:        Vladimir Ivanchenko on base of Laszlo Urban code
 35 //
 36 // Creation date: 24.06.2002
 37 //
 38 // Modifications:
 39 //
 40 // 04-12-02 Change G4DynamicParticle constructor in PostStepDoIt (V.Ivanchenko)
 41 // 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
 42 // 24-01-03 Fix for compounds (V.Ivanchenko)
 43 // 27-01-03 Make models region aware (V.Ivanchenko)
 44 // 13-02-03 Add name (V.Ivanchenko)
 45 // 10-02-04 Add lowestKinEnergy (V.Ivanchenko)
 46 // 08-04-05 Major optimisation of internal interfaces (V.Ivanchenko)
 47 // 03-08-05 Angular correlations according to PRM (V.Ivanchenko)
 48 // 13-02-06 add ComputeCrossSectionPerAtom (mma)
 49 // 21-03-06 Fix problem of initialisation in case when cuts are not defined (VI)
 50 // 07-11-07 Improve sampling of final state (A.Bogdanov)
 51 // 28-02-08 Use precomputed Z^1/3 and Log(A) (V.Ivanchenko)
 52 // 31-05-13 Use element selectors instead of local data structure (V.Ivanchenko)
 53 //
 54 // -------------------------------------------------------------------
 55 //
 56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 57 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 58 
 59 #include "G4MuBremsstrahlungModel.hh"
 60 #include "G4PhysicalConstants.hh"
 61 #include "G4SystemOfUnits.hh"
 62 #include "G4Gamma.hh"
 63 #include "G4MuonMinus.hh"
 64 #include "G4MuonPlus.hh"
 65 #include "Randomize.hh"
 66 #include "G4Material.hh"
 67 #include "G4Element.hh"
 68 #include "G4ElementVector.hh"
 69 #include "G4ProductionCutsTable.hh"
 70 #include "G4ModifiedMephi.hh"
 71 #include "G4ParticleChangeForLoss.hh"
 72 #include "G4Log.hh"
 73 
 74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 75 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 76 
 77 const G4double G4MuBremsstrahlungModel::xgi[] = 
 78   {0.03377,0.16940,0.38069,0.61931,0.83060,0.96623};
 79 const G4double G4MuBremsstrahlungModel::wgi[] = 
 80   {0.08566,0.18038,0.23396,0.23396,0.18038,0.08566};
 81 G4double G4MuBremsstrahlungModel::fDN[] = {0.0};
 82 
 83 G4MuBremsstrahlungModel::G4MuBremsstrahlungModel(const G4ParticleDefinition* p,
 84                                                  const G4String& nam)
 85   : G4VEmModel(nam),
 86     sqrte(std::sqrt(G4Exp(1.))),
 87     lowestKinEnergy(0.1*CLHEP::GeV),
 88     minThreshold(0.9*CLHEP::keV)
 89 {
 90   theGamma = G4Gamma::Gamma();
 91   nist = G4NistManager::Instance();  
 92 
 93   SetAngularDistribution(new G4ModifiedMephi());
 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 }
106 
107 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
108 
109 G4double G4MuBremsstrahlungModel::MinEnergyCut(const G4ParticleDefinition*,
110                                                const G4MaterialCutsCouple*)
111 {
112   return minThreshold;
113 }
114 
115 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
116 
117 G4double G4MuBremsstrahlungModel::MinPrimaryEnergy(const G4Material*,
118                                                    const G4ParticleDefinition*,
119                                                    G4double cut)
120 {
121   return std::max(lowestKinEnergy, cut);
122 }
123 
124 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
125 
126 void G4MuBremsstrahlungModel::SetParticle(const G4ParticleDefinition* p)
127 {
128   if(nullptr == particle) {
129     particle = p;
130     mass = particle->GetPDGMass();
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 }
136 
137 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
138 
139 void G4MuBremsstrahlungModel::Initialise(const G4ParticleDefinition* p,
140                                          const G4DataVector& cuts)
141 {
142   SetParticle(p);
143   if(nullptr == fParticleChange) { 
144     fParticleChange = GetParticleChangeForLoss(); 
145   }
146   if(IsMaster() && p == particle && lowestKinEnergy < HighEnergyLimit()) { 
147     InitialiseElementSelectors(p, cuts); 
148   }
149 }
150 
151 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
152 
153 void G4MuBremsstrahlungModel::InitialiseLocal(const G4ParticleDefinition* p,
154                                               G4VEmModel* masterModel)
155 {
156   if(p == particle && lowestKinEnergy < HighEnergyLimit()) {
157     SetElementSelectors(masterModel->GetElementSelectors());
158   }
159 }
160 
161 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
162 
163 G4double G4MuBremsstrahlungModel::ComputeDEDXPerVolume(
164                                               const G4Material* material,
165                                               const G4ParticleDefinition*,
166                                                     G4double kineticEnergy,
167                                                     G4double cutEnergy)
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 }
189 
190 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
191 
192 G4double G4MuBremsstrahlungModel::ComputMuBremLoss(G4double Z,
193                                                    G4double tkin, G4double cut)
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 }
218 
219 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
220 
221 G4double G4MuBremsstrahlungModel::ComputeMicroscopicCrossSection(
222                                            G4double tkin,
223                                            G4double Z,
224                                            G4double cut)
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 }
256 
257 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
258 
259 G4double G4MuBremsstrahlungModel::ComputeDMicroscopicCrossSection(
260                                            G4double tkin,
261                                            G4double Z,
262                                            G4double gammaEnergy)
263 //  differential cross section
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 }
308 
309 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
310 
311 G4double G4MuBremsstrahlungModel::ComputeCrossSectionPerAtom(
312                                            const G4ParticleDefinition*,
313                                                  G4double kineticEnergy,
314                                                  G4double Z, G4double,
315                                                  G4double cutEnergy,
316                                                  G4double maxEnergy)
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 }
331 
332 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
333 
334 void G4MuBremsstrahlungModel::SampleSecondaries(
335                               std::vector<G4DynamicParticle*>* vdp,
336                               const G4MaterialCutsCouple* couple,
337                               const G4DynamicParticle* dp,
338                               G4double minEnergy,
339                               G4double maxEnergy)
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()) {
391     fParticleChange->ProposeTrackStatus(fStopAndKill);
392     fParticleChange->SetProposedKineticEnergy(0.0);
393     G4DynamicParticle* newdp = new G4DynamicParticle(particle, dir, finalE);
394     vdp->push_back(newdp);
395   } else {
396     // continue tracking the primary e-/e+ otherwise
397     fParticleChange->SetProposedMomentumDirection(dir);
398     fParticleChange->SetProposedKineticEnergy(finalE);
399   }
400 }
401 
402 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
403