File indexing completed on 2025-01-18 09:58:04
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 InitialiseOnFly(const G4int Z);
0095
0096 void PrepareCache(const G4Material*);
0097
0098 G4double IsoCrossSection(const G4double kinE, const G4double loge,
0099 const G4int Z, const G4int A,
0100 const G4double temperature);
0101
0102 inline G4double GetCrossSection(const G4int Z, const G4int A);
0103
0104 inline G4bool CheckCache(const G4int Z);
0105
0106 const G4ParticleDefinition* fParticle;
0107 const G4ParticleDefinition* fNeutron;
0108 G4ParticleHPManager* fManagerHP;
0109
0110 const G4double emax;
0111 const G4double emaxT;
0112 const G4double elimit;
0113 const G4double logElimit;
0114 G4LorentzVector fLV;
0115 G4ThreeVector fBoost;
0116
0117 const G4int minZ;
0118 const G4int maxZ;
0119
0120 G4int binSearch{2};
0121 std::size_t index{0};
0122 G4bool fPrinted{false};
0123
0124
0125 const G4Material* fCurrentMat{nullptr};
0126 std::vector<std::pair<G4int, G4int> > fZA;
0127 std::vector<G4double> fIsoXS;
0128 std::vector<G4double> fTemp;
0129
0130 const G4String fDataName;
0131 const G4String fDataDirectory;
0132 G4ElementData* fData{nullptr};
0133 };
0134
0135 inline G4double
0136 G4CrossSectionHP::GetCrossSection(const G4int Z, const G4int A)
0137 {
0138 G4double res = 0.0;
0139 for (std::size_t i=0; i<fZA.size(); ++i) {
0140 if (Z == fZA[i].first && A == fZA[i].second) {
0141 res = fIsoXS[i];
0142 break;
0143 }
0144 }
0145 return res;
0146 }
0147
0148 inline G4bool G4CrossSectionHP::CheckCache(const G4int Z)
0149 {
0150 G4bool res = false;
0151 for (auto const & p : fZA) {
0152 if (Z == p.first) {
0153 res = true;
0154 break;
0155 }
0156 }
0157 return res;
0158 }
0159
0160 inline G4double G4CrossSectionHP::GetMaxHPEnergy() const
0161 {
0162 return emax;
0163 }
0164
0165 inline void G4CrossSectionHP::SetBinSearch(G4int n)
0166 {
0167 if (n > 0) { binSearch = n; }
0168 }
0169
0170 #endif