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