File indexing completed on 2025-01-18 09:59:05
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 #ifndef G4SeltzerBergerModel_h
0059 #define G4SeltzerBergerModel_h 1
0060
0061 #include "G4VEmModel.hh"
0062 #include "G4eBremsstrahlungRelModel.hh"
0063 #include "globals.hh"
0064
0065 class G4Physics2DVector;
0066 class G4SBBremTable;
0067 class G4ParticleChangeForLoss;
0068
0069 class G4SeltzerBergerModel : public G4VEmModel
0070 {
0071
0072 public:
0073
0074 explicit G4SeltzerBergerModel(const G4ParticleDefinition* p = nullptr,
0075 const G4String& nam = "eBremSB");
0076
0077 ~G4SeltzerBergerModel() override;
0078
0079 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
0080
0081 void InitialiseLocal(const G4ParticleDefinition*,
0082 G4VEmModel* masterModel) override;
0083
0084 G4double ComputeDEDXPerVolume(const G4Material*,
0085 const G4ParticleDefinition*,
0086 G4double ekin,
0087 G4double cutEnergy) override;
0088
0089 G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
0090 G4double ekin,
0091 G4double zet,
0092 G4double,
0093 G4double cutEnergy,
0094 G4double maxEnergy = DBL_MAX) override;
0095
0096 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
0097 const G4MaterialCutsCouple*,
0098 const G4DynamicParticle*,
0099 G4double cutEnergy,
0100 G4double maxEnergy) override;
0101
0102 void SetupForMaterial(const G4ParticleDefinition*,
0103 const G4Material*, G4double) override;
0104
0105 G4double MinPrimaryEnergy(const G4Material*,
0106 const G4ParticleDefinition*,
0107 G4double cutEnergy) override;
0108
0109 inline void SetBicubicInterpolationFlag(G4bool val)
0110 { fIsUseBicubicInterpolation = val; };
0111
0112
0113 G4SeltzerBergerModel & operator=(const G4SeltzerBergerModel &right) = delete;
0114 G4SeltzerBergerModel(const G4SeltzerBergerModel&) = delete;
0115
0116 private:
0117
0118 void SetParticle(const G4ParticleDefinition* p);
0119
0120 void ReadData(G4int Z);
0121
0122 G4double ComputeBremLoss(G4double cutEnergy);
0123
0124 G4double ComputeXSectionPerAtom(G4double cutEnergy);
0125
0126 G4double ComputeDXSectionPerAtom(G4double gammaEnergy);
0127
0128 G4double SampleEnergyTransfer(const G4double kineticEnergy,
0129 const G4double logKineticEnergy,
0130 const G4double cut,
0131 const G4double emax);
0132
0133 protected:
0134
0135 G4ParticleChangeForLoss* fParticleChange{nullptr};
0136
0137 private:
0138
0139 static constexpr G4int gMaxZet{101};
0140 static constexpr G4double gExpNumLimit{-12.};
0141 static G4double gYLimitData[gMaxZet];
0142 static G4Physics2DVector* gSBDCSData[gMaxZet];
0143 static G4SBBremTable* gSBSamplingTable;
0144 static const G4double gBremFactor;
0145 static const G4double gMigdalConstant;
0146
0147 G4bool fIsUseBicubicInterpolation{false};
0148 G4bool fIsUseSamplingTables{true};
0149 G4bool fIsElectron{true};
0150 G4bool fIsScatOffElectron{false};
0151 G4bool isInitializer{false};
0152
0153 G4int fCurrentIZ{0};
0154 G4int fNumWarnings{0};
0155
0156 const G4ParticleDefinition* fPrimaryParticle{nullptr};
0157 G4ParticleDefinition* fGammaParticle;
0158
0159
0160 G4double fPrimaryKinEnergy{0.};
0161 G4double fPrimaryTotalEnergy{0.};
0162 G4double fDensityFactor{0.};
0163 G4double fDensityCorr{0.};
0164 G4double fLowestKinEnergy;
0165
0166 std::size_t fIndx{0};
0167 std::size_t fIndy{0};
0168 };
0169
0170 #endif