Geant4  10.02.p03
G4PhysListFactory Class Reference

#include <G4PhysListFactory.hh>

Collaboration diagram for G4PhysListFactory:

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)
 

Private Attributes

G4String defName
 
std::vector< G4Stringlistnames_hadr
 
std::vector< G4Stringlistnames_em
 
size_t nlists_hadr
 
size_t nlists_em
 
G4int verbose
 

Detailed Description

Definition at line 44 of file G4PhysListFactory.hh.

Constructor & Destructor Documentation

◆ G4PhysListFactory()

G4PhysListFactory::G4PhysListFactory ( )

Definition at line 73 of file G4PhysListFactory.cc.

74  : defName("FTFP_BERT"),verbose(1)
75 {
76  nlists_hadr = 22;
77  G4String ss[22] = {
78  "FTFP_BERT","FTFP_BERT_TRV","FTFP_BERT_ATL","FTFP_BERT_HP","FTFP_INCLXX",
79  "FTFP_INCLXX_HP","FTF_BIC", "LBE","QBBC",
80  "QGSP_BERT","QGSP_BERT_HP","QGSP_BIC","QGSP_BIC_HP","QGSP_BIC_AllHP",
81  "QGSP_FTFP_BERT","QGSP_INCLXX","QGSP_INCLXX_HP","QGS_BIC",
82  "Shielding","ShieldingLEND","ShieldingM","NuBeam"};
83  for(size_t i=0; i<nlists_hadr; ++i) {
84  listnames_hadr.push_back(ss[i]);
85  }
86 
87  nlists_em = 8;
88  G4String s1[8] = {"","_EMV","_EMX","_EMY","_EMZ","_LIV","_PEN","__GS"};
89  for(size_t i=0; i<nlists_em; ++i) {
90  listnames_em.push_back(s1[i]);
91  }
92 }
std::vector< G4String > listnames_em
std::vector< G4String > listnames_hadr

◆ ~G4PhysListFactory()

G4PhysListFactory::~G4PhysListFactory ( )

Definition at line 94 of file G4PhysListFactory.cc.

95 {}

Member Function Documentation

◆ AvailablePhysLists()

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

Definition at line 227 of file G4PhysListFactory.cc.

228 {
229  return listnames_hadr;
230 }
std::vector< G4String > listnames_hadr

◆ AvailablePhysListsEM()

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

Definition at line 233 of file G4PhysListFactory.cc.

234 {
235  return listnames_em;
236 }
std::vector< G4String > listnames_em

◆ GetReferencePhysList()

G4VModularPhysicsList * G4PhysListFactory::GetReferencePhysList ( const G4String name)

Definition at line 118 of file G4PhysListFactory.cc.

119 {
120  // analysis on the string
121  size_t n = name.size();
122 
123  // last characters in the string
124  size_t em_opt = 0;
125  G4String em_name = "";
126 
127  // check EM options
128  if(n > 4) {
129  em_name = name.substr(n - 4, 4);
130  for(size_t i=1; i<nlists_em; ++i) {
131  if(listnames_em[i] == em_name) {
132  em_opt = i;
133  n -= 4;
134  break;
135  }
136  }
137  if(0 == em_opt) { em_name = ""; }
138  }
139 
140  // hadronic pHysics List
141  G4String had_name = name.substr(0, n);
142 
143  if(0 < verbose) {
144  G4cout << "G4PhysListFactory::GetReferencePhysList <" << had_name
145  << em_name << "> EMoption= " << em_opt << G4endl;
146  }
147  G4VModularPhysicsList* p = 0;
148  if(had_name == "FTFP_BERT") {p = new FTFP_BERT(verbose);}
149  else if(had_name == "FTFP_BERT_HP") {p = new FTFP_BERT_HP(verbose);}
150  else if(had_name == "FTFP_BERT_TRV") {p = new FTFP_BERT_TRV(verbose);}
151  else if(had_name == "FTFP_BERT_ATL") {p = new FTFP_BERT_ATL(verbose);}
152  else if(had_name == "FTFP_INCLXX") {p = new FTFP_INCLXX(verbose);}
153  else if(had_name == "FTFP_INCLXX_HP") {p = new FTFP_INCLXX_HP(verbose);}
154  else if(had_name == "FTF_BIC") {p = new FTF_BIC(verbose);}
155  else if(had_name == "LBE") {p = new LBE();}
156  else if(had_name == "QBBC") {p = new QBBC(verbose);}
157  else if(had_name == "QGSP_BERT") {p = new QGSP_BERT(verbose);}
158  else if(had_name == "QGSP_BERT_HP") {p = new QGSP_BERT_HP(verbose);}
159  else if(had_name == "QGSP_BIC") {p = new QGSP_BIC(verbose);}
160  else if(had_name == "QGSP_BIC_HP") {p = new QGSP_BIC_HP(verbose);}
161  else if(had_name == "QGSP_BIC_AllHP") {p = new QGSP_BIC_AllHP(verbose);}
162  else if(had_name == "QGSP_FTFP_BERT") {p = new QGSP_FTFP_BERT(verbose);}
163  else if(had_name == "QGSP_INCLXX") {p = new QGSP_INCLXX(verbose);}
164  else if(had_name == "QGSP_INCLXX_HP") {p = new QGSP_INCLXX_HP(verbose);}
165  else if(had_name == "QGS_BIC") {p = new QGS_BIC(verbose);}
166  else if(had_name == "Shielding") {p = new Shielding(verbose);}
167  else if(had_name == "ShieldingLEND") {p = new Shielding(verbose,"LEND");}
168  else if(had_name == "ShieldingM") {p = new Shielding(verbose,"HP","M");}
169  else if(had_name == "NuBeam") {p = new NuBeam(verbose);}
170  else {
171  G4cout << "### G4PhysListFactory WARNING: "
172  << "PhysicsList " << had_name << " is not known"
173  << G4endl;
174  }
175  if(p) {
176  G4cout << "<<< Reference Physics List " << had_name
177  << em_name << " is built" << G4endl;
178  G4int ver = p->GetVerboseLevel();
179  p->SetVerboseLevel(0);
180  if(0 < em_opt) {
181  if(1 == em_opt) {
183  } else if(2 == em_opt) {
185  } else if(3 == em_opt) {
187  } else if(4 == em_opt) {
189  } else if(5 == em_opt) {
191  } else if(6 == em_opt) {
193  } else if(7 == em_opt) {
195  }
196  }
197  p->SetVerboseLevel(ver);
198  }
199  G4cout << G4endl;
200  return p;
201 }
TQGSP_FTFP_BERT< G4VModularPhysicsList > QGSP_FTFP_BERT
TQGSP_BERT< G4VModularPhysicsList > QGSP_BERT
Definition: QGSP_BERT.hh:63
TQGSP_BIC_HP< G4VModularPhysicsList > QGSP_BIC_HP
Definition: QGSP_BIC_HP.hh:63
TFTFP_BERT_HP< G4VModularPhysicsList > FTFP_BERT_HP
Definition: FTFP_BERT_HP.hh:63
TShielding< G4VModularPhysicsList > Shielding
Definition: Shielding.hh:68
TLBE< G4VModularPhysicsList > LBE
Definition: LBE.hh:123
int G4int
Definition: G4Types.hh:78
Char_t n[5]
G4GLOB_DLL std::ostream G4cout
TQGS_BIC< G4VModularPhysicsList > QGS_BIC
Definition: QGS_BIC.hh:63
TFTFP_BERT_TRV< G4VModularPhysicsList > FTFP_BERT_TRV
void SetVerboseLevel(G4int value)
TNuBeam< G4VModularPhysicsList > NuBeam
Definition: NuBeam.hh:64
TQGSP_BIC_AllHP< G4VModularPhysicsList > QGSP_BIC_AllHP
void ReplacePhysics(G4VPhysicsConstructor *)
std::vector< G4String > listnames_em
TFTFP_BERT< G4VModularPhysicsList > FTFP_BERT
Definition: FTFP_BERT.hh:63
TINCLXXPhysicsListHelper< G4VModularPhysicsList, false, true > FTFP_INCLXX
Definition: FTFP_INCLXX.hh:41
TINCLXXPhysicsListHelper< G4VModularPhysicsList, false, false > QGSP_INCLXX
TINCLXXPhysicsListHelper< G4VModularPhysicsList, true, false > QGSP_INCLXX_HP
#define G4endl
Definition: G4ios.hh:61
TQGSP_BIC< G4VModularPhysicsList > QGSP_BIC
Definition: QGSP_BIC.hh:62
TFTF_BIC< G4VModularPhysicsList > FTF_BIC
Definition: FTF_BIC.hh:62
TFTFP_BERT_ATL< G4VModularPhysicsList > FTFP_BERT_ATL
TQGSP_BERT_HP< G4VModularPhysicsList > QGSP_BERT_HP
Definition: QGSP_BERT_HP.hh:62
Definition: QBBC.hh:44
TINCLXXPhysicsListHelper< G4VModularPhysicsList, true, true > FTFP_INCLXX_HP
Here is the call graph for this function:
Here is the caller graph for this function:

◆ IsReferencePhysList()

G4bool G4PhysListFactory::IsReferencePhysList ( const G4String name)

Definition at line 203 of file G4PhysListFactory.cc.

204 {
205  G4bool res = false;
206  size_t n = name.size();
207  if(n > 4) {
208  G4String em_name = name.substr(n - 4, 4);
209  for(size_t i=1; i<nlists_em; ++i) {
210  if(listnames_em[i] == em_name) {
211  n -= 4;
212  break;
213  }
214  }
215  }
216  G4String had_name = name.substr(0, n);
217  for(size_t i=0; i<nlists_hadr; ++i) {
218  if(had_name == listnames_hadr[i]) {
219  res = true;
220  break;
221  }
222  }
223  return res;
224 }
Char_t n[5]
bool G4bool
Definition: G4Types.hh:79
std::vector< G4String > listnames_em
std::vector< G4String > listnames_hadr
Here is the caller graph for this function:

◆ ReferencePhysList()

G4VModularPhysicsList * G4PhysListFactory::ReferencePhysList ( )

Definition at line 98 of file G4PhysListFactory.cc.

99 {
100  // instantiate PhysList by environment variable "PHYSLIST"
101  G4String name = "";
102  char* path = getenv("PHYSLIST");
103  if (path) {
104  name = G4String(path);
105  } else {
106  name = defName;
107  G4cout << "### G4PhysListFactory WARNING: "
108  << " environment variable PHYSLIST is not defined"
109  << G4endl
110  << " Default Physics Lists " << name
111  << " is instantiated"
112  << G4endl;
113  }
114  return GetReferencePhysList(name);
115 }
G4String name
Definition: TRTMaterials.hh:40
G4GLOB_DLL std::ostream G4cout
G4VModularPhysicsList * GetReferencePhysList(const G4String &)
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetVerbose()

void G4PhysListFactory::SetVerbose ( G4int  val)
inline

Definition at line 67 of file G4PhysListFactory.hh.

67 { verbose = val; }

Member Data Documentation

◆ defName

G4String G4PhysListFactory::defName
private

Definition at line 71 of file G4PhysListFactory.hh.

◆ listnames_em

std::vector<G4String> G4PhysListFactory::listnames_em
private

Definition at line 73 of file G4PhysListFactory.hh.

◆ listnames_hadr

std::vector<G4String> G4PhysListFactory::listnames_hadr
private

Definition at line 72 of file G4PhysListFactory.hh.

◆ nlists_em

size_t G4PhysListFactory::nlists_em
private

Definition at line 75 of file G4PhysListFactory.hh.

◆ nlists_hadr

size_t G4PhysListFactory::nlists_hadr
private

Definition at line 74 of file G4PhysListFactory.hh.

◆ verbose

G4int G4PhysListFactory::verbose
private

Definition at line 76 of file G4PhysListFactory.hh.


The documentation for this class was generated from the following files: