Geant4  10.01.p01
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 81579 2014-06-03 10:15:54Z 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("eCoulombScattering"),
82  cosThetaMin(1.0),
83  cosThetaMax(-1.0),
84  isCombined(combined),
85  isInitialised(false)
86 {
87  fParticleChange = 0;
91  currentMaterial = 0;
92  fixedCut = -1.0;
93 
94  pCuts = 0;
95 
96  lowEnergyThreshold = 1*keV; // particle will be killed for lower energy
97  recoilThreshold = 0.*keV; // by default does not work
98 
99  particle = 0;
100  currentCouple = 0;
101  wokvi = new G4WentzelOKandVIxSection(combined);
102 
104  mass = proton_mass_c2;
105  elecRatio = 0.0;
106 }
107 
108 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
109 
111 {
112  delete wokvi;
113 }
114 
115 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
116 
118  const G4DataVector& cuts)
119 {
120  SetupParticle(part);
121  currentCouple = 0;
122 
123  if(isCombined) {
124  cosThetaMin = 1.0;
125  G4double tet = PolarAngleLimit();
126  if(tet >= pi) { cosThetaMin = -1.0; }
127  else if(tet > 0.0) { cosThetaMin = cos(tet); }
128  }
129 
130  wokvi->Initialise(part, cosThetaMin);
131  /*
132  G4cout << "G4eCoulombScatteringModel: " << particle->GetParticleName()
133  << " 1-cos(ThetaLimit)= " << 1 - cosThetaMin
134  << " cos(thetaMax)= " << cosThetaMax
135  << G4endl;
136  */
137  pCuts = &cuts;
138  //G4ProductionCutsTable::GetProductionCutsTable()->GetEnergyCutsVector(3);
139  /*
140  G4cout << "!!! G4eCoulombScatteringModel::Initialise for "
141  << part->GetParticleName() << " cos(TetMin)= " << cosThetaMin
142  << " cos(TetMax)= " << cosThetaMax <<G4endl;
143  G4cout << "cut= " << (*pCuts)[0] << " cut1= " << (*pCuts)[1] << G4endl;
144  */
145  if(!isInitialised) {
146  isInitialised = true;
148  }
149  if(IsMaster() && mass < GeV && part->GetParticleName() != "GenericIon") {
150  InitialiseElementSelectors(part,cuts);
151  }
152 }
153 
154 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
155 
157  G4VEmModel* masterModel)
158 {
160 }
161 
162 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
163 
164 G4double
166  const G4ParticleDefinition* part,
167  G4double)
168 {
169  SetupParticle(part);
170 
171  // define cut using cuts for proton
172  G4double cut =
173  std::max(recoilThreshold, (*pCuts)[CurrentCouple()->GetIndex()]);
174 
175  // find out lightest element
176  const G4ElementVector* theElementVector = material->GetElementVector();
177  G4int nelm = material->GetNumberOfElements();
178 
179  // select lightest element
180  G4int Z = 300;
181  for (G4int j=0; j<nelm; ++j) {
182  G4int iz = (G4int)(*theElementVector)[j]->GetZ();
183  if(iz < Z) { Z = iz; }
184  }
186  G4double targetMass = G4NucleiProperties::GetNuclearMass(A, Z);
187  G4double t = std::max(cut, 0.5*(cut + sqrt(2*cut*targetMass)));
188 
189  return std::max(lowEnergyThreshold, t);
190 }
191 
192 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
193 
195  const G4ParticleDefinition* p,
196  G4double kinEnergy,
197  G4double Z, G4double,
198  G4double cutEnergy, G4double)
199 {
200  //G4cout << "### G4eCoulombScatteringModel::ComputeCrossSectionPerAtom for "
201  //<< p->GetParticleName()<<" Z= "<<Z<<" e(MeV)= "<< kinEnergy/MeV << G4endl;
202  G4double cross = 0.0;
203  elecRatio = 0.0;
204  if(p != particle) { SetupParticle(p); }
205 
206  // cross section is set to zero to avoid problems in sample secondary
207  if(kinEnergy <= 0.0) { return cross; }
209  G4double costmin = wokvi->SetupKinematic(kinEnergy, currentMaterial);
210  if(cosThetaMax < costmin) {
211  G4int iz = G4int(Z);
212  G4double cut = cutEnergy;
213  if(fixedCut > 0.0) { cut = fixedCut; }
214  costmin = wokvi->SetupTarget(iz, cut);
215  G4double costmax = cosThetaMax;
216  if(iz == 1 && costmax < 0.0 && particle == theProton) {
217  costmax = 0.0;
218  }
219  if(costmin > costmax) {
220  cross = wokvi->ComputeNuclearCrossSection(costmin, costmax)
221  + wokvi->ComputeElectronCrossSection(costmin, costmax);
222  }
223  /*
224  if(p->GetParticleName() == "mu+")
225  G4cout << "e(MeV)= " << kinEnergy/MeV << " cross(b)= " << cross/barn
226  << " 1-costmin= " << 1-costmin
227  << " 1-costmax= " << 1-costmax
228  << " 1-cosThetaMax= " << 1-cosThetaMax
229  << G4endl;
230  */
231  }
232  return cross;
233 }
234 
235 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
236 
238  std::vector<G4DynamicParticle*>* fvect,
239  const G4MaterialCutsCouple* couple,
240  const G4DynamicParticle* dp,
241  G4double cutEnergy,
242  G4double)
243 {
244  G4double kinEnergy = dp->GetKineticEnergy();
245 
246  // absorb particle below low-energy limit to avoid situation
247  // when a particle has no energy loss
248  if(kinEnergy < lowEnergyThreshold) {
252  return;
253  }
255  DefineMaterial(couple);
256  /*
257  G4cout << "G4eCoulombScatteringModel::SampleSecondaries e(MeV)= "
258  << kinEnergy << " " << particle->GetParticleName()
259  << " cut= " << cutEnergy<< G4endl;
260  */
261  // Choose nucleus
262  G4double cut = cutEnergy;
263  if(fixedCut > 0.0) { cut = fixedCut; }
264 
265  wokvi->SetupKinematic(kinEnergy, currentMaterial);
266 
267  const G4Element* currentElement =
268  SelectRandomAtom(couple,particle,kinEnergy,cut,kinEnergy);
269 
270  G4double Z = currentElement->GetZ();
271  G4int iz = G4int(Z);
272 
273  G4double costmin = wokvi->SetupTarget(iz, cut);
274  G4double costmax = cosThetaMax;
275  if(iz == 1 && costmax < 0.0 && particle == theProton) {
276  costmax = 0.0;
277  }
278 
279  if(costmin > costmax) {
280  G4double cross = wokvi->ComputeNuclearCrossSection(costmin, costmax);
281  G4double ecross = wokvi->ComputeElectronCrossSection(costmin, costmax);
282  G4double ratio = ecross/(cross + ecross);
283 
284  G4int ia = SelectIsotopeNumber(currentElement);
285  G4double targetMass = G4NucleiProperties::GetNuclearMass(ia, iz);
286  wokvi->SetTargetMass(targetMass);
287 
288  G4ThreeVector newDirection =
289  wokvi->SampleSingleScattering(costmin, costmax, ratio);
290  G4double cost = newDirection.z();
291  /*
292  G4cout << "SampleSec: e(MeV)= " << kinEnergy/MeV
293  << " 1-costmin= " << 1-costmin
294  << " 1-costmax= " << 1-costmax
295  << " 1-cost= " << 1-cost
296  << " ratio= " << ratio
297  << G4endl;
298  */
299  G4ThreeVector direction = dp->GetMomentumDirection();
300  newDirection.rotateUz(direction);
301 
303 
304  // recoil sampling assuming a small recoil
305  // and first order correction to primary 4-momentum
306  G4double mom2 = wokvi->GetMomentumSquare();
307  G4double trec = mom2*(1.0 - cost)
308  /(targetMass + (mass + kinEnergy)*(1.0 - cost));
309 
310  // the check likely not needed
311  if(trec > kinEnergy) { trec = kinEnergy; }
312  G4double finalT = kinEnergy - trec;
313  G4double edep = 0.0;
314  /*
315  G4cout<<"G4eCoulombScatteringModel: finalT= "<<finalT<<" Trec= "
316  <<trec << " Z= " << iz << " A= " << ia
317  << " tcut(keV)= " << (*pCuts)[currentMaterialIndex]/keV << G4endl;
318  */
319  G4double tcut = recoilThreshold;
320  if(pCuts) { tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]); }
321 
322  if(trec > tcut) {
324  G4ThreeVector dir = (direction*sqrt(mom2) -
325  newDirection*sqrt(finalT*(2*mass + finalT))).unit();
326  G4DynamicParticle* newdp = new G4DynamicParticle(ion, dir, trec);
327  fvect->push_back(newdp);
328  } else {
329  edep = trec;
331  }
332 
333  // finelize primary energy and energy balance
334  // this threshold may be applied only because for low-enegry
335  // e+e- msc model is applied
336  if(finalT <= lowEnergyThreshold) {
337  edep += finalT;
338  finalT = 0.0;
339  }
342  }
343 }
344 
345 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
346 
347 
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:148
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:463
const G4double pi
G4double GetZ() const
Definition: G4Element.hh:131
G4bool IsMaster() const
Definition: G4VEmModel.hh:699
G4ParticleDefinition * GetDefinition() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
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)
G4IonTable * GetIonTable() const
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:433
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
bool G4bool
Definition: G4Types.hh:79
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:566
G4eCoulombScatteringModel(G4bool combined=true)
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:783
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:791
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:643
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:186
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:525
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134
G4ThreeVector & SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio=0.0)