00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054 #ifndef G4WentzelOKandVIxSection_h
00055 #define G4WentzelOKandVIxSection_h 1
00056
00057
00058
00059 #include "globals.hh"
00060 #include "G4Material.hh"
00061 #include "G4Element.hh"
00062 #include "G4ElementVector.hh"
00063 #include "G4NistManager.hh"
00064 #include "G4ThreeVector.hh"
00065 #include "G4Pow.hh"
00066
00067 class G4ParticleDefinition;
00068
00069
00070
00071 class G4WentzelOKandVIxSection
00072 {
00073
00074 public:
00075
00076 G4WentzelOKandVIxSection();
00077
00078 virtual ~G4WentzelOKandVIxSection();
00079
00080 void Initialise(const G4ParticleDefinition*, G4double CosThetaLim);
00081
00082 void SetupParticle(const G4ParticleDefinition*);
00083
00084
00085
00086 G4double SetupTarget(G4int Z, G4double cut = DBL_MAX);
00087
00088 G4double ComputeTransportCrossSectionPerAtom(G4double CosThetaMax);
00089
00090 G4ThreeVector SampleSingleScattering(G4double CosThetaMin,
00091 G4double CosThetaMax,
00092 G4double elecRatio = 0.0);
00093
00094 inline G4double ComputeNuclearCrossSection(G4double CosThetaMin,
00095 G4double CosThetaMax);
00096
00097 inline G4double ComputeElectronCrossSection(G4double CosThetaMin,
00098 G4double CosThetaMax);
00099
00100 inline G4double SetupKinematic(G4double kinEnergy, const G4Material* mat);
00101
00102 inline void SetTargetMass(G4double value);
00103
00104 inline G4double GetMomentumSquare() const;
00105
00106 inline G4double GetCosThetaNuc() const;
00107
00108 inline G4double GetCosThetaElec() const;
00109
00110 private:
00111
00112 void ComputeMaxElectronScattering(G4double cut);
00113
00114
00115 G4WentzelOKandVIxSection & operator=(const G4WentzelOKandVIxSection &right);
00116 G4WentzelOKandVIxSection(const G4WentzelOKandVIxSection&);
00117
00118 const G4ParticleDefinition* theProton;
00119 const G4ParticleDefinition* theElectron;
00120 const G4ParticleDefinition* thePositron;
00121 const G4Material* currentMaterial;
00122
00123 G4NistManager* fNistManager;
00124 G4Pow* fG4pow;
00125
00126 G4double numlimit;
00127
00128 G4double elecXSRatio;
00129
00130
00131 G4int nwarnings;
00132 G4int nwarnlimit;
00133
00134
00135 G4double coeff;
00136 G4double cosTetMaxElec;
00137 G4double cosTetMaxNuc;
00138 G4double cosThetaMax;
00139 G4double alpha2;
00140
00141
00142 const G4ParticleDefinition* particle;
00143
00144 G4double chargeSquare;
00145 G4double charge3;
00146 G4double spin;
00147 G4double mass;
00148 G4double tkin;
00149 G4double mom2;
00150 G4double momCM2;
00151 G4double invbeta2;
00152 G4double kinFactor;
00153 G4double etag;
00154 G4double ecut;
00155 G4double lowEnergyLimit;
00156
00157
00158 G4int targetZ;
00159 G4double targetMass;
00160 G4double screenZ;
00161 G4double formfactA;
00162 G4double factorA2;
00163 G4double factB;
00164 G4double factB1;
00165 G4double factD;
00166 G4double gam0pcmp;
00167 G4double pcmp2;
00168
00169 static G4double ScreenRSquareElec[100];
00170 static G4double ScreenRSquare[100];
00171 static G4double FormFactor[100];
00172 };
00173
00174
00175
00176 inline G4double
00177 G4WentzelOKandVIxSection::SetupKinematic(G4double ekin, const G4Material* mat)
00178 {
00179 if(ekin != tkin || mat != currentMaterial) {
00180 currentMaterial = mat;
00181 tkin = ekin;
00182 mom2 = tkin*(tkin + 2.0*mass);
00183 invbeta2 = 1.0 + mass*mass/mom2;
00184 factB = spin/invbeta2;
00185 cosTetMaxNuc =
00186 std::max(cosThetaMax,1.-factorA2*mat->GetIonisation()->GetInvA23()/mom2);
00187 }
00188 return cosTetMaxNuc;
00189 }
00190
00191
00192
00193 inline void G4WentzelOKandVIxSection::SetTargetMass(G4double value)
00194 {
00195 targetMass = value;
00196 factD = std::sqrt(mom2)/value;
00197 }
00198
00199
00200
00201 inline G4double G4WentzelOKandVIxSection::GetMomentumSquare() const
00202 {
00203 return mom2;
00204 }
00205
00206
00207
00208 inline G4double G4WentzelOKandVIxSection::GetCosThetaNuc() const
00209 {
00210 return cosTetMaxNuc;
00211 }
00212
00213
00214
00215 inline G4double G4WentzelOKandVIxSection::GetCosThetaElec() const
00216 {
00217 return cosTetMaxElec;
00218 }
00219
00220
00221
00222 inline G4double
00223 G4WentzelOKandVIxSection::ComputeNuclearCrossSection(G4double cosTMin,
00224 G4double cosTMax)
00225 {
00226 G4double xsec = 0.0;
00227 if(cosTMax < cosTMin) {
00228 xsec = targetZ*kinFactor*(cosTMin - cosTMax)/
00229 ((1.0 - cosTMin + screenZ)*(1.0 - cosTMax + screenZ));
00230 }
00231 return xsec;
00232 }
00233
00234
00235
00236 inline G4double
00237 G4WentzelOKandVIxSection::ComputeElectronCrossSection(G4double cosTMin,
00238 G4double cosTMax)
00239 {
00240 G4double xsec = 0.0;
00241 G4double cost1 = std::max(cosTMin,cosTetMaxElec);
00242 G4double cost2 = std::max(cosTMax,cosTetMaxElec);
00243 if(cost1 > cost2) {
00244 xsec = kinFactor*(cost1 - cost2)/((1.0 - cost1 + screenZ)*(1.0 - cost2 + screenZ));
00245 }
00246 return xsec;
00247 }
00248
00249
00250 #endif
00251