Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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$
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
51#include "G4SystemOfUnits.hh"
52#include "Randomize.hh"
53#include "G4DataVector.hh"
54#include "G4ElementTable.hh"
56#include "G4Proton.hh"
57#include "G4ParticleTable.hh"
59#include "G4NucleiProperties.hh"
60#include "G4Pow.hh"
61#include "G4LossTableManager.hh"
62#include "G4LossTableBuilder.hh"
63#include "G4NistManager.hh"
64
65//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
66
67using namespace std;
68
70 : G4VEmModel(nam),
71 cosThetaMin(1.0),
72 cosThetaMax(-1.0),
73 isInitialised(false)
74{
79 currentMaterial = 0;
80
81 pCuts = 0;
82
83 lowEnergyThreshold = 1*keV; // particle will be killed for lower energy
84 recoilThreshold = 0.*keV; // by default does not work
85
86 particle = 0;
87 currentCouple = 0;
89
91
92 cosTetMinNuc = 1.0;
93 cosTetMaxNuc = -1.0;
94 elecRatio = 0.0;
95 mass = proton_mass_c2;
96}
97
98//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
99
101{
102 delete wokvi;
103}
104
105//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
106
108 const G4DataVector& cuts)
109{
110 SetupParticle(p);
111 currentCouple = 0;
114 /*
115 G4cout << "G4hCoulombScatteringModel: " << particle->GetParticleName()
116 << " 1-cos(ThetaLimit)= " << 1 - cosThetaMin
117 << " cos(thetaMax)= " << cosThetaMax
118 << G4endl;
119 */
121 //G4cout << "!!! G4hCoulombScatteringModel::Initialise for "
122 // << p->GetParticleName() << " cos(TetMin)= " << cosThetaMin
123 // << " cos(TetMax)= " << cosThetaMax <<G4endl;
124 // G4cout << "cut0= " << cuts[0] << " cut1= " << cuts[1] << G4endl;
125 if(!isInitialised) {
126 isInitialised = true;
128 }
129 if(mass < GeV) {
131 }
132}
133
134//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
135
137 const G4ParticleDefinition* p,
138 G4double kinEnergy,
140 G4double cutEnergy, G4double)
141{
142 //G4cout << "### G4hCoulombScatteringModel::ComputeCrossSectionPerAtom for "
143 // << p->GetParticleName()<<" Z= "<<Z<<" e(MeV)= "<< kinEnergy/MeV << G4endl;
144 G4double cross = 0.0;
145 if(p != particle) { SetupParticle(p); }
146
147 // cross section is set to zero to avoid problems in sample secondary
148 if(kinEnergy <= 0.0) { return cross; }
152 G4int iz = G4int(Z);
153 cosTetMinNuc = wokvi->SetupTarget(iz, cutEnergy);
155 if(iz == 1 && cosTetMaxNuc < 0.0 && particle == theProton) {
156 cosTetMaxNuc = 0.0;
157 }
160 cross += elecRatio;
161 if(cross > 0.0) { elecRatio /= cross; }
162 }
163 /*
164 if(p->GetParticleName() == "e-")
165 G4cout << "e(MeV)= " << kinEnergy/MeV << " cross(b)= " << cross/barn
166 << " 1-cosTetMinNuc= " << 1-cosTetMinNuc
167 << " 1-cosTetMaxNuc= " << 1-cosTetMaxNuc
168 << G4endl;
169 */
170 return cross;
171}
172
173//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
174
176 std::vector<G4DynamicParticle*>* fvect,
177 const G4MaterialCutsCouple* couple,
178 const G4DynamicParticle* dp,
179 G4double cutEnergy,
180 G4double)
181{
182 G4double kinEnergy = dp->GetKineticEnergy();
183
184 // absorb particle below low-energy limit to avoid situation
185 // when a particle has no energy loss
186 if(kinEnergy < lowEnergyThreshold) {
190 return;
191 }
193 DefineMaterial(couple);
194
195 //G4cout << "G4hCoulombScatteringModel::SampleSecondaries e(MeV)= "
196 // << kinEnergy << " " << particle->GetParticleName()
197 // << " cut= " << cutEnergy<< G4endl;
198
199 // Choose nucleus
200 const G4Element* currentElement =
201 SelectRandomAtom(couple,particle,kinEnergy,cutEnergy,kinEnergy);
202
203 G4double Z = currentElement->GetZ();
204
205 if(ComputeCrossSectionPerAtom(particle,kinEnergy, Z,
206 kinEnergy, cutEnergy, kinEnergy) == 0.0)
207 { return; }
208
209 G4int iz = G4int(Z);
210 G4int ia = SelectIsotopeNumber(currentElement);
211 G4double targetMass = G4NucleiProperties::GetNuclearMass(ia, iz);
212 wokvi->SetTargetMass(targetMass);
213
214 G4ThreeVector newDirection =
216 G4double cost = newDirection.z();
217
218 G4ThreeVector direction = dp->GetMomentumDirection();
219 newDirection.rotateUz(direction);
220
222
223 // recoil sampling assuming a small recoil
224 // and first order correction to primary 4-momentum
226 G4double trec = mom2*(1.0 - cost)/(targetMass + (mass + kinEnergy)*(1.0 - cost));
227 G4double finalT = kinEnergy - trec;
228 //G4cout<<"G4hCoulombScatteringModel: finalT= "<<finalT<<" Trec= "<<trec<<G4endl;
229 if(finalT <= lowEnergyThreshold) {
230 trec = kinEnergy;
231 finalT = 0.0;
232 }
233
236 if(pCuts) { tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]); }
237
238 if(trec > tcut) {
239 G4ParticleDefinition* ion = theParticleTable->GetIon(iz, ia, 0.0);
240 G4ThreeVector dir = (direction*sqrt(mom2) -
241 newDirection*sqrt(finalT*(2*mass + finalT))).unit();
242 G4DynamicParticle* newdp = new G4DynamicParticle(ion, dir, trec);
243 fvect->push_back(newdp);
244 } else {
247 }
248
249 return;
250}
251
252//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
253
254
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
double z() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
static G4NistManager * Instance()
static G4double GetNuclearMass(const G4double A, const G4double Z)
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
static G4ParticleTable * GetParticleTable()
G4ParticleDefinition * GetIon(G4int atomicNumber, G4int atomicMass, G4double excitationEnergy)
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:550
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:459
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.hh:478
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:377
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:123
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetTargetMass(G4double value)
G4double SetupTarget(G4int Z, G4double cut=DBL_MAX)
G4ThreeVector SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio=0.0)
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
void SetupParticle(const G4ParticleDefinition *)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax)
const G4MaterialCutsCouple * currentCouple
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
const G4ParticleDefinition * theProton
void DefineMaterial(const G4MaterialCutsCouple *)
const std::vector< G4double > * pCuts
G4hCoulombScatteringModel(const G4String &nam="eCoulombScattering")
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
const G4ParticleDefinition * particle
G4ParticleChangeForGamma * fParticleChange