File indexing completed on 2025-09-18 09:15:19
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 #ifndef G4ParticleHPChannel_h
0040 #define G4ParticleHPChannel_h 1
0041
0042 #include "G4Element.hh"
0043 #include "G4HadProjectile.hh"
0044 #include "G4Material.hh"
0045 #include "G4ParticleHPCaptureFS.hh"
0046 #include "G4ParticleHPFinalState.hh"
0047 #include "G4ParticleHPIsoData.hh"
0048 #include "G4ParticleHPManager.hh"
0049 #include "G4ParticleHPVector.hh"
0050 #include "G4StableIsotopes.hh"
0051 #include "G4WendtFissionFragmentGenerator.hh"
0052 #include "globals.hh"
0053
0054 class G4ParticleDefinition;
0055
0056 class G4ParticleHPChannel
0057 {
0058 public:
0059
0060 G4ParticleHPChannel(G4ParticleDefinition* projectile = nullptr);
0061
0062 ~G4ParticleHPChannel();
0063
0064 G4double GetXsec(G4double energy) const;
0065
0066 G4double GetWeightedXsec(G4double energy, G4int isoNumber) const;
0067
0068 G4double GetFSCrossSection(G4double energy, G4int isoNumber) const;
0069
0070 G4bool IsActive(G4int isoNumber) const
0071 {
0072 return active[isoNumber];
0073 }
0074
0075 G4bool HasFSData(G4int isoNumber) const
0076 {
0077 return theFinalStates[isoNumber]->HasFSData();
0078 }
0079
0080 G4bool HasAnyData(G4int isoNumber) const
0081 {
0082 return theFinalStates[isoNumber]->HasAnyData();
0083 }
0084
0085 G4bool Register(G4ParticleHPFinalState* theFS);
0086
0087 void Init(G4Element* theElement, const G4String& dirName);
0088
0089 void Init(G4Element* theElement, const G4String& dirName,
0090 const G4String& fsType);
0091
0092 void UpdateData(G4int A, G4int Z, G4int index, G4double abundance,
0093 G4ParticleDefinition* projectile)
0094 {
0095 UpdateData(A, Z, 0, index, abundance, projectile);
0096 }
0097
0098 void UpdateData(G4int A, G4int Z, G4int M, G4int index, G4double abundance,
0099 G4ParticleDefinition* projectile);
0100
0101 void Harmonise(G4ParticleHPVector*& theStore, G4ParticleHPVector* theNew);
0102
0103 G4HadFinalState* ApplyYourself(const G4HadProjectile& theTrack,
0104 G4int isoNumber = -1,
0105 G4bool isElastic = false);
0106
0107 G4int GetNiso() const { return niso; }
0108
0109 G4double GetN(G4int i) const { return theFinalStates[i]->GetN(); }
0110 G4double GetZ(G4int i) const { return theFinalStates[i]->GetZ(); }
0111 G4double GetM(G4int i) const { return theFinalStates[i]->GetM(); }
0112
0113 G4bool HasDataInAnyFinalState() const
0114 {
0115 G4bool result = false;
0116 for (G4int i = 0; i < niso; ++i) {
0117 if (theFinalStates[i]->HasAnyData()) {
0118 result = true;
0119 break;
0120 }
0121 }
0122 return result;
0123 }
0124
0125 void DumpInfo() const;
0126
0127 const G4String& GetFSType() { return theFSType; }
0128
0129 G4ParticleHPFinalState** GetFinalStates() const { return theFinalStates; }
0130
0131
0132 G4WendtFissionFragmentGenerator* GetWendtFissionGenerator() const;
0133
0134 G4ParticleHPChannel(G4ParticleHPChannel &) = delete;
0135 G4ParticleHPChannel & operator=(const G4ParticleHPChannel &right) = delete;
0136
0137 protected:
0138
0139 G4ParticleHPManager* fManager;
0140
0141 private:
0142
0143 G4ParticleDefinition* theProjectile;
0144
0145 G4ParticleHPVector* theChannelData;
0146 G4Element* theElement{nullptr};
0147
0148
0149 G4ParticleHPVector* theBuffer{nullptr};
0150
0151 G4ParticleHPIsoData* theIsotopeWiseData{nullptr};
0152
0153 G4ParticleHPFinalState** theFinalStates{nullptr};
0154
0155
0156 G4WendtFissionFragmentGenerator* wendtFissionGenerator{nullptr};
0157 G4bool* active{nullptr};
0158 G4int niso{-1};
0159 G4int registerCount{-1};
0160
0161 G4String theDir{""};
0162 G4String theFSType{""};
0163 };
0164
0165 #endif