File indexing completed on 2025-02-23 09:21:03
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
0054
0055
0056
0057
0058
0059 #ifndef G4ScreenedNuclearRecoil_h
0060 #define G4ScreenedNuclearRecoil_h 1
0061
0062 #include "c2_function.hh"
0063
0064 #include "G4ParticleChange.hh"
0065 #include "G4PhysicalConstants.hh"
0066 #include "G4SystemOfUnits.hh"
0067 #include "G4VDiscreteProcess.hh"
0068 #include "globals.hh"
0069
0070 #include <map>
0071 #include <vector>
0072
0073 class G4VNIELPartition;
0074
0075 typedef c2_const_ptr<G4double> G4_c2_const_ptr;
0076 typedef c2_ptr<G4double> G4_c2_ptr;
0077 typedef c2_function<G4double> G4_c2_function;
0078
0079 typedef struct G4ScreeningTables
0080 {
0081 G4double z1, z2, m1, m2, au, emin;
0082 G4_c2_const_ptr EMphiData;
0083 } G4ScreeningTables;
0084
0085
0086 class G4ScreenedCoulombCrossSectionInfo
0087 {
0088 public:
0089 G4ScreenedCoulombCrossSectionInfo() {}
0090 ~G4ScreenedCoulombCrossSectionInfo() {}
0091
0092 static const char* CVSHeaderVers()
0093 {
0094 return "G4ScreenedNuclearRecoil.hh,v 1.24 2008/05/01 19:58:59 marcus Exp GEANT4 tag ";
0095 }
0096 static const char* CVSFileVers();
0097 };
0098
0099
0100 class G4ScreenedCoulombCrossSection : public G4ScreenedCoulombCrossSectionInfo
0101 {
0102 public:
0103 G4ScreenedCoulombCrossSection() : verbosity(1) {}
0104 G4ScreenedCoulombCrossSection(const G4ScreenedCoulombCrossSection& src)
0105 : G4ScreenedCoulombCrossSectionInfo(), verbosity(src.verbosity)
0106 {}
0107 virtual ~G4ScreenedCoulombCrossSection();
0108
0109 typedef std::map<G4int, G4ScreeningTables> ScreeningMap;
0110
0111
0112 typedef std::map<G4int, class G4ParticleDefinition*> ParticleCache;
0113
0114
0115
0116 virtual void LoadData(G4String screeningKey, G4int z1, G4double m1, G4double recoilCutoff) = 0;
0117
0118
0119
0120 void BuildMFPTables(void);
0121
0122
0123 virtual G4ScreenedCoulombCrossSection* create() = 0;
0124
0125 const G4ScreeningTables* GetScreening(G4int Z) { return &(screeningData[Z]); }
0126 void SetVerbosity(G4int v) { verbosity = v; }
0127
0128
0129 G4ParticleDefinition* SelectRandomUnweightedTarget(const G4MaterialCutsCouple* couple);
0130
0131 enum
0132 {
0133 nMassMapElements = 116
0134 };
0135
0136 G4double standardmass(G4int z1) { return z1 <= nMassMapElements ? massmap[z1] : 2.5 * z1; }
0137
0138
0139 const G4_c2_function* operator[](G4int materialIndex)
0140 {
0141 return MFPTables.find(materialIndex) != MFPTables.end() ? &(MFPTables[materialIndex].get())
0142 : (G4_c2_function*)0;
0143 }
0144
0145 protected:
0146 ScreeningMap screeningData;
0147 ParticleCache targetMap;
0148 G4int verbosity;
0149 std::map<G4int, G4_c2_const_ptr> sigmaMap;
0150
0151 std::map<G4int, G4_c2_const_ptr> MFPTables;
0152
0153 private:
0154 static const G4double massmap[nMassMapElements + 1];
0155 };
0156
0157 typedef struct G4CoulombKinematicsInfo
0158 {
0159 G4double impactParameter;
0160 G4ScreenedCoulombCrossSection* crossSection;
0161 G4double a1, a2, sinTheta, cosTheta, sinZeta, cosZeta, eRecoil;
0162 G4ParticleDefinition* recoilIon;
0163 const G4Material* targetMaterial;
0164 } G4CoulombKinematicsInfo;
0165
0166 class G4ScreenedCollisionStage
0167 {
0168 public:
0169 virtual void DoCollisionStep(class G4ScreenedNuclearRecoil* master, const class G4Track& aTrack,
0170 const class G4Step& aStep) = 0;
0171 virtual ~G4ScreenedCollisionStage() {}
0172 };
0173
0174 class G4ScreenedCoulombClassicalKinematics : public G4ScreenedCoulombCrossSectionInfo,
0175 public G4ScreenedCollisionStage
0176 {
0177 public:
0178 G4ScreenedCoulombClassicalKinematics();
0179 virtual void DoCollisionStep(class G4ScreenedNuclearRecoil* master, const class G4Track& aTrack,
0180 const class G4Step& aStep);
0181
0182 G4bool DoScreeningComputation(class G4ScreenedNuclearRecoil* master,
0183 const G4ScreeningTables* screen, G4double eps, G4double beta);
0184
0185 virtual ~G4ScreenedCoulombClassicalKinematics() {}
0186
0187 protected:
0188
0189 c2_const_plugin_function_p<G4double>& phifunc;
0190 c2_linear_p<G4double>& xovereps;
0191 G4_c2_ptr diff;
0192 };
0193
0194 class G4SingleScatter : public G4ScreenedCoulombCrossSectionInfo, public G4ScreenedCollisionStage
0195 {
0196 public:
0197 G4SingleScatter() {}
0198 virtual void DoCollisionStep(class G4ScreenedNuclearRecoil* master, const class G4Track& aTrack,
0199 const class G4Step& aStep);
0200 virtual ~G4SingleScatter() {}
0201 };
0202
0203
0204
0205
0206
0207
0208 class G4ScreenedNuclearRecoil : public G4ScreenedCoulombCrossSectionInfo, public G4VDiscreteProcess
0209 {
0210 public:
0211 friend class G4ScreenedCollisionStage;
0212
0213
0214
0215
0216
0217
0218
0219
0220
0221
0222
0223
0224
0225
0226
0227
0228
0229
0230
0231
0232
0233 G4ScreenedNuclearRecoil(const G4String& processName = "ScreenedElastic",
0234 const G4String& ScreeningKey = "zbl", G4bool GenerateRecoils = 1,
0235 G4double RecoilCutoff = 100.0 * eV, G4double PhysicsCutoff = 10.0 * eV);
0236
0237
0238 virtual ~G4ScreenedNuclearRecoil();
0239
0240
0241 virtual G4double GetMeanFreePath(const G4Track&, G4double, G4ForceCondition*);
0242
0243
0244 virtual G4VParticleChange* PostStepDoIt(const G4Track& aTrack, const G4Step& aStep);
0245
0246
0247 virtual G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
0248
0249
0250 virtual void BuildPhysicsTable(const G4ParticleDefinition& aParticleType);
0251
0252
0253 virtual void DumpPhysicsTable(const G4ParticleDefinition& aParticleType);
0254
0255
0256
0257
0258
0259
0260 virtual G4bool CheckNuclearCollision(G4double A, G4double A1,
0261 G4double apsis);
0262
0263 virtual G4ScreenedCoulombCrossSection* GetNewCrossSectionHandler(void);
0264
0265
0266 G4double GetNIEL() const { return NIEL; }
0267
0268
0269 void ResetTables();
0270
0271
0272
0273
0274
0275
0276
0277
0278
0279
0280
0281
0282
0283 void SetMaxEnergyForScattering(G4double energy) { processMaxEnergy = energy; }
0284
0285
0286 std::string GetScreeningKey() const { return screeningKey; }
0287
0288
0289
0290
0291
0292 void AllowEnergyDeposition(G4bool flag) { registerDepositedEnergy = flag; }
0293
0294
0295 G4bool GetAllowEnergyDeposition() const { return registerDepositedEnergy; }
0296
0297
0298
0299
0300
0301
0302
0303
0304 void EnableRecoils(G4bool flag) { generateRecoils = flag; }
0305
0306
0307 G4bool GetEnableRecoils() const { return generateRecoils; }
0308
0309
0310
0311
0312
0313
0314
0315 void SetMFPScaling(G4double scale) { MFPScale = scale; }
0316
0317
0318 G4double GetMFPScaling() const { return MFPScale; }
0319
0320
0321
0322
0323
0324
0325
0326
0327 void AvoidNuclearReactions(G4bool flag) { avoidReactions = flag; }
0328
0329
0330 G4bool GetAvoidNuclearReactions() const { return avoidReactions; }
0331
0332
0333
0334
0335
0336
0337 void SetRecoilCutoff(G4double energy) { recoilCutoff = energy; }
0338
0339
0340 G4double GetRecoilCutoff() const { return recoilCutoff; }
0341
0342
0343
0344
0345
0346 void SetPhysicsCutoff(G4double energy)
0347 {
0348 physicsCutoff = energy;
0349 ResetTables();
0350 }
0351
0352
0353 G4double GetPhysicsCutoff() const { return physicsCutoff; }
0354
0355
0356
0357
0358 void SetNIELPartitionFunction(const G4VNIELPartition* part);
0359
0360
0361
0362
0363
0364
0365
0366
0367 void SetCrossSectionHardening(G4double fraction, G4double HardeningFactor)
0368 {
0369 hardeningFraction = fraction;
0370 hardeningFactor = HardeningFactor;
0371 }
0372
0373
0374 G4double GetHardeningFraction() const { return hardeningFraction; }
0375
0376
0377 G4double GetHardeningFactor() const { return hardeningFactor; }
0378
0379
0380 G4double GetCurrentInteractionLength() const { return currentInteractionLength; }
0381
0382
0383
0384
0385
0386 void SetExternalCrossSectionHandler(G4ScreenedCoulombCrossSection* cs)
0387 {
0388 externalCrossSectionConstructor = cs;
0389 }
0390
0391
0392 G4int GetVerboseLevel() const { return verboseLevel; }
0393
0394 std::map<G4int, G4ScreenedCoulombCrossSection*>& GetCrossSectionHandlers()
0395 {
0396 return crossSectionHandlers;
0397 }
0398
0399 void ClearStages(void);
0400
0401 void AddStage(G4ScreenedCollisionStage* stage) { collisionStages.push_back(stage); }
0402
0403 G4CoulombKinematicsInfo& GetKinematics() { return kinematics; }
0404
0405 void SetValidCollision(G4bool flag) { validCollision = flag; }
0406
0407 G4bool GetValidCollision() const { return validCollision; }
0408
0409
0410
0411 class G4ParticleChange& GetParticleChange()
0412 {
0413 return static_cast<G4ParticleChange&>(*pParticleChange);
0414 }
0415
0416
0417
0418
0419 void DepositEnergy(G4int z1, G4double a1, const G4Material* material, G4double energy);
0420
0421 protected:
0422
0423 G4double highEnergyLimit;
0424
0425
0426 G4double lowEnergyLimit;
0427
0428
0429
0430 G4double processMaxEnergy;
0431 G4String screeningKey;
0432 G4bool generateRecoils, avoidReactions;
0433 G4double recoilCutoff, physicsCutoff;
0434 G4bool registerDepositedEnergy;
0435 G4double IonizingLoss, NIEL;
0436 G4double MFPScale;
0437 G4double hardeningFraction, hardeningFactor;
0438
0439 G4ScreenedCoulombCrossSection* externalCrossSectionConstructor;
0440 std::vector<G4ScreenedCollisionStage*> collisionStages;
0441
0442 std::map<G4int, G4ScreenedCoulombCrossSection*> crossSectionHandlers;
0443
0444 G4bool validCollision;
0445 G4CoulombKinematicsInfo kinematics;
0446 const G4VNIELPartition* NIELPartitionFunction;
0447 };
0448
0449
0450
0451
0452 class G4NativeScreenedCoulombCrossSection : public G4ScreenedCoulombCrossSection
0453 {
0454 public:
0455 G4NativeScreenedCoulombCrossSection();
0456
0457 G4NativeScreenedCoulombCrossSection(const G4NativeScreenedCoulombCrossSection& src)
0458 : G4ScreenedCoulombCrossSection(src), phiMap(src.phiMap)
0459 {}
0460
0461 G4NativeScreenedCoulombCrossSection(const G4ScreenedCoulombCrossSection& src)
0462 : G4ScreenedCoulombCrossSection(src)
0463 {}
0464
0465 virtual ~G4NativeScreenedCoulombCrossSection();
0466
0467 virtual void LoadData(G4String screeningKey, G4int z1, G4double m1, G4double recoilCutoff);
0468
0469 virtual G4ScreenedCoulombCrossSection* create()
0470 {
0471 return new G4NativeScreenedCoulombCrossSection(*this);
0472 }
0473
0474
0475 std::vector<G4String> GetScreeningKeys() const;
0476
0477 typedef G4_c2_function& (*ScreeningFunc)(G4int z1, G4int z2, size_t nPoints, G4double rMax,
0478 G4double* au);
0479
0480 void AddScreeningFunction(G4String name, ScreeningFunc fn) { phiMap[name] = fn; }
0481
0482 private:
0483
0484 std::map<std::string, ScreeningFunc> phiMap;
0485 };
0486
0487 G4_c2_function& ZBLScreening(G4int z1, G4int z2, size_t npoints, G4double rMax, G4double* auval);
0488
0489 G4_c2_function& MoliereScreening(G4int z1, G4int z2, size_t npoints, G4double rMax,
0490 G4double* auval);
0491
0492 G4_c2_function& LJScreening(G4int z1, G4int z2, size_t npoints, G4double rMax, G4double* auval);
0493
0494 G4_c2_function& LJZBLScreening(G4int z1, G4int z2, size_t npoints, G4double rMax, G4double* auval);
0495
0496 #endif