Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4HadronicProcessStore.cc
Go to the documentation of this file.
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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class file
30//
31//
32// File name: G4HadronicProcessStore
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 09.05.2008
37//
38// Modifications:
39// 23.01.2009 V.Ivanchenko add destruction of processes
40// 12.05.2020 A.Ribon introduced general verbose level in hadronics
41//
42// Class Description:
43// Singleton to store hadronic processes, to provide access to processes
44// and to printout information about processes
45//
46// -------------------------------------------------------------------
47//
48//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
49//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
50
52#include "G4SystemOfUnits.hh"
53#include "G4UnitsTable.hh"
54#include "G4Element.hh"
55#include "G4ProcessManager.hh"
56#include "G4Electron.hh"
57#include "G4Proton.hh"
58#include "G4ParticleTable.hh"
64#include <algorithm>
65
66//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
67
69{
70 static thread_local auto* _instance = new G4HadronicProcessStore{};
71 return _instance;
72}
73
74//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
75
77{
78 Clean();
79 delete theEPTestMessenger;
80}
81
82//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
83
85{
86 for(auto& itr : process)
87 delete itr;
88 process.clear();
89
90 for(auto& itr : extraProcess)
91 delete itr;
92 extraProcess.clear();
93
94 m_map.clear();
95 p_map.clear();
96
97 n_extra = 0;
98 n_proc = 0;
99}
100
101//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
102
103G4HadronicProcessStore::G4HadronicProcessStore()
104{
105 n_proc = 0;
106 n_part = 0;
107 n_model= 0;
108 n_extra= 0;
109 currentProcess = nullptr;
110 currentParticle = nullptr;
111 theGenericIon =
114 buildTableStart = true;
115 buildXSTable = false;
116 theEPTestMessenger = new G4HadronicEPTestMessenger(this);
117}
118
119//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
120
122 const G4ParticleDefinition* part,
123 G4double energy,
124 const G4VProcess* proc,
125 const G4Element* element,
126 const G4Material* material)
127{
128 G4double cross = 0.;
129 G4int subType = proc->GetProcessSubType();
130 if (subType == fHadronElastic)
131 cross = GetElasticCrossSectionPerAtom(part,energy,element,material);
132 else if (subType == fHadronInelastic)
133 cross = GetInelasticCrossSectionPerAtom(part,energy,element,material);
134 else if (subType == fCapture)
135 cross = GetCaptureCrossSectionPerAtom(part,energy,element,material);
136 else if (subType == fFission)
137 cross = GetFissionCrossSectionPerAtom(part,energy,element,material);
138 else if (subType == fChargeExchange)
139 cross = GetChargeExchangeCrossSectionPerAtom(part,energy,element,material);
140 return cross;
141}
142
143//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
144
146 const G4ParticleDefinition* part,
147 G4double energy,
148 const G4VProcess* proc,
149 const G4Material* material)
150{
151 G4double cross = 0.;
152 G4int subType = proc->GetProcessSubType();
153 if (subType == fHadronElastic)
154 cross = GetElasticCrossSectionPerVolume(part,energy,material);
155 else if (subType == fHadronInelastic)
156 cross = GetInelasticCrossSectionPerVolume(part,energy,material);
157 else if (subType == fCapture)
158 cross = GetCaptureCrossSectionPerVolume(part,energy,material);
159 else if (subType == fFission)
160 cross = GetFissionCrossSectionPerVolume(part,energy,material);
161 else if (subType == fChargeExchange)
162 cross = GetChargeExchangeCrossSectionPerVolume(part,energy,material);
163 return cross;
164}
165
166//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
167
169 const G4ParticleDefinition *aParticle,
170 G4double kineticEnergy,
171 const G4Material *material)
172{
173 G4double cross = 0.0;
174 const G4ElementVector* theElementVector = material->GetElementVector();
175 const G4double* theAtomNumDensityVector =
176 material->GetVecNbOfAtomsPerVolume();
177 size_t nelm = material->GetNumberOfElements();
178 for (size_t i=0; i<nelm; ++i) {
179 const G4Element* elm = (*theElementVector)[i];
180 cross += theAtomNumDensityVector[i]*
181 GetElasticCrossSectionPerAtom(aParticle,kineticEnergy,elm,material);
182 }
183 return cross;
184}
185
186//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
187
189 const G4ParticleDefinition *aParticle,
190 G4double kineticEnergy,
191 const G4Element *anElement, const G4Material* mat)
192{
194 G4double cross = 0.0;
195 localDP.SetKineticEnergy(kineticEnergy);
196 if(hp) {
197 cross = hp->GetElementCrossSection(&localDP,anElement,mat);
198 }
199 return cross;
200}
201
202//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
203
206 G4double,
207 G4int, G4int)
208{
209 return 0.0;
210}
211
212//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
213
215 const G4ParticleDefinition *aParticle,
216 G4double kineticEnergy,
217 const G4Material *material)
218{
219 G4double cross = 0.0;
220 const G4ElementVector* theElementVector = material->GetElementVector();
221 const G4double* theAtomNumDensityVector =
222 material->GetVecNbOfAtomsPerVolume();
223 size_t nelm = material->GetNumberOfElements();
224 for (size_t i=0; i<nelm; ++i) {
225 const G4Element* elm = (*theElementVector)[i];
226 cross += theAtomNumDensityVector[i]*
227 GetInelasticCrossSectionPerAtom(aParticle,kineticEnergy,elm,material);
228 }
229 return cross;
230}
231
232//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
233
235 const G4ParticleDefinition *aParticle,
236 G4double kineticEnergy,
237 const G4Element *anElement, const G4Material* mat)
238{
240 localDP.SetKineticEnergy(kineticEnergy);
241 G4double cross = 0.0;
242 if(hp) {
243 cross = hp->GetElementCrossSection(&localDP,anElement,mat);
244 }
245 return cross;
246}
247
248//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
249
251 const G4ParticleDefinition *,
252 G4double,
253 G4int, G4int)
254{
255 return 0.0;
256}
257
258//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
259
261 const G4ParticleDefinition *aParticle,
262 G4double kineticEnergy,
263 const G4Material *material)
264{
265 G4double cross = 0.0;
266 const G4ElementVector* theElementVector = material->GetElementVector();
267 const G4double* theAtomNumDensityVector =
268 material->GetVecNbOfAtomsPerVolume();
269 size_t nelm = material->GetNumberOfElements();
270 for (size_t i=0; i<nelm; ++i) {
271 const G4Element* elm = (*theElementVector)[i];
272 cross += theAtomNumDensityVector[i]*
273 GetCaptureCrossSectionPerAtom(aParticle,kineticEnergy,elm,material);
274 }
275 return cross;
276}
277
278//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
279
281 const G4ParticleDefinition *aParticle,
282 G4double kineticEnergy,
283 const G4Element *anElement, const G4Material* mat)
284{
285 G4HadronicProcess* hp = FindProcess(aParticle, fCapture);
286 localDP.SetKineticEnergy(kineticEnergy);
287 G4double cross = 0.0;
288 if(hp) {
289 cross = hp->GetElementCrossSection(&localDP,anElement,mat);
290 }
291 return cross;
292}
293
294//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
295
297 const G4ParticleDefinition *,
298 G4double,
299 G4int, G4int)
300{
301 return 0.0;
302}
303
304//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
305
307 const G4ParticleDefinition *aParticle,
308 G4double kineticEnergy,
309 const G4Material *material)
310{
311 G4double cross = 0.0;
312 const G4ElementVector* theElementVector = material->GetElementVector();
313 const G4double* theAtomNumDensityVector =
314 material->GetVecNbOfAtomsPerVolume();
315 size_t nelm = material->GetNumberOfElements();
316 for (size_t i=0; i<nelm; i++) {
317 const G4Element* elm = (*theElementVector)[i];
318 cross += theAtomNumDensityVector[i]*
319 GetFissionCrossSectionPerAtom(aParticle,kineticEnergy,elm,material);
320 }
321 return cross;
322}
323
324//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
325
327 const G4ParticleDefinition *aParticle,
328 G4double kineticEnergy,
329 const G4Element *anElement, const G4Material* mat)
330{
331 G4HadronicProcess* hp = FindProcess(aParticle, fFission);
332 localDP.SetKineticEnergy(kineticEnergy);
333 G4double cross = 0.0;
334 if(hp) {
335 cross = hp->GetElementCrossSection(&localDP,anElement,mat);
336 }
337 return cross;
338}
339
340//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
341
343 const G4ParticleDefinition *,
344 G4double,
345 G4int, G4int)
346{
347 return 0.0;
348}
349
350//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
351
353 const G4ParticleDefinition *aParticle,
354 G4double kineticEnergy,
355 const G4Material *material)
356{
357 G4double cross = 0.0;
358 const G4ElementVector* theElementVector = material->GetElementVector();
359 const G4double* theAtomNumDensityVector =
360 material->GetVecNbOfAtomsPerVolume();
361 size_t nelm = material->GetNumberOfElements();
362 for (size_t i=0; i<nelm; ++i) {
363 const G4Element* elm = (*theElementVector)[i];
364 cross += theAtomNumDensityVector[i]*
365 GetChargeExchangeCrossSectionPerAtom(aParticle,kineticEnergy,elm,material);
366 }
367 return cross;
368}
369
370//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
371
373 const G4ParticleDefinition *aParticle,
374 G4double kineticEnergy,
375 const G4Element *anElement, const G4Material* mat)
376{
378 localDP.SetKineticEnergy(kineticEnergy);
379 G4double cross = 0.0;
380 if(hp) {
381 cross = hp->GetElementCrossSection(&localDP,anElement,mat);
382 }
383 return cross;
384}
385
386//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
387
389 const G4ParticleDefinition *,
390 G4double,
391 G4int, G4int)
392{
393 return 0.0;
394}
395
396//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
397
399{
400 for(G4int i=0; i<n_proc; ++i) {
401 if(process[i] == proc) { return; }
402 }
403 if(1 < param->GetVerboseLevel()) {
404 G4cout << "G4HadronicProcessStore::Register hadronic " << n_proc
405 << " " << proc->GetProcessName() << G4endl;
406 }
407 ++n_proc;
408 process.push_back(proc);
409}
410
411//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
412
414 const G4ParticleDefinition* part)
415{
416 G4int i=0;
417 for(; i<n_proc; ++i) {if(process[i] == proc) break;}
418 G4int j=0;
419 for(; j<n_part; ++j) {if(particle[j] == part) break;}
420
421 if(1 < param->GetVerboseLevel()) {
422 G4cout << "G4HadronicProcessStore::RegisterParticle "
423 << part->GetParticleName()
424 << " for " << proc->GetProcessName() << G4endl;
425 }
426 if(j == n_part) {
427 ++n_part;
428 particle.push_back(part);
429 wasPrinted.push_back(0);
430 }
431
432 // the pair should be added?
433 if(i < n_proc) {
434 std::multimap<PD,HP,std::less<PD> >::iterator it;
435 for(it=p_map.lower_bound(part); it!=p_map.upper_bound(part); ++it) {
436 if(it->first == part) {
437 HP process2 = (it->second);
438 if(proc == process2) { return; }
439 }
440 }
441 }
442
443 p_map.insert(std::multimap<PD,HP>::value_type(part,proc));
444}
445
446//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
447
450{
451 G4int i=0;
452 for(; i<n_proc; ++i) {if(process[i] == proc) { break; }}
453 G4int k=0;
454 for(; k<n_model; ++k) {if(model[k] == mod) { break; }}
455
456 m_map.insert(std::multimap<HP,HI>::value_type(proc,mod));
457
458 if(k == n_model) {
459 ++n_model;
460 model.push_back(mod);
461 modelName.push_back(mod->GetModelName());
462 }
463}
464
465//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
466
468{
469 for(G4int i=0; i<n_proc; ++i) {
470 if(process[i] == proc) {
471 process[i] = nullptr;
473 return;
474 }
475 }
476}
477
478//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
479
481{
482 for(G4int i=0; i<n_extra; ++i) {
483 if(extraProcess[i] == proc) { return; }
484 }
485 G4HadronicProcess* hproc = static_cast<G4HadronicProcess*>(proc);
486 if(hproc) {
487 for(G4int i=0; i<n_proc; ++i) {
488 if(process[i] == hproc) { return; }
489 }
490 }
491 if(1 < param->GetVerboseLevel()) {
492 G4cout << "Extra Process: " << n_extra
493 << " " << proc->GetProcessName() << G4endl;
494 }
495 ++n_extra;
496 extraProcess.push_back(proc);
497}
498
499//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
500
502 G4VProcess* proc,
503 const G4ParticleDefinition* part)
504{
505 G4int i=0;
506 for(; i<n_extra; ++i) { if(extraProcess[i] == proc) { break; } }
507 G4int j=0;
508 for(; j<n_part; ++j) { if(particle[j] == part) { break; } }
509
510 if(j == n_part) {
511 ++n_part;
512 particle.push_back(part);
513 wasPrinted.push_back(0);
514 }
515
516 // the pair should be added?
517 if(i < n_extra) {
518 std::multimap<PD,G4VProcess*,std::less<PD> >::iterator it;
519 for(it=ep_map.lower_bound(part); it!=ep_map.upper_bound(part); ++it) {
520 if(it->first == part) {
521 G4VProcess* process2 = (it->second);
522 if(proc == process2) { return; }
523 }
524 }
525 }
526
527 ep_map.insert(std::multimap<PD,G4VProcess*>::value_type(part,proc));
528}
529
530//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
531
533{
534 for(G4int i=0; i<n_extra; ++i) {
535 if(extraProcess[i] == proc) {
536 extraProcess[i] = nullptr;
537 if(1 < param->GetVerboseLevel()) {
538 G4cout << "Extra Process: " << i << " "
539 <<proc->GetProcessName()<< " is deregisted " << G4endl;
540 }
541 return;
542 }
543 }
544}
545
546//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
547
549{
550 buildXSTable = val;
551}
552
553//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
554
556{
557 return buildXSTable;
558}
559
560//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
561
563{
564 // Trigger particle/process/model printout only when last particle is
565 // registered
566 if(buildTableStart && part == particle[n_part - 1]) {
567 buildTableStart = false;
568 Dump(param->GetVerboseLevel());
569 if (std::getenv("G4PhysListDocDir") ) DumpHtml();
571 }
572}
573
574//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
575
577{
578 // Automatic generation of html documentation page for physics lists
579 // List processes, models and cross sections for the most important
580 // particles in descending order of importance
581
582 char* dirName = std::getenv("G4PhysListDocDir");
583 char* physListName = std::getenv("G4PhysListName");
584 if (dirName && physListName) {
585
586 // Open output file with path name
587 G4String pathName = G4String(dirName) + "/" + G4String(physListName) + ".html";
588 std::ofstream outFile;
589 outFile.open(pathName);
590
591 // Write physics list summary file
592 outFile << "<html>\n";
593 outFile << "<head>\n";
594 outFile << "<title>Physics List Summary</title>\n";
595 outFile << "</head>\n";
596 outFile << "<body>\n";
597 outFile << "<h2> Summary of Hadronic Processes, Models and Cross Sections for Physics List "
598 << G4String(physListName) << "</h2>\n";
599 outFile << "<ul>\n";
600
601 PrintHtml(G4Proton::Proton(), outFile);
602 PrintHtml(G4Neutron::Neutron(), outFile);
603 PrintHtml(G4PionPlus::PionPlus(), outFile);
605 PrintHtml(G4Gamma::Gamma(), outFile);
607// PrintHtml(G4MuonMinus::MuonMinus(), outFile);
611 PrintHtml(G4Lambda::Lambda(), outFile);
612 PrintHtml(G4Alpha::Alpha(), outFile);
614
615 outFile << "</ul>\n";
616 outFile << "</body>\n";
617 outFile << "</html>\n";
618 outFile.close();
619 }
620}
621
622//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
623
625 std::ofstream& outFile)
626{
627 // Automatic generation of html documentation page for physics lists
628 // List processes for the most important particles in descending order
629 // of importance
630
631 outFile << "<br> <li><h2><font color=\" ff0000 \">"
632 << theParticle->GetParticleName() << "</font></h2></li>\n";
633
634 typedef std::multimap<PD,HP,std::less<PD> > PDHPmap;
635 typedef std::multimap<HP,HI,std::less<HP> > HPHImap;
636
637 std::pair<PDHPmap::iterator, PDHPmap::iterator> itpart =
638 p_map.equal_range(theParticle);
639
640 // Loop over processes assigned to particle
641
642 G4HadronicProcess* theProcess;
643 for (PDHPmap::iterator it = itpart.first; it != itpart.second; ++it) {
644 theProcess = (*it).second;
645 // description is inline
646 //outFile << "<br> &nbsp;&nbsp; <b><font color=\" 0000ff \">process : <a href=\""
647 // << theProcess->GetProcessName() << ".html\"> "
648 // << theProcess->GetProcessName() << "</a></font></b>\n";
649 outFile << "<br> &nbsp;&nbsp; <b><font color=\" 0000ff \">process : "
650 << theProcess->GetProcessName() << "</font></b>\n";
651 outFile << "<ul>\n";
652 outFile << " <li>";
653 theProcess->ProcessDescription(outFile);
654 outFile << " <li><b><font color=\" 00AA00 \">models : </font></b>\n";
655 // Loop over models assigned to process
656 std::pair<HPHImap::iterator, HPHImap::iterator> itmod =
657 m_map.equal_range(theProcess);
658
659 outFile << " <ul>\n";
660 G4String physListName(std::getenv("G4PhysListName"));
661
662 for (HPHImap::iterator jt = itmod.first; jt != itmod.second; ++jt) {
663 outFile << " <li><b><a href=\"" << physListName << "_"
664 << HtmlFileName((*jt).second->GetModelName()) << "\"> "
665 << (*jt).second->GetModelName() << "</a>"
666 << " from " << (*jt).second->GetMinEnergy()/GeV
667 << " GeV to " << (*jt).second->GetMaxEnergy()/GeV
668 << " GeV </b></li>\n";
669
670 // Print ModelDescription, ignore that we overwrite files n-times.
671 PrintModelHtml((*jt).second);
672
673 }
674 outFile << " </ul>\n";
675 outFile << " </li>\n";
676
677 // List cross sections assigned to process
678 outFile << " <li><b><font color=\" 00AA00 \">cross sections : </font></b>\n";
679 outFile << " <ul>\n";
680 theProcess->GetCrossSectionDataStore()->DumpHtml(*theParticle, outFile);
681 // << " \n";
682 outFile << " </ul>\n";
683
684 outFile << " </li>\n";
685 outFile << "</ul>\n";
686
687 }
688
689 // Loop over extra (G4VProcess) processes
690
691 std::multimap<PD,G4VProcess*,std::less<PD> >::iterator itp;
692 for (itp=ep_map.lower_bound(theParticle); itp!=ep_map.upper_bound(theParticle); ++itp) {
693 if (itp->first == theParticle) {
694 G4VProcess* proc = (itp->second);
695 outFile << "<br> &nbsp;&nbsp; <b><font color=\" 0000ff \">process : "
696 << proc->GetProcessName() << "</font></b>\n";
697 outFile << "<ul>\n";
698 outFile << " <li>";
699 proc->ProcessDescription(outFile);
700 outFile << " </li>\n";
701 outFile << "</ul>\n";
702 }
703 }
704
705} // PrintHtml for particle
706
707//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
708
709void
711{
712 G4String dirName(std::getenv("G4PhysListDocDir"));
713 G4String physListName(std::getenv("G4PhysListName"));
714 G4String pathName = dirName + "/" + physListName + "_" + HtmlFileName(mod->GetModelName());
715 std::ofstream outModel;
716 outModel.open(pathName);
717 outModel << "<html>\n";
718 outModel << "<head>\n";
719 outModel << "<title>Description of " << mod->GetModelName()
720 << "</title>\n";
721 outModel << "</head>\n";
722 outModel << "<body>\n";
723
724 mod->ModelDescription(outModel);
725
726 outModel << "</body>\n";
727 outModel << "</html>\n";
728
729}
730
731//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
732//private
733G4String G4HadronicProcessStore::HtmlFileName(const G4String & in) const
734{
735 G4String str(in);
736
737 // replace blanks:
738 std::transform(str.begin(), str.end(), str.begin(), [](char ch)
739 {
740 return ch == ' ' ? '_' : ch;
741 });
742 str=str + ".html";
743 return str;
744}
745
746//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
747
749{
750 G4int level = std::max(param->GetVerboseLevel(), verb);
751 if (0 == level) return;
752
753 G4cout
754 << "\n====================================================================\n"
755 << std::setw(60) << "HADRONIC PROCESSES SUMMARY (verbose level "
756 << level << ")" << G4endl;
757
758 for (G4int i=0; i<n_part; ++i) {
759 PD part = particle[i];
760 G4String pname = part->GetParticleName();
761 G4bool yes = false;
762
763 if (level == 1 && (pname == "proton" ||
764 pname == "neutron" ||
765 pname == "deuteron" ||
766 pname == "triton" ||
767 pname == "He3" ||
768 pname == "alpha" ||
769 pname == "pi+" ||
770 pname == "pi-" ||
771 pname == "gamma" ||
772 pname == "e+" ||
773 pname == "e-" ||
774 pname == "mu+" ||
775 pname == "mu-" ||
776 pname == "kaon+" ||
777 pname == "kaon-" ||
778 pname == "lambda" ||
779 pname == "anti_lambda" ||
780 pname == "sigma-" ||
781 pname == "D-" ||
782 pname == "B-" ||
783 pname == "GenericIon" ||
784 pname == "hypertriton" ||
785 pname == "anti_neutron" ||
786 pname == "anti_proton" ||
787 pname == "anti_deuteron" ||
788 pname == "anti_triton" ||
789 pname == "anti_He3" ||
790 pname == "anti_alpha" ||
791 pname == "anti_hypertriton")) yes = true;
792 if (level > 1) yes = true;
793 if (yes) {
794 // main processes
795 std::multimap<PD,HP,std::less<PD> >::iterator it;
796
797 for (it=p_map.lower_bound(part); it!=p_map.upper_bound(part); ++it) {
798 if (it->first == part) {
799 HP proc = (it->second);
800 G4int j=0;
801 for (; j<n_proc; ++j) {
802 if (process[j] == proc) { Print(j, i); }
803 }
804 }
805 }
806
807 // extra processes
808 std::multimap<PD,G4VProcess*,std::less<PD> >::iterator itp;
809 for(itp=ep_map.lower_bound(part); itp!=ep_map.upper_bound(part); ++itp) {
810 if(itp->first == part) {
811 G4VProcess* proc = (itp->second);
812 if (wasPrinted[i] == 0) {
813 G4cout << "\n---------------------------------------------------\n"
814 << std::setw(50) << "Hadronic Processes for "
815 << part->GetParticleName() << "\n";
816 wasPrinted[i] = 1;
817 }
818 G4cout << "\n Process: " << proc->GetProcessName() << G4endl;
819 }
820 }
821 }
822 }
823
824 G4cout << "\n================================================================"
825 << G4endl;
826}
827
828//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
829
830void G4HadronicProcessStore::Print(G4int idxProc, G4int idxPart)
831{
832 G4HadronicProcess* proc = process[idxProc];
833 const G4ParticleDefinition* part = particle[idxPart];
834 if(part == nullptr || proc == nullptr) { return; }
835 if (wasPrinted[idxPart] == 0) {
836 G4cout << "\n---------------------------------------------------\n"
837 << std::setw(50) << "Hadronic Processes for "
838 << part->GetParticleName() << "\n";
839 wasPrinted[idxPart] = 1;
840 }
841
842 G4cout << "\n Process: " << proc->GetProcessName();
843
844 // Append the string "/n" (i.e. "per nucleon") on the kinetic energy of ions.
845 G4String stringEnergyPerNucleon = "";
846 if (part == G4GenericIon::Definition() ||
847 std::abs( part->GetBaryonNumber() ) > 1) {
848 stringEnergyPerNucleon = "/n";
849 }
850 // print cross section factor
851 if(param->ApplyFactorXS()) {
852 G4int pdg = part->GetPDGEncoding();
853 G4int subType = proc->GetProcessSubType();
854 G4double fact = 1.0;
855 if(subType == fHadronInelastic) {
856 if(pdg == 2212 || pdg == 2112) {
857 fact = param->XSFactorNucleonInelastic();
858 } else if(std::abs(pdg) == 211) {
859 fact = param->XSFactorPionInelastic();
860 } else {
861 fact = param->XSFactorHadronInelastic();
862 }
863 } else if(subType == fHadronElastic) {
864 if(pdg == 2212 || pdg == 2112) {
865 fact = param->XSFactorNucleonElastic();
866 } else if(std::abs(pdg) == 211) {
867 fact = param->XSFactorPionElastic();
868 } else {
869 fact = param->XSFactorHadronElastic();
870 }
871 }
872 if(std::abs(fact - 1.0) > 1.e-6) {
873 G4cout << " XSfactor= " << fact;
874 }
875 }
876
877 HI hi = 0;
878 std::multimap<HP,HI,std::less<HP> >::iterator ih;
879 for(ih=m_map.lower_bound(proc); ih!=m_map.upper_bound(proc); ++ih) {
880 if(ih->first == proc) {
881 hi = ih->second;
882 G4int i=0;
883 for(; i<n_model; ++i) {
884 if(model[i] == hi) { break; }
885 }
886 G4cout << "\n Model: " << std::setw(25) << modelName[i] << ": "
887 << G4BestUnit(hi->GetMinEnergy(), "Energy") << stringEnergyPerNucleon
888 << " ---> "
889 << G4BestUnit(hi->GetMaxEnergy(), "Energy") << stringEnergyPerNucleon;
890 }
891 }
892 G4cout << G4endl;
893
895 csds->DumpPhysicsTable(*part);
896}
897
898//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
899
901// this code is obsolete - not optimal change verbose in each thread
902{
903 G4int i;
904 for(i=0; i<n_proc; ++i) {
905 if(process[i]) { process[i]->SetVerboseLevel(val); }
906 }
907 for(i=0; i<n_model; ++i) {
908 if(model[i]) { model[i]->SetVerboseLevel(val); }
909 }
910}
911
912//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
913
915{
916 return param->GetVerboseLevel();
917}
918
919//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
920
922 const G4ParticleDefinition* part, G4HadronicProcessType subType)
923{
924 bool isNew = false;
925 G4HadronicProcess* hp = nullptr;
926 localDP.SetDefinition(part);
927
928 if(part != currentParticle) {
929 const G4ParticleDefinition* p = part;
930 if(p->GetBaryonNumber() > 4 && p->GetParticleType() == "nucleus") {
931 p = theGenericIon;
932 }
933 if(p != currentParticle) {
934 isNew = true;
935 currentParticle = p;
936 }
937 }
938 if(!isNew) {
939 if(!currentProcess) {
940 isNew = true;
941 } else if(subType == currentProcess->GetProcessSubType()) {
942 hp = currentProcess;
943 } else {
944 isNew = true;
945 }
946 }
947 if(isNew) {
948 std::multimap<PD,HP,std::less<PD> >::iterator it;
949 for(it=p_map.lower_bound(currentParticle);
950 it!=p_map.upper_bound(currentParticle); ++it) {
951 if(it->first == currentParticle &&
952 subType == (it->second)->GetProcessSubType()) {
953 hp = it->second;
954 break;
955 }
956 }
957 currentProcess = hp;
958 }
959 return hp;
960}
961
962//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
963
965{
966 G4cout << " Setting energy/momentum report level to " << level
967 << " for " << process.size() << " hadronic processes " << G4endl;
968 for (G4int i = 0; i < G4int(process.size()); ++i) {
969 process[i]->SetEpReportLevel(level);
970 }
971}
972
973//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
974
976{
977 G4cout << " Setting absolute energy/momentum test level to " << abslevel
978 << G4endl;
979 G4double rellevel = 0.0;
980 G4HadronicProcess* theProcess = 0;
981 for (G4int i = 0; i < G4int(process.size()); ++i) {
982 theProcess = process[i];
983 rellevel = theProcess->GetEnergyMomentumCheckLevels().first;
984 theProcess->SetEnergyMomentumCheckLevels(rellevel, abslevel);
985 }
986}
987
988//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
989
991{
992 G4cout << " Setting relative energy/momentum test level to " << rellevel
993 << G4endl;
994 G4double abslevel = 0.0;
995 G4HadronicProcess* theProcess = 0;
996 for (G4int i = 0; i < G4int(process.size()); ++i) {
997 theProcess = process[i];
998 abslevel = theProcess->GetEnergyMomentumCheckLevels().second;
999 theProcess->SetEnergyMomentumCheckLevels(rellevel, abslevel);
1000 }
1001}
1002
1003//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
std::vector< const G4Element * > G4ElementVector
G4HadronicProcessType
@ fChargeExchange
@ fHadronElastic
@ fHadronInelastic
#define G4BestUnit(a, b)
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4Alpha * Alpha()
Definition: G4Alpha.cc:88
void DumpHtml(const G4ParticleDefinition &, std::ofstream &) const
void DumpPhysicsTable(const G4ParticleDefinition &)
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetKineticEnergy(G4double aEnergy)
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4GenericIon * Definition()
Definition: G4GenericIon.cc:48
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:92
static G4HadronicInteractionRegistry * Instance()
virtual void ModelDescription(std::ostream &outFile) const
const G4String & GetModelName() const
G4double XSFactorPionElastic() const
static G4HadronicParameters * Instance()
G4double XSFactorNucleonElastic() const
G4double XSFactorHadronInelastic() const
G4double XSFactorPionInelastic() const
G4double XSFactorHadronElastic() const
G4double XSFactorNucleonInelastic() const
void DeRegister(G4HadronicProcess *)
G4double GetCaptureCrossSectionPerAtom(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Element *anElement, const G4Material *mat=nullptr)
G4double GetCaptureCrossSectionPerIsotope(const G4ParticleDefinition *aParticle, G4double kineticEnergy, G4int Z, G4int A)
G4double GetCrossSectionPerVolume(const G4ParticleDefinition *particle, G4double kineticEnergy, const G4VProcess *process, const G4Material *material)
G4double GetCaptureCrossSectionPerVolume(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Material *material)
G4HadronicProcess * FindProcess(const G4ParticleDefinition *, G4HadronicProcessType subType)
void RegisterParticle(G4HadronicProcess *, const G4ParticleDefinition *)
G4double GetChargeExchangeCrossSectionPerVolume(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Material *material)
void PrintHtml(const G4ParticleDefinition *, std::ofstream &)
void SetProcessAbsLevel(G4double absoluteLevel)
G4double GetChargeExchangeCrossSectionPerIsotope(const G4ParticleDefinition *aParticle, G4double kineticEnergy, G4int Z, G4int A)
G4double GetFissionCrossSectionPerIsotope(const G4ParticleDefinition *aParticle, G4double kineticEnergy, G4int Z, G4int A)
G4double GetInelasticCrossSectionPerAtom(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Element *anElement, const G4Material *mat=nullptr)
G4double GetInelasticCrossSectionPerIsotope(const G4ParticleDefinition *aParticle, G4double kineticEnergy, G4int Z, G4int A)
G4double GetInelasticCrossSectionPerVolume(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Material *material)
void SetProcessRelLevel(G4double relativeLevel)
void DeRegisterExtraProcess(G4VProcess *)
G4double GetFissionCrossSectionPerAtom(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Element *anElement, const G4Material *mat=nullptr)
G4double GetChargeExchangeCrossSectionPerAtom(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Element *anElement, const G4Material *mat=nullptr)
void RegisterExtraProcess(G4VProcess *)
G4double GetElasticCrossSectionPerVolume(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Material *material)
void RegisterParticleForExtraProcess(G4VProcess *, const G4ParticleDefinition *)
void SetEpReportLevel(G4int level)
static G4HadronicProcessStore * Instance()
G4double GetElasticCrossSectionPerAtom(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Element *anElement, const G4Material *mat=0)
void RegisterInteraction(G4HadronicProcess *, G4HadronicInteraction *)
G4double GetCrossSectionPerAtom(const G4ParticleDefinition *particle, G4double kineticEnergy, const G4VProcess *process, const G4Element *element, const G4Material *material=nullptr)
void Register(G4HadronicProcess *)
void PrintModelHtml(const G4HadronicInteraction *model) const
G4double GetElasticCrossSectionPerIsotope(const G4ParticleDefinition *aParticle, G4double kineticEnergy, G4int Z, G4int A)
void PrintInfo(const G4ParticleDefinition *)
G4double GetFissionCrossSectionPerVolume(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Material *material)
void ProcessDescription(std::ostream &outFile) const override
G4double GetElementCrossSection(const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=nullptr)
std::pair< G4double, G4double > GetEnergyMomentumCheckLevels() const
G4CrossSectionDataStore * GetCrossSectionDataStore()
void SetEnergyMomentumCheckLevels(G4double relativeLevel, G4double absoluteLevel)
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:112
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:112
static G4Lambda * Lambda()
Definition: G4Lambda.cc:107
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:185
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:201
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
const G4String & GetParticleType() const
const G4String & GetParticleName() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:97
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:97
static G4Positron * Positron()
Definition: G4Positron.cc:93
static G4Proton * Proton()
Definition: G4Proton.cc:92
virtual void ProcessDescription(std::ostream &outfile) const
Definition: G4VProcess.cc:181
G4int GetProcessSubType() const
Definition: G4VProcess.hh:404
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386