Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4eDPWACoulombScatteringModel.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 header file
30//
31//
32// File name: G4eDPWACoulombScatteringModel
33//
34// Author: Mihaly Novak
35//
36// Creation date: 02.07.2020
37//
38// Modifications:
39//
40// -------------------------------------------------------------------
41
43
44#include "G4eDPWAElasticDCS.hh"
47#include "G4DataVector.hh"
48
50#include "G4Material.hh"
51#include "G4Element.hh"
52#include "G4ElementVector.hh"
53
54#include "G4Electron.hh"
55
57#include "G4SystemOfUnits.hh"
58#include "Randomize.hh"
59#include "G4ThreeVector.hh"
60
61
63: G4VEmModel("eDPWACoulombScattering"),
64 fIsMixedModel(ismixed),
65 fIsScpCorrection(isscpcor),
66 fMuMin(mumin),
67 fTheDCS(nullptr),
68 fParticleChange(nullptr)
69{
70 SetLowEnergyLimit ( 0.0*CLHEP::eV); // ekin = 10 eV is used if (E< 10 eV)
71 SetHighEnergyLimit(100.0*CLHEP::MeV); // ekin = 100 MeV is used if (E>100 MeV)
72}
73
74
76{
77 if (IsMaster()) {
78 delete fTheDCS;
79 }
80}
81
82
84 const G4DataVector& prodcuts)
85{
86 if(!fParticleChange) {
87 fParticleChange = GetParticleChangeForGamma();
88 }
89 fMuMin = 0.5*(1.0-std::cos(PolarAngleLimit()));
90 fIsMixedModel = (fMuMin > 0.0);
91 if(IsMaster()) {
92 // clean the G4eDPWAElasticDCS object if any
93 delete fTheDCS;
94 fTheDCS = new G4eDPWAElasticDCS(pdef==G4Electron::Electron(), fIsMixedModel);
95 // init only for the elements that are used in the geometry
97 G4int numOfCouples = (G4int)theCpTable->GetTableSize();
98 for(G4int j=0; j<numOfCouples; ++j) {
99 const G4Material* mat = theCpTable->GetMaterialCutsCouple(j)->GetMaterial();
100 const G4ElementVector* elV = mat->GetElementVector();
101 std::size_t numOfElem = mat->GetNumberOfElements();
102 for (std::size_t ie = 0; ie < numOfElem; ++ie) {
103 fTheDCS->InitialiseForZ((*elV)[ie]->GetZasInt());
104 }
105 }
106 // init scattering power correction
107 if (fIsScpCorrection) {
109 }
110 // will make use of the cross sections so the above needs to be done before
111 InitialiseElementSelectors(pdef, prodcuts);
112 }
113}
114
115
117 G4VEmModel* masterModel)
118{
120 SetTheDCS(static_cast<G4eDPWACoulombScatteringModel*>(masterModel)->GetTheDCS());
121}
122
123
126 G4double ekin,
127 G4double Z,
128 G4double /*A*/,
129 G4double /*prodcut*/,
130 G4double /*emax*/)
131{
132 // Cross sections are computed by numerical integration of the pre-computed
133 // DCS data between the muMin, muMax limits where mu(theta)=0.5[1-cos(theta)].
134 // In case of single scattering model (i.e. when fMuMin=0): [muMin=0, muMax=1]
135 // In case of mixed simulation model (i.e. when fMuMin>0): [fMuMin , muMax=1]
136 // NOTE: cross sections will be zero if the kinetic enrgy is out of the
137 // [10 eV-100 MeV] range for which DCS data has been computed.
138 //
139 G4double elCS = 0.0; // elastic cross section
140 G4double tr1CS = 0.0; // first transport cross section
141 G4double tr2CS = 0.0; // second transport cross section
142 const G4double muMin = fMuMin;
143 const G4double muMax = 1.0;
144 fTheDCS->ComputeCSPerAtom((G4int)Z, ekin, elCS, tr1CS, tr2CS, muMin, muMax);
145 // scattering power correction: should be only in condensed history ioni!
146 if (fIsScpCorrection && CurrentCouple()) {
147 const G4double theScpCor = fTheDCS->ComputeScatteringPowerCorrection(CurrentCouple(), ekin);
148 elCS *= (theScpCor*(1.0+1.0/Z));
149 }
150 return std::max(0.0, elCS);
151}
152
153
154void
156 const G4MaterialCutsCouple* cp,
157 const G4DynamicParticle* dp,
159{
160 const G4double ekin = dp->GetKineticEnergy();
161 const G4double lekin = dp->GetLogKineticEnergy();
162 const G4Element* target = SelectTargetAtom(cp, dp->GetParticleDefinition(), ekin, lekin);
163 const G4int izet = target->GetZasInt();
164 // sample cosine of the polar scattering angle in (hard) elastic insteraction
165 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
166 G4double cost = 1.0;
167 if (!fIsMixedModel) {
168 G4double rndm[3];
169 rndmEngine->flatArray(3, rndm);
170 cost = fTheDCS->SampleCosineTheta(izet, lekin, rndm[0], rndm[1], rndm[2]);
171 } else {
172 //sample cost between costMax,costMin where costMax = 1-2xfMuMin;
173 const G4double costMax = 1.0-2.0*fMuMin;
174 const G4double costMin = -1.0;
175 G4double rndm[2];
176 rndmEngine->flatArray(2, rndm);
177 cost = fTheDCS->SampleCosineThetaRestricted(izet, lekin, rndm[0], rndm[1], costMin, costMax);
178 }
179 // compute the new direction in the scattering frame
180 const G4double sint = std::sqrt((1.0-cost)*(1.0+cost));
181 const G4double phi = CLHEP::twopi*rndmEngine->flat();
182 G4ThreeVector theNewDirection(sint*std::cos(phi), sint*std::sin(phi), cost);
183 // get original direction in lab frame and rotate new direction to lab frame
184 G4ThreeVector theOrgDirectionLab = dp->GetMomentumDirection();
185 theNewDirection.rotateUz(theOrgDirectionLab);
186 // set new direction
187 fParticleChange->ProposeMomentumDirection(theNewDirection);
188}
189
190
191
192
193
std::vector< const G4Element * > G4ElementVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
virtual double flat()=0
virtual void flatArray(const int size, double *vect)=0
const G4ThreeVector & GetMomentumDirection() const
G4double GetLogKineticEnergy() const
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:93
G4int GetZasInt() const
Definition: G4Element.hh:132
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:185
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:831
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:662
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:124
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:823
G4bool IsMaster() const
Definition: G4VEmModel.hh:725
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:486
const G4Element * SelectTargetAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:577
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:139
G4eDPWACoulombScatteringModel(G4bool ismixed=false, G4bool isscpcor=true, G4double mumin=0.0)
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *) override
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double ekin, G4double Z, G4double A, G4double prodcut, G4double emax) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
void SetTheDCS(G4eDPWAElasticDCS *theDCS)
void ComputeCSPerAtom(G4int iz, G4double ekin, G4double &elcs, G4double &tr1cs, G4double &tr2cs, G4double mumin=0.0, G4double mumax=1.0)
void InitialiseForZ(std::size_t iz)
void InitSCPCorrection(G4double lowEnergyLimit, G4double highEnergyLimit)
G4double ComputeScatteringPowerCorrection(const G4MaterialCutsCouple *matcut, G4double ekin)
G4double SampleCosineTheta(std::size_t iz, G4double lekin, G4double r1, G4double r2, G4double r3)
G4double SampleCosineThetaRestricted(std::size_t iz, G4double lekin, G4double r1, G4double r2, G4double costMax, G4double costMin)