Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class file
30//
31//
32// File name: G4eCoulombScatteringModel
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 22.08.2005
37//
38// Modifications: V.Ivanchenko
39//
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"
58#include "G4IonTable.hh"
60#include "G4NucleiProperties.hh"
61#include "G4Pow.hh"
62#include "G4NistManager.hh"
63
64//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
65
66using namespace std;
67
69 : G4VEmModel("eCoulombScattering"), isCombined(combined)
70{
71 fNistManager = G4NistManager::Instance();
73 theProton = G4Proton::Proton();
74
75 wokvi = new G4WentzelOKandVIxSection(isCombined);
76
77 mass = CLHEP::proton_mass_c2;
78}
79
80//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
81
86
87//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
88
90 const G4DataVector& cuts)
91{
92 SetupParticle(part);
93 currentCouple = nullptr;
94
96
97 // defined theta limit between single and multiple scattering
98 if(isCombined) {
99 if(tet >= CLHEP::pi) { cosThetaMin = -1.0; }
100 else if(tet > 0.0) { cosThetaMin = std::cos(tet); }
101
102 // single scattering without multiple
103 } else if(tet > 0.0) {
104 cosThetaMin = std::cos(std::min(tet, CLHEP::pi));
105 }
106
107 wokvi->Initialise(part, cosThetaMin);
108 pCuts = &cuts;
109 /*
110 G4cout << "G4eCoulombScatteringModel::Initialise for "
111 << part->GetParticleName() << " 1-cos(TetMin)= " << 1.0 - cosThetaMin
112 << " 1-cos(TetMax)= " << 1. - cosThetaMax << G4endl;
113 G4cout << "cut[0]= " << (*pCuts)[0] << G4endl;
114 */
115 if(nullptr == fParticleChange) {
116 fParticleChange = GetParticleChangeForGamma();
117 }
118 if(IsMaster() && mass < GeV && part->GetParticleName() != "GenericIon") {
119 InitialiseElementSelectors(part, cuts);
120 }
121}
122
123//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
124
130
131//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
132
135 const G4ParticleDefinition* part,
136 G4double)
137{
138 SetupParticle(part);
139
140 // define cut using cuts for proton
141 G4double cut =
142 std::max(recoilThreshold, (*pCuts)[CurrentCouple()->GetIndex()]);
143
144 // find out lightest element
145 const G4ElementVector* theElementVector = material->GetElementVector();
146 std::size_t nelm = material->GetNumberOfElements();
147
148 // select lightest element
149 G4int Z = 300;
150 for (std::size_t j=0; j<nelm; ++j) {
151 Z = std::min(Z,(*theElementVector)[j]->GetZasInt());
152 }
153 G4int A = G4lrint(fNistManager->GetAtomicMassAmu(Z));
155 G4double t = std::max(cut, 0.5*(cut + sqrt(2*cut*targetMass)));
156
157 return t;
158}
159
160//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
161
163 const G4ParticleDefinition* p,
164 G4double kinEnergy,
166 G4double cutEnergy, G4double)
167{
168 /*
169 G4cout << "### G4eCoulombScatteringModel::ComputeCrossSectionPerAtom for "
170 << p->GetParticleName()<<" Z= "<<Z<<" e(MeV)= "<< kinEnergy/MeV
171 << G4endl;
172 */
173 G4double cross = 0.0;
174 elecRatio = 0.0;
175 if(p != particle) { SetupParticle(p); }
176
177 // cross section is set to zero to avoid problems in sample secondary
178 if(kinEnergy <= 0.0) { return cross; }
180 G4double costmin = wokvi->SetupKinematic(kinEnergy, currentMaterial);
181
182 //G4cout << "cosThetaMax= "<<cosThetaMax<<" costmin= "<<costmin<< G4endl;
183
184 if(cosThetaMax < costmin) {
185 G4int iz = G4lrint(Z);
186 G4double cut = (0.0 < fixedCut) ? fixedCut : cutEnergy;
187 costmin = wokvi->SetupTarget(iz, cut);
188 //G4cout << "SetupTarget: Z= " << iz << " cut= " << cut << " "
189 // << costmin << G4endl;
190 G4double costmax = (1 == iz && particle == theProton && cosThetaMax < 0.0)
191 ? 0.0 : cosThetaMax;
192 if(costmin > costmax) {
193 cross = wokvi->ComputeNuclearCrossSection(costmin, costmax)
194 + wokvi->ComputeElectronCrossSection(costmin, costmax);
195 }
196 /*
197 if(p->GetParticleName() == "e-")
198 G4cout << "Z= " << Z << " e(MeV)= " << kinEnergy/MeV
199 << " cross(b)= " << cross/barn << " 1-costmin= " << 1-costmin
200 << " 1-costmax= " << 1-costmax
201 << " 1-cosThetaMax= " << 1-cosThetaMax
202 << " " << currentMaterial->GetName()
203 << G4endl;
204 */
205 }
206 //G4cout << "====== cross= " << cross << G4endl;
207 return cross;
208}
209
210//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
211
213 std::vector<G4DynamicParticle*>* fvect,
214 const G4MaterialCutsCouple* couple,
215 const G4DynamicParticle* dp,
216 G4double cutEnergy,
217 G4double)
218{
219 G4double kinEnergy = dp->GetKineticEnergy();
221 DefineMaterial(couple);
222 /*
223 G4cout << "G4eCoulombScatteringModel::SampleSecondaries e(MeV)= "
224 << kinEnergy << " " << particle->GetParticleName()
225 << " cut= " << cutEnergy<< G4endl;
226 */
227 // Choose nucleus
228 G4double cut = (0.0 < fixedCut) ? fixedCut : cutEnergy;
229
230 wokvi->SetupKinematic(kinEnergy, currentMaterial);
231
232 const G4Element* currentElement = SelectTargetAtom(couple,particle,kinEnergy,
233 dp->GetLogKineticEnergy(),cut,kinEnergy);
234 G4int iz = currentElement->GetZasInt();
235
236 G4double costmin = wokvi->SetupTarget(iz, cut);
237 G4double costmax = (1 == iz && particle == theProton && cosThetaMax < 0.0)
238 ? 0.0 : cosThetaMax;
239 if(costmin <= costmax) { return; }
240
241 G4double cross = wokvi->ComputeNuclearCrossSection(costmin, costmax);
242 G4double ecross = wokvi->ComputeElectronCrossSection(costmin, costmax);
243 G4double ratio = ecross/(cross + ecross);
244
245 G4int ia = SelectIsotopeNumber(currentElement);
246 G4double targetMass = G4NucleiProperties::GetNuclearMass(ia, iz);
247 wokvi->SetTargetMass(targetMass);
248
249 G4ThreeVector newDirection =
250 wokvi->SampleSingleScattering(costmin, costmax, ratio);
251 G4double cost = newDirection.z();
252 /*
253 G4cout << "SampleSec: e(MeV)= " << kinEnergy/MeV
254 << " 1-costmin= " << 1-costmin
255 << " 1-costmax= " << 1-costmax
256 << " 1-cost= " << 1-cost
257 << " ratio= " << ratio
258 << G4endl;
259 */
260 G4ThreeVector direction = dp->GetMomentumDirection();
261 newDirection.rotateUz(direction);
262
263 fParticleChange->ProposeMomentumDirection(newDirection);
264
265 // recoil sampling assuming a small recoil
266 // and first order correction to primary 4-momentum
267 G4double mom2 = wokvi->GetMomentumSquare();
268 G4double trec = mom2*(1.0 - cost)
269 /(targetMass + (mass + kinEnergy)*(1.0 - cost));
270
271 // the check likely not needed
272 trec = std::min(trec, kinEnergy);
273 G4double finalT = kinEnergy - trec;
274 G4double edep = 0.0;
275 /*
276 G4cout<<"G4eCoulombScatteringModel: finalT= "<<finalT<<" Trec= "
277 <<trec << " Z= " << iz << " A= " << ia
278 << " tcut(keV)= " << (*pCuts)[currentMaterialIndex]/keV << G4endl;
279 */
280 G4double tcut = recoilThreshold;
281 if(pCuts) { tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]); }
282
283 if(trec > tcut) {
284 G4ParticleDefinition* ion = theIonTable->GetIon(iz, ia, 0);
285 G4ThreeVector dir = (direction*sqrt(mom2) -
286 newDirection*sqrt(finalT*(2*mass + finalT))).unit();
287 auto newdp = new G4DynamicParticle(ion, dir, trec);
288 fvect->push_back(newdp);
289 } else {
290 edep = trec;
291 fParticleChange->ProposeNonIonizingEnergyDeposit(edep);
292 }
293
294 // finelize primary energy and energy balance
295 // this threshold may be applied only because for low-enegry
296 // e+e- msc model is applied
297 if(finalT < 0.0) {
298 edep += finalT;
299 finalT = 0.0;
300 }
301 edep = std::max(edep, 0.0);
302 fParticleChange->SetProposedKineticEnergy(finalT);
303 fParticleChange->ProposeLocalEnergyDeposit(edep);
304}
305
306//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
std::vector< const G4Element * > G4ElementVector
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 & rotateUz(const Hep3Vector &)
const G4ThreeVector & GetMomentumDirection() const
G4double GetLogKineticEnergy() 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 G4MaterialCutsCouple * CurrentCouple() const
const G4Element * SelectTargetAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
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)
virtual G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)
G4ThreeVector & SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio)
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax) override
G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double) final
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
G4eCoulombScatteringModel(G4bool combined=true)
void DefineMaterial(const G4MaterialCutsCouple *)
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void SetupParticle(const G4ParticleDefinition *)
int G4lrint(double ad)
Definition templates.hh:134