| 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 | //
|
|---|
| 26 | // $Id: G4HadronProcessStore.cc,v 1.1 2006/10/31 11:35:02 gunter Exp $
|
|---|
| 27 | // GEANT4 tag $Name: $
|
|---|
| 28 | //
|
|---|
| 29 | // -------------------------------------------------------------------
|
|---|
| 30 | //
|
|---|
| 31 | // GEANT4 Class file
|
|---|
| 32 | //
|
|---|
| 33 | //
|
|---|
| 34 | // File name: G4HadronProcessStore
|
|---|
| 35 | //
|
|---|
| 36 | // Author: Vladimir Ivanchenko
|
|---|
| 37 | //
|
|---|
| 38 | // Creation date: 03.05.2006
|
|---|
| 39 | //
|
|---|
| 40 | // Modifications:
|
|---|
| 41 | // 04.08.06 V.Ivanchenko add computation of cross sections
|
|---|
| 42 | //
|
|---|
| 43 | //
|
|---|
| 44 | // Class Description:
|
|---|
| 45 | //
|
|---|
| 46 | // -------------------------------------------------------------------
|
|---|
| 47 | //
|
|---|
| 48 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 49 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 50 |
|
|---|
| 51 | #include "G4HadronProcessStore.hh"
|
|---|
| 52 | #include "G4Element.hh"
|
|---|
| 53 | #include "G4ProcessManager.hh"
|
|---|
| 54 | #include "G4Electron.hh"
|
|---|
| 55 | #include "G4Proton.hh"
|
|---|
| 56 |
|
|---|
| 57 | G4HadronProcessStore* G4HadronProcessStore::theInstance = 0;
|
|---|
| 58 |
|
|---|
| 59 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 60 |
|
|---|
| 61 | G4HadronProcessStore* G4HadronProcessStore::Instance()
|
|---|
| 62 | {
|
|---|
| 63 | if(0 == theInstance) {
|
|---|
| 64 | static G4HadronProcessStore manager;
|
|---|
| 65 | theInstance = &manager;
|
|---|
| 66 | }
|
|---|
| 67 | return theInstance;
|
|---|
| 68 | }
|
|---|
| 69 |
|
|---|
| 70 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 71 |
|
|---|
| 72 | G4HadronProcessStore::~G4HadronProcessStore()
|
|---|
| 73 | {
|
|---|
| 74 | for (G4int i=0; i<n_proc; i++) {
|
|---|
| 75 | if( process[i] ) delete process[i];
|
|---|
| 76 | }
|
|---|
| 77 | // for (G4int j=0; j<n_model; j++) {
|
|---|
| 78 | // if( model[j] ) delete model[j];
|
|---|
| 79 | // }
|
|---|
| 80 | }
|
|---|
| 81 |
|
|---|
| 82 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 83 |
|
|---|
| 84 | G4HadronProcessStore::G4HadronProcessStore()
|
|---|
| 85 | {
|
|---|
| 86 | n_proc = 0;
|
|---|
| 87 | n_part = 0;
|
|---|
| 88 | n_model= 0;
|
|---|
| 89 | currentProcess = 0;
|
|---|
| 90 | currentParticle = 0;
|
|---|
| 91 | verbose = 1;
|
|---|
| 92 | }
|
|---|
| 93 |
|
|---|
| 94 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 95 |
|
|---|
| 96 | G4double G4HadronProcessStore::GetInelasticCrossSectionPerVolume(
|
|---|
| 97 | const G4ParticleDefinition *aParticle,
|
|---|
| 98 | G4double kineticEnergy,
|
|---|
| 99 | const G4Material *material)
|
|---|
| 100 | {
|
|---|
| 101 | G4double cross = 0.0;
|
|---|
| 102 | const G4ElementVector* theElementVector = material->GetElementVector();
|
|---|
| 103 | const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
|
|---|
| 104 | size_t nelm = material->GetNumberOfElements();
|
|---|
| 105 | for (size_t i=0; i<nelm; i++) {
|
|---|
| 106 | const G4Element* elm = (*theElementVector)[i];
|
|---|
| 107 | cross += theAtomNumDensityVector[i]*
|
|---|
| 108 | GetInelasticCrossSectionPerAtom(aParticle,kineticEnergy,elm);
|
|---|
| 109 | }
|
|---|
| 110 | return cross;
|
|---|
| 111 | }
|
|---|
| 112 |
|
|---|
| 113 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 114 |
|
|---|
| 115 | G4double G4HadronProcessStore::GetInelasticCrossSectionPerAtom(
|
|---|
| 116 | const G4ParticleDefinition *aParticle,
|
|---|
| 117 | G4double kineticEnergy,
|
|---|
| 118 | const G4Element *anElement)
|
|---|
| 119 | {
|
|---|
| 120 | G4HadronicProcess* hp = FindInelasticProcess(aParticle);
|
|---|
| 121 | localDP.SetDefinition(const_cast<G4ParticleDefinition*>(aParticle));
|
|---|
| 122 | localDP.SetKineticEnergy(kineticEnergy);
|
|---|
| 123 | G4double cross = 0.0;
|
|---|
| 124 | if(hp) cross = hp->GetMicroscopicCrossSection(&localDP,
|
|---|
| 125 | anElement,
|
|---|
| 126 | STP_Temperature);
|
|---|
| 127 | return cross;
|
|---|
| 128 | }
|
|---|
| 129 |
|
|---|
| 130 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 131 |
|
|---|
| 132 | G4double G4HadronProcessStore::GetInelasticCrossSectionPerIsotope(
|
|---|
| 133 | const G4ParticleDefinition *,
|
|---|
| 134 | G4double,
|
|---|
| 135 | G4int, G4int)
|
|---|
| 136 | {
|
|---|
| 137 | return 0.0;
|
|---|
| 138 | }
|
|---|
| 139 |
|
|---|
| 140 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 141 |
|
|---|
| 142 | G4double G4HadronProcessStore::GetElasticCrossSectionPerVolume(
|
|---|
| 143 | const G4ParticleDefinition *aParticle,
|
|---|
| 144 | G4double kineticEnergy,
|
|---|
| 145 | const G4Material *material)
|
|---|
| 146 | {
|
|---|
| 147 | G4double cross = 0.0;
|
|---|
| 148 | const G4ElementVector* theElementVector = material->GetElementVector();
|
|---|
| 149 | const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
|
|---|
| 150 | size_t nelm = material->GetNumberOfElements();
|
|---|
| 151 | for (size_t i=0; i<nelm; i++) {
|
|---|
| 152 | const G4Element* elm = (*theElementVector)[i];
|
|---|
| 153 | cross += theAtomNumDensityVector[i]*
|
|---|
| 154 | GetElasticCrossSectionPerAtom(aParticle,kineticEnergy,elm);
|
|---|
| 155 | }
|
|---|
| 156 | return cross;
|
|---|
| 157 | }
|
|---|
| 158 |
|
|---|
| 159 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 160 |
|
|---|
| 161 | G4double G4HadronProcessStore::GetElasticCrossSectionPerAtom(
|
|---|
| 162 | const G4ParticleDefinition *aParticle,
|
|---|
| 163 | G4double kineticEnergy,
|
|---|
| 164 | const G4Element *anElement)
|
|---|
| 165 | {
|
|---|
| 166 | G4HadronicProcess* hp = FindElasticProcess(aParticle);
|
|---|
| 167 | localDP.SetKineticEnergy(kineticEnergy);
|
|---|
| 168 | G4double cross = 0.0;
|
|---|
| 169 | if(hp) cross = hp->GetMicroscopicCrossSection(&localDP,
|
|---|
| 170 | anElement,
|
|---|
| 171 | STP_Temperature);
|
|---|
| 172 | return cross;
|
|---|
| 173 | }
|
|---|
| 174 |
|
|---|
| 175 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 176 |
|
|---|
| 177 | G4double G4HadronProcessStore::GetElasticCrossSectionPerIsotope(
|
|---|
| 178 | const G4ParticleDefinition*,
|
|---|
| 179 | G4double,
|
|---|
| 180 | G4int, G4int)
|
|---|
| 181 | {
|
|---|
| 182 | return 0.0;
|
|---|
| 183 | }
|
|---|
| 184 |
|
|---|
| 185 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 186 |
|
|---|
| 187 | void G4HadronProcessStore::Register(G4HadronicProcess* proc,
|
|---|
| 188 | const G4ParticleDefinition* part,
|
|---|
| 189 | G4HadronicInteraction* mod,
|
|---|
| 190 | const G4String& name)
|
|---|
| 191 | {
|
|---|
| 192 | G4int i=0;
|
|---|
| 193 | for(; i<n_proc; i++) {if(process[i] == proc) break;}
|
|---|
| 194 | G4int j=0;
|
|---|
| 195 | for(; j<n_part; j++) {if(particle[j] == part) break;}
|
|---|
| 196 | G4int k=0;
|
|---|
| 197 | for(; k<n_model; k++) {if(model[k] == mod) break;}
|
|---|
| 198 |
|
|---|
| 199 | if(i == n_proc || j == n_part)
|
|---|
| 200 | p_map.insert(std::multimap<PD,HP>::value_type(part,proc));
|
|---|
| 201 |
|
|---|
| 202 | m_map.insert(std::multimap<HP,HI>::value_type(proc,mod));
|
|---|
| 203 |
|
|---|
| 204 | if(i == n_proc) {
|
|---|
| 205 | n_proc++;
|
|---|
| 206 | process.push_back(proc);
|
|---|
| 207 | }
|
|---|
| 208 | if(j == n_part) {
|
|---|
| 209 | n_part++;
|
|---|
| 210 | particle.push_back(part);
|
|---|
| 211 | }
|
|---|
| 212 | if(k == n_model) {
|
|---|
| 213 | n_model++;
|
|---|
| 214 | model.push_back(mod);
|
|---|
| 215 | modelName.push_back(name);
|
|---|
| 216 | }
|
|---|
| 217 | }
|
|---|
| 218 |
|
|---|
| 219 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 220 |
|
|---|
| 221 | void G4HadronProcessStore::Print(const G4ParticleDefinition* part)
|
|---|
| 222 | {
|
|---|
| 223 | G4cout<<G4endl;
|
|---|
| 224 | G4cout << " Hadronic Processes for "
|
|---|
| 225 | <<part->GetParticleName() << G4endl;
|
|---|
| 226 | HP hp = 0;
|
|---|
| 227 | HI hi = 0;
|
|---|
| 228 | G4bool first;
|
|---|
| 229 | std::multimap<PD,HP,std::less<PD> >::iterator it;
|
|---|
| 230 | std::multimap<HP,HI,std::less<HP> >::iterator ih;
|
|---|
| 231 | for(it=p_map.lower_bound(part); it!=p_map.upper_bound(part); ++it) {
|
|---|
| 232 | if(it->first == part) {
|
|---|
| 233 | hp = it->second;
|
|---|
| 234 | G4cout << std::setw(10) << hp->GetProcessName()
|
|---|
| 235 | << " Models: ";
|
|---|
| 236 | first = true;
|
|---|
| 237 | for(ih=m_map.lower_bound(hp); ih!=m_map.upper_bound(hp); ++ih) {
|
|---|
| 238 | if(ih->first == hp) {
|
|---|
| 239 | hi = ih->second;
|
|---|
| 240 | G4int i=0;
|
|---|
| 241 | for(; i<n_model; i++) {
|
|---|
| 242 | if(model[i] == hi) break;
|
|---|
| 243 | }
|
|---|
| 244 | if(!first) G4cout << " ";
|
|---|
| 245 | first = false;
|
|---|
| 246 | G4cout << std::setw(10) << modelName[i]
|
|---|
| 247 | << " Emin(GeV)= "
|
|---|
| 248 | << std::setw(5) << hi->GetMinEnergy()/GeV
|
|---|
| 249 | << " Emax(GeV)= "
|
|---|
| 250 | << std::setw(7) << hi->GetMaxEnergy()/GeV
|
|---|
| 251 | << G4endl;
|
|---|
| 252 | }
|
|---|
| 253 | }
|
|---|
| 254 | }
|
|---|
| 255 | }
|
|---|
| 256 | }
|
|---|
| 257 |
|
|---|
| 258 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 259 |
|
|---|
| 260 | void G4HadronProcessStore::Dump(G4int level)
|
|---|
| 261 | {
|
|---|
| 262 | if(level > 0)
|
|---|
| 263 | G4cout << "=============================================================="
|
|---|
| 264 | << "============================="
|
|---|
| 265 | << G4endl;
|
|---|
| 266 | for(G4int i=0; i<n_part; i++) {
|
|---|
| 267 | G4String pname = (particle[i])->GetParticleName();
|
|---|
| 268 | G4bool yes = false;
|
|---|
| 269 | if(level >= 2) yes = true;
|
|---|
| 270 | else if(level == 1 && (pname == "proton" ||
|
|---|
| 271 | pname == "neutron" ||
|
|---|
| 272 | pname == "pi+" ||
|
|---|
| 273 | pname == "pi-" ||
|
|---|
| 274 | pname == "kaon+" ||
|
|---|
| 275 | pname == "kaon-" ||
|
|---|
| 276 | pname == "lambda" ||
|
|---|
| 277 | pname == "anti_neutron" ||
|
|---|
| 278 | pname == "anti_proton")) yes = true;
|
|---|
| 279 | if(yes) Print(particle[i]);
|
|---|
| 280 | }
|
|---|
| 281 | if(level > 0)
|
|---|
| 282 | G4cout << "=============================================================="
|
|---|
| 283 | << "============================="
|
|---|
| 284 | << G4endl;
|
|---|
| 285 | }
|
|---|
| 286 |
|
|---|
| 287 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 288 |
|
|---|
| 289 | void G4HadronProcessStore::SetVerbose(G4int val)
|
|---|
| 290 | {
|
|---|
| 291 | verbose = val;
|
|---|
| 292 | G4int i;
|
|---|
| 293 | for(i=0; i<n_proc; i++) {
|
|---|
| 294 | if(process[i]) process[i]->SetVerboseLevel(val);
|
|---|
| 295 | }
|
|---|
| 296 | for(i=0; i<n_model; i++) {
|
|---|
| 297 | if(model[i]) model[i]->SetVerboseLevel(val);
|
|---|
| 298 | }
|
|---|
| 299 | }
|
|---|
| 300 |
|
|---|
| 301 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 302 |
|
|---|
| 303 | G4int G4HadronProcessStore::GetVerbose()
|
|---|
| 304 | {
|
|---|
| 305 | return verbose;
|
|---|
| 306 | }
|
|---|
| 307 |
|
|---|
| 308 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 309 |
|
|---|
| 310 | G4HadronicProcess* G4HadronProcessStore::FindElasticProcess(
|
|---|
| 311 | const G4ParticleDefinition* part)
|
|---|
| 312 | {
|
|---|
| 313 | bool isNew = false;
|
|---|
| 314 | G4HadronicProcess* hp = 0;
|
|---|
| 315 |
|
|---|
| 316 | if(part != currentParticle) {
|
|---|
| 317 | isNew = true;
|
|---|
| 318 | currentParticle = part;
|
|---|
| 319 | localDP.SetDefinition(const_cast<G4ParticleDefinition*>(part));
|
|---|
| 320 | } else if(!currentProcess) {
|
|---|
| 321 | isNew = true;
|
|---|
| 322 | } else if(currentProcess->GetProcessName() == "hElastic") {
|
|---|
| 323 | hp = currentProcess;
|
|---|
| 324 | } else {
|
|---|
| 325 | isNew = true;
|
|---|
| 326 | }
|
|---|
| 327 |
|
|---|
| 328 | if(isNew) {
|
|---|
| 329 | std::multimap<PD,HP,std::less<PD> >::iterator it;
|
|---|
| 330 | for(it=p_map.lower_bound(part); it!=p_map.upper_bound(part); ++it) {
|
|---|
| 331 | if(it->first == part && (it->second)->GetProcessName() == "hElastic") {
|
|---|
| 332 | hp = it->second;
|
|---|
| 333 | break;
|
|---|
| 334 | }
|
|---|
| 335 | }
|
|---|
| 336 | currentProcess = hp;
|
|---|
| 337 | }
|
|---|
| 338 |
|
|---|
| 339 | return hp;
|
|---|
| 340 | }
|
|---|
| 341 |
|
|---|
| 342 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|
| 343 |
|
|---|
| 344 | G4HadronicProcess* G4HadronProcessStore::FindInelasticProcess(
|
|---|
| 345 | const G4ParticleDefinition* part)
|
|---|
| 346 | {
|
|---|
| 347 | bool isNew = false;
|
|---|
| 348 | G4HadronicProcess* hp = 0;
|
|---|
| 349 |
|
|---|
| 350 | if(part != currentParticle) {
|
|---|
| 351 | isNew = true;
|
|---|
| 352 | currentParticle = part;
|
|---|
| 353 | localDP.SetDefinition(const_cast<G4ParticleDefinition*>(part));
|
|---|
| 354 | } else if(!currentProcess) {
|
|---|
| 355 | isNew = true;
|
|---|
| 356 | } else if(currentProcess->GetProcessName() == "hInelastic") {
|
|---|
| 357 | hp = currentProcess;
|
|---|
| 358 | } else {
|
|---|
| 359 | isNew = true;
|
|---|
| 360 | }
|
|---|
| 361 |
|
|---|
| 362 | if(isNew) {
|
|---|
| 363 | std::multimap<PD,HP,std::less<PD> >::iterator it;
|
|---|
| 364 | for(it=p_map.lower_bound(part); it!=p_map.upper_bound(part); ++it) {
|
|---|
| 365 | if(it->first == part && (it->second)->GetProcessName() == "hInelastic") {
|
|---|
| 366 | hp = it->second;
|
|---|
| 367 | break;
|
|---|
| 368 | }
|
|---|
| 369 | }
|
|---|
| 370 | currentProcess = hp;
|
|---|
| 371 | }
|
|---|
| 372 |
|
|---|
| 373 | return hp;
|
|---|
| 374 | }
|
|---|
| 375 |
|
|---|
| 376 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
|
|---|