Geant4-11
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//
27//---------------------------------------------------------------------------
28//
29// ClassName: G4PhysListFactory
30//
31// Author: 21 April 2008 V. Ivanchenko
32//
33// Modified:
34//
35// 2014.08.05 K.L.Genser used provision for Hadronic Physics Variant M in
36// Shielding for ShieldingM
37//
38//----------------------------------------------------------------------------
39//
40
41#include "G4PhysListFactory.hh"
42#include "FTFP_BERT.hh"
43#include "FTFP_BERT_HP.hh"
44#include "FTFP_BERT_TRV.hh"
45#include "FTFP_BERT_ATL.hh"
46#include "FTFQGSP_BERT.hh"
47#include "FTFP_INCLXX.hh"
48#include "FTFP_INCLXX_HP.hh"
49#include "FTF_BIC.hh"
50#include "LBE.hh"
51#include "QBBC.hh"
52#include "QGSP_BERT.hh"
53#include "QGSP_BERT_HP.hh"
54#include "QGSP_BIC.hh"
55#include "QGSP_BIC_HP.hh"
56#include "QGSP_BIC_AllHP.hh"
57#include "QGSP_FTFP_BERT.hh"
58#include "QGS_BIC.hh"
59#include "QGSP_INCLXX.hh"
60#include "QGSP_INCLXX_HP.hh"
61#include "Shielding.hh"
62#include "ShieldingLEND.hh"
63#include "NuBeam.hh"
64
73#include "G4EmLowEPPhysics.hh"
77#include "G4UImessenger.hh"
78
80 : defName("FTFP_BERT"),verbose(ver),theMessenger(nullptr)
81{
82 nlists_hadr = 23;
83 G4String ss[23] = {
84 "FTFP_BERT","FTFP_BERT_TRV","FTFP_BERT_ATL","FTFP_BERT_HP","FTFQGSP_BERT",
85 "FTFP_INCLXX","FTFP_INCLXX_HP","FTF_BIC", "LBE","QBBC",
86 "QGSP_BERT","QGSP_BERT_HP","QGSP_BIC","QGSP_BIC_HP","QGSP_BIC_AllHP",
87 "QGSP_FTFP_BERT","QGSP_INCLXX","QGSP_INCLXX_HP","QGS_BIC",
88 "Shielding","ShieldingLEND","ShieldingM","NuBeam"};
89 for(size_t i=0; i<nlists_hadr; ++i) {
90 listnames_hadr.push_back(ss[i]);
91 }
92
93 nlists_em = 12;
94 G4String s1[12] = {"","_EMV","_EMX","_EMY","_EMZ","_LIV","_PEN",
95 "__GS","__SS","_EM0","_WVI","__LE"};
96 for(size_t i=0; i<nlists_em; ++i) {
97 listnames_em.push_back(s1[i]);
98 }
99}
100
102{
103 delete theMessenger;
104}
105
108{
109 // instantiate PhysList by environment variable "PHYSLIST"
110 G4String name = "";
111 char* path = std::getenv("PHYSLIST");
112 if (path) {
113 name = G4String(path);
114 } else {
115 name = defName;
116 G4cout << "### G4PhysListFactory WARNING: "
117 << " environment variable PHYSLIST is not defined"
118 << G4endl
119 << " Default Physics Lists " << name
120 << " is instantiated"
121 << G4endl;
122 }
124}
125
128{
129 // analysis on the string
130 size_t n = name.size();
131
132 // last characters in the string
133 size_t em_opt = 0;
134 G4String em_name = "";
135
136 // check EM options
137 if(n > 4) {
138 em_name = name.substr(n - 4, 4);
139 for(size_t i=1; i<nlists_em; ++i) {
140 if(listnames_em[i] == em_name) {
141 em_opt = i;
142 n -= 4;
143 break;
144 }
145 }
146 if(0 == em_opt) { em_name = ""; }
147 }
148
149 // hadronic pHysics List
150 G4String had_name = name.substr(0, n);
151
152 if(0 < verbose) {
153 G4cout << "G4PhysListFactory::GetReferencePhysList <" << had_name
154 << em_name << "> EMoption= " << em_opt << G4endl;
155 }
156 G4VModularPhysicsList* p = nullptr;
157 if(had_name == "FTFP_BERT") {p = new FTFP_BERT(verbose);}
158 else if(had_name == "FTFP_BERT_HP") {p = new FTFP_BERT_HP(verbose);}
159 else if(had_name == "FTFP_BERT_TRV") {p = new FTFP_BERT_TRV(verbose);}
160 else if(had_name == "FTFP_BERT_ATL") {p = new FTFP_BERT_ATL(verbose);}
161 else if(had_name == "FTFQGSP_BERT") {p = new FTFQGSP_BERT(verbose);}
162 else if(had_name == "FTFP_INCLXX") {p = new FTFP_INCLXX(verbose);}
163 else if(had_name == "FTFP_INCLXX_HP") {p = new FTFP_INCLXX_HP(verbose);}
164 else if(had_name == "FTF_BIC") {p = new FTF_BIC(verbose);}
165 else if(had_name == "LBE") {p = new LBE();}
166 else if(had_name == "QBBC") {p = new QBBC(verbose);}
167 else if(had_name == "QGSP_BERT") {p = new QGSP_BERT(verbose);}
168 else if(had_name == "QGSP_BERT_HP") {p = new QGSP_BERT_HP(verbose);}
169 else if(had_name == "QGSP_BIC") {p = new QGSP_BIC(verbose);}
170 else if(had_name == "QGSP_BIC_HP") {p = new QGSP_BIC_HP(verbose);}
171 else if(had_name == "QGSP_BIC_AllHP") {p = new QGSP_BIC_AllHP(verbose);}
172 else if(had_name == "QGSP_FTFP_BERT") {p = new QGSP_FTFP_BERT(verbose);}
173 else if(had_name == "QGSP_INCLXX") {p = new QGSP_INCLXX(verbose);}
174 else if(had_name == "QGSP_INCLXX_HP") {p = new QGSP_INCLXX_HP(verbose);}
175 else if(had_name == "QGS_BIC") {p = new QGS_BIC(verbose);}
176 else if(had_name == "Shielding") {p = new Shielding(verbose);}
177 else if(had_name == "ShieldingLEND") {p = new ShieldingLEND(verbose);}
178 else if(had_name == "ShieldingM") {p = new Shielding(verbose,"HP","M");}
179 else if(had_name == "NuBeam") {p = new NuBeam(verbose);}
180 else {
181 p = new FTFP_BERT(verbose);
183 ed << "PhysicsList " << had_name << " is not known;"
184 << " the default FTFP_BERT is created";
185 G4Exception("G4PhysListFactory: ","pl0003",JustWarning,ed,"");
186 }
187 if(nullptr != p) {
188 if(0 < em_opt && had_name != "LBE") {
189 if(1 == em_opt) {
191 } else if(2 == em_opt) {
193 } else if(3 == em_opt) {
195 } else if(4 == em_opt) {
197 } else if(5 == em_opt) {
199 } else if(6 == em_opt) {
201 } else if(7 == em_opt) {
203 } else if(8 == em_opt) {
205 } else if(9 == em_opt) {
207 } else if(10 == em_opt) {
209 } else if(11 == em_opt) {
211 }
212 }
214 }
215 if(0 < verbose) G4cout << G4endl;
216 return p;
217}
218
220{
221 G4bool res = false;
222 size_t n = name.size();
223 if(n > 4) {
224 G4String em_name = name.substr(n - 4, 4);
225 for(size_t i=1; i<nlists_em; ++i) {
226 if(listnames_em[i] == em_name) {
227 n -= 4;
228 break;
229 }
230 }
231 }
232 G4String had_name = name.substr(0, n);
233 for(size_t i=0; i<nlists_hadr; ++i) {
234 if(had_name == listnames_hadr[i]) {
235 res = true;
236 break;
237 }
238 }
239 return res;
240}
241
242const std::vector<G4String>&
244{
245 return listnames_hadr;
246}
247
248const std::vector<G4String>&
250{
251 return listnames_em;
252}
253
TINCLXXPhysicsListHelper< G4VModularPhysicsList, false, true > FTFP_INCLXX
Definition: FTFP_INCLXX.hh:40
TINCLXXPhysicsListHelper< G4VModularPhysicsList, true, true > FTFP_INCLXX_HP
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
TINCLXXPhysicsListHelper< G4VModularPhysicsList, true, false > QGSP_INCLXX_HP
TINCLXXPhysicsListHelper< G4VModularPhysicsList, false, false > QGSP_INCLXX
std::vector< G4String > listnames_hadr
const std::vector< G4String > & AvailablePhysListsEM() const
G4UImessenger * theMessenger
std::vector< G4String > listnames_em
G4bool IsReferencePhysList(const G4String &) const
G4PhysListFactory(G4int ver=1)
G4VModularPhysicsList * ReferencePhysList()
G4VModularPhysicsList * GetReferencePhysList(const G4String &)
const std::vector< G4String > & AvailablePhysLists() const
void ReplacePhysics(G4VPhysicsConstructor *)
Definition: LBE.hh:59
Definition: NuBeam.hh:42
Definition: QBBC.hh:44
const char * name(G4int ptype)