Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 "G4SystemOfUnits.hh"
58 #include "Randomize.hh"
60 #include "G4Proton.hh"
61 #include "G4ProductionCutsTable.hh"
62 #include "G4NucleiProperties.hh"
63 #include "G4NistManager.hh"
64 #include "G4ParticleTable.hh"
65 #include "G4IonTable.hh"
66 
67 #include "G4UnitsTable.hh"
68 
69 
70 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
71 
72 using namespace std;
73 
75  : G4VEmModel(nam),
76 
77  cosThetaMin(1.0),
78  isInitialised(false)
79 {
80  fNistManager = G4NistManager::Instance();
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;
93  currentMaterialIndex = -1;
94 
95  Mottcross = new G4ScreeningMottCrossSection();
96 
97 }
98 
99 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
100 
102 { delete Mottcross;}
103 
104 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
105 
107  const G4DataVector& cuts)
108 {
109  SetupParticle(p);
110  currentCouple = 0;
111  currentMaterialIndex = -1;
112  cosThetaMin = cos(PolarAngleLimit());
113  Mottcross->Initialise(p,cosThetaMin);
114 
115  pCuts = &cuts;
116  //G4ProductionCutsTable::GetProductionCutsTable()->GetEnergyCutsVector(3);
117 
118 
119  if(!isInitialised) {
120  isInitialised = true;
121  fParticleChange = GetParticleChangeForGamma();
122  }
123 }
124 
125 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
126 
128  const G4ParticleDefinition* p,
129  G4double kinEnergy,
130  G4double Z,
131  G4double ,
132  G4double,
133  G4double )
134 {
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
170  currentElement = SelectRandomAtom(couple,particle,
171  kinEnergy,cutEnergy,kinEnergy);
172 
173  G4double Z = currentElement->GetZ();
174  G4int iz = G4int(Z);
175  G4int ia = SelectIsotopeNumber(currentElement);
176 
177  //cout<<"Element "<<currentElement->GetName()<<endl;;
178 
179  G4double cross= Mottcross->GetTotalCross();
180 
181  if(cross == 0.0) { return; }
182 
183  G4ThreeVector dir = dp->GetMomentumDirection(); //old direction
184  G4ThreeVector newDirection=Mottcross->GetNewDirection();//new direction
185  newDirection.rotateUz(dir);
186 
187  fParticleChange->ProposeMomentumDirection(newDirection);
188 
189  //Recoil energy
190  G4double trec= Mottcross->GetTrec();
191 
192  //Energy after scattering
193  if(trec > kinEnergy) { trec = kinEnergy; }
194  G4double finalT = kinEnergy - trec;
195  G4double edep = 0.0;
196 
197  G4double tcut = recoilThreshold;
198  if(pCuts) {
199  tcut= std::min(tcut,(*pCuts)[currentMaterialIndex]);
200  }
201 
202  if(trec > tcut) {
203 
204  //cout<<"Trec "<<trec/eV<<endl;
205  G4ParticleDefinition* ion = theIonTable->GetIon(iz, ia, 0);
206 
207  //incident before scattering
208  G4double ptot=sqrt(Mottcross->GetMom2Lab());
209  //incident after scattering
210  G4double plab = sqrt(finalT*(finalT + 2.0*mass));
211  G4ThreeVector p2 = (ptot*dir - plab*newDirection).unit();
212  //secondary particle
213  G4DynamicParticle* newdp = new G4DynamicParticle(ion, p2, trec);
214  fvect->push_back(newdp);
215  } else if(trec > 0.0) {
216  edep = trec;
217  fParticleChange->ProposeNonIonizingEnergyDeposit(trec);
218  }
219 
220  // finelize primary energy and energy balance
221  if(finalT <= lowEnergyLimit) {
222  edep += finalT;
223  finalT = 0.0;
224  }
225  fParticleChange->SetProposedKineticEnergy(finalT);
226  fParticleChange->ProposeLocalEnergyDeposit(edep);
227 
228 }
229 
230 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
231 
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4double GetKineticEnergy() const
void Initialise(const G4ParticleDefinition *, G4double cosThetaLim)
const char * p
Definition: xmltok.h:285
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:449
G4eSingleCoulombScatteringModel(const G4String &nam="eSingleCoulombScattering")
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4IonTable * GetIonTable() const
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:426
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
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.hh:548
void DefineMaterial(const G4MaterialCutsCouple *)
void SetupParticle(const G4ParticleDefinition *)
static G4ParticleTable * GetParticleTable()
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:620
T min(const T t1, const T t2)
brief Return the smallest 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:510
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax)