00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049 #ifndef G4HadronicProcessStore_h
00050 #define G4HadronicProcessStore_h 1
00051
00052
00053 #include "globals.hh"
00054 #include "G4DynamicParticle.hh"
00055 #include "G4ThreeVector.hh"
00056 #include "G4HadronicProcess.hh"
00057 #include "G4HadronicInteraction.hh"
00058 #include "G4ParticleDefinition.hh"
00059 #include "G4HadronicProcessType.hh"
00060 #include <map>
00061 #include <vector>
00062 #include <iostream>
00063
00064 class G4Element;
00065 class G4HadronicEPTestMessenger;
00066
00067
00068 class G4HadronicProcessStore
00069 {
00070
00071 public:
00072
00073 static G4HadronicProcessStore* Instance();
00074
00075 ~G4HadronicProcessStore();
00076
00077 void Clean();
00078 G4double GetCrossSectionPerAtom(
00079 const G4ParticleDefinition* particle,
00080 G4double kineticEnergy,
00081 const G4VProcess* process,
00082 const G4Element* element);
00083
00084 G4double GetCrossSectionPerVolume(
00085 const G4ParticleDefinition* particle,
00086 G4double kineticEnergy,
00087 const G4VProcess* process,
00088 const G4Material* material);
00089
00090 G4double GetInelasticCrossSectionPerVolume(
00091 const G4ParticleDefinition *aParticle,
00092 G4double kineticEnergy,
00093 const G4Material *material);
00094
00095 G4double GetInelasticCrossSectionPerAtom(
00096 const G4ParticleDefinition *aParticle,
00097 G4double kineticEnergy,
00098 const G4Element *anElement, const G4Material* mat=0);
00099
00100 G4double GetInelasticCrossSectionPerIsotope(
00101 const G4ParticleDefinition *aParticle,
00102 G4double kineticEnergy,
00103 G4int Z, G4int A);
00104
00105 G4double GetElasticCrossSectionPerVolume(
00106 const G4ParticleDefinition *aParticle,
00107 G4double kineticEnergy,
00108 const G4Material *material);
00109
00110 G4double GetElasticCrossSectionPerAtom(
00111 const G4ParticleDefinition *aParticle,
00112 G4double kineticEnergy,
00113 const G4Element *anElement, const G4Material* mat=0);
00114
00115 G4double GetElasticCrossSectionPerIsotope(
00116 const G4ParticleDefinition *aParticle,
00117 G4double kineticEnergy,
00118 G4int Z, G4int A);
00119
00120 G4double GetCaptureCrossSectionPerVolume(
00121 const G4ParticleDefinition *aParticle,
00122 G4double kineticEnergy,
00123 const G4Material *material);
00124
00125 G4double GetCaptureCrossSectionPerAtom(
00126 const G4ParticleDefinition *aParticle,
00127 G4double kineticEnergy,
00128 const G4Element *anElement, const G4Material* mat=0);
00129
00130 G4double GetCaptureCrossSectionPerIsotope(
00131 const G4ParticleDefinition *aParticle,
00132 G4double kineticEnergy,
00133 G4int Z, G4int A);
00134
00135 G4double GetFissionCrossSectionPerVolume(
00136 const G4ParticleDefinition *aParticle,
00137 G4double kineticEnergy,
00138 const G4Material *material);
00139
00140 G4double GetFissionCrossSectionPerAtom(
00141 const G4ParticleDefinition *aParticle,
00142 G4double kineticEnergy,
00143 const G4Element *anElement, const G4Material* mat=0);
00144
00145 G4double GetFissionCrossSectionPerIsotope(
00146 const G4ParticleDefinition *aParticle,
00147 G4double kineticEnergy,
00148 G4int Z, G4int A);
00149
00150 G4double GetChargeExchangeCrossSectionPerVolume(
00151 const G4ParticleDefinition *aParticle,
00152 G4double kineticEnergy,
00153 const G4Material *material);
00154
00155 G4double GetChargeExchangeCrossSectionPerAtom(
00156 const G4ParticleDefinition *aParticle,
00157 G4double kineticEnergy,
00158 const G4Element *anElement, const G4Material* mat=0);
00159
00160 G4double GetChargeExchangeCrossSectionPerIsotope(
00161 const G4ParticleDefinition *aParticle,
00162 G4double kineticEnergy,
00163 G4int Z, G4int A);
00164
00165
00166 void Register(G4HadronicProcess*);
00167
00168 void RegisterParticle(G4HadronicProcess*,
00169 const G4ParticleDefinition*);
00170
00171 void RegisterInteraction(G4HadronicProcess*,
00172 G4HadronicInteraction*);
00173
00174 void DeRegister(G4HadronicProcess*);
00175
00176
00177 void RegisterExtraProcess(G4VProcess*);
00178
00179 void RegisterParticleForExtraProcess(G4VProcess*,
00180 const G4ParticleDefinition*);
00181
00182 void DeRegisterExtraProcess(G4VProcess*);
00183
00184 void PrintInfo(const G4ParticleDefinition*);
00185
00186 void Dump(G4int level);
00187 void DumpHtml();
00188 void PrintHtml(const G4ParticleDefinition*, std::ofstream&);
00189 void PrintModelHtml(const G4HadronicInteraction * model) const;
00190
00191 void SetVerbose(G4int val);
00192
00193 G4int GetVerbose();
00194
00195 G4HadronicProcess* FindProcess(const G4ParticleDefinition*,
00196 G4HadronicProcessType subType);
00197
00198
00199 void SetEpReportLevel(G4int level);
00200
00201 void SetProcessAbsLevel(G4double absoluteLevel);
00202
00203 void SetProcessRelLevel(G4double relativeLevel);
00204
00205 private:
00206
00207
00208 G4HadronicProcessStore();
00209
00210
00211 void Print(G4int idxProcess, G4int idxParticle);
00212
00213 static G4HadronicProcessStore* theInstance;
00214
00215 typedef const G4ParticleDefinition* PD;
00216 typedef G4HadronicProcess* HP;
00217 typedef G4HadronicInteraction* HI;
00218
00219
00220 std::vector<G4HadronicProcess*> process;
00221 std::vector<G4HadronicInteraction*> model;
00222 std::vector<G4String> modelName;
00223 std::vector<PD> particle;
00224 std::vector<G4int> wasPrinted;
00225
00226 std::multimap<PD,HP> p_map;
00227 std::multimap<HP,HI> m_map;
00228
00229
00230 std::vector<G4VProcess*> extraProcess;
00231 std::multimap<PD,G4VProcess*> ep_map;
00232
00233
00234 G4int n_proc;
00235 G4int n_model;
00236 G4int n_part;
00237 G4int n_extra;
00238
00239 G4int verbose;
00240 G4bool buildTableStart;
00241
00242
00243 HP currentProcess;
00244 PD currentParticle;
00245
00246 G4DynamicParticle localDP;
00247
00248 G4HadronicEPTestMessenger* theEPTestMessenger;
00249 };
00250
00251
00252 #endif
00253