Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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// $Id$
27//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class file
31//
32//
33// File name: G4HadronicProcessStore
34//
35// Author: Vladimir Ivanchenko
36//
37// Creation date: 09.05.2008
38//
39// Modifications:
40// 23.01.2009 V.Ivanchenko add destruction of processes
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 "G4Element.hh"
54#include "G4ProcessManager.hh"
55#include "G4Electron.hh"
56#include "G4Proton.hh"
60
61G4HadronicProcessStore* G4HadronicProcessStore::theInstance = 0;
62
63//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
64
66{
67 if(0 == theInstance) {
68 static G4HadronicProcessStore manager;
69 theInstance = &manager;
70 }
71 return theInstance;
72}
73
74//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
75
77{
78 Clean();
81 delete theEPTestMessenger;
82}
83
84//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
85
87{
88 G4int i;
89 //G4cout << "G4HadronicProcessStore::Clean() Nproc= " << n_proc
90 // << " Nextra= " << n_extra << G4endl;
91 if(n_proc > 0) {
92 for (i=0; i<n_proc; ++i) {
93 if( process[i] ) {
94 //G4cout << "G4HadronicProcessStore::Clean() delete hadronic " << i << G4endl;
95 //G4cout << process[i]->GetProcessName() << G4endl;
96 G4HadronicProcess* p = process[i];
97 process[i] = 0;
98 delete p;
99 }
100 }
101 }
102 if(n_extra > 0) {
103 for(i=0; i<n_extra; ++i) {
104 if(extraProcess[i]) {
105 //G4cout << "G4HadronicProcessStore::Clean() delete extra "
106 // << i << G4endl;
107 //G4cout << extraProcess[i]->GetProcessName() << G4endl;
108 G4VProcess* p = extraProcess[i];
109 extraProcess[i] = 0;
110 delete p;
111 }
112 }
113 }
114 //G4cout << "G4HadronicProcessStore::Clean() done" << G4endl;
115 n_extra = 0;
116 n_proc = 0;
117}
118
119//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
120
121G4HadronicProcessStore::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;
131 theEPTestMessenger = new G4HadronicEPTestMessenger(this);
132}
133//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
134
136 const G4ParticleDefinition* part,
137 G4double energy,
138 const G4VProcess* proc,
139 const G4Element* element)
140{
141 G4double cross = 0.;
142 G4int subType = proc->GetProcessSubType();
143 if (subType == fHadronElastic)
144 cross = GetElasticCrossSectionPerAtom(part,energy,element);
145 else if (subType == fHadronInelastic)
146 cross = GetInelasticCrossSectionPerAtom(part,energy,element);
147 else if (subType == fCapture)
148 cross = GetCaptureCrossSectionPerAtom(part,energy,element);
149 else if (subType == fFission)
150 cross = GetFissionCrossSectionPerAtom(part,energy,element);
151 else if (subType == fChargeExchange)
152 cross = GetChargeExchangeCrossSectionPerAtom(part,energy,element);
153 return cross;
154}
155
156//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
157
159 const G4ParticleDefinition* part,
160 G4double energy,
161 const G4VProcess* proc,
162 const G4Material* material)
163{
164 G4double cross = 0.;
165 G4int subType = proc->GetProcessSubType();
166 if (subType == fHadronElastic)
167 cross = GetElasticCrossSectionPerVolume(part,energy,material);
168 else if (subType == fHadronInelastic)
169 cross = GetInelasticCrossSectionPerVolume(part,energy,material);
170 else if (subType == fCapture)
171 cross = GetCaptureCrossSectionPerVolume(part,energy,material);
172 else if (subType == fFission)
173 cross = GetFissionCrossSectionPerVolume(part,energy,material);
174 else if (subType == fChargeExchange)
175 cross = GetChargeExchangeCrossSectionPerVolume(part,energy,material);
176 return cross;
177}
178
179//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
180
182 const G4ParticleDefinition *aParticle,
183 G4double kineticEnergy,
184 const G4Material *material)
185{
186 G4double cross = 0.0;
187 const G4ElementVector* theElementVector = material->GetElementVector();
188 const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
189 size_t nelm = material->GetNumberOfElements();
190 for (size_t i=0; i<nelm; ++i) {
191 const G4Element* elm = (*theElementVector)[i];
192 cross += theAtomNumDensityVector[i]*
193 GetElasticCrossSectionPerAtom(aParticle,kineticEnergy,elm,material);
194 }
195 return cross;
196}
197
198//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
199
201 const G4ParticleDefinition *aParticle,
202 G4double kineticEnergy,
203 const G4Element *anElement, const G4Material* mat)
204{
206 G4double cross = 0.0;
207 localDP.SetKineticEnergy(kineticEnergy);
208 if(hp) {
209 cross = hp->GetElementCrossSection(&localDP,anElement,mat);
210 }
211 return cross;
212}
213
214//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
215
218 G4double,
219 G4int, G4int)
220{
221 return 0.0;
222}
223
224//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
225
227 const G4ParticleDefinition *aParticle,
228 G4double kineticEnergy,
229 const G4Material *material)
230{
231 G4double cross = 0.0;
232 const G4ElementVector* theElementVector = material->GetElementVector();
233 const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
234 size_t nelm = material->GetNumberOfElements();
235 for (size_t i=0; i<nelm; ++i) {
236 const G4Element* elm = (*theElementVector)[i];
237 cross += theAtomNumDensityVector[i]*
238 GetInelasticCrossSectionPerAtom(aParticle,kineticEnergy,elm,material);
239 }
240 return cross;
241}
242
243//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
244
246 const G4ParticleDefinition *aParticle,
247 G4double kineticEnergy,
248 const G4Element *anElement, const G4Material* mat)
249{
251 localDP.SetKineticEnergy(kineticEnergy);
252 G4double cross = 0.0;
253 if(hp) {
254 cross = hp->GetElementCrossSection(&localDP,anElement,mat);
255 }
256 return cross;
257}
258
259//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
260
262 const G4ParticleDefinition *,
263 G4double,
264 G4int, G4int)
265{
266 return 0.0;
267}
268
269//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
270
272 const G4ParticleDefinition *aParticle,
273 G4double kineticEnergy,
274 const G4Material *material)
275{
276 G4double cross = 0.0;
277 const G4ElementVector* theElementVector = material->GetElementVector();
278 const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
279 size_t nelm = material->GetNumberOfElements();
280 for (size_t i=0; i<nelm; ++i) {
281 const G4Element* elm = (*theElementVector)[i];
282 cross += theAtomNumDensityVector[i]*
283 GetCaptureCrossSectionPerAtom(aParticle,kineticEnergy,elm,material);
284 }
285 return cross;
286}
287
288//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
289
291 const G4ParticleDefinition *aParticle,
292 G4double kineticEnergy,
293 const G4Element *anElement, const G4Material* mat)
294{
295 G4HadronicProcess* hp = FindProcess(aParticle, fCapture);
296 localDP.SetKineticEnergy(kineticEnergy);
297 G4double cross = 0.0;
298 if(hp) {
299 cross = hp->GetElementCrossSection(&localDP,anElement,mat);
300 }
301 return cross;
302}
303
304//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
305
307 const G4ParticleDefinition *,
308 G4double,
309 G4int, G4int)
310{
311 return 0.0;
312}
313
314//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
315
317 const G4ParticleDefinition *aParticle,
318 G4double kineticEnergy,
319 const G4Material *material)
320{
321 G4double cross = 0.0;
322 const G4ElementVector* theElementVector = material->GetElementVector();
323 const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
324 size_t nelm = material->GetNumberOfElements();
325 for (size_t i=0; i<nelm; i++) {
326 const G4Element* elm = (*theElementVector)[i];
327 cross += theAtomNumDensityVector[i]*
328 GetFissionCrossSectionPerAtom(aParticle,kineticEnergy,elm,material);
329 }
330 return cross;
331}
332
333//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
334
336 const G4ParticleDefinition *aParticle,
337 G4double kineticEnergy,
338 const G4Element *anElement, const G4Material* mat)
339{
340 G4HadronicProcess* hp = FindProcess(aParticle, fFission);
341 localDP.SetKineticEnergy(kineticEnergy);
342 G4double cross = 0.0;
343 if(hp) {
344 cross = hp->GetElementCrossSection(&localDP,anElement,mat);
345 }
346 return cross;
347}
348
349//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
350
352 const G4ParticleDefinition *,
353 G4double,
354 G4int, G4int)
355{
356 return 0.0;
357}
358
359//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
360
362 const G4ParticleDefinition *aParticle,
363 G4double kineticEnergy,
364 const G4Material *material)
365{
366 G4double cross = 0.0;
367 const G4ElementVector* theElementVector = material->GetElementVector();
368 const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
369 size_t nelm = material->GetNumberOfElements();
370 for (size_t i=0; i<nelm; ++i) {
371 const G4Element* elm = (*theElementVector)[i];
372 cross += theAtomNumDensityVector[i]*
373 GetChargeExchangeCrossSectionPerAtom(aParticle,kineticEnergy,elm,material);
374 }
375 return cross;
376}
377
378//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
379
381 const G4ParticleDefinition *aParticle,
382 G4double kineticEnergy,
383 const G4Element *anElement, const G4Material* mat)
384{
386 localDP.SetKineticEnergy(kineticEnergy);
387 G4double cross = 0.0;
388 if(hp) {
389 cross = hp->GetElementCrossSection(&localDP,anElement,mat);
390 }
391 return cross;
392}
393
394//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
395
397 const G4ParticleDefinition *,
398 G4double,
399 G4int, G4int)
400{
401 return 0.0;
402}
403
404//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
405
407{
408 if(0 < n_proc) {
409 for(G4int i=0; i<n_proc; ++i) {
410 if(process[i] == proc) { return; }
411 }
412 }
413 // G4cout << "G4HadronicProcessStore::Register hadronic " << n_proc
414 // << " " << proc->GetProcessName() << G4endl;
415 ++n_proc;
416 process.push_back(proc);
417}
418
419//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
420
422 const G4ParticleDefinition* part)
423{
424 G4int i=0;
425 for(; i<n_proc; ++i) {if(process[i] == proc) break;}
426 G4int j=0;
427 for(; j<n_part; ++j) {if(particle[j] == part) break;}
428
429 if(j == n_part) {
430 ++n_part;
431 particle.push_back(part);
432 wasPrinted.push_back(0);
433 }
434
435 // the pair should be added?
436 if(i < n_proc) {
437 std::multimap<PD,HP,std::less<PD> >::iterator it;
438 for(it=p_map.lower_bound(part); it!=p_map.upper_bound(part); ++it) {
439 if(it->first == part) {
440 HP process2 = (it->second);
441 if(proc == process2) { return; }
442 }
443 }
444 }
445
446 p_map.insert(std::multimap<PD,HP>::value_type(part,proc));
447}
448
449//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
450
453{
454 G4int i=0;
455 for(; i<n_proc; ++i) {if(process[i] == proc) { break; }}
456 G4int k=0;
457 for(; k<n_model; ++k) {if(model[k] == mod) { break; }}
458
459 m_map.insert(std::multimap<HP,HI>::value_type(proc,mod));
460
461 if(k == n_model) {
462 ++n_model;
463 model.push_back(mod);
464 modelName.push_back(mod->GetModelName());
465 }
466}
467
468//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
469
471{
472 if(0 == n_proc) return;
473 for(G4int i=0; i<n_proc; ++i) {
474 if(process[i] == proc) {
475 process[i] = 0;
476 return;
477 }
478 }
479}
480
481//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
482
484{
485 if(0 < n_extra) {
486 for(G4int i=0; i<n_extra; ++i) {
487 if(extraProcess[i] == proc) { return; }
488 }
489 }
490 //G4cout << "Extra Process: " << n_extra << " " << proc->GetProcessName()
491 // << " " << proc << G4endl;
492
493 n_extra++;
494 extraProcess.push_back(proc);
495}
496
497//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
498
500 G4VProcess* proc,
501 const G4ParticleDefinition* part)
502{
503 G4int i=0;
504 for(; i<n_extra; ++i) { if(extraProcess[i] == proc) { break; } }
505 G4int j=0;
506 for(; j<n_part; ++j) { if(particle[j] == part) { break; } }
507
508 if(j == n_part) {
509 ++n_part;
510 particle.push_back(part);
511 wasPrinted.push_back(0);
512 }
513
514 // the pair should be added?
515 if(i < n_extra) {
516 std::multimap<PD,G4VProcess*,std::less<PD> >::iterator it;
517 for(it=ep_map.lower_bound(part); it!=ep_map.upper_bound(part); ++it) {
518 if(it->first == part) {
519 G4VProcess* process2 = (it->second);
520 if(proc == process2) { return; }
521 }
522 }
523 }
524
525 ep_map.insert(std::multimap<PD,G4VProcess*>::value_type(part,proc));
526}
527
528//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
529
531{
532 //G4cout << "Deregister Extra Process: " << proc << " "<<proc->GetProcessName()<< G4endl;
533 if(0 == n_extra) { return; }
534 for(G4int i=0; i<n_extra; ++i) {
535 if(extraProcess[i] == proc) {
536 extraProcess[i] = 0;
537 //G4cout << "Extra Process: " << i << " is deregisted " << G4endl;
538 return;
539 }
540 }
541}
542
543
545{
546 // Trigger particle/process/model printout only when last particle is
547 // registered
548 if(buildTableStart && part == particle[n_part - 1]) {
549 buildTableStart = false;
550 Dump(verbose);
551 if (getenv("G4PhysListDocDir") ) DumpHtml();
552 }
553}
554
555
557{
558 // Automatic generation of html documentation page for physics lists
559 // List processes, models and cross sections for the most important
560 // particles in descending order of importance
561
562 char* dirName = getenv("G4PhysListDocDir");
563 char* physListName = getenv("G4PhysListName");
564 if (dirName && physListName) {
565
566 // Open output file with path name
567 G4String pathName = G4String(dirName) + "/" + G4String(physListName) + ".html";
568 std::ofstream outFile;
569 outFile.open(pathName);
570
571 // Write physics list summary file
572 outFile << "<html>\n";
573 outFile << "<head>\n";
574 outFile << "<title>Physics List Summary</title>\n";
575 outFile << "</head>\n";
576 outFile << "<body>\n";
577 outFile << "<h2> Summary of Hadronic Processes, Models and Cross Sections for Physics List "
578 << G4String(physListName) << "</h2>\n";
579 outFile << "<ul>\n";
580
581 PrintHtml(G4Proton::Proton(), outFile);
582 PrintHtml(G4Neutron::Neutron(), outFile);
583 PrintHtml(G4PionPlus::PionPlus(), outFile);
585 PrintHtml(G4Gamma::Gamma(), outFile);
587// PrintHtml(G4MuonMinus::MuonMinus(), outFile);
591 PrintHtml(G4Lambda::Lambda(), outFile);
592 PrintHtml(G4Alpha::Alpha(), outFile);
593
594 outFile << "</ul>\n";
595 outFile << "</body>\n";
596 outFile << "</html>\n";
597 outFile.close();
598 }
599}
600
601
603 std::ofstream& outFile)
604{
605 // Automatic generation of html documentation page for physics lists
606 // List processes for the most important particles in descending order
607 // of importance
608
609 outFile << "<br> <li><h2><font color=\" ff0000 \">"
610 << theParticle->GetParticleName() << "</font></h2></li>\n";
611
612 typedef std::multimap<PD,HP,std::less<PD> > PDHPmap;
613 typedef std::multimap<HP,HI,std::less<HP> > HPHImap;
614
615 std::pair<PDHPmap::iterator, PDHPmap::iterator> itpart =
616 p_map.equal_range(theParticle);
617
618 // Loop over processes assigned to particle
619
620 G4HadronicProcess* theProcess;
621 for (PDHPmap::iterator it = itpart.first; it != itpart.second; ++it) {
622 theProcess = (*it).second;
623 outFile << "<br> &nbsp;&nbsp; <b><font color=\" 0000ff \">process : <a href=\""
624 << theProcess->GetProcessName() << ".html\"> "
625 << theProcess->GetProcessName() << "</a></font></b>\n";
626 outFile << "<ul>\n";
627 outFile << " <li><b><font color=\" 00AA00 \">models : </font></b>\n";
628
629 // Loop over models assigned to process
630 std::pair<HPHImap::iterator, HPHImap::iterator> itmod =
631 m_map.equal_range(theProcess);
632
633 outFile << " <ul>\n";
634 for (HPHImap::iterator jt = itmod.first; jt != itmod.second; ++jt) {
635 outFile << " <li><b><a href=\"" << (*jt).second->GetModelName() << ".html\"> "
636 << (*jt).second->GetModelName() << "</a>"
637 << " from " << (*jt).second->GetMinEnergy()/GeV
638 << " GeV to " << (*jt).second->GetMaxEnergy()/GeV
639 << " GeV </b></li>\n";
640
641 // Print ModelDescription, ignore that we overwrite files n-times.
642 PrintModelHtml((*jt).second);
643
644 }
645 outFile << " </ul>\n";
646 outFile << " </li>\n";
647
648 // List cross sections assigned to process
649 outFile << " <li><b><font color=\" 00AA00 \">cross sections : </font></b>\n";
650 outFile << " <ul>\n";
651 theProcess->GetCrossSectionDataStore()->DumpHtml(*theParticle, outFile);
652 // << " \n";
653 outFile << " </ul>\n";
654
655 outFile << " </li>\n";
656 outFile << "</ul>\n";
657 }
658}
659
661{
662 G4String dirName(getenv("G4PhysListDocDir"));
663 G4String pathName = dirName + "/" + mod->GetModelName() + ".html";
664 std::ofstream outModel;
665 outModel.open(pathName);
666 outModel << "<html>\n";
667 outModel << "<head>\n";
668 outModel << "<title>Description of " << mod->GetModelName() << "</title>\n";
669 outModel << "</head>\n";
670 outModel << "<body>\n";
671
672 mod->ModelDescription(outModel);
673
674 outModel << "</body>\n";
675 outModel << "</html>\n";
676
677}
679{
680 if(level > 0) {
681 G4cout << "=============================================================="
682 << "=============================="
683 << G4endl;
684 G4cout << " HADRONIC PROCESSES SUMMARY (verbose level " << level
685 << ")" << G4endl;
686 }
687 for(G4int i=0; i<n_part; ++i) {
688 PD part = particle[i];
689 G4String pname = part->GetParticleName();
690 G4bool yes = false;
691 if(level >= 2) yes = true;
692
693 else if(level == 1 && (pname == "proton" ||
694 pname == "neutron" ||
695 pname == "pi+" ||
696 pname == "pi-" ||
697 pname == "gamma" ||
698 pname == "e+" ||
699 pname == "e-" ||
700 pname == "mu+" ||
701 pname == "mu-" ||
702 pname == "kaon+" ||
703 pname == "kaon-" ||
704 pname == "lambda" ||
705 pname == "GenericIon" ||
706 pname == "anti_neutron" ||
707 pname == "anti_proton")) yes = true;
708 if(yes) {
709 // main processes
710 std::multimap<PD,HP,std::less<PD> >::iterator it;
711
712 for(it=p_map.lower_bound(part); it!=p_map.upper_bound(part); ++it) {
713 if(it->first == part) {
714 HP proc = (it->second);
715 G4int j=0;
716 for(; j<n_proc; ++j) {
717 if(process[j] == proc) {
718 Print(j, i);
719 }
720 }
721 }
722 }
723 // extra processes
724 std::multimap<PD,G4VProcess*,std::less<PD> >::iterator itp;
725 for(itp=ep_map.lower_bound(part); itp!=ep_map.upper_bound(part); ++itp) {
726 if(itp->first == part) {
727 G4VProcess* proc = (itp->second);
728 if(wasPrinted[i] == 0) {
729 wasPrinted[i] = 1;
730 G4cout<<G4endl;
731 G4cout << " Hadronic Processes for <"
732 <<part->GetParticleName() << ">" << G4endl;
733 }
734 G4cout << " " << proc->GetProcessName() << G4endl;
735 }
736 }
737 }
738 }
739 if(level > 0) {
740 G4cout << "=============================================================="
741 << "=============================="
742 << G4endl;
743 }
744}
745
746
747void G4HadronicProcessStore::Print(G4int idxProc, G4int idxPart)
748{
749 G4HadronicProcess* proc = process[idxProc];
750 const G4ParticleDefinition* part = particle[idxPart];
751 if(wasPrinted[idxPart] == 0) {
752 wasPrinted[idxPart] = 1;
753 G4cout<<G4endl;
754 G4cout << " Hadronic Processes for <"
755 << part->GetParticleName() << ">" << G4endl;
756 G4cout << " ------------------------"
757 << "-----------" << G4endl;
758 }
759 HI hi = 0;
760 G4bool first;
761 std::multimap<HP,HI,std::less<HP> >::iterator ih;
762 G4cout << std::setw(20) << proc->GetProcessName()
763 << " Models: ";
764 first = true;
765 for(ih=m_map.lower_bound(proc); ih!=m_map.upper_bound(proc); ++ih) {
766 if(ih->first == proc) {
767 hi = ih->second;
768 G4int i=0;
769 for(; i<n_model; ++i) {
770 if(model[i] == hi) { break; }
771 }
772 if(!first) G4cout << " ";
773 first = false;
774 G4cout << std::setw(28) << modelName[i]
775 << ": Emin(GeV)= "
776 << std::setw(4) << hi->GetMinEnergy()/GeV
777 << " Emax(GeV)= "
778 << hi->GetMaxEnergy()/GeV
779 << G4endl;
780 }
781 }
782
783 G4cout << G4endl;
784 G4cout << std::setw(20) << proc->GetProcessName()
785 << " Crs sctns: ";
787 csds->DumpPhysicsTable(*part);
788
789}
790
791//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
792
794{
795 verbose = val;
796 G4int i;
797 for(i=0; i<n_proc; ++i) {
798 if(process[i]) { process[i]->SetVerboseLevel(val); }
799 }
800 for(i=0; i<n_model; ++i) {
801 if(model[i]) { model[i]->SetVerboseLevel(val); }
802 }
803}
804
805//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
806
808{
809 return verbose;
810}
811
812//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
813
815 const G4ParticleDefinition* part, G4HadronicProcessType subType)
816{
817 bool isNew = false;
818 G4HadronicProcess* hp = 0;
819
820 if(part != currentParticle) {
821 isNew = true;
822 currentParticle = part;
823 localDP.SetDefinition(part);
824 } else if(!currentProcess) {
825 isNew = true;
826 } else if(subType == currentProcess->GetProcessSubType()) {
827 hp = currentProcess;
828 } else {
829 isNew = true;
830 }
831
832 if(isNew) {
833 std::multimap<PD,HP,std::less<PD> >::iterator it;
834 for(it=p_map.lower_bound(part); it!=p_map.upper_bound(part); ++it) {
835 if(it->first == part && subType == (it->second)->GetProcessSubType()) {
836 hp = it->second;
837 break;
838 }
839 }
840 currentProcess = hp;
841 }
842
843 return hp;
844}
845
847{
848 G4cout << " Setting energy/momentum report level to " << level
849 << " for " << process.size() << " hadronic processes " << G4endl;
850 for (G4int i = 0; i < G4int(process.size()); ++i) {
851 process[i]->SetEpReportLevel(level);
852 }
853}
854
856{
857 G4cout << " Setting absolute energy/momentum test level to " << abslevel << G4endl;
858 G4double rellevel = 0.0;
859 G4HadronicProcess* theProcess = 0;
860 for (G4int i = 0; i < G4int(process.size()); ++i) {
861 theProcess = process[i];
862 rellevel = theProcess->GetEnergyMomentumCheckLevels().first;
863 theProcess->SetEnergyMomentumCheckLevels(rellevel, abslevel);
864 }
865}
866
868{
869 G4cout << " Setting relative energy/momentum test level to " << rellevel << G4endl;
870 G4double abslevel = 0.0;
871 G4HadronicProcess* theProcess = 0;
872 for (G4int i = 0; i < G4int(process.size()); ++i) {
873 theProcess = process[i];
874 abslevel = theProcess->GetEnergyMomentumCheckLevels().second;
875 theProcess->SetEnergyMomentumCheckLevels(rellevel, abslevel);
876 }
877}
std::vector< G4Element * > G4ElementVector
G4HadronicProcessType
@ fChargeExchange
@ fHadronElastic
@ fHadronInelastic
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
static G4CrossSectionDataSetRegistry * Instance()
void DumpPhysicsTable(const G4ParticleDefinition &)
void DumpHtml(const G4ParticleDefinition &, std::ofstream &)
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetKineticEnergy(G4double aEnergy)
static G4Electron * Electron()
Definition: G4Electron.cc:94
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4HadronicInteractionRegistry * Instance()
virtual void ModelDescription(std::ostream &outFile) const
const G4String & GetModelName() const
void DeRegister(G4HadronicProcess *)
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 GetFissionCrossSectionPerAtom(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Element *anElement, const G4Material *mat=0)
G4double GetCrossSectionPerAtom(const G4ParticleDefinition *particle, G4double kineticEnergy, const G4VProcess *process, const G4Element *element)
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 *)
void RegisterExtraProcess(G4VProcess *)
G4double GetElasticCrossSectionPerVolume(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Material *material)
G4double GetCaptureCrossSectionPerAtom(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Element *anElement, const G4Material *mat=0)
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 GetChargeExchangeCrossSectionPerAtom(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Element *anElement, const G4Material *mat=0)
G4double GetInelasticCrossSectionPerAtom(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Element *anElement, const G4Material *mat=0)
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)
std::pair< G4double, G4double > GetEnergyMomentumCheckLevels() const
G4CrossSectionDataStore * GetCrossSectionDataStore()
G4double GetElementCrossSection(const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=0)
void SetEnergyMomentumCheckLevels(G4double relativeLevel, G4double absoluteLevel)
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
static G4Lambda * Lambda()
Definition: G4Lambda.cc:108
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:189
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:205
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
const G4String & GetParticleName() const
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
static G4Positron * Positron()
Definition: G4Positron.cc:94
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4int GetProcessSubType() const
Definition: G4VProcess.hh:397
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379