File indexing completed on 2025-11-04 10:12:37
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 #ifndef G4HadronicProcessStore_h
0048 #define G4HadronicProcessStore_h 1
0049
0050
0051 #include "globals.hh"
0052 #include "G4DynamicParticle.hh"
0053 #include "G4ThreeVector.hh"
0054 #include "G4HadronicProcess.hh"
0055 #include "G4HadronicInteraction.hh"
0056 #include "G4ParticleDefinition.hh"
0057 #include "G4HadronicProcessType.hh"
0058 #include "G4ThreadLocalSingleton.hh"
0059 #include <map>
0060 #include <vector>
0061 #include <iostream>
0062
0063 class G4Element;
0064 class G4HadronicEPTestMessenger;
0065 class G4HadronicParameters;
0066
0067 class G4HadronicProcessStore
0068 {
0069
0070 friend class G4ThreadLocalSingleton<G4HadronicProcessStore>;
0071
0072 public:
0073
0074 static G4HadronicProcessStore* Instance();
0075
0076 ~G4HadronicProcessStore();
0077
0078 G4double GetCrossSectionPerAtom(
0079 const G4ParticleDefinition* particle,
0080 G4double kineticEnergy,
0081 const G4VProcess* process,
0082 const G4Element* element,
0083 const G4Material* material=nullptr);
0084
0085 G4double GetCrossSectionPerVolume(
0086 const G4ParticleDefinition* particle,
0087 G4double kineticEnergy,
0088 const G4VProcess* process,
0089 const G4Material* material);
0090
0091 G4double GetInelasticCrossSectionPerVolume(
0092 const G4ParticleDefinition *aParticle,
0093 G4double kineticEnergy,
0094 const G4Material *material);
0095
0096 G4double GetInelasticCrossSectionPerAtom(
0097 const G4ParticleDefinition *aParticle,
0098 G4double kineticEnergy,
0099 const G4Element *anElement,
0100 const G4Material* mat=nullptr);
0101
0102 G4double GetInelasticCrossSectionPerIsotope(
0103 const G4ParticleDefinition *aParticle,
0104 G4double kineticEnergy,
0105 G4int Z, G4int A);
0106
0107 G4double GetElasticCrossSectionPerVolume(
0108 const G4ParticleDefinition *aParticle,
0109 G4double kineticEnergy,
0110 const G4Material *material);
0111
0112 G4double GetElasticCrossSectionPerAtom(
0113 const G4ParticleDefinition *aParticle,
0114 G4double kineticEnergy,
0115 const G4Element *anElement, const G4Material* mat=0);
0116
0117 G4double GetElasticCrossSectionPerIsotope(
0118 const G4ParticleDefinition *aParticle,
0119 G4double kineticEnergy,
0120 G4int Z, G4int A);
0121
0122 G4double GetCaptureCrossSectionPerVolume(
0123 const G4ParticleDefinition *aParticle,
0124 G4double kineticEnergy,
0125 const G4Material *material);
0126
0127 G4double GetCaptureCrossSectionPerAtom(
0128 const G4ParticleDefinition *aParticle,
0129 G4double kineticEnergy,
0130 const G4Element *anElement,
0131 const G4Material* mat=nullptr);
0132
0133 G4double GetCaptureCrossSectionPerIsotope(
0134 const G4ParticleDefinition *aParticle,
0135 G4double kineticEnergy,
0136 G4int Z, G4int A);
0137
0138 G4double GetFissionCrossSectionPerVolume(
0139 const G4ParticleDefinition *aParticle,
0140 G4double kineticEnergy,
0141 const G4Material *material);
0142
0143 G4double GetFissionCrossSectionPerAtom(
0144 const G4ParticleDefinition *aParticle,
0145 G4double kineticEnergy,
0146 const G4Element *anElement,
0147 const G4Material* mat=nullptr);
0148
0149 G4double GetFissionCrossSectionPerIsotope(
0150 const G4ParticleDefinition *aParticle,
0151 G4double kineticEnergy,
0152 G4int Z, G4int A);
0153
0154 G4double GetChargeExchangeCrossSectionPerVolume(
0155 const G4ParticleDefinition *aParticle,
0156 G4double kineticEnergy,
0157 const G4Material *material);
0158
0159 G4double GetChargeExchangeCrossSectionPerAtom(
0160 const G4ParticleDefinition *aParticle,
0161 G4double kineticEnergy,
0162 const G4Element *anElement,
0163 const G4Material* mat=nullptr);
0164
0165 G4double GetChargeExchangeCrossSectionPerIsotope(
0166 const G4ParticleDefinition *aParticle,
0167 G4double kineticEnergy,
0168 G4int Z, G4int A);
0169
0170
0171 void Register(G4HadronicProcess*);
0172
0173 void RegisterParticle(G4HadronicProcess*,
0174 const G4ParticleDefinition*);
0175
0176 void RegisterInteraction(G4HadronicProcess*,
0177 G4HadronicInteraction*);
0178
0179 void DeRegister(G4HadronicProcess*);
0180
0181
0182 void RegisterExtraProcess(G4VProcess*);
0183
0184 void RegisterParticleForExtraProcess(G4VProcess*,
0185 const G4ParticleDefinition*);
0186
0187 void DeRegisterExtraProcess(G4VProcess*);
0188
0189 void SetBuildXSTable(G4bool val);
0190
0191 G4bool GetBuildXSTable() const;
0192
0193 void PrintInfo(const G4ParticleDefinition*);
0194
0195 void Dump(G4int level);
0196 void DumpHtml();
0197 void PrintHtml(const G4ParticleDefinition*, std::ofstream&);
0198 void PrintModelHtml(const G4HadronicInteraction * model) const;
0199
0200 void SetVerbose(G4int val);
0201 G4int GetVerbose();
0202
0203
0204 G4HadronicProcess* FindProcess(const G4ParticleDefinition*,
0205 G4HadronicProcessType subType);
0206
0207
0208 void SetEpReportLevel(G4int level);
0209
0210 void SetProcessAbsLevel(G4double absoluteLevel);
0211
0212 void SetProcessRelLevel(G4double relativeLevel);
0213
0214 private:
0215
0216
0217 G4HadronicProcessStore();
0218
0219
0220 void Print(G4int idxProcess, G4int idxParticle);
0221
0222 G4String HtmlFileName(const G4String &) const;
0223
0224 static G4ThreadLocal G4HadronicProcessStore* instance;
0225
0226 typedef const G4ParticleDefinition* PD;
0227 typedef G4HadronicProcess* HP;
0228 typedef G4HadronicInteraction* HI;
0229
0230
0231 std::vector<G4HadronicProcess*> process;
0232 std::vector<G4HadronicInteraction*> model;
0233 std::vector<G4String> modelName;
0234 std::vector<PD> particle;
0235 std::vector<G4int> wasPrinted;
0236
0237 std::multimap<PD,HP> p_map;
0238 std::multimap<HP,HI> m_map;
0239
0240
0241 std::vector<G4VProcess*> extraProcess;
0242 std::multimap<PD,G4VProcess*> ep_map;
0243
0244 G4HadronicParameters* param;
0245
0246
0247 G4int n_proc = 0;
0248 G4int n_model = 0;
0249 G4int n_part = 0;
0250 G4int n_extra = 0;
0251
0252 G4bool buildTableStart = true;
0253 G4bool buildXSTable = false;
0254
0255
0256 HP currentProcess = nullptr;
0257 PD currentParticle = nullptr;
0258 PD theGenericIon;
0259
0260 G4DynamicParticle localDP;
0261
0262 G4HadronicEPTestMessenger* theEPTestMessenger;
0263 };
0264
0265 #endif