source: trunk/source/processes/electromagnetic/standard/src/G4eCoulombScatteringModel.cc @ 1228

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

update geant4.9.3 tag

File size: 12.7 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: G4eCoulombScatteringModel.cc,v 1.78 2009/10/28 10:14:13 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name:     G4eCoulombScatteringModel
35//
36// Author:        Vladimir Ivanchenko
37//
38// Creation date: 22.08.2005
39//
40// Modifications:
41//
42// 01.08.06 V.Ivanchenko extend upper limit of table to TeV and review the
43//          logic of building - only elements from G4ElementTable
44// 08.08.06 V.Ivanchenko build internal table in ekin scale, introduce faclim
45// 19.08.06 V.Ivanchenko add inline function ScreeningParameter
46// 09.10.07 V.Ivanchenko reorganized methods, add cut dependence in scattering off e-
47// 09.06.08 V.Ivanchenko add SelectIsotope and sampling of the recoil ion
48// 16.06.09 C.Consolandi fixed computation of effective mass
49//
50//
51// Class Description:
52//
53// -------------------------------------------------------------------
54//
55//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
56//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
57
58#include "G4eCoulombScatteringModel.hh"
59#include "Randomize.hh"
60#include "G4DataVector.hh"
61#include "G4ElementTable.hh"
62#include "G4PhysicsLogVector.hh"
63#include "G4ParticleChangeForGamma.hh"
64#include "G4Electron.hh"
65#include "G4Positron.hh"
66#include "G4Proton.hh"
67#include "G4ParticleTable.hh"
68#include "G4ProductionCutsTable.hh"
69#include "G4NucleiProperties.hh"
70
71//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
72
73G4double G4eCoulombScatteringModel::ScreenRSquare[] = {0.0};
74G4double G4eCoulombScatteringModel::FormFactor[]    = {0.0};
75
76using namespace std;
77
78G4eCoulombScatteringModel::G4eCoulombScatteringModel(const G4String& nam)
79  : G4VEmModel(nam),
80    cosThetaMin(1.0),
81    cosThetaMax(-1.0),
82    q2Limit(TeV*TeV),
83    alpha2(fine_structure_const*fine_structure_const),
84    faclim(100.0),
85    isInitialised(false)
86{
87  fNistManager = G4NistManager::Instance();
88  theParticleTable = G4ParticleTable::GetParticleTable();
89  theElectron = G4Electron::Electron();
90  thePositron = G4Positron::Positron();
91  theProton   = G4Proton::Proton();
92  currentMaterial = 0; 
93  currentElement  = 0;
94  lowEnergyLimit  = 0.1*keV;
95  G4double p0 = electron_mass_c2*classic_electr_radius;
96  coeff  = twopi*p0*p0;
97  tkin = targetZ = mom2 = DBL_MIN;
98  elecXSection = nucXSection = 0.0;
99  recoilThreshold = 0.*keV;
100  ecut = DBL_MAX;
101  particle = 0;
102  currentCouple = 0;
103
104  // Thomas-Fermi screening radii
105  // Formfactors from A.V. Butkevich et al., NIM A 488 (2002) 282
106
107  if(0.0 == ScreenRSquare[0]) {
108    G4double a0 = electron_mass_c2/0.88534; 
109    G4double constn = 6.937e-6/(MeV*MeV);
110
111    ScreenRSquare[0] = alpha2*a0*a0;
112    for(G4int j=1; j<100; j++) {
113      G4double x = a0*fNistManager->GetZ13(j);
114      ScreenRSquare[j] = alpha2*x*x;
115      x = fNistManager->GetA27(j); 
116      FormFactor[j] = constn*x*x;
117    } 
118  }
119}
120
121//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
122
123G4eCoulombScatteringModel::~G4eCoulombScatteringModel()
124{}
125
126//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
127
128void G4eCoulombScatteringModel::Initialise(const G4ParticleDefinition* p,
129                                           const G4DataVector& cuts)
130{
131  SetupParticle(p);
132  currentCouple = 0;
133  elecXSection = nucXSection = 0.0;
134  tkin = targetZ = mom2 = DBL_MIN;
135  ecut = etag = DBL_MAX;
136  cosThetaMin = cos(PolarAngleLimit());
137  pCuts = G4ProductionCutsTable::GetProductionCutsTable()->GetEnergyCutsVector(3);
138  //G4cout << "!!! G4eCoulombScatteringModel::Initialise for "
139  //     << p->GetParticleName() << "  cos(TetMin)= " << cosThetaMin
140  //     << "  cos(TetMax)= " << cosThetaMax <<G4endl;
141  // G4cout << "cut0= " << cuts[0] << "  cut1= " << cuts[1] << G4endl;
142  if(!isInitialised) {
143    isInitialised = true;
144    fParticleChange = GetParticleChangeForGamma();
145  }
146  if(mass < GeV && particle->GetParticleType() != "nucleus") {
147    InitialiseElementSelectors(p,cuts);
148  }
149}
150
151//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
152
153void G4eCoulombScatteringModel::ComputeMaxElectronScattering(G4double cutEnergy)
154{
155  ecut = cutEnergy;
156  G4double tmax = tkin;
157  cosTetMaxElec = 1.0;
158  if(mass > MeV) {
159    G4double ratio = electron_mass_c2/mass;
160    G4double tau = tkin/mass;
161    tmax = 2.0*electron_mass_c2*tau*(tau + 2.)/
162      (1.0 + 2.0*ratio*(tau + 1.0) + ratio*ratio); 
163    cosTetMaxElec = 1.0 - std::min(cutEnergy, tmax)*electron_mass_c2/mom2;
164  } else {
165
166    if(particle == theElectron) tmax *= 0.5;
167    G4double t = std::min(cutEnergy, tmax);
168    G4double mom21 = t*(t + 2.0*electron_mass_c2);
169    G4double t1 = tkin - t;
170    //G4cout << "tkin= " << tkin << " t= " << t << " t1= " << t1 << G4endl;
171    if(t1 > 0.0) {
172      G4double mom22 = t1*(t1 + 2.0*mass);
173      G4double ctm = (mom2 + mom22 - mom21)*0.5/sqrt(mom2*mom22);
174      //G4cout << "ctm= " << ctm << G4endl;
175      if(ctm <  1.0) cosTetMaxElec = ctm;
176      if(ctm < -1.0) cosTetMaxElec = -1.0;
177    }
178  }
179}
180
181//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
182
183G4double G4eCoulombScatteringModel::ComputeCrossSectionPerAtom(
184                const G4ParticleDefinition* p,
185                G4double kinEnergy,
186                G4double Z, G4double,
187                G4double cutEnergy, G4double)
188{
189  //G4cout << "### G4eCoulombScatteringModel::ComputeCrossSectionPerAtom  for "
190  //  << p->GetParticleName()<<" Z= "<<Z<<" e(MeV)= "<< kinEnergy/MeV << G4endl;
191  G4double xsec = 0.0;
192  SetupParticle(p);
193  if(kinEnergy < lowEnergyLimit) return xsec;
194  SetupKinematic(kinEnergy, cutEnergy);
195  if(cosTetMaxNuc < cosTetMinNuc) {
196    SetupTarget(Z, kinEnergy);
197    xsec = CrossSectionPerAtom(); 
198  }
199  /*
200  G4cout << "e(MeV)= " << ekin/MeV << "cosTetMinNuc= " << cosTetMinNuc
201         << " cosTetMaxNuc= " << cosTetMaxNuc
202         << " cosTetMaxElec= " << cosTetMaxElec
203         << " screenZ= " << screenZ
204         << " formfactA= " << formfactA << G4endl;
205  */
206  return xsec; 
207}
208
209//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
210
211G4double G4eCoulombScatteringModel::CrossSectionPerAtom()
212{
213  // This method needs initialisation before be called
214  //G4double fac = coeff*targetZ*chargeSquare*invbeta2/mom2;
215
216  G4double meff = targetMass/(mass+targetMass);
217  G4double fac  = coeff*targetZ*chargeSquare*invbeta2/(mom2*meff*meff);
218
219  elecXSection = 0.0;
220  nucXSection  = 0.0;
221
222  G4double x  = 1.0 - cosTetMinNuc;
223  G4double x1 = x + screenZ;
224
225  if(cosTetMaxElec2 < cosTetMinNuc) {
226    elecXSection = fac*(cosTetMinNuc - cosTetMaxElec2)/
227      (x1*(1.0 - cosTetMaxElec2 + screenZ));
228    nucXSection  = elecXSection;
229  }
230
231  //G4cout << "XS tkin(MeV)= " << tkin<<" xs= " <<nucXSection
232  //     << " costmax= " << cosTetMaxNuc2
233  //     << " costmin= " << cosTetMinNuc << "  Z= " << targetZ <<G4endl;
234  if(cosTetMaxNuc2 < cosTetMinNuc) {
235    G4double s  = screenZ*formfactA;
236    G4double z1 = 1.0 - cosTetMaxNuc2 + screenZ;
237    G4double s1 = 1.0 - s;
238    G4double d  = s1/formfactA;
239    //G4cout <<"x1= "<<x1<<" z1= " <<z1<<" s= "<<s << " d= " <<d <<G4endl;
240    if(d < 0.2*x1) {
241      G4double x2 = x1*x1;
242      G4double z2 = z1*z1;
243      x = (1.0/(x1*x2) - 1.0/(z1*z2) - d*1.5*(1.0/(x2*x2) - 1.0/(z2*z2)))/
244        (3.0*formfactA*formfactA);
245    } else {
246      G4double x2 = x1 + d;
247      G4double z2 = z1 + d;
248      x = (1.0/x1 - 1.0/z1 + 1.0/x2 - 1.0/z2 - 2.0*log(z1*x2/(z2*x1))/d)/(s1*s1);
249    }
250    nucXSection += fac*targetZ*x;
251  }
252  //G4cout<<" cross(bn)= "<<nucXSection/barn<<" xsElec(bn)= "<<elecXSection/barn
253  //    << " Asc= " << screenZ << G4endl;
254 
255  return nucXSection;
256}
257
258//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
259
260void G4eCoulombScatteringModel::SampleSecondaries(
261                std::vector<G4DynamicParticle*>* fvect,
262                const G4MaterialCutsCouple* couple,
263                const G4DynamicParticle* dp,
264                G4double cutEnergy,
265                G4double)
266{
267  G4double kinEnergy = dp->GetKineticEnergy();
268  if(kinEnergy < lowEnergyLimit) return;
269  DefineMaterial(couple);
270  SetupParticle(dp->GetDefinition());
271
272  SetupKinematic(kinEnergy, cutEnergy);
273  //G4cout << "G4eCoulombScatteringModel::SampleSecondaries e(MeV)= "
274  //     << kinEnergy << "  " << particle->GetParticleName()
275  //     << " cut= " << cutEnergy<< G4endl;
276 
277  // Choose nucleus
278  currentElement = SelectRandomAtom(couple,particle,
279                                    kinEnergy,cutEnergy,kinEnergy);
280
281  SetupTarget(currentElement->GetZ(),kinEnergy);
282
283  G4int ia = SelectIsotopeNumber(currentElement);
284  targetMass = G4NucleiProperties::GetNuclearMass(ia, iz);
285 
286  G4double cost = SampleCosineTheta();
287  G4double z1   = 1.0 - cost;
288  if(z1 < 0.0) return;
289
290  G4double sint = sqrt(z1*(1.0 + cost));
291 
292  //G4cout<<"## Sampled sint= " << sint << "  Z= " << targetZ << " A= " << ia
293  //    << "  screenZ= " << screenZ << " cn= " << formfactA << G4endl;
294 
295  G4double phi  = twopi * G4UniformRand();
296
297  G4ThreeVector direction = dp->GetMomentumDirection(); 
298  G4ThreeVector newDirection(cos(phi)*sint,sin(phi)*sint,cost);
299  newDirection.rotateUz(direction);   
300
301  fParticleChange->ProposeMomentumDirection(newDirection);   
302
303  // recoil sampling assuming a small recoil
304  // and first order correction to primary 4-momentum
305  G4double q2   = 2*z1*mom2;
306  G4double trec = q2/(sqrt(targetMass*targetMass + q2) + targetMass);
307  G4double finalT = kinEnergy - trec; 
308  //G4cout<<"G4eCoulombScatteringModel: finalT= "<<finalT<<" Trec= "<<trec<<G4endl;
309  if(finalT <= lowEnergyLimit) { 
310    trec = kinEnergy; 
311    finalT = 0.0;
312  } 
313
314  fParticleChange->SetProposedKineticEnergy(finalT);
315  G4double tcut = recoilThreshold;
316  if(pCuts) { tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]); }
317
318  if(trec > tcut) {
319    G4ParticleDefinition* ion = theParticleTable->FindIon(iz, ia, 0, iz);
320    G4ThreeVector dir = (direction*sqrt(mom2) - 
321                         newDirection*sqrt(finalT*(2*mass + finalT))).unit();
322    G4DynamicParticle* newdp = new G4DynamicParticle(ion, dir, trec);
323    fvect->push_back(newdp);
324  } else {
325    fParticleChange->ProposeLocalEnergyDeposit(trec);
326    fParticleChange->ProposeNonIonizingEnergyDeposit(trec);
327  }
328 
329  return;
330}
331
332//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
333
334G4double G4eCoulombScatteringModel::SampleCosineTheta()
335{
336  G4double costm = cosTetMaxNuc2;
337  G4double formf = formfactA;
338  G4double prob  = 0.0; 
339  G4double xs = CrossSectionPerAtom();
340  if(xs > 0.0) prob = elecXSection/xs;
341
342  // scattering off e or A?
343  if(G4UniformRand() < prob) {
344    costm = cosTetMaxElec2;
345    formf = 0.0;
346  }
347
348  /* 
349  G4cout << "SampleCost: e(MeV)= " << tkin
350         << " 1-ctmaxN= " << 1. - cosTetMinNuc
351         << " 1-ctmax= " << 1. - costm
352         << " Z= " << targetZ
353         << G4endl;
354  */
355
356  if(costm >= cosTetMinNuc) return 2.0; 
357
358  G4double x1 = 1. - cosTetMinNuc + screenZ;
359  G4double x2 = 1. - costm + screenZ;
360  G4double x3 = cosTetMinNuc - costm;
361  G4double grej, z1; 
362  do {
363    z1 = x1*x2/(x1 + G4UniformRand()*x3) - screenZ;
364    grej = 1.0/(1.0 + formf*z1);
365  } while ( G4UniformRand() > grej*grej ); 
366
367  if(mass > MeV) {
368    if(G4UniformRand() > (1. - z1*0.5)/(1.0 + z1*sqrt(mom2)/targetMass)) {
369      return 2.0;
370    }
371  }
372  //G4cout << "z1= " << z1 << " cross= " << nucXSection/barn
373  //     << " crossE= " << elecXSection/barn << G4endl;
374
375  return 1.0 - z1;
376}
377
378//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
379
380
Note: See TracBrowser for help on using the repository browser.