source: trunk/source/processes/electromagnetic/standard/include/G4eCoulombScatteringModel.hh @ 1201

Last change on this file since 1201 was 1196, checked in by garnier, 15 years ago

update CVS release candidate geant4.9.3.01

File size: 8.5 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.hh,v 1.49 2009/10/10 15:16:57 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03-cand-01 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class header file
32//
33//
34// File name:     G4eCoulombScatteringModel
35//
36// Author:        Vladimir Ivanchenko
37//
38// Creation date: 19.02.2006
39//
40// Modifications:
41// 01.08.06 V.Ivanchenko extend upper limit of table to TeV and review the
42//          logic of building - only elements from G4ElementTable
43// 08.08.06 V.Ivanchenko build internal table in ekin scale, introduce faclim
44// 19.08.06 V.Ivanchenko add inline function ScreeningParameter and
45//                       make some members protected
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// 17.06.09 C.Consoalndi modified SetupTarget method - remove kinFactor
49//                                     
50//
51// Class Description:
52//
53// Implementation of eCoulombScattering of pointlike charge particle
54// on Atomic Nucleus for interval of scattering anles in Lab system
55// thetaMin - ThetaMax, nucleus recoil is neglected.
56//   The model based on analysis of J.M.Fernandez-Varea et al.
57// NIM B73(1993)447 originated from G.Wentzel Z.Phys. 40(1927)590 with
58// screening parameter from H.A.Bethe Phys. Rev. 89 (1953) 1256.
59//
60
61// -------------------------------------------------------------------
62//
63
64#ifndef G4eCoulombScatteringModel_h
65#define G4eCoulombScatteringModel_h 1
66
67#include "G4VEmModel.hh"
68#include "G4PhysicsTable.hh"
69#include "globals.hh"
70#include "G4NistManager.hh"
71#include <vector>
72
73class G4ParticleChangeForGamma;
74class G4ParticleDefinition;
75
76class G4eCoulombScatteringModel : public G4VEmModel
77{
78
79public:
80
81  G4eCoulombScatteringModel(const G4String& nam = "eCoulombScattering");
82 
83  virtual ~G4eCoulombScatteringModel();
84
85  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
86
87  virtual G4double ComputeCrossSectionPerAtom(
88                                const G4ParticleDefinition*,
89                                G4double kinEnergy, 
90                                G4double Z, 
91                                G4double A, 
92                                G4double cut,
93                                G4double emax);
94
95  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
96                                 const G4MaterialCutsCouple*,
97                                 const G4DynamicParticle*,
98                                 G4double tmin,
99                                 G4double maxEnergy);
100
101  inline void SetRecoilThreshold(G4double eth);
102
103protected:
104
105  G4double CrossSectionPerAtom();
106
107  G4double SampleCosineTheta();
108
109  inline void DefineMaterial(const G4MaterialCutsCouple*);
110
111  inline void SetupParticle(const G4ParticleDefinition*);
112
113  inline void SetupKinematic(G4double kinEnergy, G4double cut);
114 
115  inline void SetupTarget(G4double Z, G4double kinEnergy); 
116
117private:
118
119  void ComputeMaxElectronScattering(G4double cut);
120
121  // hide assignment operator
122  G4eCoulombScatteringModel & operator=(const G4eCoulombScatteringModel &right);
123  G4eCoulombScatteringModel(const  G4eCoulombScatteringModel&);
124
125protected:
126 
127  const G4ParticleDefinition* theProton;
128  const G4ParticleDefinition* theElectron;
129  const G4ParticleDefinition* thePositron;
130
131  G4ParticleTable*          theParticleTable; 
132  G4ParticleChangeForGamma* fParticleChange;
133  G4NistManager*            fNistManager;
134
135  const std::vector<G4double>* pCuts;
136
137  const G4MaterialCutsCouple* currentCouple;
138  const G4Material*           currentMaterial;
139  const G4Element*            currentElement;
140  G4int                       currentMaterialIndex;
141
142  G4double                  coeff;
143  G4double                  cosThetaMin;
144  G4double                  cosThetaMax;
145  G4double                  cosTetMinNuc;
146  G4double                  cosTetMaxNuc;
147  G4double                  cosTetMaxNuc2;
148  G4double                  cosTetMaxElec;
149  G4double                  cosTetMaxElec2;
150  G4double                  q2Limit;
151  G4double                  recoilThreshold;
152  G4double                  elecXSection;
153  G4double                  nucXSection;
154  G4double                  ecut;
155
156  // projectile
157  const G4ParticleDefinition* particle;
158
159  G4double                  chargeSquare;
160  G4double                  spin;
161  G4double                  mass;
162  G4double                  tkin;
163  G4double                  mom2;
164  G4double                  invbeta2;
165  G4double                  etag;
166  G4double                  lowEnergyLimit;
167
168  // target
169  G4double                  targetZ;
170  G4double                  targetMass;
171  G4double                  screenZ;
172  G4double                  formfactA;
173  G4int                     idxelm;
174  G4int                     iz;
175
176private:
177
178  G4double                  alpha2;
179  G4double                  faclim;
180
181  static G4double ScreenRSquare[100];
182  static G4double FormFactor[100];
183
184  G4bool                    isInitialised;             
185};
186
187//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
188
189inline
190void G4eCoulombScatteringModel::DefineMaterial(const G4MaterialCutsCouple* cup) 
191{ 
192  if(cup != currentCouple) {
193    currentCouple = cup;
194    currentMaterial = cup->GetMaterial();
195    currentMaterialIndex = currentCouple->GetIndex(); 
196  }
197}
198
199//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
200
201inline 
202void G4eCoulombScatteringModel::SetupParticle(const G4ParticleDefinition* p)
203{
204  // Initialise mass and charge
205  if(p != particle) {
206    particle = p;
207    mass = particle->GetPDGMass();
208    spin = particle->GetPDGSpin();
209    G4double q = particle->GetPDGCharge()/eplus;
210    chargeSquare = q*q;
211    tkin = 0.0;
212  }
213}
214
215//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
216
217inline void G4eCoulombScatteringModel::SetupKinematic(G4double ekin, 
218                                                      G4double cut)
219{
220  if(ekin != tkin || ecut != cut) {
221    tkin = ekin;
222    mom2 = tkin*(tkin + 2.0*mass);
223    invbeta2 = 1.0 +  mass*mass/mom2;
224    cosTetMinNuc = cosThetaMin;
225    cosTetMaxNuc = cosThetaMax;
226    if(mass < MeV && cosThetaMin < 1.0 && ekin <= 10.*cut) {
227      cosTetMinNuc = ekin*(cosThetaMin + 1.0)/(10.*cut) - 1.0;
228    }
229    ComputeMaxElectronScattering(cut);
230  }
231}
232
233//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
234 
235inline void G4eCoulombScatteringModel::SetupTarget(G4double Z, G4double e)
236{
237  if(Z != targetZ || e != etag) {
238    etag    = e; 
239    targetZ = Z;
240    iz= G4int(Z);
241    if(iz > 99) iz = 99;
242    targetMass = fNistManager->GetAtomicMassAmu(iz)*amu_c2;
243    screenZ = ScreenRSquare[iz]/mom2;
244    screenZ *=(1.13 + std::min(1.0,3.76*Z*Z*invbeta2*alpha2));
245    if(mass > MeV) { screenZ *= 2.0; }
246    formfactA = FormFactor[iz]*mom2;
247    cosTetMaxNuc2 = cosTetMaxNuc;
248    if(1 == iz && particle == theProton && cosTetMaxNuc2 < 0.0) {
249      cosTetMaxNuc2 = 0.0;
250    }
251    cosTetMaxElec2 = cosTetMaxElec;
252  } 
253} 
254
255//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
256
257inline void G4eCoulombScatteringModel::SetRecoilThreshold(G4double eth)
258{
259  recoilThreshold = eth;
260}
261
262//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
263
264#endif
Note: See TracBrowser for help on using the repository browser.