source: trunk/source/processes/electromagnetic/standard/src/G4BetheBlochModel.cc@ 1306

Last change on this file since 1306 was 1228, checked in by garnier, 16 years ago

update geant4.9.3 tag

File size: 14.8 KB
RevLine 
[819]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//
[1228]26// $Id: G4BetheBlochModel.cc,v 1.36 2009/12/03 17:26:40 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03 $
[819]28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class header file
32//
33//
34// File name: G4BetheBlochModel
35//
36// Author: Vladimir Ivanchenko on base of Laszlo Urban code
37//
38// Creation date: 03.01.2002
39//
40// Modifications:
41//
42// 04-12-02 Fix problem of G4DynamicParticle constructor (V.Ivanchenko)
43// 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
44// 27-01-03 Make models region aware (V.Ivanchenko)
45// 13-02-03 Add name (V.Ivanchenko)
46// 24-03-05 Add G4EmCorrections (V.Ivanchenko)
[961]47// 11-04-05 Major optimisation of internal interfaces (V.Ivanchenko)
[819]48// 11-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
49// 12-02-06 move G4LossTableManager::Instance()->EmCorrections()
50// in constructor (mma)
[961]51// 12-08-08 Added methods GetParticleCharge, GetChargeSquareRatio,
52// CorrectionsAlongStep needed for ions(V.Ivanchenko)
[819]53//
54// -------------------------------------------------------------------
55//
56
57
58//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
60
61#include "G4BetheBlochModel.hh"
62#include "Randomize.hh"
63#include "G4Electron.hh"
64#include "G4LossTableManager.hh"
65#include "G4EmCorrections.hh"
66#include "G4ParticleChangeForLoss.hh"
67
68//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
69
70using namespace std;
71
72G4BetheBlochModel::G4BetheBlochModel(const G4ParticleDefinition* p,
73 const G4String& nam)
74 : G4VEmModel(nam),
[961]75 particle(0),
76 tlimit(DBL_MAX),
77 twoln10(2.0*log(10.0)),
78 bg2lim(0.0169),
79 taulim(8.4146e-3),
80 isIon(false),
81 isInitialised(false)
[819]82{
[961]83 fParticleChange = 0;
[1055]84 if(p) {
85 SetGenericIon(p);
86 SetParticle(p);
87 }
[819]88 theElectron = G4Electron::Electron();
89 corr = G4LossTableManager::Instance()->EmCorrections();
[961]90 nist = G4NistManager::Instance();
91 SetLowEnergyLimit(2.0*MeV);
[819]92}
93
94//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
95
96G4BetheBlochModel::~G4BetheBlochModel()
97{}
98
99//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
100
101G4double G4BetheBlochModel::MinEnergyCut(const G4ParticleDefinition*,
102 const G4MaterialCutsCouple* couple)
103{
104 return couple->GetMaterial()->GetIonisation()->GetMeanExcitationEnergy();
105}
106
107//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
108
109void G4BetheBlochModel::Initialise(const G4ParticleDefinition* p,
110 const G4DataVector&)
111{
[1055]112 SetGenericIon(p);
113 SetParticle(p);
[819]114
[1055]115 //G4cout << "G4BetheBlochModel::Initialise for " << p->GetParticleName()
116 // << " isIon= " << isIon
117 // << G4endl;
118
[961]119 corrFactor = chargeSquare;
[1055]120 // always false before the run
121 SetDeexcitationFlag(false);
[961]122
123 if(!isInitialised) {
124 isInitialised = true;
[1055]125 fParticleChange = GetParticleChangeForLoss();
[961]126 }
[819]127}
128
129//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
130
[961]131G4double G4BetheBlochModel::GetChargeSquareRatio(const G4ParticleDefinition* p,
132 const G4Material* mat,
133 G4double kineticEnergy)
[819]134{
[961]135 // this method is called only for ions
136 G4double q2 = corr->EffectiveChargeSquareRatio(p,mat,kineticEnergy);
137 corrFactor = q2*corr->EffectiveChargeCorrection(p,mat,kineticEnergy);
138 return corrFactor;
139}
140
141//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
142
143G4double G4BetheBlochModel::GetParticleCharge(const G4ParticleDefinition* p,
144 const G4Material* mat,
145 G4double kineticEnergy)
146{
147 // this method is called only for ions
148 return corr->GetParticleCharge(p,mat,kineticEnergy);
149}
150
151//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
152
153G4double
154G4BetheBlochModel::ComputeCrossSectionPerElectron(const G4ParticleDefinition* p,
155 G4double kineticEnergy,
156 G4double cutEnergy,
157 G4double maxKinEnergy)
158{
[819]159 G4double cross = 0.0;
160 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
161 G4double maxEnergy = min(tmax,maxKinEnergy);
162 if(cutEnergy < maxEnergy) {
163
164 G4double totEnergy = kineticEnergy + mass;
165 G4double energy2 = totEnergy*totEnergy;
166 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
167
[961]168 cross = 1.0/cutEnergy - 1.0/maxEnergy
169 - beta2*log(maxEnergy/cutEnergy)/tmax;
[819]170
171 // +term for spin=1/2 particle
172 if( 0.5 == spin ) cross += 0.5*(maxEnergy - cutEnergy)/energy2;
173
[961]174 // High order correction different for hadrons and ions
175 // nevetheless they are applied to reduce high energy transfers
176 // if(!isIon)
177 //cross += corr->FiniteSizeCorrectionXS(p,currentMaterial,
178 // kineticEnergy,cutEnergy);
179
[819]180 cross *= twopi_mc2_rcl2*chargeSquare/beta2;
181 }
182
183 // G4cout << "BB: e= " << kineticEnergy << " tmin= " << cutEnergy
184 // << " tmax= " << tmax << " cross= " << cross << G4endl;
185
186 return cross;
187}
188
189//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
190
191G4double G4BetheBlochModel::ComputeCrossSectionPerAtom(
192 const G4ParticleDefinition* p,
193 G4double kineticEnergy,
194 G4double Z, G4double,
195 G4double cutEnergy,
196 G4double maxEnergy)
197{
198 G4double cross = Z*ComputeCrossSectionPerElectron
199 (p,kineticEnergy,cutEnergy,maxEnergy);
200 return cross;
201}
202
203//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
204
205G4double G4BetheBlochModel::CrossSectionPerVolume(
206 const G4Material* material,
207 const G4ParticleDefinition* p,
208 G4double kineticEnergy,
209 G4double cutEnergy,
210 G4double maxEnergy)
211{
[961]212 currentMaterial = material;
[819]213 G4double eDensity = material->GetElectronDensity();
214 G4double cross = eDensity*ComputeCrossSectionPerElectron
215 (p,kineticEnergy,cutEnergy,maxEnergy);
216 return cross;
217}
218
219//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
220
221G4double G4BetheBlochModel::ComputeDEDXPerVolume(const G4Material* material,
222 const G4ParticleDefinition* p,
223 G4double kineticEnergy,
224 G4double cut)
225{
226 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
227 G4double cutEnergy = min(cut,tmax);
228
229 G4double tau = kineticEnergy/mass;
230 G4double gam = tau + 1.0;
231 G4double bg2 = tau * (tau+2.0);
232 G4double beta2 = bg2/(gam*gam);
233
234 G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
235 G4double eexc2 = eexc*eexc;
[1196]236 //G4double cden = material->GetIonisation()->GetCdensity();
237 //G4double mden = material->GetIonisation()->GetMdensity();
238 //G4double aden = material->GetIonisation()->GetAdensity();
239 //G4double x0den = material->GetIonisation()->GetX0density();
240 //G4double x1den = material->GetIonisation()->GetX1density();
[819]241
242 G4double eDensity = material->GetElectronDensity();
243
244 G4double dedx = log(2.0*electron_mass_c2*bg2*cutEnergy/eexc2)
245 - (1.0 + cutEnergy/tmax)*beta2;
246
247 if(0.5 == spin) {
248 G4double del = 0.5*cutEnergy/(kineticEnergy + mass);
249 dedx += del*del;
250 }
251
252 // density correction
253 G4double x = log(bg2)/twoln10;
[1196]254 //if ( x >= x0den ) {
255 // dedx -= twoln10*x - cden ;
256 // if ( x < x1den ) dedx -= aden*pow((x1den-x),mden) ;
257 //}
258 dedx -= material->GetIonisation()->DensityCorrection(x);
[819]259
260 // shell correction
261 dedx -= 2.0*corr->ShellCorrection(p,material,kineticEnergy);
262
263 // now compute the total ionization loss
264
265 if (dedx < 0.0) dedx = 0.0 ;
266
267 dedx *= twopi_mc2_rcl2*chargeSquare*eDensity/beta2;
268
[961]269 //High order correction different for hadrons and ions
270 if(isIon) {
271 dedx += corr->IonBarkasCorrection(p,material,kineticEnergy);
272 } else {
273 dedx += corr->HighOrderCorrections(p,material,kineticEnergy,cutEnergy);
274 }
[819]275 return dedx;
276}
277
278//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
[1228]279/*
280void G4BetheBlochModel::CorrectionsAlongStep(const G4MaterialCutsCouple*,
281 const G4DynamicParticle*,
282 G4double&,
283 G4double&,
284 G4double)
285{}
286*/
[819]287
[961]288void G4BetheBlochModel::CorrectionsAlongStep(const G4MaterialCutsCouple* couple,
289 const G4DynamicParticle* dp,
290 G4double& eloss,
291 G4double&,
292 G4double length)
293{
[1196]294 if(isIon) {
295 const G4ParticleDefinition* p = dp->GetDefinition();
296 const G4Material* mat = couple->GetMaterial();
297 G4double preKinEnergy = dp->GetKineticEnergy();
298 G4double e = preKinEnergy - eloss*0.5;
299 if(e < 0.0) e = preKinEnergy*0.5;
[961]300
301 G4double q2 = corr->EffectiveChargeSquareRatio(p,mat,e);
302 GetModelOfFluctuations()->SetParticleAndCharge(p, q2);
[1228]303 G4double qfactor = q2*corr->EffectiveChargeCorrection(p,mat,e)/corrFactor;
304 G4double highOrder = length*corr->IonHighOrderCorrections(p,couple,e);
305 eloss *= qfactor;
306 eloss += highOrder;
307 //G4cout << "G4BetheBlochModel::CorrectionsAlongStep: e= " << preKinEnergy
308 // << " qfactor= " << qfactor
309 // << " highOrder= " << highOrder << " (" << highOrder/eloss << ")" << G4endl;
[961]310 }
311}
312
313//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
314
[819]315void G4BetheBlochModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp,
316 const G4MaterialCutsCouple*,
317 const G4DynamicParticle* dp,
318 G4double minKinEnergy,
319 G4double maxEnergy)
320{
321 G4double kineticEnergy = dp->GetKineticEnergy();
322 G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(),kineticEnergy);
323
[961]324 G4double maxKinEnergy = std::min(maxEnergy,tmax);
[819]325 if(minKinEnergy >= maxKinEnergy) return;
326
327 G4double totEnergy = kineticEnergy + mass;
328 G4double etot2 = totEnergy*totEnergy;
329 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
330
[961]331 G4double deltaKinEnergy, f;
332 G4double f1 = 0.0;
333 G4double fmax = 1.0;
334 if( 0.5 == spin ) fmax += 0.5*maxKinEnergy*maxKinEnergy/etot2;
[819]335
[961]336 // sampling without nuclear size effect
[819]337 do {
338 G4double q = G4UniformRand();
339 deltaKinEnergy = minKinEnergy*maxKinEnergy
340 /(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
341
342 f = 1.0 - beta2*deltaKinEnergy/tmax;
[961]343 if( 0.5 == spin ) {
344 f1 = 0.5*deltaKinEnergy*deltaKinEnergy/etot2;
345 f += f1;
346 }
[819]347
[961]348 } while( fmax*G4UniformRand() > f);
349
350 // projectile formfactor - suppresion of high energy
351 // delta-electron production at high energy
352
353 G4double x = formfact*deltaKinEnergy;
354 if(x > 1.e-6) {
355
356 G4double x1 = 1.0 + x;
357 G4double g = 1.0/(x1*x1);
358 if( 0.5 == spin ) {
359 G4double x2 = 0.5*electron_mass_c2*deltaKinEnergy/(mass*mass);
360 g *= (1.0 + magMoment2*(x2 - f1/f)/(1.0 + x2));
[819]361 }
[961]362 if(g > 1.0) {
363 G4cout << "### G4BetheBlochModel WARNING: g= " << g
364 << dp->GetDefinition()->GetParticleName()
365 << " Ekin(MeV)= " << kineticEnergy
366 << " delEkin(MeV)= " << deltaKinEnergy
367 << G4endl;
368 }
369 if(G4UniformRand() > g) return;
370 }
[819]371
[961]372 // delta-electron is produced
[819]373 G4double totMomentum = totEnergy*sqrt(beta2);
374 G4double deltaMomentum =
375 sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
376 G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
377 (deltaMomentum * totMomentum);
[1055]378 /*
[819]379 if(cost > 1.0) {
380 G4cout << "### G4BetheBlochModel WARNING: cost= "
381 << cost << " > 1 for "
382 << dp->GetDefinition()->GetParticleName()
383 << " Ekin(MeV)= " << kineticEnergy
384 << " p(MeV/c)= " << totMomentum
385 << " delEkin(MeV)= " << deltaKinEnergy
386 << " delMom(MeV/c)= " << deltaMomentum
387 << " tmin(MeV)= " << minKinEnergy
388 << " tmax(MeV)= " << maxKinEnergy
389 << " dir= " << dp->GetMomentumDirection()
390 << G4endl;
391 cost = 1.0;
392 }
[1055]393 */
[819]394 G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
395
396 G4double phi = twopi * G4UniformRand() ;
397
398
399 G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost);
400 G4ThreeVector direction = dp->GetMomentumDirection();
401 deltaDirection.rotateUz(direction);
402
403 // create G4DynamicParticle object for delta ray
404 G4DynamicParticle* delta = new G4DynamicParticle(theElectron,
[1055]405 deltaDirection,deltaKinEnergy);
[819]406
407 vdp->push_back(delta);
408
409 // Change kinematics of primary particle
410 kineticEnergy -= deltaKinEnergy;
411 G4ThreeVector finalP = direction*totMomentum - deltaDirection*deltaMomentum;
412 finalP = finalP.unit();
413
414 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
415 fParticleChange->SetProposedMomentumDirection(finalP);
416}
417
[1055]418//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
419
420G4double G4BetheBlochModel::MaxSecondaryEnergy(const G4ParticleDefinition* pd,
421 G4double kinEnergy)
422{
423 // here particle type is checked for any method
424 SetParticle(pd);
425 G4double tau = kinEnergy/mass;
426 G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
427 (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
428 return std::min(tmax,tlimit);
429}
430
[819]431//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracBrowser for help on using the repository browser.