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 #ifndef G4VMscModel_h
00050 #define G4VMscModel_h 1
00051
00052 #include <CLHEP/Units/SystemOfUnits.h>
00053
00054 #include "G4VEmModel.hh"
00055 #include "G4MscStepLimitType.hh"
00056 #include "globals.hh"
00057 #include "G4ThreeVector.hh"
00058 #include "G4Track.hh"
00059 #include "G4SafetyHelper.hh"
00060 #include "G4VEnergyLossProcess.hh"
00061 #include "G4LossTableManager.hh"
00062 #include "G4PhysicsTable.hh"
00063 #include "G4ThreeVector.hh"
00064 #include <vector>
00065
00066 class G4ParticleChangeForMSC;
00067
00068 class G4VMscModel : public G4VEmModel
00069 {
00070
00071 public:
00072
00073 G4VMscModel(const G4String& nam);
00074
00075 virtual ~G4VMscModel();
00076
00077 virtual G4double ComputeTruePathLengthLimit(const G4Track& track,
00078 G4double& stepLimit);
00079
00080 virtual G4double ComputeGeomPathLength(G4double truePathLength);
00081
00082 virtual G4double ComputeTrueStepLength(G4double geomPathLength);
00083
00084 virtual G4ThreeVector& SampleScattering(const G4ThreeVector&,
00085 G4double safety);
00086
00087
00088 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
00089 const G4MaterialCutsCouple*,
00090 const G4DynamicParticle*,
00091 G4double tmin,
00092 G4double tmax);
00093
00094
00095
00096
00097
00098 inline void SetStepLimitType(G4MscStepLimitType);
00099
00100 inline void SetLateralDisplasmentFlag(G4bool val);
00101
00102 inline void SetRangeFactor(G4double);
00103
00104 inline void SetGeomFactor(G4double);
00105
00106 inline void SetSkin(G4double);
00107
00108 inline void SetSampleZ(G4bool);
00109
00110
00111
00112
00113
00114 inline G4VEnergyLossProcess* GetIonisation() const;
00115
00116 inline void SetIonisation(G4VEnergyLossProcess*,
00117 const G4ParticleDefinition* part);
00118
00119
00120
00121
00122
00123 protected:
00124
00125
00126
00127 G4ParticleChangeForMSC*
00128 GetParticleChangeForMSC(const G4ParticleDefinition* p = 0);
00129
00130
00131 inline G4double ConvertTrueToGeom(G4double& tLength, G4double& gLength);
00132
00133 public:
00134
00135
00136 inline G4double ComputeSafety(const G4ThreeVector& position, G4double limit);
00137
00138
00139 inline G4double ComputeGeomLimit(const G4Track&, G4double& presafety,
00140 G4double limit);
00141
00142 inline G4double GetDEDX(const G4ParticleDefinition* part,
00143 G4double kineticEnergy,
00144 const G4MaterialCutsCouple* couple);
00145
00146 inline G4double GetRange(const G4ParticleDefinition* part,
00147 G4double kineticEnergy,
00148 const G4MaterialCutsCouple* couple);
00149
00150 inline G4double GetEnergy(const G4ParticleDefinition* part,
00151 G4double range,
00152 const G4MaterialCutsCouple* couple);
00153
00154
00155 inline
00156 G4double GetTransportMeanFreePath(const G4ParticleDefinition* part,
00157 G4double kinEnergy);
00158
00159 private:
00160
00161
00162 G4VMscModel & operator=(const G4VMscModel &right);
00163 G4VMscModel(const G4VMscModel&);
00164
00165 G4SafetyHelper* safetyHelper;
00166 G4VEnergyLossProcess* ionisation;
00167 const G4ParticleDefinition* currentPart;
00168 G4LossTableManager* man;
00169
00170 G4double dedx;
00171 G4double localtkin;
00172 G4double localrange;
00173
00174 protected:
00175
00176 G4double facrange;
00177 G4double facgeom;
00178 G4double facsafety;
00179 G4double skin;
00180 G4double dtrl;
00181 G4double lambdalimit;
00182 G4double geomMin;
00183 G4double geomMax;
00184
00185 G4ThreeVector fDisplacement;
00186 G4MscStepLimitType steppingAlgorithm;
00187
00188 G4bool samplez;
00189 G4bool latDisplasment;
00190
00191 };
00192
00193
00194
00195
00196 inline void G4VMscModel::SetLateralDisplasmentFlag(G4bool val)
00197 {
00198 latDisplasment = val;
00199 }
00200
00201
00202
00203 inline void G4VMscModel::SetSkin(G4double val)
00204 {
00205 skin = val;
00206 }
00207
00208
00209
00210 inline void G4VMscModel::SetRangeFactor(G4double val)
00211 {
00212 facrange = val;
00213 }
00214
00215
00216
00217 inline void G4VMscModel::SetGeomFactor(G4double val)
00218 {
00219 facgeom = val;
00220 }
00221
00222
00223
00224 inline void G4VMscModel::SetStepLimitType(G4MscStepLimitType val)
00225 {
00226 steppingAlgorithm = val;
00227 }
00228
00229
00230
00231 inline void G4VMscModel::SetSampleZ(G4bool val)
00232 {
00233 samplez = val;
00234 }
00235
00236
00237
00238 inline G4double G4VMscModel::ComputeSafety(const G4ThreeVector& position,
00239 G4double)
00240 {
00241 return safetyHelper->ComputeSafety(position);
00242 }
00243
00244
00245
00246 inline G4double G4VMscModel::ConvertTrueToGeom(G4double& tlength,
00247 G4double& glength)
00248 {
00249 glength = ComputeGeomPathLength(tlength);
00250
00251 return tlength;
00252 }
00253
00254
00255
00256 inline G4double G4VMscModel::ComputeGeomLimit(const G4Track& track,
00257 G4double& presafety,
00258 G4double limit)
00259 {
00260 G4double res = geomMax;
00261 if(track.GetVolume() != safetyHelper->GetWorldVolume()) {
00262 res = safetyHelper->CheckNextStep(
00263 track.GetStep()->GetPreStepPoint()->GetPosition(),
00264 track.GetMomentumDirection(),
00265 limit, presafety);
00266 }
00267 return res;
00268 }
00269
00270
00271
00272 inline G4double
00273 G4VMscModel::GetDEDX(const G4ParticleDefinition* part,
00274 G4double kinEnergy, const G4MaterialCutsCouple* couple)
00275 {
00276 G4double x;
00277 if(ionisation) { x = ionisation->GetDEDX(kinEnergy, couple); }
00278 else {
00279 G4double q = part->GetPDGCharge()/CLHEP::eplus;
00280 x = dedx*q*q;
00281 }
00282 return x;
00283 }
00284
00285
00286
00287 inline G4double
00288 G4VMscModel::GetRange(const G4ParticleDefinition* part,
00289 G4double kinEnergy, const G4MaterialCutsCouple* couple)
00290 {
00291 if(ionisation) {
00292 localrange = ionisation->GetRangeForLoss(kinEnergy, couple);
00293 } else {
00294 G4double q = part->GetPDGCharge()/CLHEP::eplus;
00295 localrange = kinEnergy/(dedx*q*q*couple->GetMaterial()->GetDensity());
00296 localtkin = kinEnergy;
00297 }
00298 return localrange;
00299 }
00300
00301
00302
00303 inline G4double
00304 G4VMscModel::GetEnergy(const G4ParticleDefinition* part,
00305 G4double range, const G4MaterialCutsCouple* couple)
00306 {
00307 G4double e;
00308 if(ionisation) { e = ionisation->GetKineticEnergy(range, couple); }
00309 else {
00310 e = localtkin;
00311 if(localrange > range) {
00312 G4double q = part->GetPDGCharge()/CLHEP::eplus;
00313 e -= (localrange - range)*dedx*q*q*couple->GetMaterial()->GetDensity();
00314 }
00315 }
00316 return e;
00317 }
00318
00319 inline G4VEnergyLossProcess* G4VMscModel::GetIonisation() const
00320 {
00321 return ionisation;
00322 }
00323
00324 inline void G4VMscModel::SetIonisation(G4VEnergyLossProcess* p,
00325 const G4ParticleDefinition* part)
00326 {
00327 ionisation = p;
00328 currentPart = part;
00329 }
00330
00331 inline G4double
00332 G4VMscModel::GetTransportMeanFreePath(const G4ParticleDefinition* part,
00333 G4double ekin)
00334 {
00335 G4double x;
00336 if(xSectionTable) {
00337 G4int idx = CurrentCouple()->GetIndex();
00338 x = (*xSectionTable)[(*theDensityIdx)[idx]]->Value(ekin)
00339 *(*theDensityFactor)[idx]/(ekin*ekin);
00340 } else {
00341 x = CrossSectionPerVolume(CurrentCouple()->GetMaterial(), part, ekin,
00342 0.0, DBL_MAX);
00343 }
00344 if(0.0 >= x) { x = DBL_MAX; }
00345 else { x = 1.0/x; }
00346 return x;
00347 }
00348
00349
00350
00351 #endif
00352