Geant4  10.02.p02
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 {
81  fParticleChange = 0;
82 
83  pCuts=0;
84  currentMaterial = 0;
85  currentElement = 0;
86  currentCouple = 0;
87 
88  lowEnergyLimit = 0*keV;
89  recoilThreshold = 0.*eV;
90  particle = 0;
91  mass=0;
93 
95 }
96 
97 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
98 
100 {
101  delete Mottcross;
102 }
103 
104 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
105 
107  const G4DataVector& cuts)
108 {
109  SetupParticle(p);
110  currentCouple = 0;
112  //cosThetaMin = cos(PolarAngleLimit());
114 
115  pCuts = &cuts;
116  //G4ProductionCutsTable::GetProductionCutsTable()->GetEnergyCutsVector(3);
117 
118  /*
119  G4cout << "!!! G4eSingleCoulombScatteringModel::Initialise for "
120  << part->GetParticleName() << " cos(TetMin)= " << cosThetaMin
121  << " cos(TetMax)= " << cosThetaMax <<G4endl;
122  G4cout << "cut= " << (*pCuts)[0] << " cut1= " << (*pCuts)[1] << G4endl;
123  */
124 
125  if(!fParticleChange) {
127  }
128 
129  if(IsMaster()) {
131  }
132 }
133 
135  G4VEmModel* masterModel)
136 {
138 }
139 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
140 
142  const G4ParticleDefinition* p,
143  G4double kinEnergy,
144  G4double Z,
145  G4double ,
146  G4double,
147  G4double )
148 {
149  SetupParticle(p);
150 
151  G4double cross =0.0;
152  if(kinEnergy < lowEnergyLimit) return cross;
153 
155 
156  //Total Cross section
157  Mottcross->SetupKinematic(kinEnergy, Z);
158  cross = Mottcross->NuclearCrossSection();
159 
160  //cout<< "Compute Cross Section....cross "<<G4BestUnit(cross,"Surface") << " cm2 "<< cross/cm2 <<" Z: "<<Z<<" kinEnergy: "<<kinEnergy<<endl;
161  return cross;
162 }
163 
164 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
165 
167  std::vector<G4DynamicParticle*>* fvect,
168  const G4MaterialCutsCouple* couple,
169  const G4DynamicParticle* dp,
170  G4double cutEnergy,
171  G4double)
172 {
173  G4double kinEnergy = dp->GetKineticEnergy();
174  //cout<<"--- kinEnergy "<<kinEnergy<<endl;
175 
176  if(kinEnergy < lowEnergyLimit) return;
177 
178  DefineMaterial(couple);
180 
181  // Choose nucleus
182  //last two :cutEnergy= min e kinEnergy=max
184  kinEnergy,cutEnergy,kinEnergy);
185 
186  G4double Z = currentElement->GetZ();
187  G4int iz = G4int(Z);
188  G4int ia = SelectIsotopeNumber(currentElement);
190 
191  //G4cout<<"..Z: "<<Z<<" ..iz: "<<iz<<" ..ia: "<<ia<<" ..mass2: "<<mass2<<G4endl;
192 
193  Mottcross->SetupKinematic(kinEnergy, Z);
194  G4double cross= Mottcross->NuclearCrossSection(); //MODIFY TO LOAD TABLE
195  if(cross == 0.0) { return; }
196  //cout<< "Energy: "<<kinEnergy/MeV<<" Z: "<<Z<<"....cross "<<G4BestUnit(cross,"Surface") << " cm2 "<< cross/cm2 <<endl;
197 
199  G4double sint = sin(z1);
200  G4double cost = sqrt(1.0 - sint*sint);
201  G4double phi = twopi* G4UniformRand();
202 
203  // kinematics in the Lab system
204  G4double ptot = dp->GetTotalMomentum();
205  G4double e1 = dp->GetTotalEnergy();
206  // Lab. system kinematics along projectile direction
207  G4LorentzVector v0 = G4LorentzVector(0, 0, ptot, e1);
208  G4double bet = ptot/(v0.e() + mass2);
209  G4double gam = 1.0/sqrt((1.0 - bet)*(1.0 + bet));
210 
211  //CM Projectile
212  G4double momCM = gam*(ptot - bet*e1);
213  G4double eCM = gam*(e1 - bet*ptot);
214  //energy & momentum after scattering of incident particle
215  G4double pxCM = momCM*sint*cos(phi);
216  G4double pyCM = momCM*sint*sin(phi);
217  G4double pzCM = momCM*cost;
218 
219  //CM--->Lab
220  G4LorentzVector v1(pxCM , pyCM, gam*(pzCM + bet*eCM), gam*(eCM + bet*pzCM));
221 
222  // Rotate to global system
223  G4ThreeVector dir = dp->GetMomentumDirection();
224  G4ThreeVector newDirection = v1.vect().unit();
225  newDirection.rotateUz(dir);
226 
228 
229  // recoil
230  v0 -= v1;
231  G4double trec = v0.e();
232  G4double edep = 0.0;
233 
234  G4double tcut = recoilThreshold;
235 
236  //G4cout<<" Energy Transfered: "<<trec/eV<<G4endl;
237 
238  if(pCuts) {
239  tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]);
240  }
241 
242  if(trec > tcut) {
244  newDirection = v0.vect().unit();
245  newDirection.rotateUz(dir);
246  G4DynamicParticle* newdp = new G4DynamicParticle(ion, newDirection, trec);
247  fvect->push_back(newdp);
248  } else if(trec > 0.0) {
249  edep = trec;
251  }
252 
253  // finelize primary energy and energy balance
254  G4double finalT = v1.e() - mass;
255  //G4cout<<"Final Energy: "<<finalT/eV<<G4endl;
256  if(finalT <= lowEnergyLimit) {
257  edep += finalT;
258  finalT = 0.0;
259  }
260  edep = std::max(edep, 0.0);
263 
264 }
265 
266 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
267 
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)
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148
G4double GetTotalEnergy() const
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:491
G4bool IsMaster() const
Definition: G4VEmModel.hh:718
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
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:802
static const G4double e1
void DefineMaterial(const G4MaterialCutsCouple *)
static const double eV
Definition: G4SIunits.hh:212
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:810
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)
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
static const double keV
Definition: G4SIunits.hh:213
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