File indexing completed on 2025-01-18 09:59:28
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 G4WentzelOKandVIxSection_h
0054 #define G4WentzelOKandVIxSection_h 1
0055
0056
0057
0058 #include "globals.hh"
0059 #include "G4Material.hh"
0060 #include "G4Element.hh"
0061 #include "G4ElementVector.hh"
0062 #include "G4NistManager.hh"
0063 #include "G4NuclearFormfactorType.hh"
0064 #include "G4ThreeVector.hh"
0065 #include "G4Pow.hh"
0066
0067 class G4ParticleDefinition;
0068 class G4ScreeningMottCrossSection;
0069
0070
0071
0072 class G4WentzelOKandVIxSection
0073 {
0074
0075 public:
0076
0077 explicit G4WentzelOKandVIxSection(G4bool combined=true);
0078
0079 virtual ~G4WentzelOKandVIxSection();
0080
0081 void Initialise(const G4ParticleDefinition*, G4double CosThetaLim);
0082
0083 void SetupParticle(const G4ParticleDefinition*);
0084
0085
0086
0087 virtual G4double SetupKinematic(G4double kinEnergy, const G4Material* mat);
0088 G4double SetupTarget(G4int Z, G4double cut);
0089
0090 G4double ComputeTransportCrossSectionPerAtom(G4double CosThetaMax);
0091
0092 G4ThreeVector& SampleSingleScattering(G4double CosThetaMin,
0093 G4double CosThetaMax,
0094 G4double elecRatio);
0095
0096 G4double ComputeSecondTransportMoment(G4double CosThetaMax);
0097
0098 inline G4double ComputeNuclearCrossSection(G4double CosThetaMin,
0099 G4double CosThetaMax);
0100
0101 inline G4double ComputeElectronCrossSection(G4double CosThetaMin,
0102 G4double CosThetaMax);
0103
0104 inline void SetTargetMass(G4double value);
0105
0106 inline G4double GetMomentumSquare() const;
0107
0108 inline G4double GetCosThetaNuc() const;
0109
0110 inline G4double GetCosThetaElec() const;
0111
0112
0113 G4WentzelOKandVIxSection & operator=
0114 (const G4WentzelOKandVIxSection &right) = delete;
0115 G4WentzelOKandVIxSection(const G4WentzelOKandVIxSection&) = delete;
0116
0117 protected:
0118
0119 void ComputeMaxElectronScattering(G4double cut);
0120
0121 void InitialiseA();
0122
0123 inline G4double FlatFormfactor(G4double x);
0124
0125 const G4ParticleDefinition* theProton;
0126 const G4ParticleDefinition* theElectron;
0127 const G4ParticleDefinition* thePositron;
0128 const G4ParticleDefinition* particle = nullptr;
0129 const G4Material* currentMaterial = nullptr;
0130
0131 G4NistManager* fNistManager;
0132 G4Pow* fG4pow;
0133
0134 G4ScreeningMottCrossSection* fMottXSection = nullptr;
0135
0136 G4ThreeVector temp;
0137
0138
0139 G4double coeff;
0140 G4double cosTetMaxElec = 1.0;
0141 G4double cosTetMaxNuc = 1.0;
0142
0143
0144
0145 G4double cosThetaMax = 1.0;
0146
0147 G4double chargeSquare = 0.0;
0148 G4double charge3 = 0.0;
0149 G4double spin = 0.0;
0150 G4double mass = 0.0;
0151 G4double tkin = 0.0;
0152 G4double mom2 = 0.0;
0153 G4double momCM2 = 0.0;
0154 G4double invbeta2 = 1.0;
0155 G4double kinFactor = 1.0;
0156 G4double etag = DBL_MAX;
0157 G4double ecut = DBL_MAX;
0158
0159
0160 G4double targetMass;
0161 G4double screenZ = 0.0;
0162 G4double formfactA = 0.0;
0163 G4double factorA2 = 0.0;
0164 G4double factB = 0.0;
0165 G4double factD = 0.0;
0166 G4double fMottFactor = 1.0;
0167 G4double gam0pcmp = 1.0;
0168 G4double pcmp2 = 1.0;
0169
0170
0171 G4int targetZ = 0;
0172 G4int nwarnings = 0;
0173
0174 G4NuclearFormfactorType fNucFormfactor = fExponentialNF;
0175
0176 G4bool isCombined;
0177
0178 static G4double ScreenRSquareElec[100];
0179 static G4double ScreenRSquare[100];
0180 static G4double FormFactor[100];
0181 };
0182
0183
0184
0185 inline void G4WentzelOKandVIxSection::SetTargetMass(G4double value)
0186 {
0187 targetMass = value;
0188 factD = std::sqrt(mom2)/value;
0189 }
0190
0191
0192
0193 inline G4double G4WentzelOKandVIxSection::GetMomentumSquare() const
0194 {
0195 return mom2;
0196 }
0197
0198
0199
0200 inline G4double G4WentzelOKandVIxSection::GetCosThetaNuc() const
0201 {
0202 return cosTetMaxNuc;
0203 }
0204
0205
0206
0207 inline G4double G4WentzelOKandVIxSection::GetCosThetaElec() const
0208 {
0209 return cosTetMaxElec;
0210 }
0211
0212
0213
0214 inline G4double
0215 G4WentzelOKandVIxSection::ComputeNuclearCrossSection(G4double cosTMin,
0216 G4double cosTMax)
0217 {
0218 return targetZ*kinFactor*fMottFactor*(cosTMin - cosTMax)/
0219 ((1.0 - cosTMin + screenZ)*(1.0 - cosTMax + screenZ));
0220 }
0221
0222
0223
0224 inline G4double
0225 G4WentzelOKandVIxSection::ComputeElectronCrossSection(G4double cosTMin,
0226 G4double cosTMax)
0227 {
0228 G4double cost1 = std::max(cosTMin,cosTetMaxElec);
0229 G4double cost2 = std::max(cosTMax,cosTetMaxElec);
0230 return (cost1 <= cost2) ? 0.0 : kinFactor*fMottFactor*(cost1 - cost2)/
0231 ((1.0 - cost1 + screenZ)*(1.0 - cost2 + screenZ));
0232 }
0233
0234 inline G4double G4WentzelOKandVIxSection::FlatFormfactor(G4double x)
0235 {
0236 return 3.0*(std::sin(x) - x*std::cos(x))/(x*x*x);
0237 }
0238
0239
0240
0241 #endif
0242