Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/management/src/G4HadronicProcessStore.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

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