Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GoudsmitSaundersonTable.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: G4GoudsmitSaundersonTable
32//
33// Author: Mihaly Novak / (Omrane Kadri)
34//
35// Creation date: 20.02.2009
36//
37// Class description:
38// Class to handle multiple scattering angular distributions precomputed by
39// using Kawrakow-Bielajew Goudsmit-Saunderson MSC model based on the screened
40// Rutherford DCS for elastic scattering of electrons/positrons [1,2]. This
41// class is used by G4GoudsmitSaundersonMscModel to sample the angular
42// deflection of electrons/positrons after travelling a given path.
43//
44// Modifications:
45// 04.03.2009 V.Ivanchenko cleanup and format according to Geant4 EM style
46// 18.05.2015 M. Novak This class has been completely replaced (only the original
47// class name was kept; class description was also inserted):
48// A new version of Kawrakow-Bielajew Goudsmit-Saunderson MSC model
49// based on the screened Rutherford DCS for elastic scattering of
50// electrons/positrons has been introduced[1,2]. The corresponding MSC
51// angular distributions over a 2D parameter grid have been recomputed
52// and the CDFs are now stored in a variable transformed (smooth) form
53// together with the corresponding rational interpolation parameters.
54// The new version is several times faster, more robust and accurate
55// compared to the earlier version (G4GoudsmitSaundersonMscModel class
56// that use these data has been also completely replaced)
57// [1] A.F.Bielajew, NIMB, 111 (1996) 195-208
58// [2] I.Kawrakow, A.F.Bielajew, NIMB 134(1998) 325-336
59// 28.04.2017 M. Novak: the GS angular distributions has been recomputed, the
60// data size has been reduced from 16 MB down to 5 MB by using a new
61// representation, the class has been modified significantly due to
62// this new data representation.
63// 23.08.2017 M. Novak: Added funtionality to handle Mott-correction to the
64// base GS angular distributions and some other factors (screening
65// parameter, first and second moments) when Mott-correction is
66// activated in the GS-MSC model.
67// 26.10.2025 M. Novak: added the related technical note as the proper reference.
68//
69// References:
70// M. Novak: https://arxiv.org/abs/2410.13361
71//
72// -----------------------------------------------------------------------------
73
74
75#ifndef G4GoudsmitSaundersonTable_h
76#define G4GoudsmitSaundersonTable_h 1
77
78#include <vector>
79
80#include "G4Types.hh"
81
84
86
87public:
90
91 void Initialise(G4double lownergylimit, G4double highenergylimit);
92
93 // structure to store one GS transformed angular distribution (for a given s/lambda_el,s/lambda_elG1)
95 G4int fNumData; // # of data points
96 G4double *fUValues; // array of transformed variables
97 G4double *fParamA; // array of interpolation parameters a
98 G4double *fParamB; // array of interpolation parameters b
99 };
100
101 void LoadMSCData();
102
103 G4bool Sampling(G4double lambdaval, G4double qval, G4double scra,
104 G4double &cost, G4double &sint, G4double lekin,
105 G4double beta2, G4int matindx, GSMSCAngularDtr **gsDtr,
106 G4int &mcekini, G4int &mcdelti, G4double &transfPar,
107 G4bool isfirst);
108
109 G4double SampleCosTheta(G4double lambdaval, G4double qval, G4double scra,
110 G4double lekin, G4double beta2, G4int matindx,
111 GSMSCAngularDtr **gsDtr, G4int &mcekini, G4int &mcdelti,
112 G4double &transfPar, G4bool isfirst);
113
114 G4double SampleGSSRCosTheta(const GSMSCAngularDtr* gsDrt, G4double transfpar);
115
116 G4double SingleScattering(G4double lambdaval, G4double scra, G4double lekin,
117 G4double beta2, G4int matindx);
118
120 G4double &qval, G4double &transfpar);
121
122 // material dependent MSC parameters (computed at initialisation) regarding
123 // Moliere's screening parameter
124 G4double GetMoliereBc(G4int matindx) { return gMoliereBc[matindx]; }
125
126 G4double GetMoliereXc2(G4int matindx) { return gMoliereXc2[matindx]; }
127
128 void GetMottCorrectionFactors(G4double logekin, G4double beta2,
129 G4int matindx, G4double &mcToScr,
130 G4double &mcToQ1, G4double &mcToG2PerG1);
131
132 // set option to activate/inactivate Mott-correction
133 void SetOptionMottCorrection(G4bool val) { fIsMottCorrection = val; }
134 // set option to activate/inactivate PWA-correction
135 void SetOptionPWACorrection(G4bool val) { fIsPWACorrection = val; }
136
137 // this method returns with the scattering power correction (to avoid double counting of sub-threshold deflections)
138 // interpolated from tables prepared at initialisation
140
141 void InitSCPCorrection();
142
143private:
144 // initialisation of material dependent Moliere's MSC parameters
145 void InitMoliereMSCParams();
146
147
148 private:
149 static G4bool gIsInitialised; // are the precomputed angular distributions already loaded in?
150 static constexpr G4int gLAMBNUM = 64; // # L=s/lambda_el in [fLAMBMIN,fLAMBMAX]
151 static constexpr G4int gQNUM1 = 15; // # Q=s/lambda_el G1 in [fQMIN1,fQMAX1] in the 1-st Q grid
152 static constexpr G4int gQNUM2 = 32; // # Q=s/lambda_el G1 in [fQMIN2,fQMAX2] in the 2-nd Q grid
153 static constexpr G4int gNUMSCR1 = 201; // # of screening parameters in the A(G1) function
154 static constexpr G4int gNUMSCR2 = 51; // # of screening parameters in the A(G1) function
155 static constexpr G4double gLAMBMIN = 1.0; // minimum s/lambda_el
156 static constexpr G4double gLAMBMAX = 100000.0; // maximum s/lambda_el
157 static constexpr G4double gQMIN1 = 0.001; // minimum s/lambda_el G1 in the 1-st Q grid
158 static constexpr G4double gQMAX1 = 0.99; // maximum s/lambda_el G1 in the 1-st Q grid
159 static constexpr G4double gQMIN2 = 0.99; // minimum s/lambda_el G1 in the 2-nd Q grid
160 static constexpr G4double gQMAX2 = 7.99; // maximum s/lambda_el G1 in the 2-nd Q grid
161 //
162 G4bool fIsElectron; // GS-table for e- (for e+ otherwise)
163 G4bool fIsMottCorrection; // flag to indicate if Mott-correction was requested to be used
164 G4bool fIsPWACorrection; // flag to indicate is PWA corrections were requested to be used
165 G4double fLogLambda0; // ln(gLAMBMIN)
166 G4double fLogDeltaLambda; // ln(gLAMBMAX/gLAMBMIN)/(gLAMBNUM-1)
167 G4double fInvLogDeltaLambda; // 1/[ln(gLAMBMAX/gLAMBMIN)/(gLAMBNUM-1)]
168 G4double fInvDeltaQ1; // 1/[(gQMAX1-gQMIN1)/(gQNUM1-1)]
169 G4double fDeltaQ2; // [(gQMAX2-gQMIN2)/(gQNUM2-1)]
170 G4double fInvDeltaQ2; // 1/[(gQMAX2-gQMIN2)/(gQNUM2-1)]
171 //
172 G4double fLowEnergyLimit;
173 G4double fHighEnergyLimit;
174 //
175 int fNumSPCEbinPerDec; // scattering power correction energy grid bins per decade
176 struct SCPCorrection {
177 bool fIsUse; //
178 double fPrCut; // sec. e- production cut energy
179 double fLEmin; // log min energy
180 double fILDel; // inverse log delta kinetic energy
181 //std::vector<double> fVEkin; // scattering power correction energies
182 std::vector<double> fVSCPC; // scattering power correction vector
183 };
184 std::vector<SCPCorrection*> fSCPCPerMatCuts;
185
186
187 // vector to store all GS transformed angular distributions (cumputed based on the Screened-Rutherford DCS)
188 static std::vector<GSMSCAngularDtr*> gGSMSCAngularDistributions1;
189 static std::vector<GSMSCAngularDtr*> gGSMSCAngularDistributions2;
190
191 //@{
192 /** Precomputed \f$ b_lambda_{c} $\f and \f$ \chi_c^{2} $\f material dependent
193 * Moliere parameters that can be used to compute the screening parameter,
194 * the elastic scattering cross section (or \f$ \lambda_{e} $\f) under the
195 * screened Rutherford cross section approximation. (These are used in
196 * G4GoudsmitSaundersonMscModel if fgIsUsePWATotalXsecData is FALSE.)
197 */
198 static std::vector<double> gMoliereBc;
199 static std::vector<double> gMoliereXc2;
200 //
201 //
202 G4GSMottCorrection *fMottCorrection;
203};
204
205#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4double SampleGSSRCosTheta(const GSMSCAngularDtr *gsDrt, G4double transfpar)
G4double SampleCosTheta(G4double lambdaval, G4double qval, G4double scra, G4double lekin, G4double beta2, G4int matindx, GSMSCAngularDtr **gsDtr, G4int &mcekini, G4int &mcdelti, G4double &transfPar, G4bool isfirst)
G4double ComputeScatteringPowerCorrection(const G4MaterialCutsCouple *matcut, G4double ekin)
void GetMottCorrectionFactors(G4double logekin, G4double beta2, G4int matindx, G4double &mcToScr, G4double &mcToQ1, G4double &mcToG2PerG1)
G4double SingleScattering(G4double lambdaval, G4double scra, G4double lekin, G4double beta2, G4int matindx)
G4bool Sampling(G4double lambdaval, G4double qval, G4double scra, G4double &cost, G4double &sint, G4double lekin, G4double beta2, G4int matindx, GSMSCAngularDtr **gsDtr, G4int &mcekini, G4int &mcdelti, G4double &transfPar, G4bool isfirst)
GSMSCAngularDtr * GetGSAngularDtr(G4double scra, G4double &lambdaval, G4double &qval, G4double &transfpar)
void Initialise(G4double lownergylimit, G4double highenergylimit)
G4double GetMoliereXc2(G4int matindx)