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

Last change on this file since 1347 was 1337, checked in by garnier, 15 years ago

tag geant4.9.4 beta 1 + modifs locales

File size: 6.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.56 2010/05/27 14:22:05 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-04-beta-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// 27.05.10 V.Ivanchenko added G4WentzelOKandVIxSection class to
50// compute cross sections and sample scattering angle
51//
52//
53// Class Description:
54//
55// Implementation of eCoulombScattering of pointlike charge particle
56// on Atomic Nucleus for interval of scattering anles in Lab system
57// thetaMin - ThetaMax, nucleus recoil is neglected.
58// The model based on analysis of J.M.Fernandez-Varea et al.
59// NIM B73(1993)447 originated from G.Wentzel Z.Phys. 40(1927)590 with
60// screening parameter from H.A.Bethe Phys. Rev. 89 (1953) 1256.
61//
62
63// -------------------------------------------------------------------
64//
65
66#ifndef G4eCoulombScatteringModel_h
67#define G4eCoulombScatteringModel_h 1
68
69#include "G4VEmModel.hh"
70#include "globals.hh"
71#include "G4MaterialCutsCouple.hh"
72#include "G4WentzelOKandVIxSection.hh"
73
74class G4ParticleChangeForGamma;
75class G4ParticleDefinition;
76class G4ParticleTable;
77class G4NistManager;
78
79class G4eCoulombScatteringModel : public G4VEmModel
80{
81
82public:
83
84 G4eCoulombScatteringModel(const G4String& nam = "eCoulombScattering");
85
86 virtual ~G4eCoulombScatteringModel();
87
88 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
89
90 virtual G4double ComputeCrossSectionPerAtom(
91 const G4ParticleDefinition*,
92 G4double kinEnergy,
93 G4double Z,
94 G4double A,
95 G4double cut,
96 G4double emax);
97
98 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
99 const G4MaterialCutsCouple*,
100 const G4DynamicParticle*,
101 G4double tmin,
102 G4double maxEnergy);
103
104 // defines low energy limit of the model
105 inline void SetLowEnergyLimit(G4double val);
106
107 // obsolete method
108 inline void SetRecoilThreshold(G4double eth);
109
110protected:
111
112 inline void DefineMaterial(const G4MaterialCutsCouple*);
113
114 inline void SetupParticle(const G4ParticleDefinition*);
115
116private:
117
118 // hide assignment operator
119 G4eCoulombScatteringModel & operator=(const G4eCoulombScatteringModel &right);
120 G4eCoulombScatteringModel(const G4eCoulombScatteringModel&);
121
122protected:
123
124 G4ParticleTable* theParticleTable;
125 G4ParticleChangeForGamma* fParticleChange;
126 G4WentzelOKandVIxSection* wokvi;
127 G4NistManager* fNistManager;
128
129 const std::vector<G4double>* pCuts;
130
131 const G4MaterialCutsCouple* currentCouple;
132 const G4Material* currentMaterial;
133 const G4Element* currentElement;
134 G4int currentMaterialIndex;
135
136 G4double cosThetaMin;
137 G4double cosThetaMax;
138 G4double cosTetMinNuc;
139 G4double cosTetMaxNuc;
140 G4double recoilThreshold;
141 G4double elecRatio;
142 G4double mass;
143
144 // projectile
145 const G4ParticleDefinition* particle;
146 const G4ParticleDefinition* theProton;
147
148 G4double lowEnergyLimit;
149
150private:
151
152 G4bool isInitialised;
153};
154
155//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
156
157inline
158void G4eCoulombScatteringModel::DefineMaterial(const G4MaterialCutsCouple* cup)
159{
160 if(cup != currentCouple) {
161 currentCouple = cup;
162 currentMaterial = cup->GetMaterial();
163 currentMaterialIndex = currentCouple->GetIndex();
164 }
165}
166
167//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
168
169inline
170void G4eCoulombScatteringModel::SetupParticle(const G4ParticleDefinition* p)
171{
172 // Initialise mass and charge
173 if(p != particle) {
174 particle = p;
175 mass = particle->GetPDGMass();
176 wokvi->SetupParticle(p);
177 }
178}
179
180//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
181
182inline void G4eCoulombScatteringModel::SetLowEnergyLimit(G4double val)
183{
184 lowEnergyLimit = val;
185}
186
187//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
188
189inline void G4eCoulombScatteringModel::SetRecoilThreshold(G4double eth)
190{
191 recoilThreshold = eth;
192}
193
194//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
195
196#endif
Note: See TracBrowser for help on using the repository browser.