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