Geant4 Cross Reference |
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 // 26 // 27 // ------------------------------------------- 27 // ------------------------------------------------------------------- 28 // 28 // 29 // GEANT4 Class header file 29 // GEANT4 Class header file 30 // 30 // 31 // 31 // 32 // File name: G4mplIonisationWithDeltaMode 32 // File name: G4mplIonisationWithDeltaModel 33 // 33 // 34 // Author: Vladimir Ivanchenko 34 // Author: Vladimir Ivanchenko 35 // 35 // 36 // Creation date: 06.09.2005 36 // Creation date: 06.09.2005 37 // 37 // 38 // Modifications: 38 // Modifications: 39 // 12.08.2007 Changing low energy approximatio 39 // 12.08.2007 Changing low energy approximation and extrapolation. 40 // Small bug fixing and refactoring 40 // Small bug fixing and refactoring (M. Vladymyrov) 41 // 13.11.2007 Use low-energy asymptotic from [ 41 // 13.11.2007 Use low-energy asymptotic from [3] (V.Ivanchenko) 42 // 42 // 43 // 43 // 44 // ------------------------------------------- 44 // ------------------------------------------------------------------- 45 // References 45 // References 46 // [1] Steven P. Ahlen: Energy loss of relativ 46 // [1] Steven P. Ahlen: Energy loss of relativistic heavy ionizing particles, 47 // S.P. Ahlen, Rev. Mod. Phys 52(1980), p1 47 // S.P. Ahlen, Rev. Mod. Phys 52(1980), p121 48 // [2] K.A. Milton arXiv:hep-ex/0602040 48 // [2] K.A. Milton arXiv:hep-ex/0602040 49 // [3] S.P. Ahlen and K. Kinoshita, Phys. Rev. 49 // [3] S.P. Ahlen and K. Kinoshita, Phys. Rev. D26 (1982) 2347 50 50 51 51 52 //....oooOO0OOooo........oooOO0OOooo........oo 52 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 53 //....oooOO0OOooo........oooOO0OOooo........oo 53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 54 54 55 #include "G4mplIonisationWithDeltaModel.hh" 55 #include "G4mplIonisationWithDeltaModel.hh" 56 #include "Randomize.hh" 56 #include "Randomize.hh" 57 #include "G4PhysicalConstants.hh" 57 #include "G4PhysicalConstants.hh" 58 #include "G4SystemOfUnits.hh" 58 #include "G4SystemOfUnits.hh" 59 #include "G4ParticleChangeForLoss.hh" 59 #include "G4ParticleChangeForLoss.hh" 60 #include "G4Electron.hh" 60 #include "G4Electron.hh" 61 #include "G4DynamicParticle.hh" 61 #include "G4DynamicParticle.hh" 62 #include "G4ProductionCutsTable.hh" 62 #include "G4ProductionCutsTable.hh" 63 #include "G4MaterialCutsCouple.hh" 63 #include "G4MaterialCutsCouple.hh" 64 #include "G4Log.hh" 64 #include "G4Log.hh" 65 #include "G4Pow.hh" 65 #include "G4Pow.hh" 66 66 67 //....oooOO0OOooo........oooOO0OOooo........oo 67 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 68 68 69 using namespace std; 69 using namespace std; 70 70 71 std::vector<G4double>* G4mplIonisationWithDelt 71 std::vector<G4double>* G4mplIonisationWithDeltaModel::dedx0 = nullptr; 72 72 73 G4mplIonisationWithDeltaModel::G4mplIonisation 73 G4mplIonisationWithDeltaModel::G4mplIonisationWithDeltaModel(G4double mCharge, 74 74 const G4String& nam) 75 : G4VEmModel(nam),G4VEmFluctuationModel(nam) 75 : G4VEmModel(nam),G4VEmFluctuationModel(nam), 76 magCharge(mCharge), 76 magCharge(mCharge), 77 twoln10(std::log(100.0)), 77 twoln10(std::log(100.0)), 78 betalow(0.01), 78 betalow(0.01), 79 betalim(0.1), 79 betalim(0.1), 80 beta2lim(betalim*betalim), 80 beta2lim(betalim*betalim), 81 bg2lim(beta2lim*(1.0 + beta2lim)) 81 bg2lim(beta2lim*(1.0 + beta2lim)) 82 { 82 { 83 nmpl = G4lrint(std::abs(magCharge) * 2 * fin 83 nmpl = G4lrint(std::abs(magCharge) * 2 * fine_structure_const); 84 if(nmpl > 6) { nmpl = 6; } 84 if(nmpl > 6) { nmpl = 6; } 85 else if(nmpl < 1) { nmpl = 1; } 85 else if(nmpl < 1) { nmpl = 1; } 86 pi_hbarc2_over_mc2 = pi * hbarc * hbarc / el 86 pi_hbarc2_over_mc2 = pi * hbarc * hbarc / electron_mass_c2; 87 chargeSquare = magCharge * magCharge; 87 chargeSquare = magCharge * magCharge; 88 dedxlim = 45.*nmpl*nmpl*GeV*cm2/g; 88 dedxlim = 45.*nmpl*nmpl*GeV*cm2/g; 89 fParticleChange = nullptr; 89 fParticleChange = nullptr; 90 theElectron = G4Electron::Electron(); 90 theElectron = G4Electron::Electron(); 91 G4cout << "### Monopole ionisation model wit 91 G4cout << "### Monopole ionisation model with d-electron production, Gmag= " 92 << magCharge/eplus << G4endl; 92 << magCharge/eplus << G4endl; 93 monopole = nullptr; 93 monopole = nullptr; 94 mass = 0.0; 94 mass = 0.0; 95 } 95 } 96 96 97 //....oooOO0OOooo........oooOO0OOooo........oo 97 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 98 98 99 G4mplIonisationWithDeltaModel::~G4mplIonisatio 99 G4mplIonisationWithDeltaModel::~G4mplIonisationWithDeltaModel() 100 { 100 { 101 if(IsMaster()) { delete dedx0; } 101 if(IsMaster()) { delete dedx0; } 102 } 102 } 103 103 104 //....oooOO0OOooo........oooOO0OOooo........oo 104 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 105 105 106 void G4mplIonisationWithDeltaModel::SetParticl 106 void G4mplIonisationWithDeltaModel::SetParticle(const G4ParticleDefinition* p) 107 { 107 { 108 monopole = p; 108 monopole = p; 109 mass = monopole->GetPDGMass(); 109 mass = monopole->GetPDGMass(); 110 G4double emin = 110 G4double emin = 111 std::min(LowEnergyLimit(),0.1*mass*(1./sqr 111 std::min(LowEnergyLimit(),0.1*mass*(1./sqrt(1. - betalow*betalow) - 1.)); 112 G4double emax = 112 G4double emax = 113 std::max(HighEnergyLimit(),10*mass*(1./sqr 113 std::max(HighEnergyLimit(),10*mass*(1./sqrt(1. - beta2lim) - 1.)); 114 SetLowEnergyLimit(emin); 114 SetLowEnergyLimit(emin); 115 SetHighEnergyLimit(emax); 115 SetHighEnergyLimit(emax); 116 } 116 } 117 117 118 //....oooOO0OOooo........oooOO0OOooo........oo 118 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 119 119 120 void 120 void 121 G4mplIonisationWithDeltaModel::Initialise(cons 121 G4mplIonisationWithDeltaModel::Initialise(const G4ParticleDefinition* p, 122 cons 122 const G4DataVector&) 123 { 123 { 124 if(!monopole) { SetParticle(p); } 124 if(!monopole) { SetParticle(p); } 125 if(!fParticleChange) { fParticleChange = Get 125 if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); } 126 if(IsMaster()) { 126 if(IsMaster()) { 127 if(!dedx0) { dedx0 = new std::vector<G4dou 127 if(!dedx0) { dedx0 = new std::vector<G4double>; } 128 G4ProductionCutsTable* theCoupleTable= 128 G4ProductionCutsTable* theCoupleTable= 129 G4ProductionCutsTable::GetProductionCuts 129 G4ProductionCutsTable::GetProductionCutsTable(); 130 G4int numOfCouples = (G4int)theCoupleTable 130 G4int numOfCouples = (G4int)theCoupleTable->GetTableSize(); 131 G4int n = (G4int)dedx0->size(); 131 G4int n = (G4int)dedx0->size(); 132 if(n < numOfCouples) { dedx0->resize(numOf 132 if(n < numOfCouples) { dedx0->resize(numOfCouples); } 133 G4Pow* g4calc = G4Pow::GetInstance(); 133 G4Pow* g4calc = G4Pow::GetInstance(); 134 134 135 // initialise vector assuming low conducti 135 // initialise vector assuming low conductivity 136 for(G4int i=0; i<numOfCouples; ++i) { 136 for(G4int i=0; i<numOfCouples; ++i) { 137 137 138 const G4Material* material = 138 const G4Material* material = 139 theCoupleTable->GetMaterialCutsCouple( 139 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial(); 140 G4double eDensity = material->GetElectro 140 G4double eDensity = material->GetElectronDensity(); 141 G4double vF2 = 2*electron_Compton_length 141 G4double vF2 = 2*electron_Compton_length*g4calc->A13(3.*pi*pi*eDensity); 142 (*dedx0)[i] = pi_hbarc2_over_mc2*eDensit 142 (*dedx0)[i] = pi_hbarc2_over_mc2*eDensity*nmpl*nmpl* 143 (G4Log(vF2/fine_structure_const) - 0.5 143 (G4Log(vF2/fine_structure_const) - 0.5)/vF2; 144 } 144 } 145 } 145 } 146 } 146 } 147 147 148 //....oooOO0OOooo........oooOO0OOooo........oo 148 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 149 149 150 G4double 150 G4double 151 G4mplIonisationWithDeltaModel::MinEnergyCut(co 151 G4mplIonisationWithDeltaModel::MinEnergyCut(const G4ParticleDefinition*, 152 co 152 const G4MaterialCutsCouple* couple) 153 { 153 { 154 return couple->GetMaterial()->GetIonisation( 154 return couple->GetMaterial()->GetIonisation()->GetMeanExcitationEnergy(); 155 } 155 } 156 156 157 //....oooOO0OOooo........oooOO0OOooo........oo 157 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 158 158 159 G4double 159 G4double 160 G4mplIonisationWithDeltaModel::ComputeDEDXPerV 160 G4mplIonisationWithDeltaModel::ComputeDEDXPerVolume(const G4Material* material, 161 161 const G4ParticleDefinition* p, 162 162 G4double kineticEnergy, 163 163 G4double maxEnergy) 164 { 164 { 165 if(!monopole) { SetParticle(p); } 165 if(!monopole) { SetParticle(p); } 166 G4double tmax = MaxSecondaryEnergy(p,kinetic 166 G4double tmax = MaxSecondaryEnergy(p,kineticEnergy); 167 G4double cutEnergy = std::min(tmax, maxEnerg 167 G4double cutEnergy = std::min(tmax, maxEnergy); 168 cutEnergy = std::max(LowEnergyLimit(), cutEn 168 cutEnergy = std::max(LowEnergyLimit(), cutEnergy); 169 G4double tau = kineticEnergy / mass; 169 G4double tau = kineticEnergy / mass; 170 G4double gam = tau + 1.0; 170 G4double gam = tau + 1.0; 171 G4double bg2 = tau * (tau + 2.0); 171 G4double bg2 = tau * (tau + 2.0); 172 G4double beta2 = bg2 / (gam * gam); 172 G4double beta2 = bg2 / (gam * gam); 173 G4double beta = sqrt(beta2); 173 G4double beta = sqrt(beta2); 174 174 175 // low-energy asymptotic formula 175 // low-energy asymptotic formula 176 G4double dedx = (*dedx0)[CurrentCouple()->Ge 176 G4double dedx = (*dedx0)[CurrentCouple()->GetIndex()]*beta; 177 177 178 // above asymptotic 178 // above asymptotic 179 if(beta > betalow) { 179 if(beta > betalow) { 180 180 181 // high energy 181 // high energy 182 if(beta >= betalim) { 182 if(beta >= betalim) { 183 dedx = ComputeDEDXAhlen(material, bg2, c 183 dedx = ComputeDEDXAhlen(material, bg2, cutEnergy); 184 184 185 } else { 185 } else { 186 G4double dedx1 = (*dedx0)[CurrentCouple( 186 G4double dedx1 = (*dedx0)[CurrentCouple()->GetIndex()]*betalow; 187 G4double dedx2 = ComputeDEDXAhlen(materi 187 G4double dedx2 = ComputeDEDXAhlen(material, bg2lim, cutEnergy); 188 188 189 // extrapolation between two formula 189 // extrapolation between two formula 190 G4double kapa2 = beta - betalow; 190 G4double kapa2 = beta - betalow; 191 G4double kapa1 = betalim - beta; 191 G4double kapa1 = betalim - beta; 192 dedx = (kapa1*dedx1 + kapa2*dedx2)/(kapa 192 dedx = (kapa1*dedx1 + kapa2*dedx2)/(kapa1 + kapa2); 193 } 193 } 194 } 194 } 195 return dedx; 195 return dedx; 196 } 196 } 197 197 198 //....oooOO0OOooo........oooOO0OOooo........oo 198 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 199 199 200 G4double 200 G4double 201 G4mplIonisationWithDeltaModel::ComputeDEDXAhle 201 G4mplIonisationWithDeltaModel::ComputeDEDXAhlen(const G4Material* material, 202 202 G4double bg2, 203 203 G4double cutEnergy) 204 { 204 { 205 G4double eDensity = material->GetElectronDen 205 G4double eDensity = material->GetElectronDensity(); 206 G4double eexc = material->GetIonisation()-> 206 G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy(); 207 207 208 // Ahlen's formula for nonconductors, [1]p15 208 // Ahlen's formula for nonconductors, [1]p157, f(5.7) 209 G4double dedx = 209 G4double dedx = 210 0.5*(G4Log(2.0*electron_mass_c2*bg2*cutEne 210 0.5*(G4Log(2.0*electron_mass_c2*bg2*cutEnergy/(eexc*eexc)) -1.0); 211 211 212 // Kazama et al. cross-section correction 212 // Kazama et al. cross-section correction 213 G4double k = 0.406; 213 G4double k = 0.406; 214 if(nmpl > 1) { k = 0.346; } 214 if(nmpl > 1) { k = 0.346; } 215 215 216 // Bloch correction 216 // Bloch correction 217 const G4double B[7] = { 0.0, 0.248, 0.672, 1 217 const G4double B[7] = { 0.0, 0.248, 0.672, 1.022, 1.243, 1.464, 1.685}; 218 218 219 dedx += 0.5 * k - B[nmpl]; 219 dedx += 0.5 * k - B[nmpl]; 220 220 221 // density effect correction 221 // density effect correction 222 G4double x = G4Log(bg2)/twoln10; 222 G4double x = G4Log(bg2)/twoln10; 223 dedx -= material->GetIonisation()->DensityCo 223 dedx -= material->GetIonisation()->DensityCorrection(x); 224 224 225 // now compute the total ionization loss 225 // now compute the total ionization loss 226 dedx *= pi_hbarc2_over_mc2 * eDensity * nmp 226 dedx *= pi_hbarc2_over_mc2 * eDensity * nmpl * nmpl; 227 227 228 dedx = std::max(dedx, 0.0); 228 dedx = std::max(dedx, 0.0); 229 return dedx; 229 return dedx; 230 } 230 } 231 231 232 //....oooOO0OOooo........oooOO0OOooo........oo 232 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 233 233 234 G4double 234 G4double 235 G4mplIonisationWithDeltaModel::ComputeCrossSec 235 G4mplIonisationWithDeltaModel::ComputeCrossSectionPerElectron( 236 con 236 const G4ParticleDefinition* p, 237 G4d 237 G4double kineticEnergy, 238 G4d 238 G4double cut, 239 G4d 239 G4double maxKinEnergy) 240 { 240 { 241 if(!monopole) { SetParticle(p); } 241 if(!monopole) { SetParticle(p); } 242 G4double tmax = MaxSecondaryEnergy(p, kineti 242 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy); 243 G4double maxEnergy = std::min(tmax, maxKinEn 243 G4double maxEnergy = std::min(tmax, maxKinEnergy); 244 G4double cutEnergy = std::max(LowEnergyLimit 244 G4double cutEnergy = std::max(LowEnergyLimit(), cut); 245 G4double cross = (cutEnergy < maxEnergy) 245 G4double cross = (cutEnergy < maxEnergy) 246 ? (0.5/cutEnergy - 0.5/maxEnergy)*pi_hbarc 246 ? (0.5/cutEnergy - 0.5/maxEnergy)*pi_hbarc2_over_mc2 * nmpl * nmpl : 0.0; 247 return cross; 247 return cross; 248 } 248 } 249 249 250 //....oooOO0OOooo........oooOO0OOooo........oo 250 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 251 251 252 G4double 252 G4double 253 G4mplIonisationWithDeltaModel::ComputeCrossSec 253 G4mplIonisationWithDeltaModel::ComputeCrossSectionPerAtom( 254 cons 254 const G4ParticleDefinition* p, 255 G4do 255 G4double kineticEnergy, 256 G4do 256 G4double Z, G4double, 257 G4do 257 G4double cutEnergy, 258 G4do 258 G4double maxEnergy) 259 { 259 { 260 G4double cross = 260 G4double cross = 261 Z*ComputeCrossSectionPerElectron(p,kinetic 261 Z*ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy); 262 return cross; 262 return cross; 263 } 263 } 264 264 265 //....oooOO0OOooo........oooOO0OOooo........oo 265 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 266 266 267 void 267 void 268 G4mplIonisationWithDeltaModel::SampleSecondari 268 G4mplIonisationWithDeltaModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp, 269 269 const G4MaterialCutsCouple*, 270 270 const G4DynamicParticle* dp, 271 271 G4double minKinEnergy, 272 272 G4double maxEnergy) 273 { 273 { 274 G4double kineticEnergy = dp->GetKineticEnerg 274 G4double kineticEnergy = dp->GetKineticEnergy(); 275 G4double tmax = MaxSecondaryEnergy(dp->GetDe 275 G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(),kineticEnergy); 276 276 277 G4double maxKinEnergy = std::min(maxEnergy,t 277 G4double maxKinEnergy = std::min(maxEnergy,tmax); 278 if(minKinEnergy >= maxKinEnergy) { return; } 278 if(minKinEnergy >= maxKinEnergy) { return; } 279 279 280 //G4cout << "G4mplIonisationWithDeltaModel:: 280 //G4cout << "G4mplIonisationWithDeltaModel::SampleSecondaries: E(GeV)= " 281 // << kineticEnergy/GeV << " M(GeV)= 281 // << kineticEnergy/GeV << " M(GeV)= " << mass/GeV 282 // << " tmin(MeV)= " << minKinEnergy 282 // << " tmin(MeV)= " << minKinEnergy/MeV << G4endl; 283 283 284 G4double totEnergy = kineticEnergy + mas 284 G4double totEnergy = kineticEnergy + mass; 285 G4double etot2 = totEnergy*totEnergy 285 G4double etot2 = totEnergy*totEnergy; 286 G4double beta2 = kineticEnergy*(kine 286 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2; 287 287 288 // sampling without nuclear size effect 288 // sampling without nuclear size effect 289 G4double q = G4UniformRand(); 289 G4double q = G4UniformRand(); 290 G4double deltaKinEnergy = minKinEnergy*maxKi 290 G4double deltaKinEnergy = minKinEnergy*maxKinEnergy 291 /(minKinEnergy*(1.0 - q) + maxKinEnergy*q) 291 /(minKinEnergy*(1.0 - q) + maxKinEnergy*q); 292 292 293 // delta-electron is produced 293 // delta-electron is produced 294 G4double totMomentum = totEnergy*sqrt(beta2) 294 G4double totMomentum = totEnergy*sqrt(beta2); 295 G4double deltaMomentum = 295 G4double deltaMomentum = 296 sqrt(deltaKinEnergy * (deltaKinEner 296 sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2)); 297 G4double cost = deltaKinEnergy * (totEnergy 297 G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) / 298 (deltaMomen 298 (deltaMomentum * totMomentum); 299 cost = std::min(cost, 1.0); 299 cost = std::min(cost, 1.0); 300 300 301 G4double sint = sqrt((1.0 - cost)*(1.0 + cos 301 G4double sint = sqrt((1.0 - cost)*(1.0 + cost)); 302 302 303 G4double phi = twopi * G4UniformRand() ; 303 G4double phi = twopi * G4UniformRand() ; 304 304 305 G4ThreeVector deltaDirection(sint*cos(phi),s 305 G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost); 306 G4ThreeVector direction = dp->GetMomentumDir 306 G4ThreeVector direction = dp->GetMomentumDirection(); 307 deltaDirection.rotateUz(direction); 307 deltaDirection.rotateUz(direction); 308 308 309 // create G4DynamicParticle object for delta 309 // create G4DynamicParticle object for delta ray 310 G4DynamicParticle* delta = 310 G4DynamicParticle* delta = 311 new G4DynamicParticle(theElectron,deltaDir 311 new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy); 312 312 313 vdp->push_back(delta); 313 vdp->push_back(delta); 314 314 315 // Change kinematics of primary particle 315 // Change kinematics of primary particle 316 kineticEnergy -= deltaKinEnergy; 316 kineticEnergy -= deltaKinEnergy; 317 G4ThreeVector finalP = direction*totMomentum 317 G4ThreeVector finalP = direction*totMomentum - deltaDirection*deltaMomentum; 318 finalP = finalP.unit(); 318 finalP = finalP.unit(); 319 319 320 fParticleChange->SetProposedKineticEnergy(ki 320 fParticleChange->SetProposedKineticEnergy(kineticEnergy); 321 fParticleChange->SetProposedMomentumDirectio 321 fParticleChange->SetProposedMomentumDirection(finalP); 322 } 322 } 323 323 324 //....oooOO0OOooo........oooOO0OOooo........oo 324 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 325 325 326 G4double G4mplIonisationWithDeltaModel::Sample 326 G4double G4mplIonisationWithDeltaModel::SampleFluctuations( 327 const G 327 const G4MaterialCutsCouple* couple, 328 const G 328 const G4DynamicParticle* dp, 329 const G 329 const G4double tcut, 330 const G 330 const G4double tmax, 331 const G 331 const G4double length, 332 const G 332 const G4double meanLoss) 333 { 333 { 334 G4double siga = Dispersion(couple->GetMateri 334 G4double siga = Dispersion(couple->GetMaterial(),dp,tcut,tmax,length); 335 G4double loss = meanLoss; 335 G4double loss = meanLoss; 336 siga = std::sqrt(siga); 336 siga = std::sqrt(siga); 337 G4double twomeanLoss = meanLoss + meanLoss; 337 G4double twomeanLoss = meanLoss + meanLoss; 338 338 339 if(twomeanLoss < siga) { 339 if(twomeanLoss < siga) { 340 G4double x; 340 G4double x; 341 do { 341 do { 342 loss = twomeanLoss*G4UniformRand(); 342 loss = twomeanLoss*G4UniformRand(); 343 x = (loss - meanLoss)/siga; 343 x = (loss - meanLoss)/siga; 344 // Loop checking, 07-Aug-2015, Vladimir 344 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko 345 } while (1.0 - 0.5*x*x < G4UniformRand()); 345 } while (1.0 - 0.5*x*x < G4UniformRand()); 346 } else { 346 } else { 347 do { 347 do { 348 loss = G4RandGauss::shoot(meanLoss,siga) 348 loss = G4RandGauss::shoot(meanLoss,siga); 349 // Loop checking, 07-Aug-2015, Vladimir 349 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko 350 } while (0.0 > loss || loss > twomeanLoss) 350 } while (0.0 > loss || loss > twomeanLoss); 351 } 351 } 352 return loss; 352 return loss; 353 } 353 } 354 354 355 //....oooOO0OOooo........oooOO0OOooo........oo 355 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 356 356 357 G4double 357 G4double 358 G4mplIonisationWithDeltaModel::Dispersion(cons 358 G4mplIonisationWithDeltaModel::Dispersion(const G4Material* material, 359 cons 359 const G4DynamicParticle* dp, 360 const G4double tcut, 360 const G4double tcut, 361 const G4double tmax, 361 const G4double tmax, 362 const G4double length) 362 const G4double length) 363 { 363 { 364 G4double siga = 0.0; 364 G4double siga = 0.0; 365 G4double tau = dp->GetKineticEnergy()/mass 365 G4double tau = dp->GetKineticEnergy()/mass; 366 if(tau > 0.0) { 366 if(tau > 0.0) { 367 const G4double beta = dp->GetBeta(); 367 const G4double beta = dp->GetBeta(); 368 siga = (tmax/(beta*beta) - 0.5*tcut) * tw 368 siga = (tmax/(beta*beta) - 0.5*tcut) * twopi_mc2_rcl2 * length 369 * material->GetElectronDensity() * charg 369 * material->GetElectronDensity() * chargeSquare; 370 } 370 } 371 return siga; 371 return siga; 372 } 372 } 373 373 374 //....oooOO0OOooo........oooOO0OOooo........oo 374 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 375 375 376 G4double 376 G4double 377 G4mplIonisationWithDeltaModel::MaxSecondaryEne 377 G4mplIonisationWithDeltaModel::MaxSecondaryEnergy(const G4ParticleDefinition*, 378 378 G4double kinEnergy) 379 { 379 { 380 G4double tau = kinEnergy/mass; 380 G4double tau = kinEnergy/mass; 381 return 2.0*electron_mass_c2*tau*(tau + 2.); 381 return 2.0*electron_mass_c2*tau*(tau + 2.); 382 } 382 } 383 383 384 //....oooOO0OOooo........oooOO0OOooo........oo 384 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 385 385