G4PhysListFactory Class Reference

#include <G4PhysListFactory.hh>


Public Member Functions

 G4PhysListFactory ()
 ~G4PhysListFactory ()
G4VModularPhysicsListGetReferencePhysList (const G4String &)
G4VModularPhysicsListReferencePhysList ()
G4bool IsReferencePhysList (const G4String &)
const std::vector< G4String > & AvailablePhysLists () const
const std::vector< G4String > & AvailablePhysListsEM () const
void SetVerbose (G4int val)


Detailed Description

Definition at line 44 of file G4PhysListFactory.hh.


Constructor & Destructor Documentation

G4PhysListFactory::G4PhysListFactory (  ) 

Definition at line 69 of file G4PhysListFactory.cc.

00070   : defName("FTFP_BERT"),verbose(1)
00071 {
00072   nlists_hadr = 19;
00073   G4String ss[19] = {
00074     "CHIPS",
00075     "FTFP_BERT","FTFP_BERT_TRV","FTFP_BERT_HP","FTF_BIC", 
00076     "LBE","LHEP","QBBC",
00077     "QGSC_BERT","QGSP","QGSP_BERT","QGSP_BERT_CHIPS","QGSP_BERT_HP",
00078     "QGSP_BIC","QGSP_BIC_HP",
00079     "QGSP_FTFP_BERT","QGS_BIC","QGSP_INCLXX",
00080     "Shielding"};
00081   for(size_t i=0; i<nlists_hadr; ++i) {
00082     listnames_hadr.push_back(ss[i]);
00083   }
00084 
00085   nlists_em = 7;
00086   G4String s1[7] = {"","_EMV","_EMX","_EMY","_EMZ","_LIV","_PEN"};
00087   for(size_t i=0; i<nlists_em; ++i) {
00088     listnames_em.push_back(s1[i]);
00089   }
00090 }

G4PhysListFactory::~G4PhysListFactory (  ) 

Definition at line 92 of file G4PhysListFactory.cc.

00093 {}


Member Function Documentation

const std::vector< G4String > & G4PhysListFactory::AvailablePhysLists (  )  const

Definition at line 221 of file G4PhysListFactory.cc.

00222 {
00223   return listnames_hadr;
00224 }

const std::vector< G4String > & G4PhysListFactory::AvailablePhysListsEM (  )  const

Definition at line 227 of file G4PhysListFactory.cc.

00228 {
00229   return listnames_em;
00230 }

G4VModularPhysicsList * G4PhysListFactory::GetReferencePhysList ( const G4String  ) 

Definition at line 116 of file G4PhysListFactory.cc.

References G4cout, G4endl, G4VModularPhysicsList::GetVerboseLevel(), CLHEP::detail::n, G4VModularPhysicsList::ReplacePhysics(), and G4VModularPhysicsList::SetVerboseLevel().

Referenced by ReferencePhysList().

00117 {
00118   // analysis on the string 
00119   size_t n = name.size();
00120 
00121   // last characters in the string
00122   size_t em_opt = 0;
00123   G4String em_name = "";
00124 
00125   // check EM options
00126   if(n > 4) {
00127     em_name = name.substr(n - 4, 4);
00128     for(size_t i=1; i<nlists_em; ++i) { 
00129       if(listnames_em[i] == em_name) { 
00130         em_opt = i;
00131         n -= 4;
00132         break; 
00133       }
00134     }
00135     if(0 == em_opt) { em_name = ""; }
00136   }
00137 
00138   // hadronic pHysics List
00139   G4String had_name = name.substr(0, n);
00140 
00141   if(0 < verbose) {
00142     G4cout << "G4PhysListFactory::GetReferencePhysList <" << had_name
00143            << em_name << ">  EMoption= " << em_opt << G4endl;
00144   }
00145   G4VModularPhysicsList* p = 0;
00146   if(had_name == "CHIPS")               {p = new CHIPS(verbose);}
00147   else if(had_name == "FTFP_BERT")      {p = new FTFP_BERT(verbose);}
00148   else if(had_name == "FTFP_BERT_TRV")  {p = new FTFP_BERT_TRV(verbose);}
00149   else if(had_name == "FTFP_BERT_HP")   {p = new FTFP_BERT_HP(verbose);}
00150   //  else if(had_name == "FTFP_BERT_DE")   {p = new FTFP_BERT_DE(verbose);}
00151   else if(had_name == "FTF_BIC")        {p = new FTF_BIC(verbose);}
00152   else if(had_name == "LBE")            {p = new LBE();}
00153   else if(had_name == "LHEP")           {p = new LHEP(verbose);}
00154   else if(had_name == "QBBC")           {p = new QBBC(verbose);}
00155   else if(had_name == "QGSC_BERT")      {p = new QGSC_BERT(verbose);}
00156   else if(had_name == "QGSP_BERT")      {p = new QGSP_BERT(verbose);}
00157   else if(had_name == "QGSP_BERT_CHIPS"){p = new QGSP_BERT_CHIPS(verbose);}
00158   else if(had_name == "QGSP_BERT_HP")   {p = new QGSP_BERT_HP(verbose);}
00159   else if(had_name == "QGSP_BIC")       {p = new QGSP_BIC(verbose);}
00160   else if(had_name == "QGSP_BIC_HP")    {p = new QGSP_BIC_HP(verbose);}
00161   else if(had_name == "QGSP_FTFP_BERT") {p = new QGSP_FTFP_BERT(verbose);}
00162   else if(had_name == "QGS_BIC")        {p = new QGS_BIC(verbose);}
00163   else if(had_name == "QGSP_INCLXX")    {p = new QGSP_INCLXX(verbose);}
00164   else if(had_name == "Shielding")      {p = new Shielding(verbose);}
00165   else if(had_name == "ShieldingLEND")  {p = new Shielding(verbose,"LEND");}
00166   else {
00167     G4cout << "### G4PhysListFactory WARNING: "
00168            << "PhysicsList " << had_name << " is not known"
00169            << G4endl;
00170   }
00171   if(p) {
00172     G4cout << "<<< Reference Physics List " << had_name
00173            << em_name << " is built" << G4endl;
00174     G4int ver = p->GetVerboseLevel();
00175     p->SetVerboseLevel(0);
00176     if(0 < em_opt) {
00177       if(1 == em_opt) { 
00178         p->ReplacePhysics(new G4EmStandardPhysics_option1(verbose)); 
00179       } else if(2 == em_opt) {
00180         p->ReplacePhysics(new G4EmStandardPhysics_option2(verbose)); 
00181       } else if(3 == em_opt) {
00182         p->ReplacePhysics(new G4EmStandardPhysics_option3(verbose)); 
00183       } else if(4 == em_opt) {
00184         p->ReplacePhysics(new G4EmStandardPhysics_option4(verbose)); 
00185       } else if(5 == em_opt) {
00186         p->ReplacePhysics(new G4EmLivermorePhysics(verbose)); 
00187       } else if(6 == em_opt) {
00188         p->ReplacePhysics(new G4EmPenelopePhysics(verbose)); 
00189       }
00190     }
00191     p->SetVerboseLevel(ver);
00192   }
00193   G4cout << G4endl;
00194   return p;
00195 }

G4bool G4PhysListFactory::IsReferencePhysList ( const G4String  ) 

Definition at line 197 of file G4PhysListFactory.cc.

References CLHEP::detail::n.

00198 {
00199   G4bool res = false;
00200   size_t n = name.size();
00201   if(n > 4) {
00202     G4String em_name = name.substr(n - 4, 4);
00203     for(size_t i=1; i<nlists_em; ++i) { 
00204       if(listnames_em[i] == em_name) { 
00205         n -= 4;
00206         break; 
00207       }
00208     }
00209   }
00210   G4String had_name = name.substr(0, n);
00211   for(size_t i=0; i<nlists_hadr; ++i) {
00212     if(had_name == listnames_hadr[i]) {
00213       res = true;
00214       break;
00215     }
00216   }
00217   return res;
00218 }

G4VModularPhysicsList * G4PhysListFactory::ReferencePhysList (  ) 

Definition at line 96 of file G4PhysListFactory.cc.

References G4cout, G4endl, and GetReferencePhysList().

00097 {
00098   // instantiate PhysList by environment variable "PHYSLIST"
00099   G4String name = "";
00100   char* path = getenv("PHYSLIST");
00101   if (path) {
00102     name = G4String(path);
00103   } else {
00104     name = defName;
00105     G4cout << "### G4PhysListFactory WARNING: "
00106            << " environment variable PHYSLIST is not defined"
00107            << G4endl
00108            << "    Default Physics Lists " << name 
00109            << " is instantiated" 
00110            << G4endl;
00111   }
00112   return GetReferencePhysList(name);
00113 }

void G4PhysListFactory::SetVerbose ( G4int  val  )  [inline]

Definition at line 67 of file G4PhysListFactory.hh.

00067 { verbose = val; }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:52:54 2013 for Geant4 by  doxygen 1.4.7