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

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

import all except CVS

File size: 7.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.20 2007/10/24 10:42:05 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-01-patch-02 $
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//
48// Class Description:
49//
50// Implementation of eCoulombScattering of pointlike charge particle
51// on Atomic Nucleus for interval of scattering anles in Lab system
52// thetaMin - ThetaMax, nucleus recoil is neglected.
53//   The model based on analysis of J.M.Fernandez-Varea et al.
54// NIM B73(1993)447 originated from G.Wentzel Z.Phys. 40(1927)590 with
55// screening parameter from H.A.Bethe Phys. Rev. 89 (1953) 1256.
56//
57
58// -------------------------------------------------------------------
59//
60
61#ifndef G4eCoulombScatteringModel_h
62#define G4eCoulombScatteringModel_h 1
63
64#include "G4VEmModel.hh"
65#include "G4PhysicsTable.hh"
66#include "G4NistManager.hh"
67#include "globals.hh"
68
69class G4ParticleChangeForGamma;
70class G4ParticleDefinition;
71
72class G4eCoulombScatteringModel : public G4VEmModel
73{
74
75public:
76
77  G4eCoulombScatteringModel(G4double thetaMin = 0.0, G4double thetaMax = pi,
78                           G4bool build = false, G4double tlim = TeV*TeV,
79                           const G4String& nam = "eCoulombScattering");
80 
81  virtual ~G4eCoulombScatteringModel();
82
83  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
84
85  virtual G4double ComputeCrossSectionPerAtom(
86                                const G4ParticleDefinition*,
87                                G4double kinEnergy, 
88                                G4double Z, 
89                                G4double A, 
90                                G4double cut,
91                                G4double emax);
92
93  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
94                                 const G4MaterialCutsCouple*,
95                                 const G4DynamicParticle*,
96                                 G4double tmin,
97                                 G4double maxEnergy);
98
99protected:
100
101  G4double ComputeElectronXSectionPerAtom(
102                                 const G4ParticleDefinition*,
103                                 G4double kinEnergy, 
104                                 G4double Z,
105                                 G4double A,
106                                 G4double cut);
107
108  virtual G4double CalculateCrossSectionPerAtom(
109                                 const G4ParticleDefinition*, 
110                                 G4double kinEnergy, 
111                                 G4double Z, G4double A);
112
113  inline void SetupParticle(const G4ParticleDefinition*);
114
115  inline void SetupKinematic(G4double kinEnergy);
116 
117  inline void SetupTarget(G4double Z, G4double A, G4double kinEnergy); 
118
119private:
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  G4ParticleChangeForGamma* fParticleChange;
132  G4NistManager*            fNistManager;
133
134  G4double                  coeff;
135  G4double                  constn;
136  G4double                  cosThetaMin;
137  G4double                  cosThetaMax;
138  G4double                  cosTetMaxNuc;
139  G4double                  cosTetMaxElec;
140  G4double                  q2Limit;
141  G4double                  elecXSection;
142  G4double                  nucXSection;
143  G4double                  ecut;
144
145  // projectile
146  const G4ParticleDefinition* particle;
147
148  G4double                  chargeSquare;
149  G4double                  spin;
150  G4double                  mass;
151  G4double                  tkin;
152  G4double                  mom2;
153  G4double                  invbeta2;
154
155  // target
156  G4double                  targetZ;
157  G4double                  targetA;
158  G4double                  screenZ;
159  G4double                  formfactA;
160
161private:
162
163  G4PhysicsTable*           theCrossSectionTable; 
164
165  G4double                  a0;
166  G4double                  lowKEnergy;
167  G4double                  highKEnergy;
168  G4double                  alpha2;
169  G4double                  faclim;
170
171  G4int                     nbins;
172  G4int                     nmax;
173  G4int                     index[100];
174
175  G4bool                    buildTable;             
176  G4bool                    isInitialised;             
177};
178
179//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
180
181inline 
182void G4eCoulombScatteringModel::SetupParticle(const G4ParticleDefinition* p)
183{
184  // Initialise mass and charge
185  if(p != particle) {
186    particle = p;
187    mass = particle->GetPDGMass();
188    spin = particle->GetPDGSpin();
189    G4double q = particle->GetPDGCharge()/eplus;
190    chargeSquare = q*q;
191    tkin = 0.0;
192  }
193}
194
195//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
196
197inline void G4eCoulombScatteringModel::SetupKinematic(G4double ekin)
198{
199  if(ekin != tkin) {
200    tkin  = ekin;
201    mom2  = tkin*(tkin + 2.0*mass);
202    invbeta2 = 1.0 +  mass*mass/mom2;
203  } 
204}
205
206//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
207 
208inline void G4eCoulombScatteringModel::SetupTarget(G4double Z, G4double A, 
209                                                   G4double e)
210{
211  if(e != tkin || Z != targetZ || A != targetA) {
212    targetZ = Z;
213    targetA = A;
214    SetupKinematic(e);
215    cosTetMaxNuc = std::max(cosThetaMax, 1.0 - 0.5*q2Limit/mom2);
216    G4double x = fNistManager->GetZ13(Z);
217    screenZ = a0*x*x*(1.13 + 3.76*invbeta2*Z*Z*chargeSquare*alpha2)/mom2;
218    if(particle == theProton && A < 1.5 && cosTetMaxNuc < 0.0) 
219      cosTetMaxNuc = 0.0;
220    // A.V. Butkevich et al., NIM A 488 (2002) 282
221    x =  fNistManager->GetLOGA(A);
222    formfactA = mom2*constn*std::exp(0.54*x);
223  } 
224} 
225
226//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
227
228#endif
Note: See TracBrowser for help on using the repository browser.