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