File indexing completed on 2025-01-18 09:58:53
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 #ifndef G4ParticleInelasticXS_h
0043 #define G4ParticleInelasticXS_h 1
0044
0045 #include "G4VCrossSectionDataSet.hh"
0046 #include "globals.hh"
0047 #include "G4ElementData.hh"
0048 #include "G4PhysicsVector.hh"
0049 #include <vector>
0050
0051 class G4DynamicParticle;
0052 class G4ParticleDefinition;
0053 class G4Element;
0054 class G4VComponentCrossSection;
0055
0056 class G4ParticleInelasticXS final : public G4VCrossSectionDataSet
0057 {
0058 public:
0059
0060 explicit G4ParticleInelasticXS(const G4ParticleDefinition*);
0061
0062 ~G4ParticleInelasticXS() override = default;
0063
0064 G4bool IsElementApplicable(const G4DynamicParticle*, G4int Z,
0065 const G4Material* mat = nullptr) final;
0066
0067 G4bool IsIsoApplicable(const G4DynamicParticle*, G4int Z, G4int A,
0068 const G4Element* elm = nullptr,
0069 const G4Material* mat = nullptr) final;
0070
0071 G4double GetElementCrossSection(const G4DynamicParticle*, G4int Z,
0072 const G4Material* mat = nullptr) final;
0073
0074 G4double ComputeCrossSectionPerElement(G4double kinEnergy, G4double loge,
0075 const G4ParticleDefinition*,
0076 const G4Element*,
0077 const G4Material*) final;
0078
0079 G4double ComputeIsoCrossSection(G4double kinEnergy, G4double loge,
0080 const G4ParticleDefinition*,
0081 G4int Z, G4int A,
0082 const G4Isotope* iso,
0083 const G4Element* elm,
0084 const G4Material* mat) final;
0085
0086 G4double GetIsoCrossSection(const G4DynamicParticle*, G4int Z, G4int A,
0087 const G4Isotope* iso = nullptr,
0088 const G4Element* elm = nullptr,
0089 const G4Material* mat = nullptr) final;
0090
0091 const G4Isotope* SelectIsotope(const G4Element*,
0092 G4double kinEnergy, G4double logE) final;
0093
0094 void BuildPhysicsTable(const G4ParticleDefinition&) final;
0095
0096 void CrossSectionDescription(std::ostream&) const final;
0097
0098 G4double ElementCrossSection(G4double kinEnergy, G4double loge, G4int Z);
0099
0100 G4double IsoCrossSection(G4double ekin, G4double logE, G4int Z, G4int A);
0101
0102 G4ParticleInelasticXS & operator=
0103 (const G4ParticleInelasticXS &right) = delete;
0104 G4ParticleInelasticXS(const G4ParticleInelasticXS&) = delete;
0105
0106 private:
0107
0108 void Initialise(G4int Z, G4int idx);
0109
0110 void InitialiseOnFly(G4int Z);
0111
0112 void FindDirectoryPath();
0113
0114 inline const G4PhysicsVector* GetPhysicsVector(G4int Z);
0115
0116 G4PhysicsVector* RetrieveVector(std::ostringstream& in, G4bool warn);
0117
0118 G4VComponentCrossSection* highEnergyXsection = nullptr;
0119 const G4ParticleDefinition* particle;
0120
0121 std::vector<G4double> temp;
0122 G4double elimit;
0123
0124 G4int index{0};
0125 G4bool isInitializer{false};
0126
0127 static const G4int MAXZINELP = 93;
0128 static G4ElementData* data[5];
0129 static G4double coeff[MAXZINELP][5];
0130 static G4String gDataDirectory[5];
0131 };
0132
0133 inline
0134 const G4PhysicsVector* G4ParticleInelasticXS::GetPhysicsVector(G4int Z)
0135 {
0136 const G4PhysicsVector* pv = data[index]->GetElementData(Z);
0137 if (pv == nullptr) {
0138 InitialiseOnFly(Z);
0139 pv = data[index]->GetElementData(Z);
0140 }
0141 return pv;
0142 }
0143
0144 #endif