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
00055 #ifndef G4UrbanMscModel93_h
00056 #define G4UrbanMscModel93_h 1
00057
00058
00059
00060 #include <CLHEP/Units/SystemOfUnits.h>
00061
00062 #include "G4VMscModel.hh"
00063 #include "G4MscStepLimitType.hh"
00064
00065 class G4ParticleChangeForMSC;
00066 class G4SafetyHelper;
00067 class G4LossTableManager;
00068
00069
00070
00071 class G4UrbanMscModel93 : public G4VMscModel
00072 {
00073
00074 public:
00075
00076 G4UrbanMscModel93(const G4String& nam = "UrbanMsc93");
00077
00078 virtual ~G4UrbanMscModel93();
00079
00080 void Initialise(const G4ParticleDefinition*, const G4DataVector&);
00081
00082 void StartTracking(G4Track*);
00083
00084 G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition* particle,
00085 G4double KineticEnergy,
00086 G4double AtomicNumber,
00087 G4double AtomicWeight=0.,
00088 G4double cut =0.,
00089 G4double emax=DBL_MAX);
00090
00091 G4ThreeVector& SampleScattering(const G4ThreeVector&, G4double safety);
00092
00093 G4double ComputeTruePathLengthLimit(const G4Track& track,
00094 G4double& currentMinimalStep);
00095
00096 G4double ComputeGeomPathLength(G4double truePathLength);
00097
00098 G4double ComputeTrueStepLength(G4double geomStepLength);
00099
00100 G4double ComputeTheta0(G4double truePathLength,
00101 G4double KineticEnergy);
00102
00103 private:
00104
00105 G4double SimpleScattering(G4double xmeanth, G4double x2meanth);
00106
00107 G4double SampleCosineTheta(G4double trueStepLength, G4double KineticEnergy);
00108
00109 G4double SampleDisplacement();
00110
00111 G4double LatCorrelation();
00112
00113 inline void SetParticle(const G4ParticleDefinition*);
00114
00115 inline void UpdateCache();
00116
00117
00118 G4UrbanMscModel93 & operator=(const G4UrbanMscModel93 &right);
00119 G4UrbanMscModel93(const G4UrbanMscModel93&);
00120
00121 const G4ParticleDefinition* particle;
00122 G4ParticleChangeForMSC* fParticleChange;
00123
00124 const G4MaterialCutsCouple* couple;
00125 G4LossTableManager* theManager;
00126
00127 G4double mass;
00128 G4double charge,ChargeSquare;
00129 G4double masslimite,lambdalimit,fr;
00130
00131 G4double taubig;
00132 G4double tausmall;
00133 G4double taulim;
00134 G4double currentTau;
00135 G4double tlimit;
00136 G4double tlimitmin;
00137 G4double tlimitminfix;
00138 G4double tgeom;
00139
00140 G4double geombig;
00141 G4double geommin;
00142 G4double geomlimit;
00143 G4double skindepth;
00144 G4double smallstep;
00145
00146 G4double presafety;
00147
00148 G4double lambda0;
00149 G4double lambdaeff;
00150 G4double tPathLength;
00151 G4double zPathLength;
00152 G4double par1,par2,par3;
00153
00154 G4double stepmin;
00155
00156 G4double currentKinEnergy;
00157 G4double currentRange;
00158 G4double rangeinit;
00159 G4double currentRadLength;
00160
00161 G4double numlim, xsi, ea, eaa;
00162
00163 G4double theta0max,rellossmax;
00164 G4double third;
00165
00166 G4int currentMaterialIndex;
00167
00168 G4double y;
00169 G4double Zold;
00170 G4double Zeff,Z2,Z23,lnZ;
00171 G4double coeffth1,coeffth2;
00172 G4double coeffc1,coeffc2;
00173 G4double scr1ini,scr2ini,scr1,scr2;
00174
00175 G4bool firstStep;
00176 G4bool inside;
00177 G4bool insideskin;
00178 };
00179
00180
00181
00182
00183 inline
00184 void G4UrbanMscModel93::SetParticle(const G4ParticleDefinition* p)
00185 {
00186 if (p != particle) {
00187 particle = p;
00188 mass = p->GetPDGMass();
00189 charge = p->GetPDGCharge()/CLHEP::eplus;
00190 ChargeSquare = charge*charge;
00191 }
00192 }
00193
00194
00195
00196 inline
00197 void G4UrbanMscModel93::UpdateCache()
00198 {
00199 lnZ = std::log(Zeff);
00200
00201 coeffth1 = (1.-8.7780e-2/Zeff)*(0.87+0.03*lnZ);
00202 coeffth2 = (4.0780e-2+1.7315e-4*Zeff)*(0.87+0.03*lnZ);
00203
00204 G4double lnZ1 = std::log(Zeff+1.);
00205 coeffc1 = 2.943-0.197*lnZ1;
00206 coeffc2 = 0.0987-0.0143*lnZ1;
00207
00208 Z2 = Zeff*Zeff;
00209 Z23 = std::exp(2.*lnZ/3.);
00210 scr1 = scr1ini*Z23;
00211 scr2 = scr2ini*Z2*ChargeSquare;
00212
00213 Zold = Zeff;
00214 }
00215
00216 #endif
00217