Geant4  10.01.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 86515 2014-11-13 09:11:45Z 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"
58 #include "G4ParticleTable.hh"
62 
63 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
64 
66 {
68  return instance.Instance();
69 }
70 
71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
72 
74 {
75  Clean();
78  delete theEPTestMessenger;
79 }
80 
81 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
82 
84 {
85  G4int i;
86  //G4cout << "G4HadronicProcessStore::Clean() Nproc= " << n_proc
87  // << " Nextra= " << n_extra << G4endl;
88  for (i=0; i<n_proc; ++i) {
89  if( process[i] ) {
90  //G4cout << "G4HadronicProcessStore::Clean() delete hadronic "
91  // << i << " " << process[i]->GetProcessName() << G4endl;
92  G4HadronicProcess* p = process[i];
93  DeRegister(p);
94  delete p;
95  }
96  }
97  for(i=0; i<n_extra; ++i) {
98  if(extraProcess[i]) {
99  // G4cout << "G4HadronicProcessStore::Clean() delete extra proc "
100  //<< i << " " << extraProcess[i]->GetProcessName() << G4endl;
101  delete extraProcess[i];
102  extraProcess[i] = 0;
103  }
104  }
105  //G4cout << "G4HadronicProcessStore::Clean() done" << G4endl;
106  n_extra = 0;
107  n_proc = 0;
108 }
109 
110 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
111 
113 {
114  n_proc = 0;
115  n_part = 0;
116  n_model= 0;
117  n_extra= 0;
118  currentProcess = 0;
119  currentParticle = 0;
120  theGenericIon =
122  verbose = 1;
123  buildTableStart = true;
125 }
126 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
127 
129  const G4ParticleDefinition* part,
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,
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 {
200  G4HadronicProcess* hp = FindProcess(aParticle, fHadronElastic);
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 
212  const G4ParticleDefinition*,
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  if(0 < n_proc) {
408  for(G4int i=0; i<n_proc; ++i) {
409  if(process[i] == proc) { return; }
410  }
411  }
412  if(1 < verbose) {
413  G4cout << "G4HadronicProcessStore::Register hadronic " << n_proc
414  << " " << proc->GetProcessName() << G4endl;
415  }
416  ++n_proc;
417  process.push_back(proc);
418 }
419 
420 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
421 
423  const G4ParticleDefinition* part)
424 {
425  G4int i=0;
426  for(; i<n_proc; ++i) {if(process[i] == proc) break;}
427  G4int j=0;
428  for(; j<n_part; ++j) {if(particle[j] == part) break;}
429 
430  if(1 < verbose) {
431  G4cout << "G4HadronicProcessStore::RegisterParticle "
432  << part->GetParticleName()
433  << " for " << proc->GetProcessName() << G4endl;
434  }
435  if(j == n_part) {
436  ++n_part;
437  particle.push_back(part);
438  wasPrinted.push_back(0);
439  }
440 
441  // the pair should be added?
442  if(i < n_proc) {
443  std::multimap<PD,HP,std::less<PD> >::iterator it;
444  for(it=p_map.lower_bound(part); it!=p_map.upper_bound(part); ++it) {
445  if(it->first == part) {
446  HP process2 = (it->second);
447  if(proc == process2) { return; }
448  }
449  }
450  }
451 
452  p_map.insert(std::multimap<PD,HP>::value_type(part,proc));
453 }
454 
455 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
456 
459 {
460  G4int i=0;
461  for(; i<n_proc; ++i) {if(process[i] == proc) { break; }}
462  G4int k=0;
463  for(; k<n_model; ++k) {if(model[k] == mod) { break; }}
464 
465  m_map.insert(std::multimap<HP,HI>::value_type(proc,mod));
466 
467  if(k == n_model) {
468  ++n_model;
469  model.push_back(mod);
470  modelName.push_back(mod->GetModelName());
471  }
472 }
473 
474 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
475 
477 {
478  for(G4int i=0; i<n_proc; ++i) {
479  if(process[i] == proc) {
480  process[i] = 0;
482  return;
483  }
484  }
485 }
486 
487 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
488 
490 {
491  if(0 < n_extra) {
492  for(G4int i=0; i<n_extra; ++i) {
493  if(extraProcess[i] == proc) { return; }
494  }
495  }
496  G4HadronicProcess* hproc = reinterpret_cast<G4HadronicProcess*>(proc);
497  if(hproc) {
498  for(G4int i=0; i<n_proc; ++i) {
499  if(process[i] == hproc) { return; }
500  }
501  }
502  if(1 < verbose) {
503  G4cout << "Extra Process: " << n_extra
504  << " " << proc->GetProcessName() << G4endl;
505  }
506  ++n_extra;
507  extraProcess.push_back(proc);
508 }
509 
510 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
511 
513  G4VProcess* proc,
514  const G4ParticleDefinition* part)
515 {
516  G4int i=0;
517  for(; i<n_extra; ++i) { if(extraProcess[i] == proc) { break; } }
518  G4int j=0;
519  for(; j<n_part; ++j) { if(particle[j] == part) { break; } }
520 
521  if(j == n_part) {
522  ++n_part;
523  particle.push_back(part);
524  wasPrinted.push_back(0);
525  }
526 
527  // the pair should be added?
528  if(i < n_extra) {
529  std::multimap<PD,G4VProcess*,std::less<PD> >::iterator it;
530  for(it=ep_map.lower_bound(part); it!=ep_map.upper_bound(part); ++it) {
531  if(it->first == part) {
532  G4VProcess* process2 = (it->second);
533  if(proc == process2) { return; }
534  }
535  }
536  }
537 
538  ep_map.insert(std::multimap<PD,G4VProcess*>::value_type(part,proc));
539 }
540 
541 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
542 
544 {
545  for(G4int i=0; i<n_extra; ++i) {
546  if(extraProcess[i] == proc) {
547  extraProcess[i] = 0;
548  if(1 < verbose) {
549  G4cout << "Extra Process: " << i << " "
550  <<proc->GetProcessName()<< " is deregisted " << G4endl;
551  }
552  return;
553  }
554  }
555 }
556 
557 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
558 
560 {
561  // Trigger particle/process/model printout only when last particle is
562  // registered
563  if(buildTableStart && part == particle[n_part - 1]) {
564  buildTableStart = false;
565  Dump(verbose);
566  if (getenv("G4PhysListDocDir") ) DumpHtml();
567  }
568 }
569 
570 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
571 
573 {
574  // Automatic generation of html documentation page for physics lists
575  // List processes, models and cross sections for the most important
576  // particles in descending order of importance
577 
578  char* dirName = getenv("G4PhysListDocDir");
579  char* physListName = getenv("G4PhysListName");
580  if (dirName && physListName) {
581 
582  // Open output file with path name
583  G4String pathName = G4String(dirName) + "/" + G4String(physListName) + ".html";
584  std::ofstream outFile;
585  outFile.open(pathName);
586 
587  // Write physics list summary file
588  outFile << "<html>\n";
589  outFile << "<head>\n";
590  outFile << "<title>Physics List Summary</title>\n";
591  outFile << "</head>\n";
592  outFile << "<body>\n";
593  outFile << "<h2> Summary of Hadronic Processes, Models and Cross Sections for Physics List "
594  << G4String(physListName) << "</h2>\n";
595  outFile << "<ul>\n";
596 
597  PrintHtml(G4Proton::Proton(), outFile);
598  PrintHtml(G4Neutron::Neutron(), outFile);
599  PrintHtml(G4PionPlus::PionPlus(), outFile);
600  PrintHtml(G4PionMinus::PionMinus(), outFile);
601  PrintHtml(G4Gamma::Gamma(), outFile);
602  PrintHtml(G4Electron::Electron(), outFile);
603 // PrintHtml(G4MuonMinus::MuonMinus(), outFile);
604  PrintHtml(G4Positron::Positron(), outFile);
605  PrintHtml(G4KaonPlus::KaonPlus(), outFile);
606  PrintHtml(G4KaonMinus::KaonMinus(), outFile);
607  PrintHtml(G4Lambda::Lambda(), outFile);
608  PrintHtml(G4Alpha::Alpha(), outFile);
609 
610  outFile << "</ul>\n";
611  outFile << "</body>\n";
612  outFile << "</html>\n";
613  outFile.close();
614  }
615 }
616 
617 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
618 
620  std::ofstream& outFile)
621 {
622  // Automatic generation of html documentation page for physics lists
623  // List processes for the most important particles in descending order
624  // of importance
625 
626  outFile << "<br> <li><h2><font color=\" ff0000 \">"
627  << theParticle->GetParticleName() << "</font></h2></li>\n";
628 
629  typedef std::multimap<PD,HP,std::less<PD> > PDHPmap;
630  typedef std::multimap<HP,HI,std::less<HP> > HPHImap;
631 
632  std::pair<PDHPmap::iterator, PDHPmap::iterator> itpart =
633  p_map.equal_range(theParticle);
634 
635  // Loop over processes assigned to particle
636 
637  G4HadronicProcess* theProcess;
638  for (PDHPmap::iterator it = itpart.first; it != itpart.second; ++it) {
639  theProcess = (*it).second;
640  // description is inline
641  //outFile << "<br> &nbsp;&nbsp; <b><font color=\" 0000ff \">process : <a href=\""
642  // << theProcess->GetProcessName() << ".html\"> "
643  // << theProcess->GetProcessName() << "</a></font></b>\n";
644  outFile << "<br> &nbsp;&nbsp; <b><font color=\" 0000ff \">process : "
645  << theProcess->GetProcessName() << "</font></b>\n";
646  outFile << "<ul>\n";
647  outFile << " <li>";
648  theProcess->ProcessDescription(outFile);
649  outFile << " <li><b><font color=\" 00AA00 \">models : </font></b>\n";
650  // Loop over models assigned to process
651  std::pair<HPHImap::iterator, HPHImap::iterator> itmod =
652  m_map.equal_range(theProcess);
653 
654  outFile << " <ul>\n";
655  for (HPHImap::iterator jt = itmod.first; jt != itmod.second; ++jt) {
656  outFile << " <li><b><a href=\"" << (*jt).second->GetModelName()
657  << ".html\"> "
658  << (*jt).second->GetModelName() << "</a>"
659  << " from " << (*jt).second->GetMinEnergy()/GeV
660  << " GeV to " << (*jt).second->GetMaxEnergy()/GeV
661  << " GeV </b></li>\n";
662 
663  // Print ModelDescription, ignore that we overwrite files n-times.
664  PrintModelHtml((*jt).second);
665 
666  }
667  outFile << " </ul>\n";
668  outFile << " </li>\n";
669 
670  // List cross sections assigned to process
671  outFile << " <li><b><font color=\" 00AA00 \">cross sections : </font></b>\n";
672  outFile << " <ul>\n";
673  theProcess->GetCrossSectionDataStore()->DumpHtml(*theParticle, outFile);
674  // << " \n";
675  outFile << " </ul>\n";
676 
677  outFile << " </li>\n";
678  outFile << "</ul>\n";
679  }
680 }
681 
682 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
683 
684 void
686 {
687  G4String dirName(getenv("G4PhysListDocDir"));
688  G4String pathName = dirName + "/" + mod->GetModelName() + ".html";
689  std::ofstream outModel;
690  outModel.open(pathName);
691  outModel << "<html>\n";
692  outModel << "<head>\n";
693  outModel << "<title>Description of " << mod->GetModelName()
694  << "</title>\n";
695  outModel << "</head>\n";
696  outModel << "<body>\n";
697 
698  mod->ModelDescription(outModel);
699 
700  outModel << "</body>\n";
701  outModel << "</html>\n";
702 
703 }
704 
705 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
706 
708 {
709  if (level == 0) return;
710 
711  G4cout
712  << "\n====================================================================\n"
713  << std::setw(60) << "HADRONIC PROCESSES SUMMARY (verbose level " << level
714  << ")" << G4endl;
715 
716  for (G4int i=0; i<n_part; ++i) {
717  PD part = particle[i];
718  G4String pname = part->GetParticleName();
719  G4bool yes = false;
720 
721  if (level == 1 && (pname == "proton" ||
722  pname == "neutron" ||
723  pname == "deuteron" ||
724  pname == "triton" ||
725  pname == "He3" ||
726  pname == "alpha" ||
727  pname == "pi+" ||
728  pname == "pi-" ||
729  pname == "gamma" ||
730  pname == "e+" ||
731  pname == "e-" ||
732  pname == "mu+" ||
733  pname == "mu-" ||
734  pname == "kaon+" ||
735  pname == "kaon-" ||
736  pname == "lambda" ||
737  pname == "GenericIon" ||
738  pname == "anti_neutron" ||
739  pname == "anti_proton" ||
740  pname == "anti_deuteron" ||
741  pname == "anti_triton" ||
742  pname == "anti_He3" ||
743  pname == "anti_alpha")) yes = true;
744  if (level > 1) yes = true;
745  if (yes) {
746  // main processes
747  std::multimap<PD,HP,std::less<PD> >::iterator it;
748 
749  for (it=p_map.lower_bound(part); it!=p_map.upper_bound(part); ++it) {
750  if (it->first == part) {
751  HP proc = (it->second);
752  G4int j=0;
753  for (; j<n_proc; ++j) {
754  if (process[j] == proc) { Print(j, i); }
755  }
756  }
757  }
758 
759  // extra processes
760  std::multimap<PD,G4VProcess*,std::less<PD> >::iterator itp;
761  for(itp=ep_map.lower_bound(part); itp!=ep_map.upper_bound(part); ++itp) {
762  if(itp->first == part) {
763  G4VProcess* proc = (itp->second);
764  if (wasPrinted[i] == 0) {
765  G4cout << "\n---------------------------------------------------\n"
766  << std::setw(50) << "Hadronic Processes for "
767  << part->GetParticleName() << "\n";
768  wasPrinted[i] = 1;
769  }
770  G4cout << "\n Process: " << proc->GetProcessName() << G4endl;
771  }
772  }
773  }
774  }
775 
776  G4cout << "\n================================================================"
777  << G4endl;
778 }
779 
780 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
781 
783 {
784  G4HadronicProcess* proc = process[idxProc];
785  const G4ParticleDefinition* part = particle[idxPart];
786  if (wasPrinted[idxPart] == 0) {
787  G4cout << "\n---------------------------------------------------\n"
788  << std::setw(50) << "Hadronic Processes for "
789  << part->GetParticleName() << "\n";
790  wasPrinted[idxPart] = 1;
791  }
792 
793  G4cout << "\n Process: " << proc->GetProcessName();
794  HI hi = 0;
795  std::multimap<HP,HI,std::less<HP> >::iterator ih;
796  for(ih=m_map.lower_bound(proc); ih!=m_map.upper_bound(proc); ++ih) {
797  if(ih->first == proc) {
798  hi = ih->second;
799  G4int i=0;
800  for(; i<n_model; ++i) {
801  if(model[i] == hi) { break; }
802  }
803  G4cout << "\n Model: " << std::setw(25) << modelName[i] << ": "
804  << G4BestUnit(hi->GetMinEnergy(), "Energy")
805  << " ---> "
806  << G4BestUnit(hi->GetMaxEnergy(), "Energy");
807  }
808  }
809  G4cout << G4endl;
810 
812  csds->DumpPhysicsTable(*part);
813 }
814 
815 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
816 
818 {
819  verbose = val;
820  G4int i;
821  for(i=0; i<n_proc; ++i) {
822  if(process[i]) { process[i]->SetVerboseLevel(val); }
823  }
824  for(i=0; i<n_model; ++i) {
825  if(model[i]) { model[i]->SetVerboseLevel(val); }
826  }
827 }
828 
829 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
830 
832 {
833  return verbose;
834 }
835 
836 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
837 
839  const G4ParticleDefinition* part, G4HadronicProcessType subType)
840 {
841  bool isNew = false;
842  G4HadronicProcess* hp = 0;
843  localDP.SetDefinition(part);
844 
845  if(part != currentParticle) {
846  const G4ParticleDefinition* p = part;
847  if(p->GetBaryonNumber() > 4 && p->GetParticleType() == "nucleus") {
848  p = theGenericIon;
849  }
850  if(p != currentParticle) {
851  isNew = true;
852  currentParticle = p;
853  }
854  }
855  if(!isNew) {
856  if(!currentProcess) {
857  isNew = true;
858  } else if(subType == currentProcess->GetProcessSubType()) {
859  hp = currentProcess;
860  } else {
861  isNew = true;
862  }
863  }
864  if(isNew) {
865  std::multimap<PD,HP,std::less<PD> >::iterator it;
866  for(it=p_map.lower_bound(currentParticle);
867  it!=p_map.upper_bound(currentParticle); ++it) {
868  if(it->first == currentParticle &&
869  subType == (it->second)->GetProcessSubType()) {
870  hp = it->second;
871  break;
872  }
873  }
874  currentProcess = hp;
875  }
876  return hp;
877 }
878 
879 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
880 
882 {
883  G4cout << " Setting energy/momentum report level to " << level
884  << " for " << process.size() << " hadronic processes " << G4endl;
885  for (G4int i = 0; i < G4int(process.size()); ++i) {
886  process[i]->SetEpReportLevel(level);
887  }
888 }
889 
890 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
891 
893 {
894  G4cout << " Setting absolute energy/momentum test level to " << abslevel
895  << G4endl;
896  G4double rellevel = 0.0;
897  G4HadronicProcess* theProcess = 0;
898  for (G4int i = 0; i < G4int(process.size()); ++i) {
899  theProcess = process[i];
900  rellevel = theProcess->GetEnergyMomentumCheckLevels().first;
901  theProcess->SetEnergyMomentumCheckLevels(rellevel, abslevel);
902  }
903 }
904 
905 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
906 
908 {
909  G4cout << " Setting relative energy/momentum test level to " << rellevel
910  << G4endl;
911  G4double abslevel = 0.0;
912  G4HadronicProcess* theProcess = 0;
913  for (G4int i = 0; i < G4int(process.size()); ++i) {
914  theProcess = process[i];
915  abslevel = theProcess->GetEnergyMomentumCheckLevels().second;
916  theProcess->SetEnergyMomentumCheckLevels(rellevel, abslevel);
917  }
918 }
919 
920 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
G4double GetMinEnergy() const
void DeRegisterExtraProcess(G4VProcess *)
void PrintModelHtml(const G4HadronicInteraction *model) const
std::vector< G4String > modelName
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
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()
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
G4double GetFissionCrossSectionPerAtom(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Element *anElement, const G4Material *mat=0)
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
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:206
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
const G4String & GetParticleType() const
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
virtual void ProcessDescription(std::ostream &outFile) const
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)
static G4ParticleTable * GetParticleTable()
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)
static MCTruthManager * instance
G4double GetCaptureCrossSectionPerVolume(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Material *material)
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
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 *)