Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4hCoulombScatteringModel.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: G4hCoulombScatteringModel.cc 96934 2016-05-18 09:10:41Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4hCoulombScatteringModel
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 08.06.2012 from G4eCoulombScatteringModel
38 //
39 // Modifications:
40 //
41 //
42 // Class Description:
43 //
44 // -------------------------------------------------------------------
45 //
46 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48 
50 #include "G4PhysicalConstants.hh"
51 #include "G4SystemOfUnits.hh"
52 #include "Randomize.hh"
53 #include "G4DataVector.hh"
54 #include "G4ElementTable.hh"
56 #include "G4Proton.hh"
57 #include "G4ParticleTable.hh"
58 #include "G4IonTable.hh"
59 #include "G4ProductionCutsTable.hh"
60 #include "G4NucleiProperties.hh"
61 #include "G4Pow.hh"
62 #include "G4LossTableManager.hh"
63 #include "G4LossTableBuilder.hh"
64 #include "G4NistManager.hh"
65 
66 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
67 
68 using namespace std;
69 
71  : G4VEmModel("hCoulombScattering"),
72  cosThetaMin(1.0),
73  cosThetaMax(-1.0),
74  isCombined(combined)
75 {
76  fParticleChange = nullptr;
77  fNistManager = G4NistManager::Instance();
79  theProton = G4Proton::Proton();
80  currentMaterial = nullptr;
81  fixedCut = -1.0;
82 
83  pCuts = 0;
84 
85  recoilThreshold = 0.*CLHEP::keV; // by default does not work
86 
87  particle = nullptr;
88  currentCouple = nullptr;
89  wokvi = new G4WentzelVIRelXSection(combined);
90 
91  currentMaterialIndex = 0;
92  mass = CLHEP::proton_mass_c2;
93  elecRatio = 0.0;
94 }
95 
96 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
97 
99 {
100  delete wokvi;
101 }
102 
103 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
104 
106  const G4DataVector& cuts)
107 {
108  SetupParticle(part);
109  currentCouple = 0;
110 
111  if(isCombined) {
112  cosThetaMin = 1.0;
114  if(tet >= pi) { cosThetaMin = -1.0; }
115  else if(tet > 0.0) { cosThetaMin = cos(tet); }
116  }
117 
118  wokvi->Initialise(part, cosThetaMin);
119  /*
120  G4cout << "G4hCoulombScatteringModel: " << particle->GetParticleName()
121  << " 1-cos(ThetaLimit)= " << 1 - cosThetaMin
122  << " cos(thetaMax)= " << cosThetaMax
123  << G4endl;
124  */
125  pCuts = &cuts;
126  //G4ProductionCutsTable::GetProductionCutsTable()->GetEnergyCutsVector(3);
127  /*
128  G4cout << "!!! G4hCoulombScatteringModel::Initialise for "
129  << part->GetParticleName() << " cos(TetMin)= " << cosThetaMin
130  << " cos(TetMax)= " << cosThetaMax <<G4endl;
131  G4cout << "cut= " << (*pCuts)[0] << " cut1= " << (*pCuts)[1] << G4endl;
132  */
133  if(!fParticleChange) {
134  fParticleChange = GetParticleChangeForGamma();
135  }
136  if(IsMaster() && mass < GeV && part->GetParticleName() != "GenericIon") {
137  InitialiseElementSelectors(part, cuts);
138  }
139 }
140 
141 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
142 
144  G4VEmModel* masterModel)
145 {
147 }
148 
149 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
150 
151 G4double
153  const G4ParticleDefinition* part,
154  G4double)
155 {
156  SetupParticle(part);
157 
158  // define cut using cuts for proton
159  G4double cut =
160  std::max(recoilThreshold, (*pCuts)[CurrentCouple()->GetIndex()]);
161 
162  // find out lightest element
163  const G4ElementVector* theElementVector = material->GetElementVector();
164  G4int nelm = material->GetNumberOfElements();
165 
166  // select lightest element
167  G4int Z = 300;
168  for (G4int j=0; j<nelm; ++j) {
169  G4int iz = G4lrint((*theElementVector)[j]->GetZ());
170  if(iz < Z) { Z = iz; }
171  }
172  G4int A = G4lrint(fNistManager->GetAtomicMassAmu(Z));
173  G4double targetMass = G4NucleiProperties::GetNuclearMass(A, Z);
174  G4double t = std::max(cut, 0.5*(cut + sqrt(2*cut*targetMass)));
175 
176  return t;
177 }
178 
179 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
180 
182  const G4ParticleDefinition* p,
183  G4double kinEnergy,
185  G4double cutEnergy, G4double)
186 {
187  //G4cout << "### G4hCoulombScatteringModel::ComputeCrossSectionPerAtom for "
188  //<< p->GetParticleName()<<" Z= "<<Z<<" e(MeV)= "<< kinEnergy/MeV << G4endl;
189  G4double cross = 0.0;
190  elecRatio = 0.0;
191  if(p != particle) { SetupParticle(p); }
192 
193  // cross section is set to zero to avoid problems in sample secondary
194  if(kinEnergy <= 0.0) { return cross; }
196 
197  G4int iz = G4int(Z);
198  G4double tmass = proton_mass_c2;
199  if(1 < iz) {
200  tmass = fNistManager->GetAtomicMassAmu(iz)*amu_c2;
201  }
202  G4double costmin =
203  wokvi->SetupKinematic(kinEnergy, currentMaterial, cutEnergy, tmass);
204  if(cosThetaMax < costmin) {
205  G4double cut = cutEnergy;
206  if(fixedCut > 0.0) { cut = fixedCut; }
207  costmin = wokvi->SetupTarget(iz, cut);
208  G4double costmax = cosThetaMax;
209  if(iz == 1 && costmax < 0.0 && particle == theProton) {
210  costmax = 0.0;
211  }
212  if(costmin > costmax) {
213  cross = wokvi->ComputeNuclearCrossSection(costmin, costmax)
214  + wokvi->ComputeElectronCrossSection(costmin, costmax);
215  }
216  /*
217  if(p->GetParticleName() == "mu+")
218  G4cout << "e(MeV)= " << kinEnergy/MeV << " cross(b)= " << cross/barn
219  << " 1-costmin= " << 1-costmin
220  << " 1-costmax= " << 1-costmax
221  << " 1-cosThetaMax= " << 1-cosThetaMax
222  << G4endl;
223  */
224  }
225  return cross;
226 }
227 
228 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
229 
231  std::vector<G4DynamicParticle*>* fvect,
232  const G4MaterialCutsCouple* couple,
233  const G4DynamicParticle* dp,
234  G4double cutEnergy,
235  G4double)
236 {
237  G4double kinEnergy = dp->GetKineticEnergy();
239  DefineMaterial(couple);
240 
241  // Choose nucleus
242  const G4Element* elm = SelectRandomAtom(couple,particle,
243  kinEnergy,cutEnergy,kinEnergy);
244 
245  G4int iz = elm->GetZasInt();
246  G4int ia = SelectIsotopeNumber(elm);
248 
249  wokvi->SetupKinematic(kinEnergy, currentMaterial, cutEnergy, mass2);
250  G4double costmin = wokvi->SetupTarget(iz, cutEnergy);
251  G4double costmax = cosThetaMax;
252  if(iz == 1 && costmax < 0.0 && particle == theProton) {
253  costmax = 0.0;
254  }
255 
256  if(costmin <= costmax) { return; }
257  G4double cross = wokvi->ComputeNuclearCrossSection(costmin, costmax);
258  G4double ecross = wokvi->ComputeElectronCrossSection(costmin, costmax);
259  G4double ratio = ecross/(cross + ecross);
260 
261  G4ThreeVector newDirection =
262  wokvi->SampleSingleScattering(costmin, costmax, ratio);
263 
264  // kinematics in the Lab system
265  G4double ptot = dp->GetTotalMomentum();
266  G4double e1 = dp->GetTotalEnergy();
267 
268  // Lab. system kinematics along projectile direction
269  G4LorentzVector v0 = G4LorentzVector(0, 0, ptot, e1 + mass2);
270  G4double bet = ptot/v0.e();
271  G4double gam = 1.0/sqrt((1.0 - bet)*(1.0 + bet));
272 
273  // CM projectile
274  G4double momCM = gam*(ptot - bet*e1);
275  G4double eCM = gam*(e1 - bet*ptot);
276  // energy & momentum after scattering of incident particle
277  G4double pxCM = momCM*newDirection.x();
278  G4double pyCM = momCM*newDirection.y();
279  G4double pzCM = momCM*newDirection.z();
280 
281  // CM--->Lab
282  G4LorentzVector v1(pxCM , pyCM, gam*(pzCM + bet*eCM), gam*(eCM + bet*pzCM));
283 
284  G4ThreeVector dir = dp->GetMomentumDirection();
285  newDirection = v1.vect().unit();
286  newDirection.rotateUz(dir);
287 
288  fParticleChange->ProposeMomentumDirection(newDirection);
289 
290  // recoil
291  v0 -= v1;
292  G4double trec = v0.e() - mass2;
293  G4double edep = 0.0;
294 
295  G4double tcut = recoilThreshold;
296  if(pCuts) {
297  tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]);
298  //G4cout<<" tcut eV "<<tcut/eV<<endl;
299  }
300 
301  if(trec > tcut) {
302  G4ParticleDefinition* ion = theIonTable->GetIon(iz, ia, 0);
303  newDirection = v0.vect().unit();
304  newDirection.rotateUz(dir);
305  G4DynamicParticle* newdp = new G4DynamicParticle(ion, newDirection, trec);
306  fvect->push_back(newdp);
307  } else if(trec > 0.0) {
308  edep = trec;
309  fParticleChange->ProposeNonIonizingEnergyDeposit(edep);
310  }
311 
312  // finelize primary energy and energy balance
313  G4double finalT = v1.e() - mass;
314  if(finalT < 0.0) {
315  edep += finalT;
316  finalT = 0.0;
317  }
318  edep = std::max(edep, 0.0);
319  fParticleChange->SetProposedKineticEnergy(finalT);
320  fParticleChange->ProposeLocalEnergyDeposit(edep);
321 }
322 
323 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)
static G4double GetNuclearMass(const G4double A, const G4double Z)
std::vector< G4Element * > G4ElementVector
static constexpr double proton_mass_c2
G4double GetKineticEnergy() const
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:146
G4double GetTotalEnergy() const
double x() const
static constexpr double keV
const char * p
Definition: xmltok.h:285
void DefineMaterial(const G4MaterialCutsCouple *)
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:503
G4bool IsMaster() const
Definition: G4VEmModel.hh:725
G4ParticleDefinition * GetDefinition() const
void SetupParticle(const G4ParticleDefinition *)
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)
double z() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax) override
G4double GetTotalMomentum() const
string material
Definition: eplot.py:19
G4IonTable * GetIonTable() const
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:454
Hep3Vector vect() const
double A(double temperature)
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
static G4double tet[DIM]
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat, G4double cut, G4double tmass)
G4double SetupTarget(G4int Z, G4double cut)
bool G4bool
Definition: G4Types.hh:79
const G4ThreeVector & GetMomentumDirection() const
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
static G4Proton * Proton()
Definition: G4Proton.cc:93
float proton_mass_c2
Definition: hepunit.py:275
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.hh:587
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:809
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:817
G4hCoulombScatteringModel(G4bool combined=true)
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
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double) final
Hep3Vector unit() const
G4int GetZasInt() const
Definition: G4Element.hh:132
double y() const
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:664
G4double GetAtomicMassAmu(const G4String &symb) const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
static constexpr double pi
Definition: G4SIunits.hh:75
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)
double G4double
Definition: G4Types.hh:76
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4ThreeVector & SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio)
float amu_c2
Definition: hepunit.py:277
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:132
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
CLHEP::HepLorentzVector G4LorentzVector