|
||||
File indexing completed on 2025-01-18 09:58:21
0001 // 0002 // ******************************************************************** 0003 // * License and Disclaimer * 0004 // * * 0005 // * The Geant4 software is copyright of the Copyright Holders of * 0006 // * the Geant4 Collaboration. It is provided under the terms and * 0007 // * conditions of the Geant4 Software License, included in the file * 0008 // * LICENSE and available at http://cern.ch/geant4/license . These * 0009 // * include a list of copyright holders. * 0010 // * * 0011 // * Neither the authors of this software system, nor their employing * 0012 // * institutes,nor the agencies providing financial support for this * 0013 // * work make any representation or warranty, express or implied, * 0014 // * regarding this software system or assume any liability for its * 0015 // * use. Please see the license in the file LICENSE and URL above * 0016 // * for the full disclaimer and the limitation of liability. * 0017 // * * 0018 // * This code implementation is the result of the scientific and * 0019 // * technical work of the GEANT4 collaboration. * 0020 // * By using, copying, modifying or distributing the software (or * 0021 // * any work based on the software) you agree to acknowledge its * 0022 // * use in resulting scientific publications, and indicate your * 0023 // * acceptance of all terms of the Geant4 Software license. * 0024 // ******************************************************************** 0025 // 0026 // 0027 // ---------------------------------------------------------------------------- 0028 // 0029 // GEANT4 Class header file 0030 // 0031 // File name: G4GoudsmitSaundersonMscModel 0032 // 0033 // Author: Mihaly Novak / (Omrane Kadri) 0034 // 0035 // Creation date: 20.02.2009 0036 // 0037 // Modifications: 0038 // 04.03.2009 V.Ivanchenko cleanup and format according to Geant4 EM style 0039 // 12.05.2010 O.Kadri: adding Qn1 and Qn12 as private doubles 0040 // 18.05.2015 M. Novak provide PLERIMINARYY version of updated class. 0041 // All algorithms of the class were revised and updated, new methods added. 0042 // A new version of Kawrakow-Bielajew Goudsmit-Saunderson MSC model 0043 // based on the screened Rutherford DCS for elastic scattering of 0044 // electrons/positrons has been introduced[1,2]. The corresponding MSC 0045 // angular distributions over a 2D parameter grid have been recomputed 0046 // and the CDFs are now stored in a variable transformed (smooth) form[2,3] 0047 // together with the corresponding rational interpolation parameters. 0048 // These angular distributions are handled by the new 0049 // G4GoudsmitSaundersonTable class that is responsible to sample if 0050 // it was no, single, few or multiple scattering case and delivers the 0051 // angular deflection (i.e. cos(theta) and sin(theta)). 0052 // Two screening options are provided: 0053 // - if fgIsUsePWATotalXsecData=TRUE i.e. SetOptionPWAScreening(TRUE) 0054 // was called before initialisation: screening parameter value A is 0055 // determined such that the first transport coefficient G1(A) 0056 // computed according to the screened Rutherford DCS for elastic 0057 // scattering will reproduce the one computed from the PWA elastic 0058 // and first transport mean free paths[4]. 0059 // - if fgIsUsePWATotalXsecData=FALSE i.e. default value or 0060 // SetOptionPWAScreening(FALSE) was called before initialisation: 0061 // screening parameter value A is computed according to Moliere's 0062 // formula (by using material dependent parameters \chi_cc2 and b_c 0063 // precomputed for each material used at initialization in 0064 // G4GoudsmitSaundersonTable) [3] 0065 // Elastic and first trasport mean free paths are used consistently. 0066 // The new version is self-consistent, several times faster, more 0067 // robust and accurate compared to the earlier version. 0068 // Spin effects as well as a more accurate energy loss correction and 0069 // computations of Lewis moments will be implemented later on. 0070 // 02.09.2015 M. Novak: first version of new step limit is provided. 0071 // fUseSafetyPlus corresponds to Urban fUseSafety (default) 0072 // fUseDistanceToBoundary corresponds to Urban fUseDistanceToBoundary 0073 // fUseSafety corresponds to EGSnrc error-free stepping algorithm 0074 // Range factor can be significantly higher at each case than in Urban. 0075 // 23.08.2017 M. Novak: added corrections to account spin effects (Mott-correction). 0076 // It can be activated by setting the fIsMottCorrection flag to be true 0077 // before initialization using the SetOptionMottCorrection() public method. 0078 // The fMottCorrection member is responsible to handle pre-computed Mott 0079 // correction (rejection) functions obtained by numerically computing 0080 // Goudsmit-Saunderson agnular distributions based on a DCS accounting spin 0081 // effects and screening corrections. The DCS used to compute the accurate 0082 // GS angular distributions is: DCS_{cor} = DCS_{SR}x[ DCS_{R}/DCS_{Mott}] where : 0083 // # DCS_{SR} is the relativistic Screened-Rutherford DCS (first Born approximate 0084 // solution of the Klein-Gordon i.e. relativistic Schrodinger equation => 0085 // scattering of spinless e- on exponentially screened Coulomb potential) 0086 // note: the default (without using Mott-correction) GS angular distributions 0087 // are based on this DCS_{SR} with Moliere's screening parameter! 0088 // # DCS_{R} is the Rutherford DCS which is the same as above but without 0089 // screening 0090 // # DCS_{Mott} is the Mott DCS i.e. solution of the Dirac equation with a bare 0091 // Coulomb potential i.e. scattering of particles with spin (e- or e+) on a 0092 // point-like unscreened Coulomb potential 0093 // # moreover, the screening parameter of the DCS_{cor} was determined such that 0094 // the DCS_{cor} with this corrected screening parameter reproduce the first 0095 // transport cross sections obtained from the corresponding most accurate DCS 0096 // (i.e. from elsepa [4]) 0097 // Unlike the default GS, the Mott-corrected angular distributions are particle type 0098 // (different for e- and e+ <= the DCS_{Mott} and the screening correction) and target 0099 // (Z and material) dependent. 0100 // 02.02.2018 M. Novak: implemented CrossSectionPerVolume interface method (used only for testing) 0101 // 0102 // Class description: 0103 // Kawrakow-Bielajew Goudsmit-Saunderson MSC model based on the screened Rutherford DCS 0104 // for elastic scattering of e-/e+. Option, to include (Mott) correction (see above), is 0105 // also available now (SetOptionMottCorrection(true)). An EGSnrc like error-free stepping 0106 // algorithm (UseSafety) is available beyond the usual Geant4 step limitation algorithms 0107 // and true to geomerty and geometry to true step length computations that were adopted 0108 // from the Urban model[5]. The most accurate setting: error-free stepping (UseSafety) 0109 // with Mott-correction (SetOptionMottCorrection(true)). 0110 // 0111 // References: 0112 // [1] A.F.Bielajew, NIMB 111 (1996) 195-208 0113 // [2] I.Kawrakow, A.F.Bielajew, NIMB 134(1998) 325-336 0114 // [3] I.Kawrakow, E.Mainegra-Hing, D.W.O.Rogers, F.Tessier,B.R.B.Walters, NRCC 0115 // Report PIRS-701 (2013) 0116 // [4] F.Salvat, A.Jablonski, C.J. Powell, CPC 165(2005) 157-190 0117 // [5] L.Urban, Preprint CERN-OPEN-2006-077 (2006) 0118 // 0119 // ----------------------------------------------------------------------------- 0120 0121 #ifndef G4GoudsmitSaundersonMscModel_h 0122 #define G4GoudsmitSaundersonMscModel_h 1 0123 0124 #include <CLHEP/Units/SystemOfUnits.h> 0125 0126 #include "G4VMscModel.hh" 0127 #include "G4PhysicsTable.hh" 0128 #include "G4MaterialCutsCouple.hh" 0129 #include "globals.hh" 0130 0131 0132 class G4DataVector; 0133 class G4ParticleChangeForMSC; 0134 class G4LossTableManager; 0135 class G4GoudsmitSaundersonTable; 0136 class G4GSPWACorrections; 0137 0138 class G4GoudsmitSaundersonMscModel : public G4VMscModel 0139 { 0140 public: 0141 0142 G4GoudsmitSaundersonMscModel(const G4String& nam = "GoudsmitSaunderson"); 0143 0144 ~G4GoudsmitSaundersonMscModel() override; 0145 0146 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override; 0147 0148 void InitialiseLocal(const G4ParticleDefinition* p, G4VEmModel* masterModel) override; 0149 0150 G4ThreeVector& SampleScattering(const G4ThreeVector&, G4double safety) override; 0151 0152 G4double ComputeTruePathLengthLimit(const G4Track& track, G4double& currentMinimalStep) override; 0153 0154 G4double ComputeGeomPathLength(G4double truePathLength) override; 0155 0156 G4double ComputeTrueStepLength(G4double geomStepLength) override; 0157 0158 // method to compute first transport cross section per Volume (i.e. macroscropic first transport cross section; this 0159 // method is used only for testing and not during a normal simulation) 0160 G4double CrossSectionPerVolume(const G4Material*, const G4ParticleDefinition*, G4double kineticEnergy, G4double cutEnergy = 0.0, G4double maxEnergy = DBL_MAX) override; 0161 0162 void StartTracking(G4Track*) override; 0163 0164 void SampleMSC(); 0165 0166 G4double GetTransportMeanFreePath(const G4ParticleDefinition*, G4double); 0167 0168 void SetOptionPWACorrection(G4bool opt) { fIsUsePWACorrection = opt; } 0169 0170 G4bool GetOptionPWACorrection() const { return fIsUsePWACorrection; } 0171 0172 void SetOptionMottCorrection(G4bool opt) { fIsUseMottCorrection = opt; } 0173 0174 G4bool GetOptionMottCorrection() const { return fIsUseMottCorrection; } 0175 0176 G4GoudsmitSaundersonTable* GetGSTable() { return fGSTable; } 0177 0178 G4GSPWACorrections* GetPWACorrection() { return fPWACorrection; } 0179 0180 // hide assignment operator 0181 G4GoudsmitSaundersonMscModel & operator=(const G4GoudsmitSaundersonMscModel &right) = delete; 0182 G4GoudsmitSaundersonMscModel(const G4GoudsmitSaundersonMscModel&) = delete; 0183 0184 private: 0185 inline void SetParticle(const G4ParticleDefinition* p); 0186 0187 inline G4double GetLambda(G4double); 0188 0189 G4double GetTransportMeanFreePathOnly(const G4ParticleDefinition*,G4double); 0190 0191 inline G4double Randomizetlimit(); 0192 0193 private: 0194 CLHEP::HepRandomEngine* rndmEngineMod; 0195 // 0196 G4double currentKinEnergy; 0197 G4double currentRange; 0198 // 0199 G4double fr; 0200 G4double rangeinit; 0201 G4double geombig; 0202 G4double geomlimit; 0203 G4double tlimit; 0204 G4double tgeom; 0205 // 0206 G4double par1; 0207 G4double par2; 0208 G4double par3; 0209 G4double tlimitminfix2; 0210 G4double tausmall; 0211 G4double mass; 0212 G4double taulim; 0213 // 0214 // 0215 G4double presafety; 0216 G4double fZeff; 0217 // 0218 G4int charge; 0219 G4int currentMaterialIndex; 0220 // 0221 G4bool firstStep; 0222 // 0223 G4LossTableManager* theManager; 0224 const G4ParticleDefinition* particle; 0225 G4ParticleChangeForMSC* fParticleChange; 0226 const G4MaterialCutsCouple* currentCouple; 0227 0228 G4GoudsmitSaundersonTable* fGSTable; 0229 G4GSPWACorrections* fPWACorrection; 0230 0231 G4bool fIsUsePWACorrection; 0232 G4bool fIsUseMottCorrection; 0233 // 0234 G4double fLambda0; // elastic mean free path 0235 G4double fLambda1; // first transport mean free path 0236 G4double fScrA; // screening parameter 0237 G4double fG1; // first transport coef. 0238 // in case of Mott-correction 0239 G4double fMCtoScrA; 0240 G4double fMCtoQ1; 0241 G4double fMCtoG2PerG1; 0242 // 0243 G4double fTheTrueStepLenght; 0244 G4double fTheTransportDistance; 0245 G4double fTheZPathLenght; 0246 // 0247 G4ThreeVector fTheDisplacementVector; 0248 G4ThreeVector fTheNewDirection; 0249 // 0250 G4bool fIsEndedUpOnBoundary; // step ended up on boundary i.e. transportation is the winer 0251 G4bool fIsMultipleSacettring; 0252 G4bool fIsSingleScattering; 0253 G4bool fIsEverythingWasDone; 0254 G4bool fIsNoScatteringInMSC; 0255 G4bool fIsNoDisplace; 0256 G4bool fIsInsideSkin; 0257 G4bool fIsWasOnBoundary; 0258 G4bool fIsFirstRealStep; 0259 // 0260 static G4bool gIsUseAccurate; 0261 static G4bool gIsOptimizationOn; 0262 }; 0263 0264 //////////////////////////////////////////////////////////////////////////////// 0265 inline 0266 void G4GoudsmitSaundersonMscModel::SetParticle(const G4ParticleDefinition* p) 0267 { 0268 if (p != particle) { 0269 particle = p; 0270 charge = (G4int)(p->GetPDGCharge()/CLHEP::eplus); 0271 mass = p->GetPDGMass(); 0272 } 0273 } 0274 0275 0276 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0277 inline 0278 G4double G4GoudsmitSaundersonMscModel::Randomizetlimit() 0279 { 0280 G4double temptlimit; 0281 do { 0282 temptlimit = G4RandGauss::shoot(rndmEngineMod,tlimit,0.1*tlimit); 0283 } while ( (temptlimit<0.) || (temptlimit>2.*tlimit)); 0284 0285 return temptlimit; 0286 } 0287 0288 0289 0290 #endif
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |