File indexing completed on 2025-09-17 08:58:13
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 #ifndef G4CrossSectionHP_h
0034 #define G4CrossSectionHP_h 1
0035
0036 #include "G4VCrossSectionDataSet.hh"
0037 #include "globals.hh"
0038 #include "G4ElementData.hh"
0039 #include "G4PhysicsVector.hh"
0040 #include "G4LorentzVector.hh"
0041 #include "G4ThreeVector.hh"
0042 #include <vector>
0043
0044 class G4DynamicParticle;
0045 class G4ParticleDefinition;
0046 class G4ParticleHPManager;
0047 class G4Element;
0048 class G4Material;
0049
0050 class G4CrossSectionHP : public G4VCrossSectionDataSet
0051 {
0052 public:
0053 explicit G4CrossSectionHP(const G4ParticleDefinition*,
0054 const G4String& nameData,
0055 const G4String& nameDir, G4double emaxHP,
0056 G4int zmin, G4int zmax);
0057
0058 ~G4CrossSectionHP() override = default;
0059
0060 G4bool IsIsoApplicable(const G4DynamicParticle*, G4int Z, G4int A,
0061 const G4Element*, const G4Material*) override;
0062
0063 G4double GetIsoCrossSection(const G4DynamicParticle*, G4int Z, G4int A,
0064 const G4Isotope*, const G4Element*,
0065 const G4Material*) override;
0066
0067 G4double ComputeIsoCrossSection(G4double kinEnergy, G4double loge,
0068 const G4ParticleDefinition*,
0069 G4int Z, G4int A,
0070 const G4Isotope* iso,
0071 const G4Element* elm,
0072 const G4Material* mat) override;
0073
0074 const G4Isotope* SelectIsotope(const G4Element*, G4double kinEnergy,
0075 G4double logE) override;
0076
0077 void BuildPhysicsTable(const G4ParticleDefinition&) override;
0078
0079 void DumpPhysicsTable(const G4ParticleDefinition&) override;
0080
0081 G4CrossSectionHP & operator=(const G4CrossSectionHP &right) = delete;
0082 G4CrossSectionHP(const G4CrossSectionHP&) = delete;
0083
0084 protected:
0085
0086 inline G4double GetMaxHPEnergy() const;
0087
0088 inline void SetBinSearch(G4int n);
0089
0090 private:
0091
0092 void Initialise(const G4int Z);
0093
0094 void PrepareCache(const G4Material*);
0095
0096 G4double IsoCrossSection(const G4double kinE, const G4double loge,
0097 const G4int Z, const G4int A,
0098 const G4double temperature);
0099
0100 inline G4double GetCrossSection(const G4int Z, const G4int A);
0101
0102 inline G4bool CheckCache(const G4int Z);
0103
0104 const G4ParticleDefinition* fParticle;
0105 const G4ParticleDefinition* fNeutron;
0106 G4ParticleHPManager* fManagerHP;
0107
0108 const G4double emax;
0109 const G4double emaxT;
0110 const G4double elimit;
0111 const G4double logElimit;
0112 G4LorentzVector fLV;
0113 G4ThreeVector fBoost;
0114
0115 const G4int minZ;
0116 const G4int maxZ;
0117
0118 G4int binSearch{2};
0119 std::size_t index{0};
0120 G4bool fPrinted{false};
0121
0122
0123 const G4Material* fCurrentMat{nullptr};
0124 std::vector<std::pair<G4int, G4int> > fZA;
0125 std::vector<G4double> fIsoXS;
0126 std::vector<G4double> fTemp;
0127
0128 const G4String fDataName;
0129 const G4String fDataDirectory;
0130 G4ElementData* fData{nullptr};
0131 };
0132
0133 inline G4double
0134 G4CrossSectionHP::GetCrossSection(const G4int Z, const G4int A)
0135 {
0136 G4double res = 0.0;
0137 for (std::size_t i=0; i<fZA.size(); ++i) {
0138 if (Z == fZA[i].first && A == fZA[i].second) {
0139 res = fIsoXS[i];
0140 break;
0141 }
0142 }
0143 return res;
0144 }
0145
0146 inline G4bool G4CrossSectionHP::CheckCache(const G4int Z)
0147 {
0148 G4bool res = false;
0149 for (auto const & p : fZA) {
0150 if (Z == p.first) {
0151 res = true;
0152 break;
0153 }
0154 }
0155 return res;
0156 }
0157
0158 inline G4double G4CrossSectionHP::GetMaxHPEnergy() const
0159 {
0160 return emax;
0161 }
0162
0163 inline void G4CrossSectionHP::SetBinSearch(G4int n)
0164 {
0165 if (n > 0) { binSearch = n; }
0166 }
0167
0168 #endif