Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 96934 2016-05-18 09:10:41Z 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 {
86  fParticleChange = nullptr;
87  fNistManager = G4NistManager::Instance();
89  theProton = G4Proton::Proton();
90  currentMaterial = 0;
91  fixedCut = -1.0;
92 
93  pCuts = nullptr;
94 
95  recoilThreshold = 0.*keV; // by default does not work
96 
97  particle = nullptr;
98  currentCouple = nullptr;
99  wokvi = new G4WentzelOKandVIxSection(combined);
100 
101  currentMaterialIndex = 0;
102  mass = proton_mass_c2;
103  elecRatio = 0.0;
104 }
105 
106 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
107 
109 {
110  delete wokvi;
111 }
112 
113 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
114 
116  const G4DataVector& cuts)
117 {
118  SetupParticle(part);
119  currentCouple = 0;
120 
121  if(isCombined) {
122  cosThetaMin = 1.0;
124  if(tet >= pi) { cosThetaMin = -1.0; }
125  else if(tet > 0.0) { cosThetaMin = cos(tet); }
126  }
127 
128  wokvi->Initialise(part, cosThetaMin);
129  /*
130  G4cout << "G4eCoulombScatteringModel: " << particle->GetParticleName()
131  << " 1-cos(ThetaLimit)= " << 1 - cosThetaMin
132  << " cos(thetaMax)= " << cosThetaMax
133  << G4endl;
134  */
135  pCuts = &cuts;
136  //G4ProductionCutsTable::GetProductionCutsTable()->GetEnergyCutsVector(3);
137  /*
138  G4cout << "!!! G4eCoulombScatteringModel::Initialise for "
139  << part->GetParticleName() << " cos(TetMin)= " << cosThetaMin
140  << " cos(TetMax)= " << cosThetaMax <<G4endl;
141  G4cout << "cut= " << (*pCuts)[0] << " cut1= " << (*pCuts)[1] << G4endl;
142  */
143  if(!fParticleChange) {
144  fParticleChange = GetParticleChangeForGamma();
145  }
146  if(IsMaster() && mass < GeV && part->GetParticleName() != "GenericIon") {
147  InitialiseElementSelectors(part,cuts);
148  }
149 }
150 
151 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
152 
154  G4VEmModel* masterModel)
155 {
157 }
158 
159 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
160 
161 G4double
163  const G4ParticleDefinition* part,
164  G4double)
165 {
166  SetupParticle(part);
167 
168  // define cut using cuts for proton
169  G4double cut =
170  std::max(recoilThreshold, (*pCuts)[CurrentCouple()->GetIndex()]);
171 
172  // find out lightest element
173  const G4ElementVector* theElementVector = material->GetElementVector();
174  G4int nelm = material->GetNumberOfElements();
175 
176  // select lightest element
177  G4int Z = 300;
178  for (G4int j=0; j<nelm; ++j) {
179  Z = std::min(Z,(*theElementVector)[j]->GetZasInt());
180  }
181  G4int A = G4lrint(fNistManager->GetAtomicMassAmu(Z));
182  G4double targetMass = G4NucleiProperties::GetNuclearMass(A, Z);
183  G4double t = std::max(cut, 0.5*(cut + sqrt(2*cut*targetMass)));
184 
185  return t;
186 }
187 
188 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
189 
191  const G4ParticleDefinition* p,
192  G4double kinEnergy,
194  G4double cutEnergy, G4double)
195 {
196  //G4cout << "### G4eCoulombScatteringModel::ComputeCrossSectionPerAtom for "
197  //<< p->GetParticleName()<<" Z= "<<Z<<" e(MeV)= "<< kinEnergy/MeV << G4endl;
198  G4double cross = 0.0;
199  elecRatio = 0.0;
200  if(p != particle) { SetupParticle(p); }
201 
202  // cross section is set to zero to avoid problems in sample secondary
203  if(kinEnergy <= 0.0) { return cross; }
205  G4double costmin = wokvi->SetupKinematic(kinEnergy, currentMaterial);
206  if(cosThetaMax < costmin) {
207  G4int iz = G4lrint(Z);
208  G4double cut = cutEnergy;
209  if(fixedCut > 0.0) { cut = fixedCut; }
210  costmin = wokvi->SetupTarget(iz, cut);
211  G4double costmax = cosThetaMax;
212  if(iz == 1 && costmax < 0.0 && particle == theProton) {
213  costmax = 0.0;
214  }
215  if(costmin > costmax) {
216  cross = wokvi->ComputeNuclearCrossSection(costmin, costmax)
217  + wokvi->ComputeElectronCrossSection(costmin, costmax);
218  }
219  /*
220  if(p->GetParticleName() == "mu+")
221  G4cout << "e(MeV)= " << kinEnergy/MeV << " cross(b)= " << cross/barn
222  << " 1-costmin= " << 1-costmin
223  << " 1-costmax= " << 1-costmax
224  << " 1-cosThetaMax= " << 1-cosThetaMax
225  << G4endl;
226  */
227  }
228  return cross;
229 }
230 
231 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
232 
234  std::vector<G4DynamicParticle*>* fvect,
235  const G4MaterialCutsCouple* couple,
236  const G4DynamicParticle* dp,
237  G4double cutEnergy,
238  G4double)
239 {
240  G4double kinEnergy = dp->GetKineticEnergy();
242  DefineMaterial(couple);
243  /*
244  G4cout << "G4eCoulombScatteringModel::SampleSecondaries e(MeV)= "
245  << kinEnergy << " " << particle->GetParticleName()
246  << " cut= " << cutEnergy<< G4endl;
247  */
248  // Choose nucleus
249  G4double cut = cutEnergy;
250  if(fixedCut > 0.0) { cut = fixedCut; }
251 
252  wokvi->SetupKinematic(kinEnergy, currentMaterial);
253 
254  const G4Element* currentElement =
255  SelectRandomAtom(couple,particle,kinEnergy,cut,kinEnergy);
256 
257  G4int iz = currentElement->GetZasInt();
258 
259  G4double costmin = wokvi->SetupTarget(iz, cut);
260  G4double costmax = cosThetaMax;
261  if(iz == 1 && costmax < 0.0 && particle == theProton) {
262  costmax = 0.0;
263  }
264 
265  if(costmin > costmax) {
266  G4double cross = wokvi->ComputeNuclearCrossSection(costmin, costmax);
267  G4double ecross = wokvi->ComputeElectronCrossSection(costmin, costmax);
268  G4double ratio = ecross/(cross + ecross);
269 
270  G4int ia = SelectIsotopeNumber(currentElement);
271  G4double targetMass = G4NucleiProperties::GetNuclearMass(ia, iz);
272  wokvi->SetTargetMass(targetMass);
273 
274  G4ThreeVector newDirection =
275  wokvi->SampleSingleScattering(costmin, costmax, ratio);
276  G4double cost = newDirection.z();
277  /*
278  G4cout << "SampleSec: e(MeV)= " << kinEnergy/MeV
279  << " 1-costmin= " << 1-costmin
280  << " 1-costmax= " << 1-costmax
281  << " 1-cost= " << 1-cost
282  << " ratio= " << ratio
283  << G4endl;
284  */
285  G4ThreeVector direction = dp->GetMomentumDirection();
286  newDirection.rotateUz(direction);
287 
288  fParticleChange->ProposeMomentumDirection(newDirection);
289 
290  // recoil sampling assuming a small recoil
291  // and first order correction to primary 4-momentum
292  G4double mom2 = wokvi->GetMomentumSquare();
293  G4double trec = mom2*(1.0 - cost)
294  /(targetMass + (mass + kinEnergy)*(1.0 - cost));
295 
296  // the check likely not needed
297  if(trec > kinEnergy) { trec = kinEnergy; }
298  G4double finalT = kinEnergy - trec;
299  G4double edep = 0.0;
300  /*
301  G4cout<<"G4eCoulombScatteringModel: finalT= "<<finalT<<" Trec= "
302  <<trec << " Z= " << iz << " A= " << ia
303  << " tcut(keV)= " << (*pCuts)[currentMaterialIndex]/keV << G4endl;
304  */
305  G4double tcut = recoilThreshold;
306  if(pCuts) { tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]); }
307 
308  if(trec > tcut) {
309  G4ParticleDefinition* ion = theIonTable->GetIon(iz, ia, 0);
310  G4ThreeVector dir = (direction*sqrt(mom2) -
311  newDirection*sqrt(finalT*(2*mass + finalT))).unit();
312  G4DynamicParticle* newdp = new G4DynamicParticle(ion, dir, trec);
313  fvect->push_back(newdp);
314  } else {
315  edep = trec;
316  fParticleChange->ProposeNonIonizingEnergyDeposit(edep);
317  }
318 
319  // finelize primary energy and energy balance
320  // this threshold may be applied only because for low-enegry
321  // e+e- msc model is applied
322  if(finalT < 0.0) {
323  edep += finalT;
324  finalT = 0.0;
325  if(edep < 0.0) { edep = 0.0; }
326  }
327  fParticleChange->SetProposedKineticEnergy(finalT);
328  fParticleChange->ProposeLocalEnergyDeposit(edep);
329  }
330 }
331 
332 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
333 
334 
static G4double GetNuclearMass(const G4double A, const G4double Z)
std::vector< G4Element * > G4ElementVector
G4double GetKineticEnergy() const
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:146
const char * p
Definition: xmltok.h:285
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:503
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4bool IsMaster() const
Definition: G4VEmModel.hh:725
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)
double z() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
string material
Definition: eplot.py:19
G4IonTable * GetIonTable() const
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:454
double A(double temperature)
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
static G4double tet[DIM]
bool G4bool
Definition: G4Types.hh:79
const G4ThreeVector & GetMomentumDirection() const
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double) final
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
static G4Proton * Proton()
Definition: G4Proton.cc:93
float proton_mass_c2
Definition: hepunit.py:275
void DefineMaterial(const G4MaterialCutsCouple *)
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.hh:587
G4eCoulombScatteringModel(G4bool combined=true)
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:809
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:817
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)
G4int GetZasInt() const
Definition: G4Element.hh:132
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:664
G4double GetAtomicMassAmu(const G4String &symb) const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
void SetProposedKineticEnergy(G4double proposedKinEnergy)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax) override
static constexpr double pi
Definition: G4SIunits.hh:75
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
double G4double
Definition: G4Types.hh:76
void SetupParticle(const G4ParticleDefinition *)
static constexpr double keV
Definition: G4SIunits.hh:216
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:544
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132
G4ThreeVector & SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio=0.0)