Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PhysListFactory Class Reference

#include <G4PhysListFactory.hh>

Public Member Functions

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

Detailed Description

Definition at line 45 of file G4PhysListFactory.hh.

Constructor & Destructor Documentation

◆ G4PhysListFactory()

G4PhysListFactory::G4PhysListFactory ( G4int ver = 1)

Definition at line 82 of file G4PhysListFactory.cc.

83 : defName("FTFP_BERT"),verbose(ver),theMessenger(nullptr)
84{
85 nlists_hadr = 36;
86 G4String ss[36] = {
87 "FTFP_BERT","FTFP_BERT_TRV","FTFP_BERT_ATL","FTFP_BERT_HP","FTFQGSP_BERT",
88 "FTFP_INCLXX","FTFP_INCLXX_HP","FTF_BIC","LBE","QBBC",
89 "QGSP_BERT","QGSP_BERT_HP","QGSP_BIC","QGSP_BIC_HP","QGSP_BIC_AllHP",
90 "QGSP_FTFP_BERT","QGSP_INCLXX","QGSP_INCLXX_HP","QGS_BIC",
91 "Shielding","ShieldingLEND","ShieldingLIQMD","ShieldingM","NuBeam",
92 "Shielding_HP","ShieldingLIQMD_HP","ShieldingM_HP",
93 "FTFP_BERT_HPT","FTFP_INCLXX_HPT","QGSP_BERT_HPT","QGSP_BIC_HPT",
94 "QGSP_BIC_AllHPT","QGSP_INCLXX_HPT","Shielding_HPT","ShieldingLIQMD_HPT",
95 "ShieldingM_HPT"};
96 for(size_t i=0; i<nlists_hadr; ++i) {
97 listnames_hadr.push_back(ss[i]);
98 }
99
100 nlists_em = 12;
101 G4String s1[12] = {"","_EMV","_EMX","_EMY","_EMZ","_LIV","_PEN",
102 "__GS","__SS","_EM0","_WVI","__LE"};
103 for(size_t i=0; i<nlists_em; ++i) {
104 listnames_em.push_back(s1[i]);
105 }
106}

◆ ~G4PhysListFactory()

G4PhysListFactory::~G4PhysListFactory ( )

Definition at line 108 of file G4PhysListFactory.cc.

109{
110 delete theMessenger;
111}

Member Function Documentation

◆ AvailablePhysLists()

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

Definition at line 271 of file G4PhysListFactory.cc.

272{
273 return listnames_hadr;
274}

◆ AvailablePhysListsEM()

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

Definition at line 277 of file G4PhysListFactory.cc.

278{
279 return listnames_em;
280}

◆ GetReferencePhysList()

G4VModularPhysicsList * G4PhysListFactory::GetReferencePhysList ( const G4String & name)

Definition at line 134 of file G4PhysListFactory.cc.

135{
136 // analysis on the string
137 size_t n = name.size();
138
139 // last characters in the string
140 size_t em_opt = 0;
141 G4String em_name = "";
142
143 // check EM options
144 if(n > 4) {
145 em_name = name.substr(n - 4, 4);
146 for(size_t i=1; i<nlists_em; ++i) {
147 if(listnames_em[i] == em_name) {
148 em_opt = i;
149 n -= 4;
150 break;
151 }
152 }
153 if(0 == em_opt) { em_name = ""; }
154 }
155
156 // hadronic pHysics List
157 G4String had_name = name.substr(0, n);
158
159 if(0 < verbose) {
160 G4cout << "G4PhysListFactory::GetReferencePhysList <" << had_name
161 << em_name << "> EMoption= " << em_opt << G4endl;
162 }
163 G4VModularPhysicsList* p = nullptr;
164 if(had_name == "FTFP_BERT") {p = new FTFP_BERT(verbose);}
165 else if(had_name == "FTFP_BERT_HP") {p = new FTFP_BERT_HP(verbose);}
166 else if(had_name == "FTFP_BERT_TRV") {p = new FTFP_BERT_TRV(verbose);}
167 else if(had_name == "FTFP_BERT_ATL") {p = new FTFP_BERT_ATL(verbose);}
168 else if(had_name == "FTFQGSP_BERT") {p = new FTFQGSP_BERT(verbose);}
169 else if(had_name == "FTFP_INCLXX") {p = new FTFP_INCLXX(verbose);}
170 else if(had_name == "FTFP_INCLXX_HP") {p = new FTFP_INCLXX_HP(verbose);}
171 else if(had_name == "FTF_BIC") {p = new FTF_BIC(verbose);}
172 else if(had_name == "LBE") {p = new LBE();}
173 else if(had_name == "QBBC") {p = new QBBC(verbose);}
174 else if(had_name == "QGSP_BERT") {p = new QGSP_BERT(verbose);}
175 else if(had_name == "QGSP_BERT_HP") {p = new QGSP_BERT_HP(verbose);}
176 else if(had_name == "QGSP_BIC") {p = new QGSP_BIC(verbose);}
177 else if(had_name == "QGSP_BIC_HP") {p = new QGSP_BIC_HP(verbose);}
178 else if(had_name == "QGSP_BIC_AllHP") {p = new QGSP_BIC_AllHP(verbose);}
179 else if(had_name == "QGSP_FTFP_BERT") {p = new QGSP_FTFP_BERT(verbose);}
180 else if(had_name == "QGSP_INCLXX") {p = new QGSP_INCLXX(verbose);}
181 else if(had_name == "QGSP_INCLXX_HP") {p = new QGSP_INCLXX_HP(verbose);}
182 else if(had_name == "QGS_BIC") {p = new QGS_BIC(verbose);}
183 else if(had_name == "Shielding") {p = new Shielding(verbose);}
184 else if(had_name == "ShieldingLEND") {p = new ShieldingLEND(verbose);}
185 else if(had_name == "ShieldingLIQMD") {p = new Shielding(verbose,"HP","",true);}
186 else if(had_name == "ShieldingM") {p = new Shielding(verbose,"HP","M");}
187 else if(had_name == "NuBeam") {p = new NuBeam(verbose);}
188 else if(had_name == "Shielding_HP") {p = new Shielding(verbose);}
189 else if(had_name == "ShieldingLIQMD_HP") {p = new Shielding(verbose,"HP","",true);}
190 else if(had_name == "ShieldingM_HP") {p = new Shielding(verbose,"HP","M");}
191 else if(had_name == "FTFP_BERT_HPT") {p = new FTFP_BERT_HP(verbose);
193 else if(had_name == "FTFP_INCLXX_HPT") {p = new FTFP_INCLXX_HP(verbose);
195 else if(had_name == "QGSP_BERT_HPT") {p = new QGSP_BERT_HP(verbose);
197 else if(had_name == "QGSP_BIC_HPT") {p = new QGSP_BIC_HPT(verbose);}
198 else if(had_name == "QGSP_BIC_AllHPT") {p = new QGSP_BIC_AllHP(verbose);
200 else if(had_name == "QGSP_INCLXX_HPT") {p = new QGSP_INCLXX_HP(verbose);
202 else if(had_name == "Shielding_HPT") {p = new Shielding(verbose);
204 else if(had_name == "ShieldingLIQMD_HPT") {p = new Shielding(verbose,"HP","",true);
206 else if(had_name == "ShieldingM_HPT") {p = new Shielding(verbose,"HP","M");
208 else {
209 p = new FTFP_BERT(verbose);
211 ed << "PhysicsList " << had_name << " is not known;"
212 << " the default FTFP_BERT is created";
213 G4Exception("G4PhysListFactory: ","pl0003",JustWarning,ed,"");
214 }
215 if(nullptr != p) {
216 if(0 < em_opt && had_name != "LBE") {
217 if(1 == em_opt) {
219 } else if(2 == em_opt) {
221 } else if(3 == em_opt) {
223 } else if(4 == em_opt) {
225 } else if(5 == em_opt) {
226 p->ReplacePhysics(new G4EmLivermorePhysics(verbose));
227 } else if(6 == em_opt) {
228 p->ReplacePhysics(new G4EmPenelopePhysics(verbose));
229 } else if(7 == em_opt) {
230 p->ReplacePhysics(new G4EmStandardPhysicsGS(verbose));
231 } else if(8 == em_opt) {
232 p->ReplacePhysics(new G4EmStandardPhysicsSS(verbose));
233 } else if(9 == em_opt) {
234 p->ReplacePhysics(new G4EmStandardPhysics(verbose));
235 } else if(10 == em_opt) {
236 p->ReplacePhysics(new G4EmStandardPhysicsWVI(verbose));
237 } else if(11 == em_opt) {
238 p->ReplacePhysics(new G4EmLowEPPhysics(verbose));
239 }
240 }
241 theMessenger = new G4PhysListFactoryMessenger(p);
242 }
243 if(0 < verbose) G4cout << G4endl;
244 return p;
245}
TINCLXXPhysicsListHelper< G4VModularPhysicsList, false, true > FTFP_INCLXX
TINCLXXPhysicsListHelper< G4VModularPhysicsList, true, true > FTFP_INCLXX_HP
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
TINCLXXPhysicsListHelper< G4VModularPhysicsList, true, false > QGSP_INCLXX_HP
TINCLXXPhysicsListHelper< G4VModularPhysicsList, false, false > QGSP_INCLXX
void RegisterPhysics(G4VPhysicsConstructor *)
void ReplacePhysics(G4VPhysicsConstructor *)
Definition LBE.hh:59
Definition QBBC.hh:44
const char * name(G4int ptype)

Referenced by ReferencePhysList().

◆ IsReferencePhysList()

G4bool G4PhysListFactory::IsReferencePhysList ( const G4String & name) const

Definition at line 247 of file G4PhysListFactory.cc.

248{
249 G4bool res = false;
250 size_t n = name.size();
251 if(n > 4) {
252 G4String em_name = name.substr(n - 4, 4);
253 for(size_t i=1; i<nlists_em; ++i) {
254 if(listnames_em[i] == em_name) {
255 n -= 4;
256 break;
257 }
258 }
259 }
260 G4String had_name = name.substr(0, n);
261 for(size_t i=0; i<nlists_hadr; ++i) {
262 if(had_name == listnames_hadr[i]) {
263 res = true;
264 break;
265 }
266 }
267 return res;
268}
bool G4bool
Definition G4Types.hh:86

◆ ReferencePhysList()

G4VModularPhysicsList * G4PhysListFactory::ReferencePhysList ( )

Definition at line 114 of file G4PhysListFactory.cc.

115{
116 // instantiate PhysList by environment variable "PHYSLIST"
117 G4String name = "";
118 char* path = std::getenv("PHYSLIST");
119 if (path) {
120 name = G4String(path);
121 } else {
122 name = defName;
123 G4cout << "### G4PhysListFactory WARNING: "
124 << " environment variable PHYSLIST is not defined"
125 << G4endl
126 << " Default Physics Lists " << name
127 << " is instantiated"
128 << G4endl;
129 }
130 return GetReferencePhysList(name);
131}
G4VModularPhysicsList * GetReferencePhysList(const G4String &)

◆ SetVerbose()

void G4PhysListFactory::SetVerbose ( G4int val)
inline

Definition at line 68 of file G4PhysListFactory.hh.

68{ verbose = val; }

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