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