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

Last change on this file since 846 was 819, checked in by garnier, 16 years ago

import all except CVS

File size: 11.3 KB
Line 
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// $Id: G4BetheBlochModel.cc,v 1.13 2007/05/22 17:34:36 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-01-patch-02 $
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)
47// 11-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
48// 11-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
49// 12-02-06 move G4LossTableManager::Instance()->EmCorrections()
50//          in constructor (mma)
51//
52// -------------------------------------------------------------------
53//
54
55
56//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
58
59#include "G4BetheBlochModel.hh"
60#include "Randomize.hh"
61#include "G4Electron.hh"
62#include "G4LossTableManager.hh"
63#include "G4EmCorrections.hh"
64#include "G4ParticleChangeForLoss.hh"
65
66//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
67
68using namespace std;
69
70G4BetheBlochModel::G4BetheBlochModel(const G4ParticleDefinition* p, 
71                                     const G4String& nam)
72  : G4VEmModel(nam),
73  particle(0),
74  tlimit(DBL_MAX),
75  twoln10(2.0*log(10.0)),
76  bg2lim(0.0169),
77  taulim(8.4146e-3),
78  isIon(false)
79{
80  if(p) SetParticle(p);
81  theElectron = G4Electron::Electron();
82  corr = G4LossTableManager::Instance()->EmCorrections(); 
83}
84
85//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
86
87G4BetheBlochModel::~G4BetheBlochModel()
88{}
89
90//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
91
92G4double G4BetheBlochModel::MinEnergyCut(const G4ParticleDefinition*,
93                                         const G4MaterialCutsCouple* couple)
94{
95  return couple->GetMaterial()->GetIonisation()->GetMeanExcitationEnergy();
96}
97
98//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
99
100void G4BetheBlochModel::Initialise(const G4ParticleDefinition* p,
101                                   const G4DataVector&)
102{
103  if (!particle) SetParticle(p);
104  G4String pname = particle->GetParticleName();
105  if (particle->GetParticleType() == "nucleus" &&
106     pname != "deuteron" && pname != "triton") isIon = true;
107
108  if (pParticleChange) 
109    fParticleChange = reinterpret_cast<G4ParticleChangeForLoss*>
110                                                              (pParticleChange);
111  else 
112    fParticleChange = new G4ParticleChangeForLoss();
113}
114
115//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
116
117G4double G4BetheBlochModel::ComputeCrossSectionPerElectron(
118                                           const G4ParticleDefinition* p,
119                                                 G4double kineticEnergy,
120                                                 G4double cutEnergy,
121                                                 G4double maxKinEnergy)                                           
122{
123  G4double cross = 0.0;
124  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
125  G4double maxEnergy = min(tmax,maxKinEnergy);
126  if(cutEnergy < maxEnergy) {
127
128    G4double totEnergy = kineticEnergy + mass;
129    G4double energy2   = totEnergy*totEnergy;
130    G4double beta2     = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
131
132    cross = 1.0/cutEnergy - 1.0/maxEnergy - beta2*log(maxEnergy/cutEnergy)/tmax;
133
134    // +term for spin=1/2 particle
135    if( 0.5 == spin ) cross += 0.5*(maxEnergy - cutEnergy)/energy2;
136
137    cross *= twopi_mc2_rcl2*chargeSquare/beta2;
138  }
139 
140   // G4cout << "BB: e= " << kineticEnergy << " tmin= " << cutEnergy
141   //        << " tmax= " << tmax << " cross= " << cross << G4endl;
142 
143  return cross;
144}
145
146//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
147
148G4double G4BetheBlochModel::ComputeCrossSectionPerAtom(
149                                           const G4ParticleDefinition* p,
150                                                 G4double kineticEnergy,
151                                                 G4double Z, G4double,
152                                                 G4double cutEnergy,
153                                                 G4double maxEnergy)
154{
155  G4double cross = Z*ComputeCrossSectionPerElectron
156                                         (p,kineticEnergy,cutEnergy,maxEnergy);
157  return cross;
158}
159
160//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
161
162G4double G4BetheBlochModel::CrossSectionPerVolume(
163                                           const G4Material* material,
164                                           const G4ParticleDefinition* p,
165                                                 G4double kineticEnergy,
166                                                 G4double cutEnergy,
167                                                 G4double maxEnergy)
168{
169  G4double eDensity = material->GetElectronDensity();
170  G4double cross = eDensity*ComputeCrossSectionPerElectron
171                                         (p,kineticEnergy,cutEnergy,maxEnergy);
172  return cross;
173}
174
175//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
176
177G4double G4BetheBlochModel::ComputeDEDXPerVolume(const G4Material* material,
178                                                 const G4ParticleDefinition* p,
179                                                 G4double kineticEnergy,
180                                                 G4double cut)
181{
182  G4double tmax      = MaxSecondaryEnergy(p, kineticEnergy);
183  G4double cutEnergy = min(cut,tmax);
184
185  G4double tau   = kineticEnergy/mass;
186  G4double gam   = tau + 1.0;
187  G4double bg2   = tau * (tau+2.0);
188  G4double beta2 = bg2/(gam*gam);
189
190  G4double eexc  = material->GetIonisation()->GetMeanExcitationEnergy();
191  G4double eexc2 = eexc*eexc;
192  G4double cden  = material->GetIonisation()->GetCdensity();
193  G4double mden  = material->GetIonisation()->GetMdensity();
194  G4double aden  = material->GetIonisation()->GetAdensity();
195  G4double x0den = material->GetIonisation()->GetX0density();
196  G4double x1den = material->GetIonisation()->GetX1density();
197
198  G4double eDensity = material->GetElectronDensity();
199
200  G4double dedx = log(2.0*electron_mass_c2*bg2*cutEnergy/eexc2)
201                - (1.0 + cutEnergy/tmax)*beta2;
202
203  if(0.5 == spin) {
204    G4double del = 0.5*cutEnergy/(kineticEnergy + mass);
205    dedx += del*del;
206  }
207
208  // density correction
209  G4double x = log(bg2)/twoln10;
210  if ( x >= x0den ) {
211    dedx -= twoln10*x - cden ;
212    if ( x < x1den ) dedx -= aden*pow((x1den-x),mden) ;
213  }
214
215  // shell correction
216  dedx -= 2.0*corr->ShellCorrection(p,material,kineticEnergy);
217
218  // now compute the total ionization loss
219
220  if (dedx < 0.0) dedx = 0.0 ;
221
222  dedx *= twopi_mc2_rcl2*chargeSquare*eDensity/beta2;
223
224  //High order correction only for hadrons
225  if(!isIon) dedx += corr->HighOrderCorrections(p,material,kineticEnergy);
226
227  return dedx;
228}
229
230//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
231
232void G4BetheBlochModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp,
233                                          const G4MaterialCutsCouple*,
234                                          const G4DynamicParticle* dp,
235                                          G4double minKinEnergy,
236                                          G4double maxEnergy)
237{
238  G4double kineticEnergy = dp->GetKineticEnergy();
239  G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(),kineticEnergy);
240
241  G4double maxKinEnergy = min(maxEnergy,tmax);
242  if(minKinEnergy >= maxKinEnergy) return;
243
244  G4double totEnergy     = kineticEnergy + mass;
245  G4double etot2         = totEnergy*totEnergy;
246  G4double beta2         = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
247
248  G4double deltaKinEnergy, f;
249
250  // sampling follows ...
251  do {
252    G4double q = G4UniformRand();
253    deltaKinEnergy = minKinEnergy*maxKinEnergy
254                    /(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
255
256    f = 1.0 - beta2*deltaKinEnergy/tmax;
257    if( 0.5 == spin ) f += 0.5*deltaKinEnergy*deltaKinEnergy/etot2;
258
259    if(f > 1.0) {
260        G4cout << "G4BetheBlochModel::SampleSecondary Warning! "
261               << "Majorant 1.0 < "
262               << f << " for Edelta= " << deltaKinEnergy
263               << G4endl;
264    }
265
266  } while( G4UniformRand() > f );
267
268  G4double totMomentum = totEnergy*sqrt(beta2);
269  G4double deltaMomentum =
270           sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
271  G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
272                                   (deltaMomentum * totMomentum);
273  if(cost > 1.0) {
274    G4cout << "### G4BetheBlochModel WARNING: cost= " 
275           << cost << " > 1 for "
276           << dp->GetDefinition()->GetParticleName()
277           << " Ekin(MeV)= " <<  kineticEnergy
278           << " p(MeV/c)= " <<  totMomentum
279           << " delEkin(MeV)= " << deltaKinEnergy
280           << " delMom(MeV/c)= " << deltaMomentum
281           << " tmin(MeV)= " << minKinEnergy
282           << " tmax(MeV)= " << maxKinEnergy
283           << " dir= " << dp->GetMomentumDirection()
284           << G4endl;
285    cost = 1.0;
286  }
287  G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
288
289  G4double phi = twopi * G4UniformRand() ;
290
291
292  G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost);
293  G4ThreeVector direction = dp->GetMomentumDirection();
294  deltaDirection.rotateUz(direction);
295
296  // create G4DynamicParticle object for delta ray
297  G4DynamicParticle* delta = new G4DynamicParticle(theElectron,
298                                                 deltaDirection,deltaKinEnergy);
299
300  vdp->push_back(delta);
301
302  // Change kinematics of primary particle
303  kineticEnergy       -= deltaKinEnergy;
304  G4ThreeVector finalP = direction*totMomentum - deltaDirection*deltaMomentum;
305  finalP               = finalP.unit();
306 
307  fParticleChange->SetProposedKineticEnergy(kineticEnergy);
308  fParticleChange->SetProposedMomentumDirection(finalP);
309}
310
311//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracBrowser for help on using the repository browser.