Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4eCoulombScatteringModel.cc
Go to the documentation of this file.
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$
27//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class file
31//
32//
33// File name: G4eCoulombScatteringModel
34//
35// Author: Vladimir Ivanchenko
36//
37// Creation date: 22.08.2005
38//
39// Modifications:
40//
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
45// 09.10.07 V.Ivanchenko reorganized methods, add cut dependence in scattering off e-
46// 09.06.08 V.Ivanchenko add SelectIsotope and sampling of the recoil ion
47// 16.06.09 C.Consolandi fixed computation of effective mass
48// 27.05.10 V.Ivanchenko added G4WentzelOKandVIxSection class to
49// compute cross sections and sample scattering angle
50//
51//
52// Class Description:
53//
54// -------------------------------------------------------------------
55//
56//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
57//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
58
61#include "G4SystemOfUnits.hh"
62#include "Randomize.hh"
63#include "G4DataVector.hh"
64#include "G4ElementTable.hh"
66#include "G4Proton.hh"
67#include "G4ParticleTable.hh"
69#include "G4NucleiProperties.hh"
70#include "G4Pow.hh"
71#include "G4LossTableManager.hh"
72#include "G4LossTableBuilder.hh"
73#include "G4NistManager.hh"
74
75//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
76
77using namespace std;
78
80 : G4VEmModel(nam),
81 cosThetaMin(1.0),
82 cosThetaMax(-1.0),
83 isInitialised(false)
84{
89 currentMaterial = 0;
90
91 pCuts = 0;
92
93 lowEnergyThreshold = 1*keV; // particle will be killed for lower energy
94 recoilThreshold = 0.*keV; // by default does not work
95
96 particle = 0;
97 currentCouple = 0;
99
101
102 cosTetMinNuc = 1.0;
103 cosTetMaxNuc = -1.0;
104 elecRatio = 0.0;
105 mass = proton_mass_c2;
106}
107
108//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
109
111{
112 delete wokvi;
113}
114
115//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
116
118 const G4DataVector& cuts)
119{
120 SetupParticle(p);
121 currentCouple = 0;
124 /*
125 G4cout << "G4eCoulombScatteringModel: " << particle->GetParticleName()
126 << " 1-cos(ThetaLimit)= " << 1 - cosThetaMin
127 << " cos(thetaMax)= " << cosThetaMax
128 << G4endl;
129 */
131 /*
132 G4cout << "!!! G4eCoulombScatteringModel::Initialise for "
133 << p->GetParticleName() << " cos(TetMin)= " << cosThetaMin
134 << " cos(TetMax)= " << cosThetaMax <<G4endl;
135 G4cout << "cut0= " << cuts[0] << " cut1= " << cuts[1] << G4endl;
136 */
137 if(!isInitialised) {
138 isInitialised = true;
140 }
141 if(mass < GeV) {
143 }
144}
145
146//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
147
149 const G4ParticleDefinition* p,
150 G4double kinEnergy,
152 G4double cutEnergy, G4double)
153{
154 //G4cout << "### G4eCoulombScatteringModel::ComputeCrossSectionPerAtom for "
155 // << p->GetParticleName()<<" Z= "<<Z<<" e(MeV)= "<< kinEnergy/MeV << G4endl;
156 G4double cross = 0.0;
157 if(p != particle) { SetupParticle(p); }
158
159 // cross section is set to zero to avoid problems in sample secondary
160 if(kinEnergy <= 0.0) { return cross; }
164 G4int iz = G4int(Z);
165 cosTetMinNuc = wokvi->SetupTarget(iz, cutEnergy);
167 if(iz == 1 && cosTetMaxNuc < 0.0 && particle == theProton) {
168 cosTetMaxNuc = 0.0;
169 }
172 cross += elecRatio;
173 if(cross > 0.0) { elecRatio /= cross; }
174 }
175 /*
176 if(p->GetParticleName() == "e-")
177 G4cout << "e(MeV)= " << kinEnergy/MeV << " cross(b)= " << cross/barn
178 << " 1-cosTetMinNuc= " << 1-cosTetMinNuc
179 << " 1-cosTetMaxNuc= " << 1-cosTetMaxNuc
180 << G4endl;
181 */
182 return cross;
183}
184
185//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
186
188 std::vector<G4DynamicParticle*>* fvect,
189 const G4MaterialCutsCouple* couple,
190 const G4DynamicParticle* dp,
191 G4double cutEnergy,
192 G4double)
193{
194 G4double kinEnergy = dp->GetKineticEnergy();
195
196 // absorb particle below low-energy limit to avoid situation
197 // when a particle has no energy loss
198 if(kinEnergy < lowEnergyThreshold) {
202 return;
203 }
205 DefineMaterial(couple);
206 /*
207 G4cout << "G4eCoulombScatteringModel::SampleSecondaries e(MeV)= "
208 << kinEnergy << " " << particle->GetParticleName()
209 << " cut= " << cutEnergy<< G4endl;
210 */
211 // Choose nucleus
212 const G4Element* currentElement =
213 SelectRandomAtom(couple,particle,kinEnergy,cutEnergy,kinEnergy);
214
215 G4double Z = currentElement->GetZ();
216
217 if(ComputeCrossSectionPerAtom(particle,kinEnergy, Z,
218 kinEnergy, cutEnergy, kinEnergy) == 0.0)
219 { return; }
220
221 G4int iz = G4int(Z);
222 G4int ia = SelectIsotopeNumber(currentElement);
223 G4double targetMass = G4NucleiProperties::GetNuclearMass(ia, iz);
224 wokvi->SetTargetMass(targetMass);
225
226 G4ThreeVector newDirection =
228 G4double cost = newDirection.z();
229
230 G4ThreeVector direction = dp->GetMomentumDirection();
231 newDirection.rotateUz(direction);
232
234
235 // recoil sampling assuming a small recoil
236 // and first order correction to primary 4-momentum
238 G4double trec = mom2*(1.0 - cost)/(targetMass + (mass + kinEnergy)*(1.0 - cost));
239 G4double finalT = kinEnergy - trec;
240 //G4cout<<"G4eCoulombScatteringModel: finalT= "<<finalT<<" Trec= "<<trec<<G4endl;
241 if(finalT <= lowEnergyThreshold) {
242 trec = kinEnergy;
243 finalT = 0.0;
244 }
245
248 if(pCuts) { tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]); }
249
250 if(trec > tcut) {
251 G4ParticleDefinition* ion = theParticleTable->GetIon(iz, ia, 0.0);
252 G4ThreeVector dir = (direction*sqrt(mom2) -
253 newDirection*sqrt(finalT*(2*mass + finalT))).unit();
254 G4DynamicParticle* newdp = new G4DynamicParticle(ion, dir, trec);
255 fvect->push_back(newdp);
256 } else {
259 }
260
261 return;
262}
263
264//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
265
266
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
double z() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
static G4NistManager * Instance()
static G4double GetNuclearMass(const G4double A, const G4double Z)
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
static G4ParticleTable * GetParticleTable()
G4ParticleDefinition * GetIon(G4int atomicNumber, G4int atomicMass, G4double excitationEnergy)
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:550
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:459
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.hh:478
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:377
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:123
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
G4ThreeVector SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio=0.0)
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)
G4double SetupTarget(G4int Z, G4double cut=DBL_MAX)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4eCoulombScatteringModel(const G4String &nam="eCoulombScattering")
G4WentzelOKandVIxSection * wokvi
const G4ParticleDefinition * particle
const G4MaterialCutsCouple * currentCouple
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4ParticleChangeForGamma * fParticleChange
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax)
void DefineMaterial(const G4MaterialCutsCouple *)
const std::vector< G4double > * pCuts
const G4ParticleDefinition * theProton
void SetupParticle(const G4ParticleDefinition *)