G4HadronicProcessStore.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 // $Id$
00027 //
00028 // -------------------------------------------------------------------
00029 //
00030 // GEANT4 Class file
00031 //
00032 //
00033 // File name:     G4HadronicProcessStore
00034 //
00035 // Author:        Vladimir Ivanchenko
00036 //
00037 // Creation date: 09.05.2008
00038 //
00039 // Modifications:
00040 // 23.01.2009 V.Ivanchenko add destruction of processes
00041 //
00042 // Class Description:
00043 // Singleton to store hadronic processes, to provide access to processes
00044 // and to printout information about processes
00045 //
00046 // -------------------------------------------------------------------
00047 //
00048 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00049 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00050 
00051 #include "G4HadronicProcessStore.hh"
00052 #include "G4SystemOfUnits.hh"
00053 #include "G4Element.hh"
00054 #include "G4ProcessManager.hh"
00055 #include "G4Electron.hh"
00056 #include "G4Proton.hh"
00057 #include "G4HadronicInteractionRegistry.hh"
00058 #include "G4CrossSectionDataSetRegistry.hh"
00059 #include "G4HadronicEPTestMessenger.hh"
00060 
00061 G4HadronicProcessStore* G4HadronicProcessStore::theInstance = 0;
00062 
00063 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00064 
00065 G4HadronicProcessStore* G4HadronicProcessStore::Instance()
00066 {
00067   if(0 == theInstance) {
00068     static G4HadronicProcessStore manager;
00069     theInstance = &manager;
00070   }
00071   return theInstance;
00072 }
00073 
00074 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00075 
00076 G4HadronicProcessStore::~G4HadronicProcessStore()
00077 {
00078   Clean();
00079   G4HadronicInteractionRegistry::Instance()->Clean();
00080   G4CrossSectionDataSetRegistry::Instance()->Clean();
00081   delete theEPTestMessenger;
00082 }
00083 
00084 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00085 
00086 void G4HadronicProcessStore::Clean()
00087 {
00088   G4int i;
00089   //G4cout << "G4HadronicProcessStore::Clean() Nproc= " << n_proc
00090   //     << "  Nextra= " << n_extra << G4endl;
00091   if(n_proc > 0) {
00092     for (i=0; i<n_proc; ++i) {
00093       if( process[i] ) {
00094         //G4cout << "G4HadronicProcessStore::Clean() delete hadronic " << i << G4endl;
00095         //G4cout <<  process[i]->GetProcessName() << G4endl;
00096         G4HadronicProcess* p = process[i]; 
00097         process[i] = 0;
00098         delete p;
00099       }
00100     }
00101   }
00102   if(n_extra > 0) {
00103     for(i=0; i<n_extra; ++i) {
00104       if(extraProcess[i]) {
00105         //G4cout << "G4HadronicProcessStore::Clean() delete extra "  
00106         //       << i << G4endl;
00107         //G4cout << extraProcess[i]->GetProcessName() << G4endl;
00108         G4VProcess* p = extraProcess[i]; 
00109         extraProcess[i] = 0;
00110         delete p;
00111       }
00112     }
00113   }
00114   //G4cout << "G4HadronicProcessStore::Clean() done" << G4endl; 
00115   n_extra = 0;
00116   n_proc = 0;
00117 }
00118 
00119 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00120 
00121 G4HadronicProcessStore::G4HadronicProcessStore()
00122 {
00123   n_proc = 0;
00124   n_part = 0;
00125   n_model= 0;
00126   n_extra= 0;
00127   currentProcess  = 0;
00128   currentParticle = 0;
00129   verbose = 1;
00130   buildTableStart = true;
00131   theEPTestMessenger = new G4HadronicEPTestMessenger(this);
00132 }
00133 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00134 
00135 G4double G4HadronicProcessStore::GetCrossSectionPerAtom(
00136                                  const G4ParticleDefinition* part,
00137                                  G4double energy,
00138                                  const G4VProcess* proc,
00139                                  const G4Element*  element)
00140 {
00141   G4double cross = 0.;    
00142   G4int subType = proc->GetProcessSubType();      
00143   if (subType == fHadronElastic)   
00144     cross = GetElasticCrossSectionPerAtom(part,energy,element);
00145   else if (subType == fHadronInelastic)   
00146     cross = GetInelasticCrossSectionPerAtom(part,energy,element);
00147   else if (subType == fCapture)   
00148     cross = GetCaptureCrossSectionPerAtom(part,energy,element);      
00149   else if (subType == fFission)   
00150     cross = GetFissionCrossSectionPerAtom(part,energy,element); 
00151   else if (subType == fChargeExchange)   
00152     cross = GetChargeExchangeCrossSectionPerAtom(part,energy,element);
00153   return cross;
00154 }          
00155 
00156 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00157 
00158 G4double G4HadronicProcessStore::GetCrossSectionPerVolume(
00159                                  const G4ParticleDefinition* part,
00160                                  G4double energy,
00161                                  const G4VProcess* proc,
00162                                  const G4Material* material)
00163 {
00164   G4double cross = 0.;    
00165   G4int subType = proc->GetProcessSubType();      
00166   if (subType == fHadronElastic)   
00167     cross = GetElasticCrossSectionPerVolume(part,energy,material);
00168   else if (subType == fHadronInelastic)   
00169     cross = GetInelasticCrossSectionPerVolume(part,energy,material);
00170   else if (subType == fCapture)   
00171     cross = GetCaptureCrossSectionPerVolume(part,energy,material);      
00172   else if (subType == fFission)   
00173     cross = GetFissionCrossSectionPerVolume(part,energy,material); 
00174   else if (subType == fChargeExchange)   
00175     cross = GetChargeExchangeCrossSectionPerVolume(part,energy,material);
00176   return cross;
00177 }          
00178 
00179 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00180 
00181 G4double G4HadronicProcessStore::GetElasticCrossSectionPerVolume(
00182     const G4ParticleDefinition *aParticle,
00183     G4double kineticEnergy,
00184     const G4Material *material)
00185 {
00186   G4double cross = 0.0;
00187   const G4ElementVector* theElementVector = material->GetElementVector();
00188   const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
00189   size_t nelm = material->GetNumberOfElements();
00190   for (size_t i=0; i<nelm; ++i) {
00191     const G4Element* elm = (*theElementVector)[i];
00192     cross += theAtomNumDensityVector[i]*
00193       GetElasticCrossSectionPerAtom(aParticle,kineticEnergy,elm,material);
00194   }
00195   return cross;
00196 }
00197 
00198 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00199 
00200 G4double G4HadronicProcessStore::GetElasticCrossSectionPerAtom(
00201     const G4ParticleDefinition *aParticle,
00202     G4double kineticEnergy,
00203     const G4Element *anElement, const G4Material* mat)
00204 {
00205   G4HadronicProcess* hp = FindProcess(aParticle, fHadronElastic);
00206   G4double cross = 0.0;
00207   localDP.SetKineticEnergy(kineticEnergy);
00208   if(hp) {
00209     cross = hp->GetElementCrossSection(&localDP,anElement,mat);
00210   }
00211   return cross;
00212 }
00213 
00214 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00215 
00216 G4double G4HadronicProcessStore::GetElasticCrossSectionPerIsotope(
00217     const G4ParticleDefinition*,
00218     G4double,
00219     G4int, G4int)
00220 {
00221   return 0.0;
00222 }
00223 
00224 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00225 
00226 G4double G4HadronicProcessStore::GetInelasticCrossSectionPerVolume(
00227     const G4ParticleDefinition *aParticle,
00228     G4double kineticEnergy,
00229     const G4Material *material)
00230 {
00231   G4double cross = 0.0;
00232   const G4ElementVector* theElementVector = material->GetElementVector();
00233   const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
00234   size_t nelm = material->GetNumberOfElements();
00235   for (size_t i=0; i<nelm; ++i) {
00236     const G4Element* elm = (*theElementVector)[i];
00237     cross += theAtomNumDensityVector[i]*
00238       GetInelasticCrossSectionPerAtom(aParticle,kineticEnergy,elm,material);
00239   }
00240   return cross;
00241 }
00242 
00243 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00244 
00245 G4double G4HadronicProcessStore::GetInelasticCrossSectionPerAtom(
00246     const G4ParticleDefinition *aParticle,
00247     G4double kineticEnergy,
00248     const G4Element *anElement, const G4Material* mat)
00249 {
00250   G4HadronicProcess* hp = FindProcess(aParticle, fHadronInelastic);
00251   localDP.SetKineticEnergy(kineticEnergy);
00252   G4double cross = 0.0;
00253   if(hp) { 
00254     cross = hp->GetElementCrossSection(&localDP,anElement,mat);
00255   }
00256   return cross;
00257 }
00258 
00259 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00260 
00261 G4double G4HadronicProcessStore::GetInelasticCrossSectionPerIsotope(
00262     const G4ParticleDefinition *,
00263     G4double,
00264     G4int, G4int)
00265 {
00266   return 0.0;
00267 }
00268 
00269 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00270 
00271 G4double G4HadronicProcessStore::GetCaptureCrossSectionPerVolume(
00272     const G4ParticleDefinition *aParticle,
00273     G4double kineticEnergy,
00274     const G4Material *material)
00275 {
00276   G4double cross = 0.0;
00277   const G4ElementVector* theElementVector = material->GetElementVector();
00278   const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
00279   size_t nelm = material->GetNumberOfElements();
00280   for (size_t i=0; i<nelm; ++i) {
00281     const G4Element* elm = (*theElementVector)[i];
00282     cross += theAtomNumDensityVector[i]*
00283       GetCaptureCrossSectionPerAtom(aParticle,kineticEnergy,elm,material);
00284   }
00285   return cross;
00286 }
00287 
00288 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00289 
00290 G4double G4HadronicProcessStore::GetCaptureCrossSectionPerAtom(
00291     const G4ParticleDefinition *aParticle,
00292     G4double kineticEnergy,
00293     const G4Element *anElement, const G4Material* mat)
00294 {
00295   G4HadronicProcess* hp = FindProcess(aParticle, fCapture);
00296   localDP.SetKineticEnergy(kineticEnergy);
00297   G4double cross = 0.0;
00298   if(hp) {
00299     cross = hp->GetElementCrossSection(&localDP,anElement,mat);
00300   }
00301   return cross;
00302 }
00303 
00304 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00305 
00306 G4double G4HadronicProcessStore::GetCaptureCrossSectionPerIsotope(
00307     const G4ParticleDefinition *,
00308     G4double,
00309     G4int, G4int)
00310 {
00311   return 0.0;
00312 }
00313 
00314 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00315 
00316 G4double G4HadronicProcessStore::GetFissionCrossSectionPerVolume(
00317     const G4ParticleDefinition *aParticle,
00318     G4double kineticEnergy,
00319     const G4Material *material)
00320 {
00321   G4double cross = 0.0;
00322   const G4ElementVector* theElementVector = material->GetElementVector();
00323   const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
00324   size_t nelm = material->GetNumberOfElements();
00325   for (size_t i=0; i<nelm; i++) {
00326     const G4Element* elm = (*theElementVector)[i];
00327     cross += theAtomNumDensityVector[i]*
00328       GetFissionCrossSectionPerAtom(aParticle,kineticEnergy,elm,material);
00329   }
00330   return cross;
00331 }
00332 
00333 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00334 
00335 G4double G4HadronicProcessStore::GetFissionCrossSectionPerAtom(
00336     const G4ParticleDefinition *aParticle,
00337     G4double kineticEnergy,
00338     const G4Element *anElement, const G4Material* mat)
00339 {
00340   G4HadronicProcess* hp = FindProcess(aParticle, fFission);
00341   localDP.SetKineticEnergy(kineticEnergy);
00342   G4double cross = 0.0;
00343   if(hp) {
00344     cross = hp->GetElementCrossSection(&localDP,anElement,mat);
00345   }
00346   return cross;
00347 }
00348 
00349 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00350 
00351 G4double G4HadronicProcessStore::GetFissionCrossSectionPerIsotope(
00352     const G4ParticleDefinition *,
00353     G4double,
00354     G4int, G4int)
00355 {
00356   return 0.0;
00357 }
00358 
00359 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00360 
00361 G4double G4HadronicProcessStore::GetChargeExchangeCrossSectionPerVolume(
00362     const G4ParticleDefinition *aParticle,
00363     G4double kineticEnergy,
00364     const G4Material *material)
00365 {
00366   G4double cross = 0.0;
00367   const G4ElementVector* theElementVector = material->GetElementVector();
00368   const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
00369   size_t nelm = material->GetNumberOfElements();
00370   for (size_t i=0; i<nelm; ++i) {
00371     const G4Element* elm = (*theElementVector)[i];
00372     cross += theAtomNumDensityVector[i]*
00373       GetChargeExchangeCrossSectionPerAtom(aParticle,kineticEnergy,elm,material);
00374   }
00375   return cross;
00376 }
00377 
00378 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00379 
00380 G4double G4HadronicProcessStore::GetChargeExchangeCrossSectionPerAtom(
00381     const G4ParticleDefinition *aParticle,
00382     G4double kineticEnergy,
00383     const G4Element *anElement, const G4Material* mat)
00384 {
00385   G4HadronicProcess* hp = FindProcess(aParticle, fChargeExchange);
00386   localDP.SetKineticEnergy(kineticEnergy);
00387   G4double cross = 0.0;
00388   if(hp) {
00389     cross = hp->GetElementCrossSection(&localDP,anElement,mat);
00390   }
00391   return cross;
00392 }
00393 
00394 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00395 
00396 G4double G4HadronicProcessStore::GetChargeExchangeCrossSectionPerIsotope(
00397     const G4ParticleDefinition *,
00398     G4double,
00399     G4int, G4int)
00400 {
00401   return 0.0;
00402 }
00403 
00404 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00405 
00406 void G4HadronicProcessStore::Register(G4HadronicProcess* proc) 
00407 { 
00408   if(0 < n_proc) {
00409     for(G4int i=0; i<n_proc; ++i) {
00410       if(process[i] == proc) { return; }
00411     }
00412   }
00413   //  G4cout << "G4HadronicProcessStore::Register hadronic " << n_proc
00414   //     << "  " << proc->GetProcessName() << G4endl;
00415   ++n_proc;
00416   process.push_back(proc);
00417 }
00418 
00419 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00420 
00421 void G4HadronicProcessStore::RegisterParticle(G4HadronicProcess* proc, 
00422                                               const G4ParticleDefinition* part) 
00423 { 
00424   G4int i=0;
00425   for(; i<n_proc; ++i) {if(process[i] == proc) break;}
00426   G4int j=0;
00427   for(; j<n_part; ++j) {if(particle[j] == part) break;}
00428 
00429   if(j == n_part) {
00430     ++n_part;
00431     particle.push_back(part);
00432     wasPrinted.push_back(0);
00433   }
00434   
00435   // the pair should be added?
00436   if(i < n_proc) {
00437     std::multimap<PD,HP,std::less<PD> >::iterator it;
00438     for(it=p_map.lower_bound(part); it!=p_map.upper_bound(part); ++it) {
00439       if(it->first == part) {
00440         HP process2 = (it->second);
00441         if(proc == process2) { return; }
00442       }
00443     }
00444   }
00445   
00446   p_map.insert(std::multimap<PD,HP>::value_type(part,proc));
00447 }
00448 
00449 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00450 
00451 void G4HadronicProcessStore::RegisterInteraction(G4HadronicProcess* proc,
00452                                                  G4HadronicInteraction* mod)
00453 {
00454   G4int i=0;
00455   for(; i<n_proc; ++i) {if(process[i] == proc) { break; }}
00456   G4int k=0;
00457   for(; k<n_model; ++k) {if(model[k] == mod) { break; }}
00458    
00459   m_map.insert(std::multimap<HP,HI>::value_type(proc,mod));
00460     
00461   if(k == n_model) {
00462     ++n_model;
00463     model.push_back(mod);
00464     modelName.push_back(mod->GetModelName());
00465   }
00466 }
00467 
00468 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00469 
00470 void G4HadronicProcessStore::DeRegister(G4HadronicProcess* proc)
00471 {
00472   if(0 == n_proc) return;
00473   for(G4int i=0; i<n_proc; ++i) {
00474     if(process[i] == proc) {
00475       process[i] = 0;
00476       return;
00477     }
00478   }
00479 } 
00480 
00481 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00482 
00483 void G4HadronicProcessStore::RegisterExtraProcess(G4VProcess* proc)
00484 {
00485   if(0 < n_extra) {
00486     for(G4int i=0; i<n_extra; ++i) {
00487       if(extraProcess[i] == proc) { return; }
00488     }
00489   }
00490   //G4cout << "Extra Process: " << n_extra << "  " <<  proc->GetProcessName() 
00491   //     << "  " << proc << G4endl;
00492     
00493   n_extra++;
00494   extraProcess.push_back(proc);
00495 }
00496 
00497 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00498 
00499 void G4HadronicProcessStore::RegisterParticleForExtraProcess(
00500                              G4VProcess* proc,
00501                              const G4ParticleDefinition* part)
00502 {
00503   G4int i=0;
00504   for(; i<n_extra; ++i) { if(extraProcess[i] == proc) { break; } }
00505   G4int j=0;
00506   for(; j<n_part; ++j) { if(particle[j] == part) { break; } }
00507 
00508   if(j == n_part) {
00509     ++n_part;
00510     particle.push_back(part);
00511     wasPrinted.push_back(0);
00512   }
00513   
00514   // the pair should be added?
00515   if(i < n_extra) {
00516     std::multimap<PD,G4VProcess*,std::less<PD> >::iterator it;
00517     for(it=ep_map.lower_bound(part); it!=ep_map.upper_bound(part); ++it) {
00518       if(it->first == part) {
00519         G4VProcess* process2 = (it->second);
00520         if(proc == process2) { return; }
00521       }
00522     }
00523   }
00524   
00525   ep_map.insert(std::multimap<PD,G4VProcess*>::value_type(part,proc));
00526 } 
00527 
00528 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00529 
00530 void G4HadronicProcessStore::DeRegisterExtraProcess(G4VProcess* proc)
00531 {
00532   //G4cout << "Deregister Extra Process: " << proc << "  "<<proc->GetProcessName()<< G4endl;
00533   if(0 == n_extra) { return; }
00534   for(G4int i=0; i<n_extra; ++i) {
00535     if(extraProcess[i] == proc) {
00536       extraProcess[i] = 0;
00537       //G4cout << "Extra Process: " << i << " is deregisted " << G4endl;
00538       return;
00539     }
00540   }
00541 }
00542 
00543 
00544 void G4HadronicProcessStore::PrintInfo(const G4ParticleDefinition* part) 
00545 {
00546   // Trigger particle/process/model printout only when last particle is 
00547   // registered
00548   if(buildTableStart && part == particle[n_part - 1]) {
00549     buildTableStart = false;
00550     Dump(verbose);
00551     if (getenv("G4PhysListDocDir") ) DumpHtml();
00552   }
00553 }
00554 
00555 
00556 void G4HadronicProcessStore::DumpHtml()
00557 {
00558   // Automatic generation of html documentation page for physics lists
00559   // List processes, models and cross sections for the most important
00560   // particles in descending order of importance
00561 
00562   char* dirName = getenv("G4PhysListDocDir");
00563   char* physListName = getenv("G4PhysListName");
00564   if (dirName && physListName) {
00565 
00566     // Open output file with path name
00567     G4String pathName = G4String(dirName) + "/" + G4String(physListName) + ".html";
00568     std::ofstream outFile;
00569     outFile.open(pathName);
00570 
00571     // Write physics list summary file
00572     outFile << "<html>\n";
00573     outFile << "<head>\n";
00574     outFile << "<title>Physics List Summary</title>\n";
00575     outFile << "</head>\n";
00576     outFile << "<body>\n";
00577     outFile << "<h2> Summary of Hadronic Processes, Models and Cross Sections for Physics List "
00578             << G4String(physListName) << "</h2>\n";
00579     outFile << "<ul>\n";
00580 
00581     PrintHtml(G4Proton::Proton(), outFile);
00582     PrintHtml(G4Neutron::Neutron(), outFile);
00583     PrintHtml(G4PionPlus::PionPlus(), outFile); 
00584     PrintHtml(G4PionMinus::PionMinus(), outFile);
00585     PrintHtml(G4Gamma::Gamma(), outFile);
00586     PrintHtml(G4Electron::Electron(), outFile);
00587 //    PrintHtml(G4MuonMinus::MuonMinus(), outFile);
00588     PrintHtml(G4Positron::Positron(), outFile);
00589     PrintHtml(G4KaonPlus::KaonPlus(), outFile);
00590     PrintHtml(G4KaonMinus::KaonMinus(), outFile);
00591     PrintHtml(G4Lambda::Lambda(), outFile);
00592     PrintHtml(G4Alpha::Alpha(), outFile);
00593 
00594     outFile << "</ul>\n";
00595     outFile << "</body>\n";
00596     outFile << "</html>\n";
00597     outFile.close();
00598   }
00599 }
00600 
00601 
00602 void G4HadronicProcessStore::PrintHtml(const G4ParticleDefinition* theParticle,
00603                                        std::ofstream& outFile)
00604 {
00605   // Automatic generation of html documentation page for physics lists
00606   // List processes for the most important particles in descending order
00607   // of importance
00608  
00609   outFile << "<br> <li><h2><font color=\" ff0000 \">" 
00610           << theParticle->GetParticleName() << "</font></h2></li>\n";
00611 
00612   typedef std::multimap<PD,HP,std::less<PD> > PDHPmap;
00613   typedef std::multimap<HP,HI,std::less<HP> > HPHImap;
00614 
00615   std::pair<PDHPmap::iterator, PDHPmap::iterator> itpart =
00616                         p_map.equal_range(theParticle);
00617 
00618   // Loop over processes assigned to particle
00619 
00620   G4HadronicProcess* theProcess;
00621   for (PDHPmap::iterator it = itpart.first; it != itpart.second; ++it) {
00622     theProcess = (*it).second;
00623     outFile << "<br> &nbsp;&nbsp; <b><font color=\" 0000ff \">process : <a href=\"" 
00624             << theProcess->GetProcessName() << ".html\"> "
00625             << theProcess->GetProcessName() << "</a></font></b>\n";
00626     outFile << "<ul>\n";
00627     outFile << "  <li><b><font color=\" 00AA00 \">models : </font></b>\n";
00628 
00629     // Loop over models assigned to process
00630     std::pair<HPHImap::iterator, HPHImap::iterator> itmod =
00631                         m_map.equal_range(theProcess);
00632 
00633     outFile << "    <ul>\n"; 
00634     for (HPHImap::iterator jt = itmod.first; jt != itmod.second; ++jt) {
00635       outFile << "    <li><b><a href=\"" << (*jt).second->GetModelName() << ".html\"> "
00636               << (*jt).second->GetModelName() << "</a>" 
00637               << " from " << (*jt).second->GetMinEnergy()/GeV
00638               << " GeV to " << (*jt).second->GetMaxEnergy()/GeV
00639               << " GeV </b></li>\n";
00640 
00641       // Print ModelDescription, ignore that we overwrite files n-times.
00642       PrintModelHtml((*jt).second);
00643 
00644     }
00645     outFile << "    </ul>\n";
00646     outFile << "  </li>\n";
00647 
00648     // List cross sections assigned to process
00649     outFile << "  <li><b><font color=\" 00AA00 \">cross sections : </font></b>\n";
00650     outFile << "    <ul>\n";
00651     theProcess->GetCrossSectionDataStore()->DumpHtml(*theParticle, outFile);
00652     //        << " \n";
00653     outFile << "    </ul>\n";
00654 
00655     outFile << "  </li>\n";
00656     outFile << "</ul>\n";
00657   }
00658 }
00659 
00660 void G4HadronicProcessStore::PrintModelHtml(const G4HadronicInteraction * mod) const
00661 {
00662         G4String dirName(getenv("G4PhysListDocDir"));
00663         G4String pathName = dirName + "/" + mod->GetModelName() + ".html";
00664         std::ofstream outModel;
00665         outModel.open(pathName);
00666         outModel << "<html>\n";
00667         outModel << "<head>\n";
00668         outModel << "<title>Description of " << mod->GetModelName() << "</title>\n";
00669         outModel << "</head>\n";
00670         outModel << "<body>\n";
00671 
00672         mod->ModelDescription(outModel);
00673 
00674         outModel << "</body>\n";
00675         outModel << "</html>\n";
00676 
00677 }
00678 void G4HadronicProcessStore::Dump(G4int level)
00679 {
00680   if(level > 0) {
00681     G4cout << "=============================================================="
00682            << "=============================="
00683            << G4endl;
00684       G4cout << "             HADRONIC PROCESSES SUMMARY (verbose level " << level
00685              << ")" << G4endl;
00686   }
00687   for(G4int i=0; i<n_part; ++i) {
00688     PD part = particle[i];
00689     G4String pname = part->GetParticleName();
00690     G4bool yes = false;
00691     if(level >= 2) yes = true;
00692 
00693     else if(level == 1 && (pname == "proton" || 
00694                            pname == "neutron" ||
00695                            pname == "pi+" ||
00696                            pname == "pi-" ||
00697                            pname == "gamma" ||
00698                            pname == "e+" ||
00699                            pname == "e-" ||
00700                            pname == "mu+" ||
00701                            pname == "mu-" ||
00702                            pname == "kaon+" ||
00703                            pname == "kaon-" ||
00704                            pname == "lambda" ||
00705                            pname == "GenericIon" ||
00706                            pname == "anti_neutron" ||
00707                            pname == "anti_proton")) yes = true;
00708     if(yes) {
00709       // main processes
00710       std::multimap<PD,HP,std::less<PD> >::iterator it;
00711 
00712       for(it=p_map.lower_bound(part); it!=p_map.upper_bound(part); ++it) {
00713         if(it->first == part) {
00714           HP proc = (it->second);
00715           G4int j=0;
00716           for(; j<n_proc; ++j) {
00717             if(process[j] == proc) {
00718               Print(j, i);
00719             }
00720           }
00721         }
00722       }
00723       // extra processes
00724       std::multimap<PD,G4VProcess*,std::less<PD> >::iterator itp;
00725       for(itp=ep_map.lower_bound(part); itp!=ep_map.upper_bound(part); ++itp) {
00726         if(itp->first == part) {
00727           G4VProcess* proc = (itp->second);
00728           if(wasPrinted[i] == 0) {
00729             wasPrinted[i] = 1;
00730             G4cout<<G4endl;
00731             G4cout << "                     Hadronic Processes for <" 
00732                    <<part->GetParticleName() << ">" << G4endl; 
00733           }
00734           G4cout << "          " << proc->GetProcessName() << G4endl;
00735         }
00736       }
00737     }
00738   }
00739   if(level > 0) {
00740     G4cout << "=============================================================="
00741            << "=============================="
00742            << G4endl;
00743   }
00744 }
00745 
00746 
00747 void G4HadronicProcessStore::Print(G4int idxProc, G4int idxPart)
00748 {
00749   G4HadronicProcess* proc = process[idxProc];
00750   const G4ParticleDefinition* part = particle[idxPart];
00751   if(wasPrinted[idxPart] == 0) {
00752     wasPrinted[idxPart] = 1;
00753     G4cout<<G4endl;
00754     G4cout << "                                  Hadronic Processes for <" 
00755            << part->GetParticleName() << ">" << G4endl;
00756     G4cout << "                                  ------------------------"
00757            << "-----------" << G4endl;
00758   }
00759   HI hi = 0;
00760   G4bool first;
00761   std::multimap<HP,HI,std::less<HP> >::iterator ih;
00762   G4cout << std::setw(20) << proc->GetProcessName()  
00763          << "  Models: ";
00764   first = true;
00765   for(ih=m_map.lower_bound(proc); ih!=m_map.upper_bound(proc); ++ih) {
00766     if(ih->first == proc) {
00767       hi = ih->second;
00768       G4int i=0;
00769       for(; i<n_model; ++i) {
00770         if(model[i] == hi) { break; }
00771       }
00772       if(!first) G4cout << "                              ";
00773       first = false;
00774       G4cout << std::setw(28) << modelName[i] 
00775              << ": Emin(GeV)= "  
00776              << std::setw(4) << hi->GetMinEnergy()/GeV
00777              << "  Emax(GeV)= " 
00778              << hi->GetMaxEnergy()/GeV
00779              << G4endl;
00780     }
00781   }
00782 
00783   G4cout << G4endl;
00784   G4cout << std::setw(20) << proc->GetProcessName()
00785          << "  Crs sctns: ";
00786   G4CrossSectionDataStore* csds = proc->GetCrossSectionDataStore();
00787   csds->DumpPhysicsTable(*part);
00788 
00789 }
00790 
00791 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00792 
00793 void G4HadronicProcessStore::SetVerbose(G4int val)
00794 {
00795   verbose = val;
00796   G4int i;
00797   for(i=0; i<n_proc; ++i) {
00798     if(process[i]) { process[i]->SetVerboseLevel(val); }
00799   }
00800   for(i=0; i<n_model; ++i) {
00801     if(model[i]) { model[i]->SetVerboseLevel(val); }
00802   }
00803 }
00804 
00805 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00806 
00807 G4int G4HadronicProcessStore::GetVerbose()
00808 {
00809   return verbose;
00810 }
00811 
00812 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
00813 
00814 G4HadronicProcess* G4HadronicProcessStore::FindProcess(
00815    const G4ParticleDefinition* part, G4HadronicProcessType subType)
00816 {
00817   bool isNew = false;
00818   G4HadronicProcess* hp = 0;
00819 
00820   if(part != currentParticle) {
00821     isNew = true;
00822     currentParticle = part;
00823     localDP.SetDefinition(part);
00824   } else if(!currentProcess) {
00825     isNew = true;
00826   } else if(subType == currentProcess->GetProcessSubType()) {
00827     hp = currentProcess;
00828   } else {
00829     isNew = true;
00830   }
00831 
00832   if(isNew) {
00833     std::multimap<PD,HP,std::less<PD> >::iterator it;
00834     for(it=p_map.lower_bound(part); it!=p_map.upper_bound(part); ++it) {
00835       if(it->first == part && subType == (it->second)->GetProcessSubType()) {
00836         hp = it->second;
00837         break;
00838       }
00839     }  
00840     currentProcess = hp;
00841   }
00842 
00843   return hp;
00844 }
00845 
00846 void G4HadronicProcessStore::SetEpReportLevel(G4int level)
00847 {
00848   G4cout << " Setting energy/momentum report level to " << level 
00849          << " for " << process.size() << " hadronic processes " << G4endl;
00850   for (G4int i = 0; i < G4int(process.size()); ++i) {
00851     process[i]->SetEpReportLevel(level);
00852   }
00853 }
00854 
00855 void G4HadronicProcessStore::SetProcessAbsLevel(G4double abslevel)
00856 {
00857   G4cout << " Setting absolute energy/momentum test level to " << abslevel << G4endl;
00858   G4double rellevel = 0.0;
00859   G4HadronicProcess* theProcess = 0;
00860   for (G4int i = 0; i < G4int(process.size()); ++i) {
00861     theProcess = process[i];
00862     rellevel = theProcess->GetEnergyMomentumCheckLevels().first;
00863     theProcess->SetEnergyMomentumCheckLevels(rellevel, abslevel);
00864   }
00865 }
00866 
00867 void G4HadronicProcessStore::SetProcessRelLevel(G4double rellevel)
00868 {
00869   G4cout << " Setting relative energy/momentum test level to " << rellevel << G4endl;
00870   G4double abslevel = 0.0;
00871   G4HadronicProcess* theProcess = 0;
00872   for (G4int i = 0; i < G4int(process.size()); ++i) {
00873     theProcess = process[i];
00874     abslevel = theProcess->GetEnergyMomentumCheckLevels().second;
00875     theProcess->SetEnergyMomentumCheckLevels(rellevel, abslevel);
00876   }
00877 }

Generated on Mon May 27 17:48:26 2013 for Geant4 by  doxygen 1.4.7