Geant4  10.00.p02
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: G4eCoulombScatteringModel.cc 79188 2014-02-20 09:22:48Z gcosmo $
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 
60 #include "G4PhysicalConstants.hh"
61 #include "G4SystemOfUnits.hh"
62 #include "Randomize.hh"
63 #include "G4DataVector.hh"
64 #include "G4ElementTable.hh"
66 #include "G4Proton.hh"
67 #include "G4ParticleTable.hh"
68 #include "G4IonTable.hh"
69 #include "G4ProductionCutsTable.hh"
70 #include "G4NucleiProperties.hh"
71 #include "G4Pow.hh"
72 #include "G4LossTableManager.hh"
73 #include "G4LossTableBuilder.hh"
74 #include "G4NistManager.hh"
75 
76 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
77 
78 using namespace std;
79 
81  : G4VEmModel(nam),
82  cosThetaMin(1.0),
83  cosThetaMax(-1.0),
84  isInitialised(false)
85 {
86  fParticleChange = 0;
90  currentMaterial = 0;
91  fixedCut = -1.0;
92 
93  pCuts = 0;
94 
95  lowEnergyThreshold = 1*keV; // particle will be killed for lower energy
96  recoilThreshold = 0.*keV; // by default does not work
97 
98  particle = 0;
99  currentCouple = 0;
101 
103 
104  cosTetMinNuc = 1.0;
105  cosTetMaxNuc = -1.0;
106  elecRatio = 0.0;
107  mass = proton_mass_c2;
108 }
109 
110 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
111 
113 {
114  delete wokvi;
115 }
116 
117 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
118 
120  const G4DataVector& cuts)
121 {
122  SetupParticle(part);
123  currentCouple = 0;
124  cosThetaMin = cos(PolarAngleLimit());
125  wokvi->Initialise(part, cosThetaMin);
126  /*
127  G4cout << "G4eCoulombScatteringModel: " << particle->GetParticleName()
128  << " 1-cos(ThetaLimit)= " << 1 - cosThetaMin
129  << " cos(thetaMax)= " << cosThetaMax
130  << G4endl;
131  */
132  pCuts = &cuts;
133  // G4ProductionCutsTable::GetProductionCutsTable()->GetEnergyCutsVector(3);
134  /*
135  G4cout << "!!! G4eCoulombScatteringModel::Initialise for "
136  << part->GetParticleName() << " cos(TetMin)= " << cosThetaMin
137  << " cos(TetMax)= " << cosThetaMax <<G4endl;
138  G4cout << "cut= " << pCuts[0] << " cut1= " << pCuts[1] << G4endl;
139  */
140  if(!isInitialised) {
141  isInitialised = true;
143  }
144  if(IsMaster() && mass < GeV && part->GetParticleName() != "GenericIon") {
145  InitialiseElementSelectors(part,cuts);
146  }
147 }
148 
149 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
150 
152  G4VEmModel* masterModel)
153 {
155 }
156 
157 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
158 
159 G4double
161  const G4ParticleDefinition* part,
162  G4double)
163 {
164  SetupParticle(part);
165 
166  // define cut using cuts for proton
167  G4double cut =
168  std::max(recoilThreshold, (*pCuts)[CurrentCouple()->GetIndex()]);
169 
170  // find out lightest element
171  const G4ElementVector* theElementVector = material->GetElementVector();
172  G4int nelm = material->GetNumberOfElements();
173  G4int Z = 300;
174  for (G4int j=0; j<nelm; ++j) {
175  G4int iz = (G4int)(*theElementVector)[j]->GetZ();
176  if(iz < Z) { Z = iz; }
177  }
179  G4double targetMass = G4NucleiProperties::GetNuclearMass(A, Z);
180  G4double t = std::max(cut, 0.5*(cut + sqrt(2*cut*targetMass)));
181 
182  return std::max(lowEnergyThreshold, t);
183 }
184 
185 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
186 
188  const G4ParticleDefinition* p,
189  G4double kinEnergy,
190  G4double Z, G4double,
191  G4double cutEnergy, G4double)
192 {
193  //G4cout << "### G4eCoulombScatteringModel::ComputeCrossSectionPerAtom for "
194  //<< p->GetParticleName()<<" Z= "<<Z<<" e(MeV)= "<< kinEnergy/MeV << G4endl;
195  G4double cross = 0.0;
196  if(p != particle) { SetupParticle(p); }
197 
198  // cross section is set to zero to avoid problems in sample secondary
199  if(kinEnergy <= 0.0) { return cross; }
202  if(cosThetaMax < cosTetMinNuc) {
203  G4int iz = G4int(Z);
204  G4double cut = cutEnergy;
205  if(fixedCut > 0.0) { cut = fixedCut; }
206  cosTetMinNuc = wokvi->SetupTarget(iz, cut);
208  if(iz == 1 && cosTetMaxNuc < 0.0 && particle == theProton) {
209  cosTetMaxNuc = 0.0;
210  }
213  cross += elecRatio;
214  if(cross > 0.0) { elecRatio /= cross; }
215  }
216  /*
217  if(p->GetParticleName() == "e-")
218  G4cout << "e(MeV)= " << kinEnergy/MeV << " cross(b)= " << cross/barn
219  << " 1-cosTetMinNuc= " << 1-cosTetMinNuc
220  << " 1-cosTetMaxNuc= " << 1-cosTetMaxNuc
221  << G4endl;
222  */
223  return cross;
224 }
225 
226 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
227 
229  std::vector<G4DynamicParticle*>* fvect,
230  const G4MaterialCutsCouple* couple,
231  const G4DynamicParticle* dp,
232  G4double cutEnergy,
233  G4double)
234 {
235  G4double kinEnergy = dp->GetKineticEnergy();
236 
237  // absorb particle below low-energy limit to avoid situation
238  // when a particle has no energy loss
239  if(kinEnergy < lowEnergyThreshold) {
243  return;
244  }
246  DefineMaterial(couple);
247  /*
248  G4cout << "G4eCoulombScatteringModel::SampleSecondaries e(MeV)= "
249  << kinEnergy << " " << particle->GetParticleName()
250  << " cut= " << cutEnergy<< G4endl;
251  */
252  // Choose nucleus
253  G4double cut = cutEnergy;
254  if(fixedCut > 0.0) { cut = fixedCut; }
255 
256  const G4Element* currentElement =
257  SelectRandomAtom(couple,particle,kinEnergy,cut,kinEnergy);
258 
259  G4double Z = currentElement->GetZ();
260 
261  if(ComputeCrossSectionPerAtom(particle,kinEnergy, Z,
262  kinEnergy, cut, kinEnergy) == 0.0)
263  { return; }
264 
265  G4int iz = G4int(Z);
266  G4int ia = SelectIsotopeNumber(currentElement);
267  G4double targetMass = G4NucleiProperties::GetNuclearMass(ia, iz);
268  wokvi->SetTargetMass(targetMass);
269 
270  G4ThreeVector newDirection =
272  G4double cost = newDirection.z();
273 
274  G4ThreeVector direction = dp->GetMomentumDirection();
275  newDirection.rotateUz(direction);
276 
278 
279  // recoil sampling assuming a small recoil
280  // and first order correction to primary 4-momentum
281  G4double mom2 = wokvi->GetMomentumSquare();
282  G4double trec = mom2*(1.0 - cost)
283  /(targetMass + (mass + kinEnergy)*(1.0 - cost));
284 
285  // the check likely not needed
286  if(trec > kinEnergy) { trec = kinEnergy; }
287  G4double finalT = kinEnergy - trec;
288  G4double edep = 0.0;
289  //G4cout<<"G4eCoulombScatteringModel: finalT= "<<finalT<<" Trec= "
290  // <<trec << " Z= " << iz << " A= " << ia<<G4endl;
291 
292  G4double tcut = recoilThreshold;
293  if(pCuts) { tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]); }
294 
295  if(trec > tcut) {
297  G4ThreeVector dir = (direction*sqrt(mom2) -
298  newDirection*sqrt(finalT*(2*mass + finalT))).unit();
299  G4DynamicParticle* newdp = new G4DynamicParticle(ion, dir, trec);
300  fvect->push_back(newdp);
301  } else {
302  edep = trec;
304  }
305 
306  // finelize primary energy and energy balance
307  // this threshold may be applied only because for low-enegry
308  // e+e- msc model is applied
309  if(finalT <= lowEnergyThreshold) {
310  edep += finalT;
311  finalT = 0.0;
312  }
315 }
316 
317 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
318 
319 
const G4ParticleDefinition * particle
static G4double GetNuclearMass(const G4double A, const G4double Z)
std::vector< G4Element * > G4ElementVector
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:135
G4ThreeVector SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio=0.0)
const std::vector< G4double > * pCuts
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
const G4MaterialCutsCouple * currentCouple
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:449
G4double GetZ() const
Definition: G4Element.hh:131
G4bool IsMaster() const
Definition: G4VEmModel.hh:676
G4ParticleDefinition * GetDefinition() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4double SetupTarget(G4int Z, G4double cut=DBL_MAX)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4eCoulombScatteringModel(const G4String &nam="eCoulombScattering")
G4IonTable * GetIonTable() const
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:426
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
const G4ThreeVector & GetMomentumDirection() const
G4double iz
Definition: TRTMaterials.hh:39
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
static G4Proton * Proton()
Definition: G4Proton.cc:93
void DefineMaterial(const G4MaterialCutsCouple *)
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.hh:548
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:760
static const G4double A[nN]
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double)
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:768
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
static G4ParticleTable * GetParticleTable()
int G4lrint(double ad)
Definition: templates.hh:163
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:620
G4double GetAtomicMassAmu(const G4String &symb) const
const G4ParticleDefinition * theProton
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
void SetProposedKineticEnergy(G4double proposedKinEnergy)
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
static const double keV
Definition: G4SIunits.hh:195
double G4double
Definition: G4Types.hh:76
G4WentzelOKandVIxSection * wokvi
void SetupParticle(const G4ParticleDefinition *)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax)
G4ParticleChangeForGamma * fParticleChange
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:510
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121