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
G4eSingleCoulombScatteringModel.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// G4eSingleCoulombScatteringModel.cc
27// -------------------------------------------------------------------
28//
29// GEANT4 Class header file
30//
31// File name: G4eSingleCoulombScatteringModel
32//
33// Author: Cristina Consolandi
34//
35// Creation date: 20.10.2012
36//
37// Class Description:
38// Single Scattering model for electron-nuclei interaction.
39// Suitable for high energy electrons and low scattering angles.
40//
41//
42// Reference:
43// M.J. Boschini et al.
44// "Non Ionizing Energy Loss induced by Electrons in the Space Environment"
45// Proc. of the 13th International Conference on Particle Physics and Advanced Technology
46// (13th ICPPAT, Como 3-7/10/2011), World Scientific (Singapore).
47// Available at: http://arxiv.org/abs/1111.4042v4
48//
49//
50// -------------------------------------------------------------------
51//
52//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
53
54
56#include "G4SystemOfUnits.hh"
57#include "Randomize.hh"
59#include "G4Proton.hh"
61#include "G4NucleiProperties.hh"
62
63#include "G4UnitsTable.hh"
64
65
66//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
67
68using namespace std;
69
71 : G4VEmModel(nam),
72
73 cosThetaMin(1.0),
74 isInitialised(false)
75{
79
80 pCuts=0;
83 currentCouple = 0;
84
85 lowEnergyLimit = 0*eV;
86 recoilThreshold = 0.*eV;
87 particle = 0;
88 mass=0;
90
92
93}
94
95
96//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
97
99{ delete Mottcross;}
100
101//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
102
104 const G4DataVector& )
105{
106 SetupParticle(p);
107 currentCouple = 0;
111
113
114
115 if(!isInitialised) {
116 isInitialised = true;
118 }
119}
120
121//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
122
124 const G4ParticleDefinition* p,
125 G4double kinEnergy,
126 G4double Z,
127 G4double ,
128 G4double,
129 G4double )
130{
131
132 SetupParticle(p);
133
134 G4double cross =0.0;
135 if(kinEnergy < lowEnergyLimit) return cross;
136
138
139 //Total Cross section
140 Mottcross->SetupKinematic(kinEnergy, Z);
142
143 //cout<< "....cross "<<G4BestUnit(cross,"Surface") << " cm2 "<< cross/cm2 <<endl;
144 return cross;
145}
146
147//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
148
150 std::vector<G4DynamicParticle*>* fvect,
151 const G4MaterialCutsCouple* couple,
152 const G4DynamicParticle* dp,
153 G4double cutEnergy,
154 G4double)
155{
156 G4double kinEnergy = dp->GetKineticEnergy();
157 //cout<<"--- kinEnergy "<<kinEnergy<<endl;
158
159
160 if(kinEnergy < lowEnergyLimit) return;
161
162 DefineMaterial(couple);
164
165 // Choose nucleus
167 kinEnergy,cutEnergy,kinEnergy);//last two :cutEnergy= min e kinEnergy=max
168
170 G4int iz = G4int(Z);
172
173 //cout<<"Element "<<currentElement->GetName()<<endl;;
174
176
177
178 if(cross == 0.0)return;
179
180 G4ThreeVector dir = dp->GetMomentumDirection(); //old direction
181 G4ThreeVector newDirection=Mottcross->GetNewDirection();//new direction
182 newDirection.rotateUz(dir);
183
185
186 //Recoil energy
187 G4double trec= Mottcross->GetTrec();
188 //Energy after scattering
189 G4double finalT = kinEnergy - trec;
190
191
192 if(finalT <= lowEnergyLimit) {
193 trec = kinEnergy;
194 finalT = 0.0;
195 }
196
198
200 if(pCuts) { tcut= std::min(tcut,(*pCuts)[currentMaterialIndex]);
201 }
202
203
204 if(trec > tcut) {
205
206 //cout<<"Trec "<<trec/eV<<endl;
207 G4ParticleDefinition* ion = theParticleTable->GetIon(iz, ia, 0.0);
208
209 //incident before scattering
210 G4double ptot=sqrt(Mottcross->GetMom2Lab());
211 //incident after scattering
212 G4double plab = sqrt(finalT*(finalT + 2.0*mass));
213 G4ThreeVector p2 = (ptot*dir - plab*newDirection).unit();
214 //secondary particle
215 G4DynamicParticle* newdp = new G4DynamicParticle(ion, p2, trec);
216 fvect->push_back(newdp);
217 }
218
219 else if(trec > 0.0) {
221 if(trec< tcut) fParticleChange->ProposeLocalEnergyDeposit(trec);
222 }
223
224
225}
226
227//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
228
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
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()
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()
void SetupKinematic(G4double kinEnergy, G4double Z)
void Initialise(const G4ParticleDefinition *, G4double cosThetaLim)
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 ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetupParticle(const G4ParticleDefinition *)
G4eSingleCoulombScatteringModel(const G4String &nam="eSingleCoulombScattering")
void DefineMaterial(const G4MaterialCutsCouple *)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax)