File indexing completed on 2025-01-18 09:58:12
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042
0043
0044
0045
0046
0047
0048
0049
0050
0051
0052
0053 #ifndef G4EmCorrections_h
0054 #define G4EmCorrections_h 1
0055
0056 #include "globals.hh"
0057 #include "G4ionEffectiveCharge.hh"
0058 #include "G4Material.hh"
0059 #include "G4ParticleDefinition.hh"
0060
0061 class G4VEmModel;
0062 class G4PhysicsVector;
0063 class G4IonTable;
0064 class G4MaterialCutsCouple;
0065 class G4PhysicsFreeVector;
0066 class G4Pow;
0067
0068 class G4EmCorrections
0069 {
0070
0071 public:
0072
0073 explicit G4EmCorrections(G4int verb);
0074
0075 ~G4EmCorrections();
0076
0077 G4double HighOrderCorrections(const G4ParticleDefinition*,
0078 const G4Material*,
0079 const G4double kineticEnergy,
0080 const G4double cutEnergy);
0081
0082 G4double IonHighOrderCorrections(const G4ParticleDefinition*,
0083 const G4MaterialCutsCouple*,
0084 const G4double kineticEnergy);
0085
0086 G4double ComputeIonCorrections(const G4ParticleDefinition*,
0087 const G4Material*,
0088 const G4double kineticEnergy);
0089
0090 G4double IonBarkasCorrection(const G4ParticleDefinition*,
0091 const G4Material*,
0092 const G4double kineticEnergy);
0093
0094 G4double Bethe(const G4ParticleDefinition*,
0095 const G4Material*,
0096 const G4double kineticEnergy);
0097
0098 G4double SpinCorrection(const G4ParticleDefinition*,
0099 const G4Material*,
0100 const G4double kineticEnergy);
0101
0102 G4double KShellCorrection(const G4ParticleDefinition*,
0103 const G4Material*,
0104 const G4double kineticEnergy);
0105
0106 G4double LShellCorrection(const G4ParticleDefinition*,
0107 const G4Material*,
0108 const G4double kineticEnergy);
0109
0110 G4double ShellCorrection(const G4ParticleDefinition*,
0111 const G4Material*,
0112 const G4double kineticEnergy);
0113
0114 G4double ShellCorrectionSTD(const G4ParticleDefinition*,
0115 const G4Material*,
0116 const G4double kineticEnergy);
0117
0118 G4double DensityCorrection(const G4ParticleDefinition*,
0119 const G4Material*,
0120 const G4double kineticEnergy);
0121
0122 G4double BarkasCorrection(const G4ParticleDefinition*,
0123 const G4Material*,
0124 const G4double kineticEnergy,
0125 const G4bool isInitialized = false);
0126
0127 G4double BlochCorrection(const G4ParticleDefinition*,
0128 const G4Material*,
0129 const G4double kineticEnergy,
0130 const G4bool isInitialized = false);
0131
0132 G4double MottCorrection(const G4ParticleDefinition*,
0133 const G4Material*,
0134 const G4double kineticEnergy,
0135 const G4bool isInitialized = false);
0136
0137 void AddStoppingData(const G4int Z, const G4int A,
0138 const G4String& materialName,
0139 G4PhysicsVector* dVector);
0140
0141 void InitialiseForNewRun();
0142
0143
0144 G4double EffectiveChargeCorrection(const G4ParticleDefinition*,
0145 const G4Material*,
0146 const G4double kineticEnergy);
0147
0148
0149 inline G4double GetParticleCharge(const G4ParticleDefinition*,
0150 const G4Material*,
0151 const G4double kineticEnergy);
0152
0153 inline
0154 G4double EffectiveChargeSquareRatio(const G4ParticleDefinition*,
0155 const G4Material*,
0156 const G4double kineticEnergy);
0157
0158
0159 inline void SetIonisationModels(G4VEmModel* mod1 = nullptr,
0160 G4VEmModel* mod2 = nullptr);
0161
0162 inline G4int GetNumberOfStoppingVectors() const;
0163
0164 inline void SetVerbose(G4int verb);
0165
0166
0167 G4EmCorrections & operator=(const G4EmCorrections &right) = delete;
0168 G4EmCorrections(const G4EmCorrections&) = delete;
0169
0170 private:
0171
0172 void Initialise();
0173
0174 void BuildCorrectionVector();
0175
0176 void SetupKinematics(const G4ParticleDefinition*,
0177 const G4Material*,
0178 const G4double kineticEnergy);
0179
0180 G4double KShell(const G4double theta, const G4double eta);
0181
0182 G4double LShell(const G4double theta, const G4double eta);
0183
0184 G4int Index(const G4double x, const G4double* y, const G4int n) const;
0185
0186 G4double Value(const G4double xv, const G4double x1, const G4double x2,
0187 const G4double y1, const G4double y2) const;
0188
0189 G4double Value2(const G4double xv, const G4double yv,
0190 const G4double x1, const G4double x2,
0191 const G4double y1, const G4double y2,
0192 const G4double z11, const G4double z21,
0193 const G4double z12, const G4double z22) const;
0194
0195 G4Pow* g4calc;
0196 G4IonTable* ionTable;
0197
0198 const G4ParticleDefinition* particle = nullptr;
0199 const G4ParticleDefinition* curParticle = nullptr;
0200 const G4Material* material = nullptr;
0201 const G4Material* curMaterial = nullptr;
0202 const G4ElementVector* theElementVector = nullptr;
0203 const G4double* atomDensity = nullptr;
0204
0205 G4PhysicsVector* curVector = nullptr;
0206
0207 G4VEmModel* ionLEModel = nullptr;
0208 G4VEmModel* ionHEModel = nullptr;
0209
0210 G4double kinEnergy = 0.0;
0211 G4double mass = 0.0;
0212 G4double massFactor = 1.0;
0213 G4double tau = 0.0;
0214 G4double gamma = 1.0;
0215 G4double bg2 = 0.0;
0216 G4double beta2 = 0.0;
0217 G4double beta = 0.0;
0218 G4double ba2 = 0.0;
0219 G4double tmax = 0.0;
0220 G4double charge = 0.0;
0221 G4double q2 = 0.0;
0222 G4double eth;
0223 G4double eCorrMin;
0224 G4double eCorrMax;
0225
0226 std::size_t ncouples = 0;
0227 std::size_t idxBarkas = 0;
0228 G4int nK = 20;
0229 G4int nL = 26;
0230 G4int nEtaK = 29;
0231 G4int nEtaL = 28;
0232 G4int nbinCorr = 52;
0233 G4int numberOfElements = 0;
0234
0235
0236 G4int nIons = 0;
0237 G4int idx = 0;
0238 G4int currentZ = 0;
0239
0240 G4int verbose;
0241 G4bool isInitializer = false;
0242
0243 std::vector<G4int> Zion;
0244 std::vector<G4int> Aion;
0245 std::vector<G4String> materialName;
0246 std::vector<const G4ParticleDefinition*> ionList;
0247
0248 std::map< G4int, std::vector<G4double> > thcorr;
0249
0250 std::vector<const G4Material*> currmat;
0251 std::vector<const G4Material*> materialList;
0252 std::vector<G4PhysicsVector*> stopData;
0253
0254 G4ionEffectiveCharge effCharge;
0255
0256 static const G4double ZD[11];
0257 static const G4double UK[20];
0258 static const G4double VK[20];
0259 static G4double ZK[20];
0260 static const G4double Eta[29];
0261 static G4double CK[20][29];
0262 static G4double CL[26][28];
0263 static const G4double UL[26];
0264 static G4double VL[26];
0265
0266 static G4double sWmaxBarkas;
0267 static G4PhysicsFreeVector* sBarkasCorr;
0268 static G4PhysicsFreeVector* sThetaK;
0269 static G4PhysicsFreeVector* sThetaL;
0270 };
0271
0272 inline G4int
0273 G4EmCorrections::Index(const G4double x, const G4double* y, const G4int n) const
0274 {
0275 G4int iddd = n-1;
0276
0277 do {--iddd;} while (iddd>0 && x<y[iddd]);
0278 return iddd;
0279 }
0280
0281 inline G4double G4EmCorrections::Value(const G4double xv, const G4double x1,
0282 const G4double x2,
0283 const G4double y1, const G4double y2) const
0284 {
0285 return y1 + (y2 - y1)*(xv - x1)/(x2 - x1);
0286 }
0287
0288 inline G4double G4EmCorrections::Value2(const G4double xv, const G4double yv,
0289 const G4double x1, const G4double x2,
0290 const G4double y1, const G4double y2,
0291 const G4double z11, const G4double z21,
0292 const G4double z12, const G4double z22) const
0293 {
0294 return ( z11*(x2-xv)*(y2-yv) + z22*(xv-x1)*(yv-y1) +
0295 z12*(x2-xv)*(yv-y1) + z21*(xv-x1)*(y2-yv) )
0296 / ((x2-x1)*(y2-y1));
0297 }
0298
0299 inline void
0300 G4EmCorrections::SetIonisationModels(G4VEmModel* mod1, G4VEmModel* mod2)
0301 {
0302 if(nullptr != mod1) { ionLEModel = mod1; }
0303 if(nullptr != mod2) { ionHEModel = mod2; }
0304 }
0305
0306 inline G4int G4EmCorrections::GetNumberOfStoppingVectors() const
0307 {
0308 return nIons;
0309 }
0310
0311 inline G4double
0312 G4EmCorrections::GetParticleCharge(const G4ParticleDefinition* p,
0313 const G4Material* mat,
0314 const G4double kineticEnergy)
0315 {
0316 return effCharge.EffectiveCharge(p,mat,kineticEnergy);
0317 }
0318
0319 inline G4double
0320 G4EmCorrections::EffectiveChargeSquareRatio(const G4ParticleDefinition* p,
0321 const G4Material* mat,
0322 const G4double kineticEnergy)
0323 {
0324 return effCharge.EffectiveChargeSquareRatio(p,mat,kineticEnergy);
0325 }
0326
0327 inline void G4EmCorrections::SetVerbose(G4int verb)
0328 {
0329 verbose = verb;
0330 }
0331
0332
0333
0334 #endif