Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/muons/src/G4MuBetheBlochModel.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 ]

Diff markup

Differences between /processes/electromagnetic/muons/src/G4MuBetheBlochModel.cc (Version 11.3.0) and /processes/electromagnetic/muons/src/G4MuBetheBlochModel.cc (Version 10.1)


  1 //                                                  1 //
  2 // *******************************************      2 // ********************************************************************
  3 // * License and Disclaimer                         3 // * License and Disclaimer                                           *
  4 // *                                                4 // *                                                                  *
  5 // * The  Geant4 software  is  copyright of th      5 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
  6 // * the Geant4 Collaboration.  It is provided      6 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
  7 // * conditions of the Geant4 Software License      7 // * conditions of the Geant4 Software License,  included in the file *
  8 // * LICENSE and available at  http://cern.ch/      8 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
  9 // * include a list of copyright holders.           9 // * include a list of copyright holders.                             *
 10 // *                                               10 // *                                                                  *
 11 // * Neither the authors of this software syst     11 // * Neither the authors of this software system, nor their employing *
 12 // * institutes,nor the agencies providing fin     12 // * institutes,nor the agencies providing financial support for this *
 13 // * work  make  any representation or  warran     13 // * work  make  any representation or  warranty, express or implied, *
 14 // * regarding  this  software system or assum     14 // * regarding  this  software system or assume any liability for its *
 15 // * use.  Please see the license in the file      15 // * use.  Please see the license in the file  LICENSE  and URL above *
 16 // * for the full disclaimer and the limitatio     16 // * for the full disclaimer and the limitation of liability.         *
 17 // *                                               17 // *                                                                  *
 18 // * This  code  implementation is the result      18 // * This  code  implementation is the result of  the  scientific and *
 19 // * technical work of the GEANT4 collaboratio     19 // * technical work of the GEANT4 collaboration.                      *
 20 // * By using,  copying,  modifying or  distri     20 // * By using,  copying,  modifying or  distributing the software (or *
 21 // * any work based  on the software)  you  ag     21 // * any work based  on the software)  you  agree  to acknowledge its *
 22 // * use  in  resulting  scientific  publicati     22 // * use  in  resulting  scientific  publications,  and indicate your *
 23 // * acceptance of all terms of the Geant4 Sof     23 // * acceptance of all terms of the Geant4 Software license.          *
 24 // *******************************************     24 // ********************************************************************
 25 //                                                 25 //
                                                   >>  26 // $Id: G4MuBetheBlochModel.cc 85023 2014-10-23 09:56:39Z gcosmo $
 26 //                                                 27 //
 27 // -------------------------------------------     28 // -------------------------------------------------------------------
 28 //                                                 29 //
 29 // GEANT4 Class header file                        30 // GEANT4 Class header file
 30 //                                                 31 //
 31 //                                                 32 //
 32 // File name:     G4MuBetheBlochModel              33 // File name:     G4MuBetheBlochModel
 33 //                                                 34 //
 34 // Author:        Vladimir Ivanchenko on base      35 // Author:        Vladimir Ivanchenko on base of Laszlo Urban code
 35 //                                                 36 //
 36 // Creation date: 09.08.2002                       37 // Creation date: 09.08.2002
 37 //                                                 38 //
 38 // Modifications:                                  39 // Modifications:
 39 //                                                 40 //
 40 // 04-12-02 Fix problem of G4DynamicParticle c     41 // 04-12-02 Fix problem of G4DynamicParticle constructor (V.Ivanchenko)
 41 // 23-12-02 Change interface in order to move      42 // 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
 42 // 27-01-03 Make models region aware (V.Ivanch     43 // 27-01-03 Make models region aware (V.Ivanchenko)
 43 // 13-02-03 Add name (V.Ivanchenko)                44 // 13-02-03 Add name (V.Ivanchenko)
 44 // 10-02-04 Calculation of radiative correctio     45 // 10-02-04 Calculation of radiative corrections using R.Kokoulin model (V.I)
 45 // 08-04-05 Major optimisation of internal int     46 // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
 46 // 12-04-05 Add usage of G4EmCorrections (V.Iv     47 // 12-04-05 Add usage of G4EmCorrections (V.Ivanchenko)
 47 // 13-02-06 ComputeCrossSectionPerElectron, Co     48 // 13-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
 48 //                                                 49 //
 49                                                    50 
 50 //                                                 51 //
 51 // -------------------------------------------     52 // -------------------------------------------------------------------
 52 //                                                 53 //
 53                                                    54 
 54 //....oooOO0OOooo........oooOO0OOooo........oo     55 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 55 //....oooOO0OOooo........oooOO0OOooo........oo     56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 56                                                    57 
 57 #include "G4MuBetheBlochModel.hh"                  58 #include "G4MuBetheBlochModel.hh"
 58 #include "G4PhysicalConstants.hh"                  59 #include "G4PhysicalConstants.hh"
 59 #include "G4SystemOfUnits.hh"                      60 #include "G4SystemOfUnits.hh"
 60 #include "Randomize.hh"                            61 #include "Randomize.hh"
 61 #include "G4Electron.hh"                           62 #include "G4Electron.hh"
 62 #include "G4LossTableManager.hh"                   63 #include "G4LossTableManager.hh"
 63 #include "G4EmCorrections.hh"                      64 #include "G4EmCorrections.hh"
 64 #include "G4ParticleChangeForLoss.hh"              65 #include "G4ParticleChangeForLoss.hh"
 65 #include "G4Log.hh"                                66 #include "G4Log.hh"
 66 #include "G4Exp.hh"                                67 #include "G4Exp.hh"
 67 #include "G4DeltaAngle.hh"                     << 
 68                                                    68 
 69 G4double G4MuBetheBlochModel::xgi[]={ 0.0199,      69 G4double G4MuBetheBlochModel::xgi[]={ 0.0199, 0.1017, 0.2372, 0.4083, 0.5917,
 70                                       0.7628,      70                                       0.7628, 0.8983, 0.9801 };
 71                                                <<  71               
 72 G4double G4MuBetheBlochModel::wgi[]={ 0.0506,      72 G4double G4MuBetheBlochModel::wgi[]={ 0.0506, 0.1112, 0.1569, 0.1813, 0.1813,
 73                                       0.1569,      73                                       0.1569, 0.1112, 0.0506 };
 74                                                    74 
 75 //....oooOO0OOooo........oooOO0OOooo........oo     75 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 76                                                    76 
                                                   >>  77 using namespace std;
                                                   >>  78 
 77 G4MuBetheBlochModel::G4MuBetheBlochModel(const     79 G4MuBetheBlochModel::G4MuBetheBlochModel(const G4ParticleDefinition* p,
 78                                          const     80                                          const G4String& nam)
 79   : G4VEmModel(nam),                               81   : G4VEmModel(nam),
 80     limitRadCorrection(250.*CLHEP::MeV),       <<  82   particle(0),
 81     limitKinEnergy(100.*CLHEP::keV),           <<  83   limitKinEnergy(100.*keV),
 82     logLimitKinEnergy(G4Log(limitKinEnergy)),  <<  84   logLimitKinEnergy(G4Log(limitKinEnergy)),
 83     twoln10(2.0*G4Log(10.0)),                  <<  85   twoln10(2.0*G4Log(10.0)),
 84     alphaprime(CLHEP::fine_structure_const/CLH <<  86   bg2lim(0.0169),
                                                   >>  87   taulim(8.4146e-3),
                                                   >>  88   alphaprime(fine_structure_const/twopi)
 85 {                                                  89 {
 86   theElectron = G4Electron::Electron();            90   theElectron = G4Electron::Electron();
 87   corr = G4LossTableManager::Instance()->EmCor     91   corr = G4LossTableManager::Instance()->EmCorrections();
 88   if(nullptr != p) { SetParticle(p); }         <<  92   fParticleChange = 0;
                                                   >>  93 
                                                   >>  94   // initial initialisation of memeber should be overwritten
                                                   >>  95   // by SetParticle
                                                   >>  96   mass = massSquare = ratio = 1.0;
                                                   >>  97 
                                                   >>  98   if(p) { SetParticle(p); }
 89 }                                                  99 }
 90                                                   100 
 91 //....oooOO0OOooo........oooOO0OOooo........oo    101 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 92                                                   102 
                                                   >> 103 G4MuBetheBlochModel::~G4MuBetheBlochModel()
                                                   >> 104 {}
                                                   >> 105 
                                                   >> 106 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
                                                   >> 107 
 93 G4double G4MuBetheBlochModel::MinEnergyCut(con    108 G4double G4MuBetheBlochModel::MinEnergyCut(const G4ParticleDefinition*,
 94                                            con    109                                            const G4MaterialCutsCouple* couple)
 95 {                                                 110 {
 96   return couple->GetMaterial()->GetIonisation(    111   return couple->GetMaterial()->GetIonisation()->GetMeanExcitationEnergy();
 97 }                                                 112 }
 98                                                   113 
 99 //....oooOO0OOooo........oooOO0OOooo........oo    114 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
100                                                   115 
101 G4double G4MuBetheBlochModel::MaxSecondaryEner    116 G4double G4MuBetheBlochModel::MaxSecondaryEnergy(const G4ParticleDefinition*,
102                                                << 117              G4double kinEnergy) 
103 {                                                 118 {
104   G4double tau  = kinEnergy/mass;                 119   G4double tau  = kinEnergy/mass;
105   G4double tmax = 2.0*CLHEP::electron_mass_c2* << 120   G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
106                   (1. + 2.0*(tau + 1.)*ratio +    121                   (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
107   return tmax;                                    122   return tmax;
108 }                                                 123 }
109                                                   124 
110 //....oooOO0OOooo........oooOO0OOooo........oo    125 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
111                                                   126 
112 void G4MuBetheBlochModel::SetParticle(const G4 << 
113 {                                              << 
114   if(nullptr == particle) {                    << 
115     particle = p;                              << 
116     mass = particle->GetPDGMass();             << 
117     massSquare = mass*mass;                    << 
118     ratio = CLHEP::electron_mass_c2/mass;      << 
119   }                                            << 
120 }                                              << 
121                                                << 
122 //....oooOO0OOooo........oooOO0OOooo........oo << 
123                                                << 
124 void G4MuBetheBlochModel::Initialise(const G4P    127 void G4MuBetheBlochModel::Initialise(const G4ParticleDefinition* p,
125                                      const G4D    128                                      const G4DataVector&)
126 {                                                 129 {
127   SetParticle(p);                              << 130   if(p) { SetParticle(p); }
128   if(nullptr == fParticleChange) {             << 131   if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
129     fParticleChange = GetParticleChangeForLoss << 
130     if(UseAngularGeneratorFlag() && nullptr == << 
131       SetAngularDistribution(new G4DeltaAngle( << 
132     }                                          << 
133   }                                            << 
134 }                                                 132 }
135                                                   133 
136 //....oooOO0OOooo........oooOO0OOooo........oo    134 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
137                                                   135 
138 G4double G4MuBetheBlochModel::ComputeCrossSect    136 G4double G4MuBetheBlochModel::ComputeCrossSectionPerElectron(
139                                            con    137                                            const G4ParticleDefinition* p,
140                                                   138                                                  G4double kineticEnergy,
141                                                   139                                                  G4double cutEnergy,
142                                                   140                                                  G4double maxKinEnergy)
143 {                                                 141 {
144   G4double cross = 0.0;                           142   G4double cross = 0.0;
145   G4double tmax = MaxSecondaryEnergy(p, kineti    143   G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
146   G4double maxEnergy = std::min(tmax, maxKinEn << 144   G4double maxEnergy = min(tmax,maxKinEnergy);
147   if(cutEnergy < maxEnergy) {                     145   if(cutEnergy < maxEnergy) {
148                                                   146 
149     G4double totEnergy = kineticEnergy + mass;    147     G4double totEnergy = kineticEnergy + mass;
150     G4double energy2 = totEnergy*totEnergy;    << 148     G4double energy2   = totEnergy*totEnergy;
151     G4double beta2 = kineticEnergy*(kineticEne << 149     G4double beta2     = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
152                                                   150 
153     cross = 1.0/cutEnergy - 1.0/maxEnergy -    << 151     cross = 1.0/cutEnergy - 1.0/maxEnergy - beta2*G4Log(maxEnergy/cutEnergy)/tmax
154       beta2*G4Log(maxEnergy/cutEnergy)/tmax +  << 152           + 0.5*(maxEnergy - cutEnergy)/energy2;
155       0.5*(maxEnergy - cutEnergy)/energy2;     << 
156                                                   153 
157     // radiative corrections of R. Kokoulin       154     // radiative corrections of R. Kokoulin
158     if (maxEnergy > limitKinEnergy && kineticE << 155     if (maxEnergy > limitKinEnergy) {
159                                                   156 
160       G4double logtmax = G4Log(maxEnergy);        157       G4double logtmax = G4Log(maxEnergy);
161       G4double logtmin = G4Log(std::max(cutEne << 158       G4double logtmin = G4Log(max(cutEnergy,limitKinEnergy));
162       G4double logstep = logtmax - logtmin;       159       G4double logstep = logtmax - logtmin;
163       G4double dcross  = 0.0;                     160       G4double dcross  = 0.0;
164                                                   161 
165       for (G4int ll=0; ll<8; ++ll) {           << 162       for (G4int ll=0; ll<8; ll++)
                                                   >> 163       {
166         G4double ep = G4Exp(logtmin + xgi[ll]*    164         G4double ep = G4Exp(logtmin + xgi[ll]*logstep);
167         G4double a1 = G4Log(1.0 + 2.0*ep/CLHEP << 165         G4double a1 = G4Log(1.0 + 2.0*ep/electron_mass_c2);
168         G4double a3 = G4Log(4.0*totEnergy*(tot    166         G4double a3 = G4Log(4.0*totEnergy*(totEnergy - ep)/massSquare);
169         dcross += wgi[ll]*(1.0/ep - beta2/tmax    167         dcross += wgi[ll]*(1.0/ep - beta2/tmax + 0.5*ep/energy2)*a1*(a3 - a1);
170       }                                           168       }
                                                   >> 169 
171       cross += dcross*logstep*alphaprime;         170       cross += dcross*logstep*alphaprime;
172     }                                             171     }
173     cross *= CLHEP::twopi_mc2_rcl2/beta2;      << 172 
                                                   >> 173     cross *= twopi_mc2_rcl2/beta2;
                                                   >> 174 
174   }                                               175   }
                                                   >> 176 
175   //  G4cout << "tmin= " << cutEnergy << " tma    177   //  G4cout << "tmin= " << cutEnergy << " tmax= " << tmax
176   //         << " cross= " << cross << G4endl; << 178   //         << " cross= " << cross << G4endl;
                                                   >> 179   
177   return cross;                                   180   return cross;
178 }                                                 181 }
179                                                   182 
180 //....oooOO0OOooo........oooOO0OOooo........oo    183 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
181                                                   184 
182 G4double G4MuBetheBlochModel::ComputeCrossSect    185 G4double G4MuBetheBlochModel::ComputeCrossSectionPerAtom(
183                                            con    186                                            const G4ParticleDefinition* p,
184                                                   187                                                  G4double kineticEnergy,
185                                                << 188              G4double Z, G4double,
186                                                   189                                                  G4double cutEnergy,
187                                                   190                                                  G4double maxEnergy)
188 {                                                 191 {
189   G4double cross = Z*ComputeCrossSectionPerEle    192   G4double cross = Z*ComputeCrossSectionPerElectron
190                                          (p,ki    193                                          (p,kineticEnergy,cutEnergy,maxEnergy);
191   return cross;                                   194   return cross;
192 }                                                 195 }
193                                                   196 
194 //....oooOO0OOooo........oooOO0OOooo........oo    197 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
195                                                   198 
196 G4double G4MuBetheBlochModel::CrossSectionPerV    199 G4double G4MuBetheBlochModel::CrossSectionPerVolume(
197                                            con << 200              const G4Material* material,
198                                            con    201                                            const G4ParticleDefinition* p,
199                                                   202                                                  G4double kineticEnergy,
200                                                   203                                                  G4double cutEnergy,
201                                                   204                                                  G4double maxEnergy)
202 {                                                 205 {
203   G4double eDensity = material->GetElectronDen    206   G4double eDensity = material->GetElectronDensity();
204   G4double cross = eDensity*ComputeCrossSectio    207   G4double cross = eDensity*ComputeCrossSectionPerElectron
205                                          (p,ki    208                                          (p,kineticEnergy,cutEnergy,maxEnergy);
206   return cross;                                   209   return cross;
207 }                                                 210 }
208                                                   211 
209 //....oooOO0OOooo........oooOO0OOooo........oo    212 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
210                                                   213 
211 G4double G4MuBetheBlochModel::ComputeDEDXPerVo    214 G4double G4MuBetheBlochModel::ComputeDEDXPerVolume(const G4Material* material,
212                                                << 215               const G4ParticleDefinition* p,
213                                                << 216               G4double kineticEnergy,
214                                                << 217               G4double cut)
215 {                                                 218 {
216   G4double tmax  = MaxSecondaryEnergy(p, kinet    219   G4double tmax  = MaxSecondaryEnergy(p, kineticEnergy);
217   G4double tau   = kineticEnergy/mass;            220   G4double tau   = kineticEnergy/mass;
218   G4double cutEnergy = std::min(cut, tmax);    << 221   G4double cutEnergy = min(cut,tmax);
219   G4double gam   = tau + 1.0;                     222   G4double gam   = tau + 1.0;
220   G4double bg2   = tau * (tau+2.0);               223   G4double bg2   = tau * (tau+2.0);
221   G4double beta2 = bg2/(gam*gam);                 224   G4double beta2 = bg2/(gam*gam);
222                                                   225 
223   G4double eexc  = material->GetIonisation()->    226   G4double eexc  = material->GetIonisation()->GetMeanExcitationEnergy();
224   G4double eexc2 = eexc*eexc;                     227   G4double eexc2 = eexc*eexc;
225                                                   228 
226   G4double eDensity = material->GetElectronDen    229   G4double eDensity = material->GetElectronDensity();
227                                                   230 
228   G4double dedx = G4Log(2.0*CLHEP::electron_ma << 231   G4double dedx = G4Log(2.0*electron_mass_c2*bg2*cutEnergy/eexc2)
229                  -(1.0 + cutEnergy/tmax)*beta2    232                  -(1.0 + cutEnergy/tmax)*beta2;
230                                                   233 
231   G4double totEnergy = kineticEnergy + mass;      234   G4double totEnergy = kineticEnergy + mass;
232   G4double del = 0.5*cutEnergy/totEnergy;         235   G4double del = 0.5*cutEnergy/totEnergy;
233   dedx += del*del;                                236   dedx += del*del;
234                                                   237 
235   // density correction                           238   // density correction
236   G4double x = G4Log(bg2)/twoln10;                239   G4double x = G4Log(bg2)/twoln10;
237   dedx -= material->GetIonisation()->DensityCo    240   dedx -= material->GetIonisation()->DensityCorrection(x);
238                                                   241 
239   // shell and high order corrections          << 242   // shell correction
240   dedx -= 2.0*corr->ShellCorrection(p,material    243   dedx -= 2.0*corr->ShellCorrection(p,material,kineticEnergy);
241                                                   244 
                                                   >> 245   // now compute the total ionization loss
                                                   >> 246 
                                                   >> 247   if (dedx < 0.0) dedx = 0.0 ;
                                                   >> 248 
242   // radiative corrections of R. Kokoulin         249   // radiative corrections of R. Kokoulin
243   if (cutEnergy > limitKinEnergy && kineticEne << 250   if (cutEnergy > limitKinEnergy) {
244                                                   251 
245     G4double logtmax = G4Log(cutEnergy);          252     G4double logtmax = G4Log(cutEnergy);
246     G4double logstep = logtmax - logLimitKinEn    253     G4double logstep = logtmax - logLimitKinEnergy;
247     G4double dloss = 0.0;                         254     G4double dloss = 0.0;
248     G4double ftot2= 0.5/(totEnergy*totEnergy);    255     G4double ftot2= 0.5/(totEnergy*totEnergy);
249                                                   256 
250     for (G4int ll=0; ll<8; ++ll) {             << 257     for (G4int ll=0; ll<8; ll++)
                                                   >> 258     {
251       G4double ep = G4Exp(logLimitKinEnergy +     259       G4double ep = G4Exp(logLimitKinEnergy + xgi[ll]*logstep);
252       G4double a1 = G4Log(1.0 + 2.0*ep/CLHEP:: << 260       G4double a1 = G4Log(1.0 + 2.0*ep/electron_mass_c2);
253       G4double a3 = G4Log(4.0*totEnergy*(totEn    261       G4double a3 = G4Log(4.0*totEnergy*(totEnergy - ep)/massSquare);
254       dloss += wgi[ll]*(1.0 - beta2*ep/tmax +     262       dloss += wgi[ll]*(1.0 - beta2*ep/tmax + ep*ep*ftot2)*a1*(a3 - a1);
255     }                                             263     }
256     dedx += dloss*logstep*alphaprime;             264     dedx += dloss*logstep*alphaprime;
257   }                                               265   }
258   dedx *= CLHEP::twopi_mc2_rcl2*eDensity/beta2 << 266 
                                                   >> 267   dedx *= twopi_mc2_rcl2*eDensity/beta2;
259                                                   268 
260   //High order corrections                        269   //High order corrections
261   dedx += corr->HighOrderCorrections(p,materia    270   dedx += corr->HighOrderCorrections(p,material,kineticEnergy,cutEnergy);
262   dedx = std::max(dedx, 0.);                   << 271 
263   return dedx;                                    272   return dedx;
264 }                                                 273 }
265                                                   274 
266 //....oooOO0OOooo........oooOO0OOooo........oo    275 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
267                                                   276 
268 void G4MuBetheBlochModel::SampleSecondaries(   << 277 void G4MuBetheBlochModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp,
269                           std::vector<G4Dynami << 278               const G4MaterialCutsCouple*,
270         const G4MaterialCutsCouple* couple,    << 279               const G4DynamicParticle* dp,
271         const G4DynamicParticle* dp,           << 280               G4double minKinEnergy,
272         G4double minKinEnergy,                 << 281               G4double maxEnergy)
273         G4double maxEnergy)                    << 
274 {                                                 282 {
275   G4double kineticEnergy = dp->GetKineticEnerg << 283   G4double tmax = MaxSecondaryKinEnergy(dp);
276   G4double tmax = MaxSecondaryEnergy(dp->GetDe << 284   G4double maxKinEnergy = min(maxEnergy,tmax);
277   G4double maxKinEnergy = std::min(maxEnergy,  << 
278   if(minKinEnergy >= maxKinEnergy) { return; }    285   if(minKinEnergy >= maxKinEnergy) { return; }
279                                                   286 
280   G4double totEnergy = kineticEnergy + mass;   << 287   G4double kineticEnergy = dp->GetKineticEnergy();
281   G4double etot2 = totEnergy*totEnergy;        << 288   G4double totEnergy     = kineticEnergy + mass;
282   G4double beta2 = kineticEnergy*(kineticEnerg << 289   G4double etot2         = totEnergy*totEnergy;
                                                   >> 290   G4double beta2         = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
283                                                   291  
284   G4double grej  = 1.;                            292   G4double grej  = 1.;
285   G4bool radC = (tmax > limitKinEnergy && kine << 293   if(tmax > limitKinEnergy) {
286   if(radC) {                                   << 294     G4double a0    = G4Log(2.*totEnergy/mass);
287     G4double a0 = G4Log(2.*totEnergy/mass);    << 295     grej  += alphaprime*a0*a0;
288     grej += alphaprime*a0*a0;                  << 
289   }                                               296   }
290                                                   297 
291   G4double tkin, f;                            << 298   G4double deltaKinEnergy, f;
292                                                   299 
293   // sampling follows ...                         300   // sampling follows ...
294   do {                                            301   do {
295     G4double q = G4UniformRand();                 302     G4double q = G4UniformRand();
296     tkin = minKinEnergy*maxKinEnergy/(minKinEn << 303     deltaKinEnergy = minKinEnergy*maxKinEnergy
297     f = 1.0 - beta2*tkin/tmax + 0.5*tkin*tkin/ << 304                     /(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
                                                   >> 305 
                                                   >> 306 
                                                   >> 307     f = 1.0 - beta2*deltaKinEnergy/tmax 
                                                   >> 308             + 0.5*deltaKinEnergy*deltaKinEnergy/etot2;
298                                                   309 
299     if(radC && tkin > limitKinEnergy) {        << 310     if(deltaKinEnergy > limitKinEnergy) {
300       G4double a1 = G4Log(1.0 + 2.0*tkin/CLHEP << 311       G4double a1 = G4Log(1.0 + 2.0*deltaKinEnergy/electron_mass_c2);
301       G4double a3 = G4Log(4.0*totEnergy*(totEn << 312       G4double a3 = G4Log(4.0*totEnergy*(totEnergy - deltaKinEnergy)/massSquare);
302       f *= (1. + alphaprime*a1*(a3 - a1));        313       f *= (1. + alphaprime*a1*(a3 - a1));
303     }                                             314     }
304                                                   315 
305     if(f > grej) {                                316     if(f > grej) {
306         G4cout << "G4MuBetheBlochModel::Sample    317         G4cout << "G4MuBetheBlochModel::SampleSecondary Warning! "
307                << "Majorant " << grej << " < "    318                << "Majorant " << grej << " < "
308                << f << " for edelta= " << tkin << 319                << f << " for edelta= " << deltaKinEnergy
309                << " tmin= " << minKinEnergy <<    320                << " tmin= " << minKinEnergy << " max= " << maxKinEnergy
310                << G4endl;                         321                << G4endl;
311     }                                             322     }
312     // Loop checking, 03-Aug-2015, Vladimir Iv << 
313   } while( grej*G4UniformRand() > f );         << 
314                                                   323 
315   G4ThreeVector deltaDirection;                << 
316                                                   324 
317   if(UseAngularGeneratorFlag()) {              << 325   } while( grej*G4UniformRand() > f );
318     const G4Material* mat = couple->GetMateria << 
319     deltaDirection = GetAngularDistribution()- << 
320                      SelectRandomAtomNumber(ma << 
321   } else {                                     << 
322                                                << 
323     G4double deltaMom = std::sqrt(tkin * (tkin << 
324     G4double totalMom = totEnergy*std::sqrt(be << 
325     G4double cost = tkin * (totEnergy + CLHEP: << 
326       (deltaMom * totalMom);                   << 
327     cost = std::min(cost, 1.0);                << 
328     const G4double sint = std::sqrt((1.0 - cos << 
329     const G4double phi = twopi*G4UniformRand() << 
330                                                   326 
331     deltaDirection.set(sint*std::cos(phi),sint << 327   G4double deltaMomentum =
332     deltaDirection.rotateUz(dp->GetMomentumDir << 328            sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
333   }                                            << 329   G4double totalMomentum = totEnergy*sqrt(beta2);
334   // create G4DynamicParticle object for delta << 330   G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
335   auto delta = new G4DynamicParticle(theElectr << 331                                    (deltaMomentum * totalMomentum);
336   vdp->push_back(delta);                       << 332 
                                                   >> 333   G4double sint = sqrt(1.0 - cost*cost);
                                                   >> 334 
                                                   >> 335   G4double phi = twopi * G4UniformRand() ;
                                                   >> 336 
                                                   >> 337   G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost) ;
                                                   >> 338   G4ThreeVector direction = dp->GetMomentumDirection();
                                                   >> 339   deltaDirection.rotateUz(direction);
337                                                   340 
338   // primary change                               341   // primary change
339   kineticEnergy -= tkin;                       << 342   kineticEnergy -= deltaKinEnergy;
340   G4ThreeVector dir = dp->GetMomentum() - delt << 343   G4ThreeVector dir = totalMomentum*direction - deltaMomentum*deltaDirection;
341   dir = dir.unit();                            << 344   direction = dir.unit();
342   fParticleChange->SetProposedKineticEnergy(ki    345   fParticleChange->SetProposedKineticEnergy(kineticEnergy);
343   fParticleChange->SetProposedMomentumDirectio << 346   fParticleChange->SetProposedMomentumDirection(direction);
                                                   >> 347 
                                                   >> 348   // create G4DynamicParticle object for delta ray
                                                   >> 349   G4DynamicParticle* delta = new G4DynamicParticle(theElectron,
                                                   >> 350                                                  deltaDirection,deltaKinEnergy);
                                                   >> 351   vdp->push_back(delta);
344 }                                                 352 }
345                                                   353 
346 //....oooOO0OOooo........oooOO0OOooo........oo    354 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
347                                                   355