Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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// -------------------------------------------------------------------
27//
28// GEANT4 Class file
29//
30//
31// File name: G4hCoulombScatteringModel
32//
33// Author: Vladimir Ivanchenko
34//
35// Creation date: 08.06.2012 from G4eCoulombScatteringModel
36//
37// Modifications:
38//
39//
40// Class Description:
41//
42// -------------------------------------------------------------------
43//
44//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
45//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
46
49#include "G4SystemOfUnits.hh"
50#include "Randomize.hh"
51#include "G4DataVector.hh"
52#include "G4ElementTable.hh"
54#include "G4Proton.hh"
55#include "G4ParticleTable.hh"
56#include "G4IonTable.hh"
58#include "G4NucleiProperties.hh"
59#include "G4Pow.hh"
60#include "G4NistManager.hh"
61
62//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
63
65 : G4VEmModel("hCoulombScattering"),
66 cosThetaMin(1.0),
67 cosThetaMax(-1.0),
68 isCombined(combined)
69{
70 fParticleChange = nullptr;
71 fNistManager = G4NistManager::Instance();
73 theProton = G4Proton::Proton();
74 currentMaterial = nullptr;
75 fixedCut = -1.0;
76
77 pCuts = nullptr;
78
79 recoilThreshold = 0.0; // by default does not work
80
81 particle = nullptr;
82 currentCouple = nullptr;
83 wokvi = new G4WentzelVIRelXSection();
84
85 currentMaterialIndex = 0;
86 mass = CLHEP::proton_mass_c2;
87 elecRatio = 0.0;
88}
89
90//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
91
96
97//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
98
100 const G4DataVector& cuts)
101{
102 SetupParticle(part);
103 currentCouple = nullptr;
104
105 // defined theta limit between single and multiple scattering
106 isCombined = true;
108
109 if(tet <= 0.0) {
110 cosThetaMin = 1.0;
111 isCombined = false;
112 } else if(tet >= CLHEP::pi) {
113 cosThetaMin = -1.0;
114 } else {
115 cosThetaMin = std::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 < CLHEP::GeV && part->GetParticleName() != "GenericIon") {
137 InitialiseElementSelectors(part, cuts);
138 }
139}
140
141//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
142
148
149//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
150
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 std::size_t nelm = material->GetNumberOfElements();
165
166 // select lightest element
167 G4int Z = 300;
168 for (std::size_t j=0; j<nelm; ++j) {
169 Z = std::min(Z,(*theElementVector)[j]->GetZasInt());
170 }
171 G4int A = G4lrint(fNistManager->GetAtomicMassAmu(Z));
173 G4double t = std::max(cut, 0.5*(cut + std::sqrt(2*cut*targetMass)));
174
175 return t;
176}
177
178//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
179
181 const G4ParticleDefinition* p,
182 G4double kinEnergy,
184 G4double cutEnergy, G4double)
185{
186 //G4cout << "### G4hCoulombScatteringModel::ComputeCrossSectionPerAtom for "
187 //<< p->GetParticleName()<<" Z= "<<Z<<" e(MeV)= "<< kinEnergy/MeV << G4endl;
188 G4double cross = 0.0;
189 elecRatio = 0.0;
190 if(p != particle) { SetupParticle(p); }
191
192 // cross section is set to zero to avoid problems in sample secondary
193 if(kinEnergy <= 0.0) { return cross; }
195
196 G4int iz = G4lrint(Z);
197 G4double tmass = (1 == iz) ? proton_mass_c2 :
198 fNistManager->GetAtomicMassAmu(iz)*amu_c2;
199 wokvi->SetTargetMass(tmass);
200
201 G4double costmin =
202 wokvi->SetupKinematic(kinEnergy, currentMaterial);
203
204 if(cosThetaMax < costmin) {
205 G4double cut = (0.0 < fixedCut) ? fixedCut : cutEnergy;
206 costmin = wokvi->SetupTarget(iz, cut);
207 G4double costmax =
208 (1 == iz && particle == theProton && cosThetaMax < 0.0)
209 ? 0.0 : cosThetaMax;
210 if(costmin > costmax) {
211 cross = wokvi->ComputeNuclearCrossSection(costmin, costmax)
212 + wokvi->ComputeElectronCrossSection(costmin, costmax);
213 }
214 /*
215 if(p->GetParticleName() == "mu+")
216 G4cout << "e(MeV)= " << kinEnergy/MeV << " cross(b)= " << cross/barn
217 << " 1-costmin= " << 1-costmin
218 << " 1-costmax= " << 1-costmax
219 << " 1-cosThetaMax= " << 1-cosThetaMax
220 << G4endl;
221 */
222 }
223 return cross;
224}
225
226//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
227
229 std::vector<G4DynamicParticle*>* fvect,
230 const G4MaterialCutsCouple* couple,
231 const G4DynamicParticle* dp,
232 G4double cutEnergy,
233 G4double)
234{
235 G4double kinEnergy = dp->GetKineticEnergy();
237 DefineMaterial(couple);
238
239 // Choose nucleus
240 G4double cut = (0.0 < fixedCut) ? fixedCut : cutEnergy;
241
242 const G4Element* elm = SelectRandomAtom(couple,particle,
243 kinEnergy,cut,kinEnergy);
244
245 G4int iz = elm->GetZasInt();
246 G4int ia = SelectIsotopeNumber(elm);
248
249 wokvi->SetTargetMass(mass2);
250 wokvi->SetupKinematic(kinEnergy, currentMaterial);
251 G4double costmin = wokvi->SetupTarget(iz, cut);
252 G4double costmax = (1 == iz && particle == theProton && cosThetaMax < 0.0)
253 ? 0.0 : cosThetaMax;
254 if(costmin <= costmax) { return; }
255
256 G4double cross = wokvi->ComputeNuclearCrossSection(costmin, costmax);
257 G4double ecross = wokvi->ComputeElectronCrossSection(costmin, costmax);
258 G4double ratio = ecross/(cross + ecross);
259
260 G4ThreeVector newDirection =
261 wokvi->SampleSingleScattering(costmin, costmax, ratio);
262
263 // kinematics in the Lab system
264 G4double ptot = std::sqrt(kinEnergy*(kinEnergy + 2.0*mass));
265 G4double e1 = mass + kinEnergy;
266
267 // Lab. system kinematics along projectile direction
268 G4LorentzVector v0 = G4LorentzVector(0, 0, ptot, e1+mass2);
269 G4LorentzVector v1 = G4LorentzVector(0, 0, ptot, e1);
270 G4ThreeVector bst = v0.boostVector();
271 v1.boost(-bst);
272 // CM projectile
273 G4double momCM = v1.pz();
274
275 // Momentum after scattering of incident particle
276 v1.setX(momCM*newDirection.x());
277 v1.setY(momCM*newDirection.y());
278 v1.setZ(momCM*newDirection.z());
279
280 // CM--->Lab
281 v1.boost(bst);
282
284 newDirection = v1.vect().unit();
285 newDirection.rotateUz(dir);
286
287 fParticleChange->ProposeMomentumDirection(newDirection);
288
289 // recoil
290 v0 -= v1;
291 G4double trec = std::max(v0.e() - mass2, 0.0);
292 G4double edep = 0.0;
293
294 G4double tcut = recoilThreshold;
295 if(pCuts) { tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]); }
296
297 if(trec > tcut) {
298 G4ParticleDefinition* ion = theIonTable->GetIon(iz, ia, 0);
299 newDirection = v0.vect().unit();
300 newDirection.rotateUz(dir);
301 auto newdp = new G4DynamicParticle(ion, newDirection, trec);
302 fvect->push_back(newdp);
303 } else if(trec > 0.0) {
304 edep = trec;
305 fParticleChange->ProposeNonIonizingEnergyDeposit(edep);
306 }
307
308 // finelize primary energy and energy balance
309 G4double finalT = v1.e() - mass;
310 if(finalT < 0.0) {
311 edep += finalT;
312 finalT = 0.0;
313 }
314 edep = std::max(edep, 0.0);
315 fParticleChange->SetProposedKineticEnergy(finalT);
316 fParticleChange->ProposeLocalEnergyDeposit(edep);
317}
318
319//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
std::vector< const G4Element * > G4ElementVector
CLHEP::HepLorentzVector G4LorentzVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
double z() const
Hep3Vector unit() const
double x() const
double y() const
Hep3Vector & rotateUz(const Hep3Vector &)
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
Hep3Vector vect() const
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4int GetZasInt() const
Definition G4Element.hh:120
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
const G4ElementVector * GetElementVector() const
std::size_t GetNumberOfElements() const
static G4NistManager * Instance()
G4double GetAtomicMassAmu(const G4String &symb) const
static G4double GetNuclearMass(const G4double A, const G4double Z)
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
static G4Proton * Proton()
Definition G4Proton.cc:90
G4int SelectIsotopeNumber(const G4Element *) const
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
G4double PolarAngleLimit() const
G4ParticleChangeForGamma * GetParticleChangeForGamma()
std::vector< G4EmElementSelector * > * GetElementSelectors()
G4bool IsMaster() const
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
const G4MaterialCutsCouple * CurrentCouple() const
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double SetupTarget(G4int Z, G4double cut)
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)
G4ThreeVector & SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio)
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat) override
void SetupParticle(const G4ParticleDefinition *)
G4hCoulombScatteringModel(G4bool combined=true)
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax) override
void DefineMaterial(const G4MaterialCutsCouple *)
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double) final
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
int G4lrint(double ad)
Definition templates.hh:134