Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GoudsmitSaundersonMscModel.hh
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// File name: G4GoudsmitSaundersonMscModel
32//
33// Author: Mihaly Novak / (Omrane Kadri)
34//
35// Creation date: 20.02.2009
36//
37// Modifications:
38// 04.03.2009 V.Ivanchenko cleanup and format according to Geant4 EM style
39// 12.05.2010 O.Kadri: adding Qn1 and Qn12 as private doubles
40// 18.05.2015 M. Novak provide PLERIMINARYY version of updated class.
41// All algorithms of the class were revised and updated, new methods added.
42// A new version of Kawrakow-Bielajew Goudsmit-Saunderson MSC model
43// based on the screened Rutherford DCS for elastic scattering of
44// electrons/positrons has been introduced[1,2]. The corresponding MSC
45// angular distributions over a 2D parameter grid have been recomputed
46// and the CDFs are now stored in a variable transformed (smooth) form[2,3]
47// together with the corresponding rational interpolation parameters.
48// These angular distributions are handled by the new
49// G4GoudsmitSaundersonTable class that is responsible to sample if
50// it was no, single, few or multiple scattering case and delivers the
51// angular deflection (i.e. cos(theta) and sin(theta)).
52// Two screening options are provided:
53// - if fgIsUsePWATotalXsecData=TRUE i.e. SetOptionPWAScreening(TRUE)
54// was called before initialisation: screening parameter value A is
55// determined such that the first transport coefficient G1(A)
56// computed according to the screened Rutherford DCS for elastic
57// scattering will reproduce the one computed from the PWA elastic
58// and first transport mean free paths[4].
59// - if fgIsUsePWATotalXsecData=FALSE i.e. default value or
60// SetOptionPWAScreening(FALSE) was called before initialisation:
61// screening parameter value A is computed according to Moliere's
62// formula (by using material dependent parameters \chi_cc2 and b_c
63// precomputed for each material used at initialization in
64// G4GoudsmitSaundersonTable) [3]
65// Elastic and first trasport mean free paths are used consistently.
66// The new version is self-consistent, several times faster, more
67// robust and accurate compared to the earlier version.
68// Spin effects as well as a more accurate energy loss correction and
69// computations of Lewis moments will be implemented later on.
70// [1] A.F.Bielajew, NIMB 111 (1996) 195-208
71// [2] I.Kawrakow, A.F.Bielajew, NIMB 134(1998) 325-336
72// [3] I.Kawrakow, E.Mainegra-Hing, D.W.O.Rogers, F.Tessier,B.R.B.Walters,
73// NRCC Report PIRS-701 (2013)
74// [4] F.Salvat, A.Jablonski, C.J. Powell, CPC 165(2005) 157-190
75// 02.09.2015 M. Novak: first version of new step limit is provided.
76// fUseSafetyPlus corresponds to Urban fUseSafety (default)
77// fUseDistanceToBoundary corresponds to Urban fUseDistanceToBoundary
78// fUseSafety corresponds to EGSnrc error-free stepping algorithm
79// Range factor can be significantly higher at each case than in Urban.
80// 23.08.2017 M. Novak: added corrections to account spin effects (Mott-correction).
81// It can be activated by setting the fIsMottCorrection flag to be true
82// before initialization using the SetOptionMottCorrection() public method.
83// The fMottCorrection member is responsible to handle pre-computed Mott
84// correction (rejection) functions obtained by numerically computing
85// Goudsmit-Saunderson agnular distributions based on a DCS accounting spin
86// effects and screening corrections. The DCS used to compute the accurate
87// GS angular distributions is: DCS_{cor} = DCS_{SR}x[ DCS_{R}/DCS_{Mott}] where :
88// # DCS_{SR} is the relativistic Screened-Rutherford DCS (first Born approximate
89// solution of the Klein-Gordon i.e. relativistic Schrodinger equation =>
90// scattering of spinless e- on exponentially screened Coulomb potential)
91// note: the default (without using Mott-correction) GS angular distributions
92// are based on this DCS_{SR} with Moliere's screening parameter!
93// # DCS_{R} is the Rutherford DCS which is the same as above but without
94// screening
95// # DCS_{Mott} is the Mott DCS i.e. solution of the Dirac equation with a bare
96// Coulomb potential i.e. scattering of particles with spin (e- or e+) on a
97// point-like unscreened Coulomb potential
98// # moreover, the screening parameter of the DCS_{cor} was determined such that
99// the DCS_{cor} with this corrected screening parameter reproduce the first
100// transport cross sections obtained from the corresponding most accurate DCS
101// (i.e. from elsepa [4])
102// Unlike the default GS, the Mott-corrected angular distributions are particle type
103// (different for e- and e+ <= the DCS_{Mott} and the screening correction) and target
104// (Z and material) dependent.
105// 02.02.2018 M. Novak: implemented CrossSectionPerVolume interface method (used only for testing)
106// 26.10.2025 M. Novak: the model has only its accurate stepping and boundary crossing algorithms
107// left as the only option that ensures the expected precision, especially when activating
108// its Mott correction option (that also activates the screeing and scattering power
109// corrections). The model has been used for describing e-/e+ MSC (below 100 MeV kinetic)
110// energy in the option4, Penelope and Livermore EM physics constructors since Geant4 10.6.
111//
112//
113// Class description:
114// Kawrakow-Bielajew Goudsmit-Saunderson MSC model based on the screened Rutherford DCS
115// for elastic scattering of e-/e+. Option, to include Mott correction, is also available
116// that also activates the screening and scattering power corrections leading to the most
117// precise settings of the model. With the accurate electron stepping and boundary crossing
118// algorithm the model provides very precise e-/e+ simulation and tracking independently
119// from the target material and geometrical configurations similarly to EGSnrc. All details
120// are available in the corresponding technical note (M. Novak: https://arxiv.org/abs/2410.13361).
121//
122// References:
123// M. Novak: https://arxiv.org/abs/2410.13361
124//
125// -----------------------------------------------------------------------------
126
127#ifndef G4GoudsmitSaundersonMscModel_h
128#define G4GoudsmitSaundersonMscModel_h 1
129
131
132#include "G4VMscModel.hh"
133#include "G4PhysicsTable.hh"
135#include "globals.hh"
136
137
138class G4DataVector;
143
145{
146public:
147
148 G4GoudsmitSaundersonMscModel(const G4String& nam = "GoudsmitSaunderson");
149
151
152 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
153
154 void InitialiseLocal(const G4ParticleDefinition* p, G4VEmModel* masterModel) override;
155
156 G4ThreeVector& SampleScattering(const G4ThreeVector&, G4double safety) override;
157
158 G4double ComputeTruePathLengthLimit(const G4Track& track, G4double& currentMinimalStep) override;
159
160 G4double ComputeGeomPathLength(G4double truePathLength) override;
161
162 G4double ComputeTrueStepLength(G4double geomStepLength) override;
163
164 // method to compute first transport cross section per Volume (i.e. macroscropic first transport cross section; this
165 // method is used only for testing and not during a normal simulation)
166 G4double CrossSectionPerVolume(const G4Material*, const G4ParticleDefinition*, G4double kineticEnergy, G4double cutEnergy = 0.0, G4double maxEnergy = DBL_MAX) override;
167
168 void StartTracking(G4Track*) override;
169
170 void SampleMSC();
171
173
174 void SetOptionPWACorrection(G4bool opt) { fIsUsePWACorrection = opt; }
175
176 G4bool GetOptionPWACorrection() const { return fIsUsePWACorrection; }
177
178 void SetOptionMottCorrection(G4bool opt) { fIsUseMottCorrection = opt; }
179
180 G4bool GetOptionMottCorrection() const { return fIsUseMottCorrection; }
181
182 void SetOptionOptimisation(G4bool opt) { fIsUseOptimisation = opt; }
183
184 G4bool GetOptionOptimisation() const { return fIsUseOptimisation; }
185
186 G4GoudsmitSaundersonTable* GetGSTable() { return fGSTable; }
187
188 G4GSPWACorrections* GetPWACorrection() { return fPWACorrection; }
189
190 // hide assignment operator
193
194private:
195 inline void SetParticle(const G4ParticleDefinition* p);
196
197 inline G4double GetLambda(G4double);
198
199 G4double GetTransportMeanFreePathOnly(const G4ParticleDefinition*,G4double);
200
201private:
202
203 G4double currentKinEnergy;
204 G4double currentRange;
205 G4double presafety;
206 G4int currentMaterialIndex;
207 //
208 const G4ParticleDefinition* particle;
209 G4ParticleChangeForMSC* fParticleChange;
210 const G4MaterialCutsCouple* currentCouple;
211
213 G4GSPWACorrections* fPWACorrection;
214
215 G4bool fIsUsePWACorrection;
216 G4bool fIsUseMottCorrection;
217 G4bool fIsUseOptimisation;
218 //
219 G4double fLambda0; // elastic mean free path
220 G4double fLambda1; // first transport mean free path
221 G4double fScrA; // screening parameter
222 G4double fG1; // first transport coef.
223 // in case of Mott-correction
224 G4double fMCtoScrA;
225 G4double fMCtoQ1;
226 G4double fMCtoG2PerG1;
227 //
228 G4double fTheTrueStepLenght;
229 G4double fTheZPathLenght;
230 //
231 G4ThreeVector fTheDisplacementVector;
232 G4ThreeVector fTheNewDirection;
233 //
234 G4bool fIsEndedUpOnBoundary;
235 G4bool fIsMultipleScattering;
236 G4bool fIsSingleScattering;
237 G4bool fIsNoScatteringInMSC;
238 G4bool fIsSimplified;
239};
240
241////////////////////////////////////////////////////////////////////////////////
242inline
243void G4GoudsmitSaundersonMscModel::SetParticle(const G4ParticleDefinition* p)
244{
245 if (p != particle) {
246 particle = p;
247 }
248}
249
250#endif
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &currentMinimalStep) override
G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety) override
G4double ComputeGeomPathLength(G4double truePathLength) override
G4GoudsmitSaundersonMscModel(const G4String &nam="GoudsmitSaunderson")
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4GoudsmitSaundersonMscModel(const G4GoudsmitSaundersonMscModel &)=delete
G4GoudsmitSaundersonTable * GetGSTable()
G4double ComputeTrueStepLength(G4double geomStepLength) override
void InitialiseLocal(const G4ParticleDefinition *p, G4VEmModel *masterModel) override
G4GoudsmitSaundersonMscModel & operator=(const G4GoudsmitSaundersonMscModel &right)=delete
G4double GetTransportMeanFreePath(const G4ParticleDefinition *, G4double)
G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) override
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67
G4VMscModel(const G4String &nam)
#define DBL_MAX
Definition templates.hh:62