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
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070 #ifndef G4VMultipleScattering_h
00071 #define G4VMultipleScattering_h 1
00072
00073 #include "G4VContinuousDiscreteProcess.hh"
00074 #include "G4LossTableManager.hh"
00075 #include "globals.hh"
00076 #include "G4Material.hh"
00077 #include "G4MaterialCutsCouple.hh"
00078 #include "G4ParticleChangeForMSC.hh"
00079 #include "G4Track.hh"
00080 #include "G4Step.hh"
00081 #include "G4EmModelManager.hh"
00082 #include "G4VEmModel.hh"
00083 #include "G4MscStepLimitType.hh"
00084
00085 class G4ParticleDefinition;
00086 class G4DataVector;
00087 class G4PhysicsTable;
00088 class G4PhysicsVector;
00089
00090
00091
00092 class G4VMultipleScattering : public G4VContinuousDiscreteProcess
00093 {
00094 public:
00095
00096 G4VMultipleScattering(const G4String& name = "msc",
00097 G4ProcessType type = fElectromagnetic);
00098
00099 virtual ~G4VMultipleScattering();
00100
00101
00102
00103
00104
00105 virtual G4bool IsApplicable(const G4ParticleDefinition& p) = 0;
00106
00107
00108 virtual void PrintInfo() = 0;
00109
00110 protected:
00111
00112 virtual void InitialiseProcess(const G4ParticleDefinition*) = 0;
00113
00114
00115
00116
00117 public:
00118
00119
00120
00121
00122
00123
00124 void PreparePhysicsTable(const G4ParticleDefinition&);
00125
00126
00127 void BuildPhysicsTable(const G4ParticleDefinition&);
00128
00129
00130 void PrintInfoDefinition();
00131
00132 G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&);
00133
00134 G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&);
00135
00136
00137
00138 G4bool StorePhysicsTable(const G4ParticleDefinition*,
00139 const G4String& directory,
00140 G4bool ascii = false);
00141
00142
00143
00144
00145
00146
00147 G4bool RetrievePhysicsTable(const G4ParticleDefinition*,
00148 const G4String& directory,
00149 G4bool ascii);
00150
00151
00152
00153
00154
00155
00156
00157
00158 G4double AlongStepGetPhysicalInteractionLength(
00159 const G4Track&,
00160 G4double previousStepSize,
00161 G4double currentMinimalStep,
00162 G4double& currentSafety,
00163 G4GPILSelection* selection);
00164
00165
00166
00167 G4double PostStepGetPhysicalInteractionLength(
00168 const G4Track&,
00169 G4double previousStepSize,
00170 G4ForceCondition* condition);
00171
00172
00173 inline G4double ContinuousStepLimit(const G4Track& track,
00174 G4double previousStepSize,
00175 G4double currentMinimalStep,
00176 G4double& currentSafety);
00177
00178
00179
00180
00181
00182
00183 G4PhysicsVector* PhysicsVector(const G4MaterialCutsCouple*);
00184
00185 inline void SetBinning(G4int nbins);
00186 inline G4int Binning() const;
00187
00188 inline void SetMinKinEnergy(G4double e);
00189 inline G4double MinKinEnergy() const;
00190
00191
00192 inline void SetMaxKinEnergy(G4double e);
00193 inline G4double MaxKinEnergy() const;
00194
00195 inline void SetBuildLambdaTable(G4bool val);
00196
00197 inline G4PhysicsTable* LambdaTable() const;
00198
00199
00200
00201
00202
00203 inline const G4ParticleDefinition* Particle() const;
00204 inline void SetParticle(const G4ParticleDefinition*);
00205
00206
00207
00208
00209
00210 inline void AddEmModel(G4int, G4VEmModel*, const G4Region* region = 0);
00211
00212 inline G4VEmModel* SelectModelForMaterial(G4double kinEnergy,
00213 size_t& idxRegion) const;
00214
00215
00216 inline G4VEmModel* GetModelByIndex(G4int idx = 0);
00217
00218
00219
00220
00221
00222 inline void SetLateralDisplasmentFlag(G4bool val);
00223
00224
00225 inline void SetSkin(G4double val);
00226
00227
00228 inline void SetRangeFactor(G4double val);
00229
00230
00231 inline void SetGeomFactor(G4double val);
00232
00233
00234 inline void SetStepLimitType(G4MscStepLimitType val);
00235
00236
00237 protected:
00238
00239
00240 G4double GetMeanFreePath(const G4Track& track,
00241 G4double,
00242 G4ForceCondition* condition);
00243
00244
00245
00246
00247
00248
00249 G4double GetContinuousStepLimit(const G4Track& track,
00250 G4double previousStepSize,
00251 G4double currentMinimalStep,
00252 G4double& currentSafety);
00253
00254 inline G4double GetLambda(const G4ParticleDefinition* p, G4double& kineticEnergy);
00255
00256
00257 inline G4double GetMscContinuousStepLimit(const G4Track& track,
00258 G4double previousStepSize,
00259 G4double currentMinimalStep,
00260 G4double& currentSafety);
00261
00262 inline G4VEmModel* SelectModel(G4double kinEnergy);
00263
00264
00265 inline const G4MaterialCutsCouple* CurrentMaterialCutsCouple() const;
00266
00267
00268 inline void DefineMaterial(const G4MaterialCutsCouple* couple);
00269
00270
00271
00272
00273
00274
00275 inline G4double Skin() const;
00276
00277 inline G4double RangeFactor() const;
00278
00279 inline G4double GeomFactor() const;
00280
00281 inline G4MscStepLimitType StepLimitType() const;
00282
00283 inline G4bool LateralDisplasmentFlag() const;
00284
00285 private:
00286
00287
00288
00289 G4VMultipleScattering(G4VMultipleScattering &);
00290 G4VMultipleScattering & operator=(const G4VMultipleScattering &right);
00291
00292
00293
00294 protected:
00295
00296 G4GPILSelection valueGPILSelectionMSC;
00297 G4ParticleChangeForMSC fParticleChange;
00298
00299 private:
00300
00301 G4EmModelManager* modelManager;
00302 G4VEmModel* currentModel;
00303 G4PhysicsTable* theLambdaTable;
00304
00305
00306 const G4ParticleDefinition* firstParticle;
00307 const G4ParticleDefinition* currentParticle;
00308 const G4MaterialCutsCouple* currentCouple;
00309 size_t currentMaterialIndex;
00310
00311 G4int nBins;
00312
00313 G4MscStepLimitType stepLimit;
00314
00315 G4double minKinEnergy;
00316 G4double maxKinEnergy;
00317 G4double skin;
00318 G4double facrange;
00319 G4double facgeom;
00320
00321 G4bool latDisplasment;
00322 G4bool buildLambdaTable;
00323 };
00324
00325
00326
00327
00328 inline void G4VMultipleScattering::DefineMaterial(const G4MaterialCutsCouple* couple)
00329 {
00330 if(couple != currentCouple) {
00331 currentCouple = couple;
00332 currentMaterialIndex = couple->GetIndex();
00333 }
00334 }
00335
00336
00337
00338 inline G4double G4VMultipleScattering::GetMscContinuousStepLimit(
00339 const G4Track& track,
00340 G4double,
00341 G4double currentMinimalStep,
00342 G4double&)
00343 {
00344 G4double x = currentMinimalStep;
00345 G4double e = track.GetKineticEnergy();
00346 DefineMaterial(track.GetMaterialCutsCouple());
00347 currentModel = SelectModel(e);
00348 if(x > 0.0 && e > 0.0) {
00349 G4double tPathLength =
00350 currentModel->ComputeTruePathLengthLimit(track, theLambdaTable, x);
00351 if (tPathLength < x) valueGPILSelectionMSC = CandidateForSelection;
00352 x = currentModel->ComputeGeomPathLength(tPathLength);
00353
00354
00355
00356 }
00357 return x;
00358 }
00359
00360
00361
00362 inline G4double G4VMultipleScattering::ContinuousStepLimit(
00363 const G4Track& track,
00364 G4double previousStepSize,
00365 G4double currentMinimalStep,
00366 G4double& currentSafety)
00367 {
00368 return GetMscContinuousStepLimit(track,previousStepSize,currentMinimalStep,
00369 currentSafety);
00370 }
00371
00372
00373
00374 inline G4double G4VMultipleScattering::GetLambda(const G4ParticleDefinition* p, G4double& e)
00375 {
00376 G4double x;
00377 if(theLambdaTable) {
00378 G4bool b;
00379 x = ((*theLambdaTable)[currentMaterialIndex])->GetValue(e, b);
00380 } else {
00381 x = currentModel->CrossSection(currentCouple,p,e);
00382 }
00383 if(x > DBL_MIN) x = 1./x;
00384 else x = DBL_MAX;
00385 return x;
00386 }
00387
00388
00389
00390 inline G4VEmModel* G4VMultipleScattering::SelectModel(G4double kinEnergy)
00391 {
00392 return modelManager->SelectModel(kinEnergy, currentMaterialIndex);
00393 }
00394
00395
00396
00397 inline G4VEmModel* G4VMultipleScattering::SelectModelForMaterial(
00398 G4double kinEnergy, size_t& idxRegion) const
00399 {
00400 return modelManager->SelectModel(kinEnergy, idxRegion);
00401 }
00402
00403
00404
00405 inline void G4VMultipleScattering::SetBinning(G4int nbins)
00406 {
00407 nBins = nbins;
00408 }
00409
00410
00411
00412 inline G4int G4VMultipleScattering::Binning() const
00413 {
00414 return nBins;
00415 }
00416
00417
00418
00419 inline void G4VMultipleScattering::SetMinKinEnergy(G4double e)
00420 {
00421 minKinEnergy = e;
00422 }
00423
00424
00425
00426 inline G4double G4VMultipleScattering::MinKinEnergy() const
00427 {
00428 return minKinEnergy;
00429 }
00430
00431
00432
00433 inline void G4VMultipleScattering::SetMaxKinEnergy(G4double e)
00434 {
00435 maxKinEnergy = e;
00436 }
00437
00438
00439
00440 inline G4double G4VMultipleScattering::MaxKinEnergy() const
00441 {
00442 return maxKinEnergy;
00443 }
00444
00445
00446
00447 inline G4bool G4VMultipleScattering::LateralDisplasmentFlag() const
00448 {
00449 return latDisplasment;
00450 }
00451
00452
00453
00454 inline void G4VMultipleScattering::SetLateralDisplasmentFlag(G4bool val)
00455 {
00456 latDisplasment = val;
00457 }
00458
00459
00460
00461 inline G4double G4VMultipleScattering::Skin() const
00462 {
00463 return skin;
00464 }
00465
00466
00467
00468 inline void G4VMultipleScattering::SetSkin(G4double val)
00469 {
00470 if(val <= 0.99999) {
00471 skin = 0.0;
00472 stepLimit = fUseSafety;
00473 } else {
00474 skin = val;
00475 stepLimit = fUseDistanceToBoundary;
00476 }
00477 }
00478
00479
00480
00481 inline G4double G4VMultipleScattering::RangeFactor() const
00482 {
00483 return facrange;
00484 }
00485
00486
00487
00488 inline void G4VMultipleScattering::SetRangeFactor(G4double val)
00489 {
00490 if(val > 0.0) facrange = val;
00491 }
00492
00493
00494
00495 inline G4double G4VMultipleScattering::GeomFactor() const
00496 {
00497 return facgeom;
00498 }
00499
00500
00501
00502 inline void G4VMultipleScattering::SetGeomFactor(G4double val)
00503 {
00504 if(val > 0.0) facgeom = val;
00505 }
00506
00507
00508
00509 inline G4MscStepLimitType G4VMultipleScattering::StepLimitType() const
00510 {
00511 return stepLimit;
00512 }
00513
00514
00515
00516 inline void G4VMultipleScattering::SetStepLimitType(G4MscStepLimitType val)
00517 {
00518 stepLimit = val;
00519 if(val == fMinimal) {
00520 skin = 0;
00521 facrange = 0.2;
00522 } else if(val == fUseSafety) {
00523 skin = 0;
00524 }
00525 }
00526
00527
00528
00529 inline void G4VMultipleScattering::SetBuildLambdaTable(G4bool val)
00530 {
00531 buildLambdaTable = val;
00532 }
00533
00534
00535
00536 inline const G4ParticleDefinition* G4VMultipleScattering::Particle() const
00537 {
00538 return currentParticle;
00539 }
00540
00541
00542
00543 inline G4PhysicsTable* G4VMultipleScattering::LambdaTable() const
00544 {
00545 return theLambdaTable;
00546 }
00547
00548
00549
00550 inline
00551 const G4MaterialCutsCouple* G4VMultipleScattering::CurrentMaterialCutsCouple() const
00552 {
00553 return currentCouple;
00554 }
00555
00556
00557
00558 inline void G4VMultipleScattering::AddEmModel(G4int order, G4VEmModel* p,
00559 const G4Region* region)
00560 {
00561 G4VEmFluctuationModel* fm = 0;
00562 modelManager->AddEmModel(order, p, fm, region);
00563 if(p)p->SetParticleChange(pParticleChange);
00564 }
00565
00566
00567
00568 inline G4VEmModel* G4VMultipleScattering::GetModelByIndex(G4int idx)
00569 {
00570 return modelManager->GetModel(idx);
00571 }
00572
00573
00574
00575 #endif