| [968] | 1 | //
|
|---|
| 2 | // ********************************************************************
|
|---|
| 3 | // * License and Disclaimer *
|
|---|
| 4 | // * *
|
|---|
| 5 | // * The Geant4 software is copyright of the Copyright Holders of *
|
|---|
| 6 | // * the Geant4 Collaboration. It is provided under the terms and *
|
|---|
| 7 | // * conditions of the Geant4 Software License, included in the file *
|
|---|
| 8 | // * LICENSE and available at http://cern.ch/geant4/license . These *
|
|---|
| 9 | // * include a list of copyright holders. *
|
|---|
| 10 | // * *
|
|---|
| 11 | // * Neither the authors of this software system, nor their employing *
|
|---|
| 12 | // * institutes,nor the agencies providing financial support for this *
|
|---|
| 13 | // * work make any representation or warranty, express or implied, *
|
|---|
| 14 | // * regarding this software system or assume any liability for its *
|
|---|
| 15 | // * use. Please see the license in the file LICENSE and URL above *
|
|---|
| 16 | // * for the full disclaimer and the limitation of liability. *
|
|---|
| 17 | // * *
|
|---|
| 18 | // * This code implementation is the result of the scientific and *
|
|---|
| 19 | // * technical work of the GEANT4 collaboration. *
|
|---|
| 20 | // * By using, copying, modifying or distributing the software (or *
|
|---|
| 21 | // * any work based on the software) you agree to acknowledge its *
|
|---|
| 22 | // * use in resulting scientific publications, and indicate your *
|
|---|
| 23 | // * acceptance of all terms of the Geant4 Software license. *
|
|---|
| 24 | // ********************************************************************
|
|---|
| 25 | //
|
|---|
| [1340] | 26 | // $Id: G4HadronicProcessStore.cc,v 1.16 2010/08/17 09:47:54 vnivanch Exp $
|
|---|
| 27 | // GEANT4 tag $Name: hadr-man-V09-03-04 $
|
|---|
| [968] | 28 | //
|
|---|
| 29 | // -------------------------------------------------------------------
|
|---|
| 30 | //
|
|---|
| 31 | // GEANT4 Class file
|
|---|
| 32 | //
|
|---|
| 33 | //
|
|---|
| 34 | // File name: G4HadronicProcessStore
|
|---|
| 35 | //
|
|---|
| 36 | // Author: Vladimir Ivanchenko
|
|---|
| 37 | //
|
|---|
| 38 | // Creation date: 09.05.2008
|
|---|
| 39 | //
|
|---|
| 40 | // Modifications:
|
|---|
| [1055] | 41 | // 23.01.2009 V.Ivanchenko add destruction of processes
|
|---|
| [968] | 42 | //
|
|---|
| 43 | // Class Description:
|
|---|
| [1055] | 44 | // Singleton to store hadronic processes, to provide access to processes
|
|---|
| 45 | // and to printout information about processes
|
|---|
| [968] | 46 | //
|
|---|
| 47 | // -------------------------------------------------------------------
|
|---|
| 48 | //
|
|---|
| 49 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 50 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 51 |
|
|---|
| 52 | #include "G4HadronicProcessStore.hh"
|
|---|
| 53 | #include "G4Element.hh"
|
|---|
| 54 | #include "G4ProcessManager.hh"
|
|---|
| 55 | #include "G4Electron.hh"
|
|---|
| 56 | #include "G4Proton.hh"
|
|---|
| [1055] | 57 | #include "G4HadronicInteractionRegistry.hh"
|
|---|
| 58 | #include "G4CrossSectionDataSetRegistry.hh"
|
|---|
| [1315] | 59 | #include "G4HadronicEPTestMessenger.hh"
|
|---|
| [968] | 60 |
|
|---|
| 61 | G4HadronicProcessStore* G4HadronicProcessStore::theInstance = 0;
|
|---|
| 62 |
|
|---|
| 63 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 64 |
|
|---|
| 65 | G4HadronicProcessStore* G4HadronicProcessStore::Instance()
|
|---|
| 66 | {
|
|---|
| 67 | if(0 == theInstance) {
|
|---|
| 68 | static G4HadronicProcessStore manager;
|
|---|
| 69 | theInstance = &manager;
|
|---|
| 70 | }
|
|---|
| 71 | return theInstance;
|
|---|
| 72 | }
|
|---|
| 73 |
|
|---|
| 74 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 75 |
|
|---|
| 76 | G4HadronicProcessStore::~G4HadronicProcessStore()
|
|---|
| 77 | {
|
|---|
| [1055] | 78 | Clean();
|
|---|
| 79 | G4HadronicInteractionRegistry::Instance()->Clean();
|
|---|
| 80 | G4CrossSectionDataSetRegistry::Instance()->Clean();
|
|---|
| [1315] | 81 | delete theEPTestMessenger;
|
|---|
| [1055] | 82 | }
|
|---|
| 83 |
|
|---|
| 84 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 85 |
|
|---|
| 86 | void G4HadronicProcessStore::Clean()
|
|---|
| 87 | {
|
|---|
| 88 | G4int i;
|
|---|
| 89 | //G4cout << "G4HadronicProcessStore::Clean() Nproc= " << n_proc
|
|---|
| 90 | // << " Nextra= " << n_extra << G4endl;
|
|---|
| 91 | if(n_proc > 0) {
|
|---|
| [1340] | 92 | for (i=0; i<n_proc; ++i) {
|
|---|
| [1055] | 93 | if( process[i] ) {
|
|---|
| 94 | //G4cout << "G4HadronicProcessStore::Clean() delete hadronic " << i << G4endl;
|
|---|
| 95 | //G4cout << process[i]->GetProcessName() << G4endl;
|
|---|
| [1228] | 96 | G4HadronicProcess* p = process[i];
|
|---|
| [1055] | 97 | process[i] = 0;
|
|---|
| [1228] | 98 | delete p;
|
|---|
| [1055] | 99 | }
|
|---|
| 100 | }
|
|---|
| [968] | 101 | }
|
|---|
| [1055] | 102 | if(n_extra > 0) {
|
|---|
| [1340] | 103 | for(i=0; i<n_extra; ++i) {
|
|---|
| [1055] | 104 | if(extraProcess[i]) {
|
|---|
| 105 | //G4cout << "G4HadronicProcessStore::Clean() delete extra "
|
|---|
| 106 | // << i << G4endl;
|
|---|
| 107 | //G4cout << extraProcess[i]->GetProcessName() << G4endl;
|
|---|
| [1228] | 108 | G4VProcess* p = extraProcess[i];
|
|---|
| [1055] | 109 | extraProcess[i] = 0;
|
|---|
| [1228] | 110 | delete p;
|
|---|
| [1055] | 111 | }
|
|---|
| 112 | }
|
|---|
| 113 | }
|
|---|
| 114 | //G4cout << "G4HadronicProcessStore::Clean() done" << G4endl;
|
|---|
| 115 | n_extra = 0;
|
|---|
| 116 | n_proc = 0;
|
|---|
| [968] | 117 | }
|
|---|
| 118 |
|
|---|
| 119 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 120 |
|
|---|
| 121 | G4HadronicProcessStore::G4HadronicProcessStore()
|
|---|
| 122 | {
|
|---|
| 123 | n_proc = 0;
|
|---|
| 124 | n_part = 0;
|
|---|
| 125 | n_model= 0;
|
|---|
| 126 | n_extra= 0;
|
|---|
| 127 | currentProcess = 0;
|
|---|
| 128 | currentParticle = 0;
|
|---|
| 129 | verbose = 1;
|
|---|
| 130 | buildTableStart = true;
|
|---|
| [1315] | 131 | theEPTestMessenger = new G4HadronicEPTestMessenger(this);
|
|---|
| [968] | 132 | }
|
|---|
| 133 |
|
|---|
| 134 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 135 |
|
|---|
| 136 | G4double G4HadronicProcessStore::GetElasticCrossSectionPerVolume(
|
|---|
| 137 | const G4ParticleDefinition *aParticle,
|
|---|
| 138 | G4double kineticEnergy,
|
|---|
| 139 | const G4Material *material)
|
|---|
| 140 | {
|
|---|
| 141 | G4double cross = 0.0;
|
|---|
| 142 | const G4ElementVector* theElementVector = material->GetElementVector();
|
|---|
| 143 | const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
|
|---|
| 144 | size_t nelm = material->GetNumberOfElements();
|
|---|
| [1340] | 145 | for (size_t i=0; i<nelm; ++i) {
|
|---|
| [968] | 146 | const G4Element* elm = (*theElementVector)[i];
|
|---|
| 147 | cross += theAtomNumDensityVector[i]*
|
|---|
| 148 | GetElasticCrossSectionPerAtom(aParticle,kineticEnergy,elm);
|
|---|
| 149 | }
|
|---|
| 150 | return cross;
|
|---|
| 151 | }
|
|---|
| 152 |
|
|---|
| 153 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 154 |
|
|---|
| 155 | G4double G4HadronicProcessStore::GetElasticCrossSectionPerAtom(
|
|---|
| 156 | const G4ParticleDefinition *aParticle,
|
|---|
| 157 | G4double kineticEnergy,
|
|---|
| 158 | const G4Element *anElement)
|
|---|
| 159 | {
|
|---|
| 160 | G4HadronicProcess* hp = FindProcess(aParticle, fHadronElastic);
|
|---|
| [1340] | 161 | G4double cross = 0.0;
|
|---|
| [968] | 162 | localDP.SetKineticEnergy(kineticEnergy);
|
|---|
| [1340] | 163 | if(hp) {
|
|---|
| 164 | cross = hp->GetMicroscopicCrossSection(&localDP,
|
|---|
| 165 | anElement,
|
|---|
| 166 | STP_Temperature);
|
|---|
| 167 | }
|
|---|
| [968] | 168 | return cross;
|
|---|
| 169 | }
|
|---|
| 170 |
|
|---|
| 171 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 172 |
|
|---|
| 173 | G4double G4HadronicProcessStore::GetElasticCrossSectionPerIsotope(
|
|---|
| 174 | const G4ParticleDefinition*,
|
|---|
| 175 | G4double,
|
|---|
| 176 | G4int, G4int)
|
|---|
| 177 | {
|
|---|
| 178 | return 0.0;
|
|---|
| 179 | }
|
|---|
| 180 |
|
|---|
| 181 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 182 |
|
|---|
| 183 | G4double G4HadronicProcessStore::GetInelasticCrossSectionPerVolume(
|
|---|
| 184 | const G4ParticleDefinition *aParticle,
|
|---|
| 185 | G4double kineticEnergy,
|
|---|
| 186 | const G4Material *material)
|
|---|
| 187 | {
|
|---|
| 188 | G4double cross = 0.0;
|
|---|
| 189 | const G4ElementVector* theElementVector = material->GetElementVector();
|
|---|
| 190 | const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
|
|---|
| 191 | size_t nelm = material->GetNumberOfElements();
|
|---|
| [1340] | 192 | for (size_t i=0; i<nelm; ++i) {
|
|---|
| [968] | 193 | const G4Element* elm = (*theElementVector)[i];
|
|---|
| 194 | cross += theAtomNumDensityVector[i]*
|
|---|
| 195 | GetInelasticCrossSectionPerAtom(aParticle,kineticEnergy,elm);
|
|---|
| 196 | }
|
|---|
| 197 | return cross;
|
|---|
| 198 | }
|
|---|
| 199 |
|
|---|
| 200 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 201 |
|
|---|
| 202 | G4double G4HadronicProcessStore::GetInelasticCrossSectionPerAtom(
|
|---|
| 203 | const G4ParticleDefinition *aParticle,
|
|---|
| 204 | G4double kineticEnergy,
|
|---|
| 205 | const G4Element *anElement)
|
|---|
| 206 | {
|
|---|
| 207 | G4HadronicProcess* hp = FindProcess(aParticle, fHadronInelastic);
|
|---|
| 208 | localDP.SetKineticEnergy(kineticEnergy);
|
|---|
| 209 | G4double cross = 0.0;
|
|---|
| [1340] | 210 | if(hp) {
|
|---|
| 211 | cross = hp->GetMicroscopicCrossSection(&localDP,
|
|---|
| 212 | anElement,
|
|---|
| 213 | STP_Temperature);
|
|---|
| 214 | }
|
|---|
| [968] | 215 | return cross;
|
|---|
| 216 | }
|
|---|
| 217 |
|
|---|
| 218 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 219 |
|
|---|
| 220 | G4double G4HadronicProcessStore::GetInelasticCrossSectionPerIsotope(
|
|---|
| 221 | const G4ParticleDefinition *,
|
|---|
| 222 | G4double,
|
|---|
| 223 | G4int, G4int)
|
|---|
| 224 | {
|
|---|
| 225 | return 0.0;
|
|---|
| 226 | }
|
|---|
| 227 |
|
|---|
| 228 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 229 |
|
|---|
| 230 | G4double G4HadronicProcessStore::GetCaptureCrossSectionPerVolume(
|
|---|
| 231 | const G4ParticleDefinition *aParticle,
|
|---|
| 232 | G4double kineticEnergy,
|
|---|
| 233 | const G4Material *material)
|
|---|
| 234 | {
|
|---|
| 235 | G4double cross = 0.0;
|
|---|
| 236 | const G4ElementVector* theElementVector = material->GetElementVector();
|
|---|
| 237 | const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
|
|---|
| 238 | size_t nelm = material->GetNumberOfElements();
|
|---|
| [1340] | 239 | for (size_t i=0; i<nelm; ++i) {
|
|---|
| [968] | 240 | const G4Element* elm = (*theElementVector)[i];
|
|---|
| 241 | cross += theAtomNumDensityVector[i]*
|
|---|
| 242 | GetCaptureCrossSectionPerAtom(aParticle,kineticEnergy,elm);
|
|---|
| 243 | }
|
|---|
| 244 | return cross;
|
|---|
| 245 | }
|
|---|
| 246 |
|
|---|
| 247 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 248 |
|
|---|
| 249 | G4double G4HadronicProcessStore::GetCaptureCrossSectionPerAtom(
|
|---|
| 250 | const G4ParticleDefinition *aParticle,
|
|---|
| 251 | G4double kineticEnergy,
|
|---|
| 252 | const G4Element *anElement)
|
|---|
| 253 | {
|
|---|
| 254 | G4HadronicProcess* hp = FindProcess(aParticle, fCapture);
|
|---|
| 255 | localDP.SetKineticEnergy(kineticEnergy);
|
|---|
| 256 | G4double cross = 0.0;
|
|---|
| [1340] | 257 | if(hp) {
|
|---|
| 258 | cross = hp->GetMicroscopicCrossSection(&localDP,
|
|---|
| 259 | anElement,
|
|---|
| 260 | STP_Temperature);
|
|---|
| 261 | }
|
|---|
| [968] | 262 | return cross;
|
|---|
| 263 | }
|
|---|
| 264 |
|
|---|
| 265 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 266 |
|
|---|
| 267 | G4double G4HadronicProcessStore::GetCaptureCrossSectionPerIsotope(
|
|---|
| 268 | const G4ParticleDefinition *,
|
|---|
| 269 | G4double,
|
|---|
| 270 | G4int, G4int)
|
|---|
| 271 | {
|
|---|
| 272 | return 0.0;
|
|---|
| 273 | }
|
|---|
| 274 |
|
|---|
| 275 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 276 |
|
|---|
| 277 | G4double G4HadronicProcessStore::GetFissionCrossSectionPerVolume(
|
|---|
| 278 | const G4ParticleDefinition *aParticle,
|
|---|
| 279 | G4double kineticEnergy,
|
|---|
| 280 | const G4Material *material)
|
|---|
| 281 | {
|
|---|
| 282 | G4double cross = 0.0;
|
|---|
| 283 | const G4ElementVector* theElementVector = material->GetElementVector();
|
|---|
| 284 | const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
|
|---|
| 285 | size_t nelm = material->GetNumberOfElements();
|
|---|
| 286 | for (size_t i=0; i<nelm; i++) {
|
|---|
| 287 | const G4Element* elm = (*theElementVector)[i];
|
|---|
| 288 | cross += theAtomNumDensityVector[i]*
|
|---|
| 289 | GetFissionCrossSectionPerAtom(aParticle,kineticEnergy,elm);
|
|---|
| 290 | }
|
|---|
| 291 | return cross;
|
|---|
| 292 | }
|
|---|
| 293 |
|
|---|
| 294 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 295 |
|
|---|
| 296 | G4double G4HadronicProcessStore::GetFissionCrossSectionPerAtom(
|
|---|
| 297 | const G4ParticleDefinition *aParticle,
|
|---|
| 298 | G4double kineticEnergy,
|
|---|
| 299 | const G4Element *anElement)
|
|---|
| 300 | {
|
|---|
| 301 | G4HadronicProcess* hp = FindProcess(aParticle, fFission);
|
|---|
| 302 | localDP.SetKineticEnergy(kineticEnergy);
|
|---|
| 303 | G4double cross = 0.0;
|
|---|
| [1340] | 304 | if(hp) {
|
|---|
| 305 | cross = hp->GetMicroscopicCrossSection(&localDP,
|
|---|
| 306 | anElement,
|
|---|
| 307 | STP_Temperature);
|
|---|
| 308 | }
|
|---|
| [968] | 309 | return cross;
|
|---|
| 310 | }
|
|---|
| 311 |
|
|---|
| 312 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 313 |
|
|---|
| 314 | G4double G4HadronicProcessStore::GetFissionCrossSectionPerIsotope(
|
|---|
| 315 | const G4ParticleDefinition *,
|
|---|
| 316 | G4double,
|
|---|
| 317 | G4int, G4int)
|
|---|
| 318 | {
|
|---|
| 319 | return 0.0;
|
|---|
| 320 | }
|
|---|
| 321 |
|
|---|
| 322 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 323 |
|
|---|
| 324 | G4double G4HadronicProcessStore::GetChargeExchangeCrossSectionPerVolume(
|
|---|
| 325 | const G4ParticleDefinition *aParticle,
|
|---|
| 326 | G4double kineticEnergy,
|
|---|
| 327 | const G4Material *material)
|
|---|
| 328 | {
|
|---|
| 329 | G4double cross = 0.0;
|
|---|
| 330 | const G4ElementVector* theElementVector = material->GetElementVector();
|
|---|
| 331 | const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
|
|---|
| 332 | size_t nelm = material->GetNumberOfElements();
|
|---|
| [1340] | 333 | for (size_t i=0; i<nelm; ++i) {
|
|---|
| [968] | 334 | const G4Element* elm = (*theElementVector)[i];
|
|---|
| 335 | cross += theAtomNumDensityVector[i]*
|
|---|
| 336 | GetChargeExchangeCrossSectionPerAtom(aParticle,kineticEnergy,elm);
|
|---|
| 337 | }
|
|---|
| 338 | return cross;
|
|---|
| 339 | }
|
|---|
| 340 |
|
|---|
| 341 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 342 |
|
|---|
| 343 | G4double G4HadronicProcessStore::GetChargeExchangeCrossSectionPerAtom(
|
|---|
| 344 | const G4ParticleDefinition *aParticle,
|
|---|
| 345 | G4double kineticEnergy,
|
|---|
| 346 | const G4Element *anElement)
|
|---|
| 347 | {
|
|---|
| 348 | G4HadronicProcess* hp = FindProcess(aParticle, fChargeExchange);
|
|---|
| 349 | localDP.SetKineticEnergy(kineticEnergy);
|
|---|
| 350 | G4double cross = 0.0;
|
|---|
| [1340] | 351 | if(hp) {
|
|---|
| 352 | cross = hp->GetMicroscopicCrossSection(&localDP,
|
|---|
| 353 | anElement,
|
|---|
| 354 | STP_Temperature);
|
|---|
| 355 | }
|
|---|
| [968] | 356 | return cross;
|
|---|
| 357 | }
|
|---|
| 358 |
|
|---|
| 359 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 360 |
|
|---|
| 361 | G4double G4HadronicProcessStore::GetChargeExchangeCrossSectionPerIsotope(
|
|---|
| 362 | const G4ParticleDefinition *,
|
|---|
| 363 | G4double,
|
|---|
| 364 | G4int, G4int)
|
|---|
| 365 | {
|
|---|
| 366 | return 0.0;
|
|---|
| 367 | }
|
|---|
| 368 |
|
|---|
| 369 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 370 |
|
|---|
| 371 | void G4HadronicProcessStore::Register(G4HadronicProcess* proc)
|
|---|
| 372 | {
|
|---|
| [1055] | 373 | if(0 < n_proc) {
|
|---|
| [1340] | 374 | for(G4int i=0; i<n_proc; ++i) {
|
|---|
| [1228] | 375 | if(process[i] == proc) { return; }
|
|---|
| [1055] | 376 | }
|
|---|
| 377 | }
|
|---|
| [1315] | 378 | // G4cout << "G4HadronicProcessStore::Register hadronic " << n_proc
|
|---|
| [1055] | 379 | // << " " << proc->GetProcessName() << G4endl;
|
|---|
| [1340] | 380 | ++n_proc;
|
|---|
| [968] | 381 | process.push_back(proc);
|
|---|
| 382 | }
|
|---|
| 383 |
|
|---|
| 384 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 385 |
|
|---|
| 386 | void G4HadronicProcessStore::RegisterParticle(G4HadronicProcess* proc,
|
|---|
| 387 | const G4ParticleDefinition* part)
|
|---|
| 388 | {
|
|---|
| 389 | G4int i=0;
|
|---|
| [1340] | 390 | for(; i<n_proc; ++i) {if(process[i] == proc) break;}
|
|---|
| [968] | 391 | G4int j=0;
|
|---|
| [1340] | 392 | for(; j<n_part; ++j) {if(particle[j] == part) break;}
|
|---|
| [968] | 393 |
|
|---|
| 394 | if(j == n_part) {
|
|---|
| [1340] | 395 | ++n_part;
|
|---|
| [968] | 396 | particle.push_back(part);
|
|---|
| 397 | wasPrinted.push_back(0);
|
|---|
| 398 | }
|
|---|
| 399 |
|
|---|
| 400 | // the pair should be added?
|
|---|
| 401 | if(i < n_proc) {
|
|---|
| 402 | std::multimap<PD,HP,std::less<PD> >::iterator it;
|
|---|
| 403 | for(it=p_map.lower_bound(part); it!=p_map.upper_bound(part); ++it) {
|
|---|
| 404 | if(it->first == part) {
|
|---|
| 405 | HP process = (it->second);
|
|---|
| [1340] | 406 | if(proc == process) { return; }
|
|---|
| [968] | 407 | }
|
|---|
| 408 | }
|
|---|
| 409 | }
|
|---|
| 410 |
|
|---|
| 411 | p_map.insert(std::multimap<PD,HP>::value_type(part,proc));
|
|---|
| 412 | }
|
|---|
| 413 |
|
|---|
| 414 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 415 |
|
|---|
| 416 | void G4HadronicProcessStore::RegisterInteraction(G4HadronicProcess* proc,
|
|---|
| 417 | G4HadronicInteraction* mod)
|
|---|
| 418 | {
|
|---|
| 419 | G4int i=0;
|
|---|
| [1340] | 420 | for(; i<n_proc; ++i) {if(process[i] == proc) { break; }}
|
|---|
| [968] | 421 | G4int k=0;
|
|---|
| [1340] | 422 | for(; k<n_model; ++k) {if(model[k] == mod) { break; }}
|
|---|
| [968] | 423 |
|
|---|
| 424 | m_map.insert(std::multimap<HP,HI>::value_type(proc,mod));
|
|---|
| 425 |
|
|---|
| 426 | if(k == n_model) {
|
|---|
| [1340] | 427 | ++n_model;
|
|---|
| [968] | 428 | model.push_back(mod);
|
|---|
| 429 | modelName.push_back(mod->GetModelName());
|
|---|
| 430 | }
|
|---|
| 431 | }
|
|---|
| 432 |
|
|---|
| 433 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 434 |
|
|---|
| 435 | void G4HadronicProcessStore::DeRegister(G4HadronicProcess* proc)
|
|---|
| 436 | {
|
|---|
| [1055] | 437 | if(0 == n_proc) return;
|
|---|
| [1340] | 438 | for(G4int i=0; i<n_proc; ++i) {
|
|---|
| [968] | 439 | if(process[i] == proc) {
|
|---|
| 440 | process[i] = 0;
|
|---|
| [1055] | 441 | return;
|
|---|
| [968] | 442 | }
|
|---|
| 443 | }
|
|---|
| 444 | }
|
|---|
| 445 |
|
|---|
| 446 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 447 |
|
|---|
| 448 | void G4HadronicProcessStore::RegisterExtraProcess(G4VProcess* proc)
|
|---|
| 449 | {
|
|---|
| [1055] | 450 | if(0 < n_extra) {
|
|---|
| [1340] | 451 | for(G4int i=0; i<n_extra; ++i) {
|
|---|
| 452 | if(extraProcess[i] == proc) { return; }
|
|---|
| [1055] | 453 | }
|
|---|
| 454 | }
|
|---|
| 455 | //G4cout << "Extra Process: " << n_extra << " " << proc->GetProcessName()
|
|---|
| 456 | // << " " << proc << G4endl;
|
|---|
| [968] | 457 |
|
|---|
| 458 | n_extra++;
|
|---|
| 459 | extraProcess.push_back(proc);
|
|---|
| 460 | }
|
|---|
| 461 |
|
|---|
| 462 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 463 |
|
|---|
| 464 | void G4HadronicProcessStore::RegisterParticleForExtraProcess(
|
|---|
| 465 | G4VProcess* proc,
|
|---|
| 466 | const G4ParticleDefinition* part)
|
|---|
| 467 | {
|
|---|
| 468 | G4int i=0;
|
|---|
| [1340] | 469 | for(; i<n_extra; ++i) { if(extraProcess[i] == proc) { break; } }
|
|---|
| [968] | 470 | G4int j=0;
|
|---|
| [1340] | 471 | for(; j<n_part; ++j) { if(particle[j] == part) { break; } }
|
|---|
| [968] | 472 |
|
|---|
| 473 | if(j == n_part) {
|
|---|
| [1340] | 474 | ++n_part;
|
|---|
| [968] | 475 | particle.push_back(part);
|
|---|
| 476 | wasPrinted.push_back(0);
|
|---|
| 477 | }
|
|---|
| 478 |
|
|---|
| 479 | // the pair should be added?
|
|---|
| 480 | if(i < n_extra) {
|
|---|
| 481 | std::multimap<PD,G4VProcess*,std::less<PD> >::iterator it;
|
|---|
| 482 | for(it=ep_map.lower_bound(part); it!=ep_map.upper_bound(part); ++it) {
|
|---|
| 483 | if(it->first == part) {
|
|---|
| 484 | G4VProcess* process = (it->second);
|
|---|
| [1340] | 485 | if(proc == process) { return; }
|
|---|
| [968] | 486 | }
|
|---|
| 487 | }
|
|---|
| 488 | }
|
|---|
| 489 |
|
|---|
| 490 | ep_map.insert(std::multimap<PD,G4VProcess*>::value_type(part,proc));
|
|---|
| 491 | }
|
|---|
| 492 |
|
|---|
| 493 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 494 |
|
|---|
| 495 | void G4HadronicProcessStore::DeRegisterExtraProcess(G4VProcess* proc)
|
|---|
| 496 | {
|
|---|
| [1340] | 497 | //G4cout << "Deregister Extra Process: " << proc << " "<<proc->GetProcessName()<< G4endl;
|
|---|
| 498 | if(0 == n_extra) { return; }
|
|---|
| 499 | for(G4int i=0; i<n_extra; ++i) {
|
|---|
| [968] | 500 | if(extraProcess[i] == proc) {
|
|---|
| 501 | extraProcess[i] = 0;
|
|---|
| [1055] | 502 | //G4cout << "Extra Process: " << i << " is deregisted " << G4endl;
|
|---|
| 503 | return;
|
|---|
| [968] | 504 | }
|
|---|
| 505 | }
|
|---|
| 506 | }
|
|---|
| 507 |
|
|---|
| 508 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 509 |
|
|---|
| 510 | void G4HadronicProcessStore::PrintInfo(const G4ParticleDefinition* part)
|
|---|
| 511 | {
|
|---|
| 512 | if(buildTableStart && part == particle[n_part - 1]) {
|
|---|
| 513 | buildTableStart = false;
|
|---|
| 514 | Dump(verbose);
|
|---|
| 515 | }
|
|---|
| 516 | }
|
|---|
| 517 |
|
|---|
| 518 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 519 |
|
|---|
| 520 | void G4HadronicProcessStore::Dump(G4int level)
|
|---|
| 521 | {
|
|---|
| 522 | if(level > 0) {
|
|---|
| 523 | G4cout << "=============================================================="
|
|---|
| 524 | << "=============================="
|
|---|
| 525 | << G4endl;
|
|---|
| 526 | G4cout << " HADRONIC PROCESSES SUMMARY (verbose level " << level
|
|---|
| 527 | << ")" << G4endl;
|
|---|
| 528 | }
|
|---|
| [1340] | 529 | for(G4int i=0; i<n_part; ++i) {
|
|---|
| [968] | 530 | PD part = particle[i];
|
|---|
| 531 | G4String pname = part->GetParticleName();
|
|---|
| 532 | G4bool yes = false;
|
|---|
| 533 | if(level >= 2) yes = true;
|
|---|
| 534 | else if(level == 1 && (pname == "proton" ||
|
|---|
| 535 | pname == "neutron" ||
|
|---|
| 536 | pname == "pi+" ||
|
|---|
| 537 | pname == "pi-" ||
|
|---|
| 538 | pname == "gamma" ||
|
|---|
| 539 | pname == "e-" ||
|
|---|
| 540 | pname == "mu-" ||
|
|---|
| 541 | pname == "kaon+" ||
|
|---|
| 542 | pname == "kaon-" ||
|
|---|
| 543 | pname == "lambda" ||
|
|---|
| 544 | pname == "anti_neutron" ||
|
|---|
| 545 | pname == "anti_proton")) yes = true;
|
|---|
| 546 | if(yes) {
|
|---|
| 547 | // main processes
|
|---|
| 548 | std::multimap<PD,HP,std::less<PD> >::iterator it;
|
|---|
| 549 | for(it=p_map.lower_bound(part); it!=p_map.upper_bound(part); ++it) {
|
|---|
| 550 | if(it->first == part) {
|
|---|
| 551 | HP proc = (it->second);
|
|---|
| 552 | G4int j=0;
|
|---|
| [1340] | 553 | for(; j<n_proc; ++j) {
|
|---|
| [968] | 554 | if(process[j] == proc) {
|
|---|
| 555 | Print(j, i);
|
|---|
| 556 | }
|
|---|
| 557 | }
|
|---|
| 558 | }
|
|---|
| 559 | }
|
|---|
| 560 | // extra processes
|
|---|
| 561 | std::multimap<PD,G4VProcess*,std::less<PD> >::iterator itp;
|
|---|
| 562 | for(itp=ep_map.lower_bound(part); itp!=ep_map.upper_bound(part); ++itp) {
|
|---|
| 563 | if(itp->first == part) {
|
|---|
| 564 | G4VProcess* proc = (itp->second);
|
|---|
| 565 | if(wasPrinted[i] == 0) {
|
|---|
| 566 | wasPrinted[i] = 1;
|
|---|
| 567 | G4cout<<G4endl;
|
|---|
| 568 | G4cout << " Hadronic Processes for <"
|
|---|
| 569 | <<part->GetParticleName() << ">" << G4endl;
|
|---|
| 570 | }
|
|---|
| 571 | G4cout << " " << proc->GetProcessName() << G4endl;
|
|---|
| 572 | }
|
|---|
| 573 | }
|
|---|
| 574 | }
|
|---|
| 575 | }
|
|---|
| 576 | if(level > 0) {
|
|---|
| 577 | G4cout << "=============================================================="
|
|---|
| 578 | << "=============================="
|
|---|
| 579 | << G4endl;
|
|---|
| 580 | }
|
|---|
| 581 | }
|
|---|
| 582 |
|
|---|
| 583 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 584 |
|
|---|
| 585 | void G4HadronicProcessStore::Print(G4int idxProc, G4int idxPart)
|
|---|
| 586 | {
|
|---|
| 587 | G4HadronicProcess* proc = process[idxProc];
|
|---|
| 588 | const G4ParticleDefinition* part = particle[idxPart];
|
|---|
| 589 | if(wasPrinted[idxPart] == 0) {
|
|---|
| 590 | wasPrinted[idxPart] = 1;
|
|---|
| 591 | G4cout<<G4endl;
|
|---|
| 592 | G4cout << " Hadronic Processes for <"
|
|---|
| 593 | <<part->GetParticleName() << ">" << G4endl;
|
|---|
| 594 | }
|
|---|
| 595 | HI hi = 0;
|
|---|
| 596 | G4bool first;
|
|---|
| 597 | std::multimap<HP,HI,std::less<HP> >::iterator ih;
|
|---|
| 598 | G4cout << std::setw(20) << proc->GetProcessName()
|
|---|
| 599 | << " Models: ";
|
|---|
| 600 | first = true;
|
|---|
| 601 | for(ih=m_map.lower_bound(proc); ih!=m_map.upper_bound(proc); ++ih) {
|
|---|
| 602 | if(ih->first == proc) {
|
|---|
| 603 | hi = ih->second;
|
|---|
| 604 | G4int i=0;
|
|---|
| [1340] | 605 | for(; i<n_model; ++i) {
|
|---|
| 606 | if(model[i] == hi) { break; }
|
|---|
| [968] | 607 | }
|
|---|
| 608 | if(!first) G4cout << " ";
|
|---|
| 609 | first = false;
|
|---|
| 610 | G4cout << std::setw(25) << modelName[i]
|
|---|
| 611 | << ": Emin(GeV)= "
|
|---|
| 612 | << std::setw(5) << hi->GetMinEnergy()/GeV
|
|---|
| 613 | << " Emax(GeV)= "
|
|---|
| 614 | << hi->GetMaxEnergy()/GeV
|
|---|
| 615 | << G4endl;
|
|---|
| 616 | }
|
|---|
| 617 | }
|
|---|
| 618 | }
|
|---|
| 619 |
|
|---|
| 620 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 621 |
|
|---|
| 622 | void G4HadronicProcessStore::SetVerbose(G4int val)
|
|---|
| 623 | {
|
|---|
| 624 | verbose = val;
|
|---|
| 625 | G4int i;
|
|---|
| [1340] | 626 | for(i=0; i<n_proc; ++i) {
|
|---|
| 627 | if(process[i]) { process[i]->SetVerboseLevel(val); }
|
|---|
| [968] | 628 | }
|
|---|
| [1340] | 629 | for(i=0; i<n_model; ++i) {
|
|---|
| 630 | if(model[i]) { model[i]->SetVerboseLevel(val); }
|
|---|
| [968] | 631 | }
|
|---|
| 632 | }
|
|---|
| 633 |
|
|---|
| 634 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 635 |
|
|---|
| 636 | G4int G4HadronicProcessStore::GetVerbose()
|
|---|
| 637 | {
|
|---|
| 638 | return verbose;
|
|---|
| 639 | }
|
|---|
| 640 |
|
|---|
| 641 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 642 |
|
|---|
| 643 | G4HadronicProcess* G4HadronicProcessStore::FindProcess(
|
|---|
| 644 | const G4ParticleDefinition* part, G4HadronicProcessType subType)
|
|---|
| 645 | {
|
|---|
| 646 | bool isNew = false;
|
|---|
| 647 | G4HadronicProcess* hp = 0;
|
|---|
| 648 |
|
|---|
| 649 | if(part != currentParticle) {
|
|---|
| 650 | isNew = true;
|
|---|
| 651 | currentParticle = part;
|
|---|
| [1340] | 652 | localDP.SetDefinition(part);
|
|---|
| [968] | 653 | } else if(!currentProcess) {
|
|---|
| 654 | isNew = true;
|
|---|
| 655 | } else if(subType == currentProcess->GetProcessSubType()) {
|
|---|
| 656 | hp = currentProcess;
|
|---|
| 657 | } else {
|
|---|
| 658 | isNew = true;
|
|---|
| 659 | }
|
|---|
| 660 |
|
|---|
| 661 | if(isNew) {
|
|---|
| 662 | std::multimap<PD,HP,std::less<PD> >::iterator it;
|
|---|
| 663 | for(it=p_map.lower_bound(part); it!=p_map.upper_bound(part); ++it) {
|
|---|
| 664 | if(it->first == part && subType == (it->second)->GetProcessSubType()) {
|
|---|
| 665 | hp = it->second;
|
|---|
| 666 | break;
|
|---|
| 667 | }
|
|---|
| 668 | }
|
|---|
| 669 | currentProcess = hp;
|
|---|
| 670 | }
|
|---|
| 671 |
|
|---|
| 672 | return hp;
|
|---|
| 673 | }
|
|---|
| 674 |
|
|---|
| [1315] | 675 | void G4HadronicProcessStore::SetEpReportLevel(G4int level)
|
|---|
| 676 | {
|
|---|
| 677 | G4cout << " Setting energy/momentum report level to " << level
|
|---|
| 678 | << " for " << process.size() << " hadronic processes " << G4endl;
|
|---|
| [1340] | 679 | for (G4int i = 0; i < G4int(process.size()); ++i) {
|
|---|
| [1315] | 680 | process[i]->SetEpReportLevel(level);
|
|---|
| 681 | }
|
|---|
| 682 | }
|
|---|
| 683 |
|
|---|
| 684 | void G4HadronicProcessStore::SetProcessAbsLevel(G4double abslevel)
|
|---|
| 685 | {
|
|---|
| 686 | G4cout << " Setting absolute energy/momentum test level to " << abslevel << G4endl;
|
|---|
| 687 | G4double rellevel = 0.0;
|
|---|
| 688 | G4HadronicProcess* theProcess = 0;
|
|---|
| [1340] | 689 | for (G4int i = 0; i < G4int(process.size()); ++i) {
|
|---|
| [1315] | 690 | theProcess = process[i];
|
|---|
| 691 | rellevel = theProcess->GetEnergyMomentumCheckLevels().first;
|
|---|
| 692 | theProcess->SetEnergyMomentumCheckLevels(rellevel, abslevel);
|
|---|
| 693 | }
|
|---|
| 694 | }
|
|---|
| 695 |
|
|---|
| 696 | void G4HadronicProcessStore::SetProcessRelLevel(G4double rellevel)
|
|---|
| 697 | {
|
|---|
| 698 | G4cout << " Setting relative energy/momentum test level to " << rellevel << G4endl;
|
|---|
| 699 | G4double abslevel = 0.0;
|
|---|
| 700 | G4HadronicProcess* theProcess = 0;
|
|---|
| [1340] | 701 | for (G4int i = 0; i < G4int(process.size()); ++i) {
|
|---|
| [1315] | 702 | theProcess = process[i];
|
|---|
| 703 | abslevel = theProcess->GetEnergyMomentumCheckLevels().second;
|
|---|
| 704 | theProcess->SetEnergyMomentumCheckLevels(rellevel, abslevel);
|
|---|
| 705 | }
|
|---|
| 706 | }
|
|---|