Geant4  10.01.p03
G4PhysListFactory.cc
Go to the documentation of this file.
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 // $Id: G4PhysListFactory.cc 83413 2014-08-21 15:21:43Z gcosmo $
27 //
28 //---------------------------------------------------------------------------
29 //
30 // ClassName: G4PhysListFactory
31 //
32 // Author: 21 April 2008 V. Ivanchenko
33 //
34 // Modified:
35 //
36 // 2014.08.05 K.L.Genser used provision for Hadronic Physics Variant M in
37 // Shielding for ShieldingM
38 //
39 //----------------------------------------------------------------------------
40 //
41 
42 #include "G4PhysListFactory.hh"
43 #include "FTFP_BERT.hh"
44 #include "FTFP_BERT_HP.hh"
45 #include "FTFP_BERT_TRV.hh"
46 #include "FTFP_INCLXX.hh"
47 #include "FTFP_INCLXX_HP.hh"
48 #include "FTF_BIC.hh"
49 #include "LBE.hh"
50 #include "QBBC.hh"
51 #include "QGSP_BERT.hh"
52 #include "QGSP_BERT_HP.hh"
53 #include "QGSP_BIC.hh"
54 #include "QGSP_BIC_HP.hh"
55 #include "QGSP_FTFP_BERT.hh"
56 #include "QGS_BIC.hh"
57 #include "QGSP_INCLXX.hh"
58 #include "QGSP_INCLXX_HP.hh"
59 #include "Shielding.hh"
60 #include "NuBeam.hh"
61 
62 #include "G4EmStandardPhysics.hh"
67 #include "G4EmLivermorePhysics.hh"
68 #include "G4EmPenelopePhysics.hh"
69 
71  : defName("FTFP_BERT"),verbose(1)
72 {
73  nlists_hadr = 20;
74  G4String ss[20] = {
75  "FTFP_BERT","FTFP_BERT_TRV","FTFP_BERT_HP","FTFP_INCLXX",
76  "FTFP_INCLXX_HP","FTF_BIC", "LBE","QBBC",
77  "QGSP_BERT","QGSP_BERT_HP","QGSP_BIC","QGSP_BIC_HP",
78  "QGSP_FTFP_BERT","QGSP_INCLXX","QGSP_INCLXX_HP","QGS_BIC",
79  "Shielding","ShieldingLEND","ShieldingM","NuBeam"};
80  for(size_t i=0; i<nlists_hadr; ++i) {
81  listnames_hadr.push_back(ss[i]);
82  }
83 
84  nlists_em = 7;
85  G4String s1[7] = {"","_EMV","_EMX","_EMY","_EMZ","_LIV","_PEN"};
86  for(size_t i=0; i<nlists_em; ++i) {
87  listnames_em.push_back(s1[i]);
88  }
89 }
90 
92 {}
93 
96 {
97  // instantiate PhysList by environment variable "PHYSLIST"
98  G4String name = "";
99  char* path = getenv("PHYSLIST");
100  if (path) {
101  name = G4String(path);
102  } else {
103  name = defName;
104  G4cout << "### G4PhysListFactory WARNING: "
105  << " environment variable PHYSLIST is not defined"
106  << G4endl
107  << " Default Physics Lists " << name
108  << " is instantiated"
109  << G4endl;
110  }
111  return GetReferencePhysList(name);
112 }
113 
116 {
117  // analysis on the string
118  size_t n = name.size();
119 
120  // last characters in the string
121  size_t em_opt = 0;
122  G4String em_name = "";
123 
124  // check EM options
125  if(n > 4) {
126  em_name = name.substr(n - 4, 4);
127  for(size_t i=1; i<nlists_em; ++i) {
128  if(listnames_em[i] == em_name) {
129  em_opt = i;
130  n -= 4;
131  break;
132  }
133  }
134  if(0 == em_opt) { em_name = ""; }
135  }
136 
137  // hadronic pHysics List
138  G4String had_name = name.substr(0, n);
139 
140  if(0 < verbose) {
141  G4cout << "G4PhysListFactory::GetReferencePhysList <" << had_name
142  << em_name << "> EMoption= " << em_opt << G4endl;
143  }
144  G4VModularPhysicsList* p = 0;
145  if(had_name == "FTFP_BERT") {p = new FTFP_BERT(verbose);}
146  else if(had_name == "FTFP_BERT_HP") {p = new FTFP_BERT_HP(verbose);}
147  else if(had_name == "FTFP_BERT_TRV") {p = new FTFP_BERT_TRV(verbose);}
148  else if(had_name == "FTFP_INCLXX") {p = new FTFP_INCLXX(verbose);}
149  else if(had_name == "FTFP_INCLXX_HP") {p = new FTFP_INCLXX_HP(verbose);}
150  else if(had_name == "FTF_BIC") {p = new FTF_BIC(verbose);}
151  else if(had_name == "LBE") {p = new LBE();}
152  else if(had_name == "QBBC") {p = new QBBC(verbose);}
153  else if(had_name == "QGSP_BERT") {p = new QGSP_BERT(verbose);}
154  else if(had_name == "QGSP_BERT_HP") {p = new QGSP_BERT_HP(verbose);}
155  else if(had_name == "QGSP_BIC") {p = new QGSP_BIC(verbose);}
156  else if(had_name == "QGSP_BIC_HP") {p = new QGSP_BIC_HP(verbose);}
157  else if(had_name == "QGSP_FTFP_BERT") {p = new QGSP_FTFP_BERT(verbose);}
158  else if(had_name == "QGSP_INCLXX") {p = new QGSP_INCLXX(verbose);}
159  else if(had_name == "QGSP_INCLXX_HP") {p = new QGSP_INCLXX_HP(verbose);}
160  else if(had_name == "QGS_BIC") {p = new QGS_BIC(verbose);}
161  else if(had_name == "Shielding") {p = new Shielding(verbose);}
162  else if(had_name == "ShieldingLEND") {p = new Shielding(verbose,"LEND");}
163  else if(had_name == "ShieldingM") {p = new Shielding(verbose,"HP","M");}
164  else if(had_name == "NuBeam") {p = new NuBeam(verbose);}
165  else {
166  G4cout << "### G4PhysListFactory WARNING: "
167  << "PhysicsList " << had_name << " is not known"
168  << G4endl;
169  }
170  if(p) {
171  G4cout << "<<< Reference Physics List " << had_name
172  << em_name << " is built" << G4endl;
173  G4int ver = p->GetVerboseLevel();
174  p->SetVerboseLevel(0);
175  if(0 < em_opt) {
176  if(1 == em_opt) {
178  } else if(2 == em_opt) {
180  } else if(3 == em_opt) {
182  } else if(4 == em_opt) {
184  } else if(5 == em_opt) {
186  } else if(6 == em_opt) {
188  }
189  }
190  p->SetVerboseLevel(ver);
191  }
192  G4cout << G4endl;
193  return p;
194 }
195 
197 {
198  G4bool res = false;
199  size_t n = name.size();
200  if(n > 4) {
201  G4String em_name = name.substr(n - 4, 4);
202  for(size_t i=1; i<nlists_em; ++i) {
203  if(listnames_em[i] == em_name) {
204  n -= 4;
205  break;
206  }
207  }
208  }
209  G4String had_name = name.substr(0, n);
210  for(size_t i=0; i<nlists_hadr; ++i) {
211  if(had_name == listnames_hadr[i]) {
212  res = true;
213  break;
214  }
215  }
216  return res;
217 }
218 
219 const std::vector<G4String>&
221 {
222  return listnames_hadr;
223 }
224 
225 const std::vector<G4String>&
227 {
228  return listnames_em;
229 }
230 
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
G4String name
Definition: TRTMaterials.hh:40
TLBE< G4VModularPhysicsList > LBE
Definition: LBE.hh:123
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
G4VModularPhysicsList * GetReferencePhysList(const G4String &)
G4bool IsReferencePhysList(const G4String &)
bool G4bool
Definition: G4Types.hh:79
TQGS_BIC< G4VModularPhysicsList > QGS_BIC
Definition: QGS_BIC.hh:63
const std::vector< G4String > & AvailablePhysListsEM() const
TFTFP_BERT_TRV< G4VModularPhysicsList > FTFP_BERT_TRV
const G4int n
void SetVerboseLevel(G4int value)
TNuBeam< G4VModularPhysicsList > NuBeam
Definition: NuBeam.hh:64
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
std::vector< G4String > listnames_hadr
TINCLXXPhysicsListHelper< G4VModularPhysicsList, false, false > QGSP_INCLXX
TINCLXXPhysicsListHelper< G4VModularPhysicsList, true, false > QGSP_INCLXX_HP
#define G4endl
Definition: G4ios.hh:61
G4VModularPhysicsList * ReferencePhysList()
TQGSP_BIC< G4VModularPhysicsList > QGSP_BIC
Definition: QGSP_BIC.hh:62
TFTF_BIC< G4VModularPhysicsList > FTF_BIC
Definition: FTF_BIC.hh:62
TQGSP_BERT_HP< G4VModularPhysicsList > QGSP_BERT_HP
Definition: QGSP_BERT_HP.hh:62
Definition: QBBC.hh:44
TINCLXXPhysicsListHelper< G4VModularPhysicsList, true, true > FTFP_INCLXX_HP
const std::vector< G4String > & AvailablePhysLists() const