Geant4  10.00.p03
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 76536 2013-11-12 15:17: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(nam),
72  cosThetaMin(1.0),
73  cosThetaMax(-1.0),
74  isInitialised(false)
75 {
76  fParticleChange = 0;
80  currentMaterial = 0;
81 
82  pCuts = 0;
83 
84  lowEnergyThreshold = 1*keV; // particle will be killed for lower energy
85  recoilThreshold = 0.*keV; // by default does not work
86 
87  particle = 0;
88  currentCouple = 0;
90 
92 
93  cosTetMinNuc = 1.0;
94  cosTetMaxNuc = -1.0;
95  elecRatio = 0.0;
96  mass = proton_mass_c2;
97 }
98 
99 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
100 
102 {
103  delete wokvi;
104 }
105 
106 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
107 
109  const G4DataVector& cuts)
110 {
111  SetupParticle(p);
112  currentCouple = 0;
113  cosThetaMin = cos(PolarAngleLimit());
115  /*
116  G4cout << "G4hCoulombScatteringModel: " << particle->GetParticleName()
117  << " 1-cos(ThetaLimit)= " << 1 - cosThetaMin
118  << " cos(thetaMax)= " << cosThetaMax
119  << G4endl;
120  */
122  //G4cout << "!!! G4hCoulombScatteringModel::Initialise for "
123  // << p->GetParticleName() << " cos(TetMin)= " << cosThetaMin
124  // << " cos(TetMax)= " << cosThetaMax <<G4endl;
125  // G4cout << "cut0= " << cuts[0] << " cut1= " << cuts[1] << G4endl;
126  if(!isInitialised) {
127  isInitialised = true;
129  }
130  if(mass < GeV) {
132  }
133 }
134 
135 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
136 
138  const G4ParticleDefinition* p,
139  G4double kinEnergy,
140  G4double Z, G4double,
141  G4double cutEnergy, G4double)
142 {
143  //G4cout << "### G4hCoulombScatteringModel::ComputeCrossSectionPerAtom for "
144  // << p->GetParticleName()<<" Z= "<<Z<<" e(MeV)= "<< kinEnergy/MeV << G4endl;
145  G4double cross = 0.0;
146  if(p != particle) { SetupParticle(p); }
147 
148  // cross section is set to zero to avoid problems in sample secondary
149  if(kinEnergy <= 0.0) { return cross; }
152  if(cosThetaMax < cosTetMinNuc) {
153  G4int iz = G4int(Z);
154  cosTetMinNuc = wokvi->SetupTarget(iz, cutEnergy);
156  if(iz == 1 && cosTetMaxNuc < 0.0 && particle == theProton) {
157  cosTetMaxNuc = 0.0;
158  }
161  cross += elecRatio;
162  if(cross > 0.0) { elecRatio /= cross; }
163  }
164  /*
165  if(p->GetParticleName() == "e-")
166  G4cout << "e(MeV)= " << kinEnergy/MeV << " cross(b)= " << cross/barn
167  << " 1-cosTetMinNuc= " << 1-cosTetMinNuc
168  << " 1-cosTetMaxNuc= " << 1-cosTetMaxNuc
169  << G4endl;
170  */
171  return cross;
172 }
173 
174 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
175 
177  std::vector<G4DynamicParticle*>* fvect,
178  const G4MaterialCutsCouple* couple,
179  const G4DynamicParticle* dp,
180  G4double cutEnergy,
181  G4double)
182 {
183  G4double kinEnergy = dp->GetKineticEnergy();
184 
185  // absorb particle below low-energy limit to avoid situation
186  // when a particle has no energy loss
187  if(kinEnergy < lowEnergyThreshold) {
191  return;
192  }
194  DefineMaterial(couple);
195 
196  //G4cout << "G4hCoulombScatteringModel::SampleSecondaries e(MeV)= "
197  // << kinEnergy << " " << particle->GetParticleName()
198  // << " cut= " << cutEnergy<< G4endl;
199 
200  // Choose nucleus
201  const G4Element* currentElement =
202  SelectRandomAtom(couple,particle,kinEnergy,cutEnergy,kinEnergy);
203 
204  G4double Z = currentElement->GetZ();
205 
206  if(ComputeCrossSectionPerAtom(particle,kinEnergy, Z,
207  kinEnergy, cutEnergy, kinEnergy) == 0.0)
208  { return; }
209 
210  G4int iz = G4int(Z);
211  G4int ia = SelectIsotopeNumber(currentElement);
212  G4double targetMass = G4NucleiProperties::GetNuclearMass(ia, iz);
213  wokvi->SetTargetMass(targetMass);
214 
215  G4ThreeVector newDirection =
217  G4double cost = newDirection.z();
218 
219  G4ThreeVector direction = dp->GetMomentumDirection();
220  newDirection.rotateUz(direction);
221 
223 
224  // recoil sampling assuming a small recoil
225  // and first order correction to primary 4-momentum
226  G4double mom2 = wokvi->GetMomentumSquare();
227  G4double trec = mom2*(1.0 - cost)
228  /(targetMass + (mass + kinEnergy)*(1.0 - cost));
229 
230  // the check likely not needed
231  if(trec > kinEnergy) { trec = kinEnergy; }
232  G4double finalT = kinEnergy - trec;
233  G4double edep = 0.0;
234  //G4cout<<"G4eCoulombScatteringModel: finalT= "<<finalT<<" Trec= "
235  // <<trec << " Z= " << iz << " A= " << ia<<G4endl;
236 
237  G4double tcut = recoilThreshold;
238  if(pCuts) { tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]); }
239 
240  if(trec > tcut) {
242  G4ThreeVector dir = (direction*sqrt(mom2) -
243  newDirection*sqrt(finalT*(2*mass + finalT))).unit();
244  G4DynamicParticle* newdp = new G4DynamicParticle(ion, dir, trec);
245  fvect->push_back(newdp);
246  } else {
247  edep = trec;
249  }
250 
251  // finelize primary energy and energy balance
252  // this threshold may be applied only because for low-enegry
253  // e+e- msc model is applied
254  if(finalT <= lowEnergyThreshold) {
255  edep += finalT;
256  finalT = 0.0;
257  }
260 }
261 
262 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
263 
264 
const G4ParticleDefinition * theProton
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)
G4hCoulombScatteringModel(const G4String &nam="eCoulombScattering")
static G4double GetNuclearMass(const G4double A, const G4double Z)
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
G4ThreeVector SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio=0.0)
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:135
void DefineMaterial(const G4MaterialCutsCouple *)
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:462
G4double GetZ() const
Definition: G4Element.hh:131
const std::vector< G4double > * pCuts
G4ParticleDefinition * GetDefinition() const
void SetupParticle(const G4ParticleDefinition *)
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
G4double SetupTarget(G4int Z, G4double cut=DBL_MAX)
const G4MaterialCutsCouple * currentCouple
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double GetMomentumSquare() const
G4IonTable * GetIonTable() const
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:426
const G4ParticleDefinition * particle
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
const G4ThreeVector & GetMomentumDirection() const
G4double iz
Definition: TRTMaterials.hh:39
static G4Proton * Proton()
Definition: G4Proton.cc:93
static const double GeV
Definition: G4SIunits.hh:196
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.hh:548
G4ParticleChangeForGamma * fParticleChange
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
static G4ProductionCutsTable * GetProductionCutsTable()
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax)
static G4ParticleTable * GetParticleTable()
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:620
void SetProposedKineticEnergy(G4double proposedKinEnergy)
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
static const double keV
Definition: G4SIunits.hh:195
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
double G4double
Definition: G4Types.hh:76
void SetTargetMass(G4double value)
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