Geant4  10.02.p01
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. "Non Ionizing Energy Loss induced by Electrons
44 // in the Space Environment" Proc. of the 13th International Conference
45 // on Particle Physics and Advanced Technology
46 //
47 // (13th ICPPAT, Como 3-7/10/2011), World Scientific (Singapore).
48 // Available at: http://arxiv.org/abs/1111.4042v4
49 //
50 //
51 // -------------------------------------------------------------------
52 //
53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
54 
55 
57 #include "G4PhysicalConstants.hh"
58 #include "G4SystemOfUnits.hh"
59 #include "Randomize.hh"
61 #include "G4Proton.hh"
62 #include "G4ProductionCutsTable.hh"
63 #include "G4NucleiProperties.hh"
64 #include "G4NistManager.hh"
65 #include "G4ParticleTable.hh"
66 #include "G4IonTable.hh"
67 
68 #include "G4UnitsTable.hh"
69 
70 
71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
72 
73 using namespace std;
74 
76  : G4VEmModel(nam),
77  cosThetaMin(1.0),
78  isInitialised(false)
79 {
82  fParticleChange = 0;
83 
84  pCuts=0;
85  currentMaterial = 0;
86  currentElement = 0;
87  currentCouple = 0;
88 
89  lowEnergyLimit = 0*eV;
90  recoilThreshold = 0.*eV;
91  particle = 0;
92  mass=0;
94 
96 }
97 
98 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
99 
101 {
102  delete Mottcross;
103 }
104 
105 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
106 
108  const G4DataVector& cuts)
109 {
110  SetupParticle(p);
111  currentCouple = 0;
113  //cosThetaMin = cos(PolarAngleLimit());
115 
116  pCuts = &cuts;
117  //G4ProductionCutsTable::GetProductionCutsTable()->GetEnergyCutsVector(3);
118 
119 
120  if(!isInitialised) {
121  isInitialised = true;
123  }
124 }
125 
126 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
127 
129  const G4ParticleDefinition* p,
130  G4double kinEnergy,
131  G4double Z,
132  G4double ,
133  G4double,
134  G4double )
135 {
136  SetupParticle(p);
137 
138  G4double cross =0.0;
139  if(kinEnergy < lowEnergyLimit) return cross;
140 
142 
143  //Total Cross section
144  Mottcross->SetupKinematic(kinEnergy, Z);
145  cross = Mottcross->NuclearCrossSection();
146 
147  //cout<< "....cross "<<G4BestUnit(cross,"Surface") << " cm2 "<< cross/cm2 <<endl;
148  return cross;
149 }
150 
151 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
152 
154  std::vector<G4DynamicParticle*>* fvect,
155  const G4MaterialCutsCouple* couple,
156  const G4DynamicParticle* dp,
157  G4double cutEnergy,
158  G4double)
159 {
160  G4double kinEnergy = dp->GetKineticEnergy();
161  //cout<<"--- kinEnergy "<<kinEnergy<<endl;
162 
163  if(kinEnergy < lowEnergyLimit) return;
164 
165  DefineMaterial(couple);
167 
168  // Choose nucleus
169  //last two :cutEnergy= min e kinEnergy=max
171  kinEnergy,cutEnergy,kinEnergy);
172 
173  G4double Z = currentElement->GetZ();
174  G4int iz = G4int(Z);
175  G4int ia = SelectIsotopeNumber(currentElement);
177 
178  //G4cout<<"..Z: "<<Z<<" ..iz: "<<iz<<" ..ia: "<<ia<<" ..mass2: "<<mass2<<G4endl;
179 
180  Mottcross->SetupKinematic(kinEnergy, Z);
181 
183  G4double sint = sin(z1);
184  G4double cost = sqrt(1.0 - sint*sint);
185  G4double phi = twopi* G4UniformRand();
186 
187  // kinematics in the Lab system
188  G4double ptot = dp->GetTotalMomentum();
189  G4double e1 = dp->GetTotalEnergy();
190  // Lab. system kinematics along projectile direction
191  G4LorentzVector v0 = G4LorentzVector(0, 0, ptot, e1);
192  G4double bet = ptot/(v0.e() + mass2);
193  G4double gam = 1.0/sqrt((1.0 - bet)*(1.0 + bet));
194 
195  //CM Projectile
196  G4double momCM = gam*(ptot - bet*e1);
197  G4double eCM = gam*(e1 - bet*ptot);
198  //energy & momentum after scattering of incident particle
199  G4double pxCM = momCM*sint*cos(phi);
200  G4double pyCM = momCM*sint*sin(phi);
201  G4double pzCM = momCM*cost;
202 
203  //CM--->Lab
204  G4LorentzVector v1(pxCM , pyCM, gam*(pzCM + bet*eCM), gam*(eCM + bet*pzCM));
205 
206  // Rotate to global system
207  G4ThreeVector dir = dp->GetMomentumDirection();
208  G4ThreeVector newDirection = v1.vect().unit();
209  newDirection.rotateUz(dir);
210 
212 
213  // recoil
214  v0 -= v1;
215  G4double trec = v0.e();
216  G4double edep = 0.0;
217 
218  G4double tcut = recoilThreshold;
219 
220  //G4cout<<" Energy Transfered: "<<trec/eV<<G4endl;
221 
222  if(pCuts) {
223  tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]);
224  }
225 
226  if(trec > tcut) {
228  newDirection = v0.vect().unit();
229  newDirection.rotateUz(dir);
230  G4DynamicParticle* newdp = new G4DynamicParticle(ion, newDirection, trec);
231  fvect->push_back(newdp);
232  } else if(trec > 0.0) {
233  edep = trec;
235  }
236 
237  // finelize primary energy and energy balance
238  G4double finalT = v1.e() - mass;
239  //G4cout<<"Final Energy: "<<finalT/eV<<G4endl;
240  if(finalT <= lowEnergyLimit) {
241  edep += finalT;
242  finalT = 0.0;
243  }
244  edep = std::max(edep, 0.0);
247 
248 }
249 
250 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
251 
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
void Initialise(const G4ParticleDefinition *, G4double cosThetaLim)
G4double GetTotalEnergy() const
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:491
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double GetTotalMomentum() const
G4IonTable * GetIonTable() const
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:452
#define G4UniformRand()
Definition: Randomize.hh:97
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
const G4ThreeVector & GetMomentumDirection() const
G4double iz
Definition: TRTMaterials.hh:39
static const double twopi
Definition: G4SIunits.hh:75
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.hh:585
static const G4double e1
void DefineMaterial(const G4MaterialCutsCouple *)
static const double eV
Definition: G4SIunits.hh:212
void SetupParticle(const G4ParticleDefinition *)
G4eSingleCoulombScatteringModel(const G4String &nam="eSingleCoulombScat")
static G4ParticleTable * GetParticleTable()
T max(const T t1, const T t2)
brief Return the largest of the two arguments
void SetProposedKineticEnergy(G4double proposedKinEnergy)
double G4double
Definition: G4Types.hh:76
void SetupKinematic(G4double kinEnergy, G4double Z)
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:544
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax)
CLHEP::HepLorentzVector G4LorentzVector