Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
PhysicsList.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 // Authors: Susanna Guatelli, susanna@uow.edu.au,
27 // Authors: Jeremy Davis, jad028@uowmail.edu.au
28 //
29 // Code based on the hadrontherapy advanced example
30 
31 #include "PhysicsList.hh"
32 #include "PhysicsListMessenger.hh"
33 #include "G4PhysListFactory.hh"
34 #include "G4VPhysicsConstructor.hh"
35 
36 // Physic lists (contained inside the Geant4 distribution)
38 #include "G4EmLivermorePhysics.hh"
39 #include "G4EmPenelopePhysics.hh"
40 #include "G4DecayPhysics.hh"
45 #include "G4Decay.hh"
46 #include "G4StepLimiter.hh"
47 #include "G4LossTableManager.hh"
48 #include "G4UnitsTable.hh"
49 #include "G4SystemOfUnits.hh"
50 #include "G4ProcessManager.hh"
51 
52 #include "G4IonFluctuations.hh"
54 #include "G4EmProcessOptions.hh"
57 
60 {
63  cutForGamma = defaultCutValue;
64  cutForElectron = defaultCutValue;
65  cutForPositron = defaultCutValue;
66 
67  helIsRegisted = false;
68  bicIsRegisted = false;
69  biciIsRegisted = false;
70  locIonIonInelasticIsRegistered = false;
71  radioactiveDecayIsRegisted = false;
72 
73  pMessenger = new PhysicsListMessenger(this);
74 
75  SetVerboseLevel(1);
76 
77  // EM physics
78  emPhysicsList = new G4EmStandardPhysics_option3(1);
79  emName = G4String("emstandard_opt3");
80 
81  // Decay physics and all particles
82  decPhysicsList = new G4DecayPhysics();
83 }
84 
87 {
88  delete pMessenger;
89  delete emPhysicsList;
90  delete decPhysicsList;
91  for(size_t i=0; i<hadronPhys.size(); i++) {delete hadronPhys[i];}
92 }
93 
96 {
97  G4PhysListFactory factory;
98  G4VModularPhysicsList* phys =factory.GetReferencePhysList(name);
99  G4int i=0;
100  const G4VPhysicsConstructor* elem= phys->GetPhysics(i);
101  G4VPhysicsConstructor* tmp = const_cast<G4VPhysicsConstructor*> (elem);
102  while (elem !=0)
103  {
104  RegisterPhysics(tmp);
105  elem= phys->GetPhysics(++i) ;
106  tmp = const_cast<G4VPhysicsConstructor*> (elem);
107  }
108 }
109 
112 {
113  decPhysicsList->ConstructParticle();
114 }
115 
118 {
119  // transportation
120  //
122 
123  // electromagnetic physics list
124  //
125  emPhysicsList->ConstructProcess();
126  em_config.AddModels();
127 
128  // decay physics list
129  //
130  decPhysicsList->ConstructProcess();
131 
132  // hadronic physics lists
133  for(size_t i=0; i<hadronPhys.size(); i++) {
134  hadronPhys[i]->ConstructProcess();
135  }
136 
137  // step limitation (as a full process)
138  //
139  // AddStepMax();
140 }
141 
144 {
145 
146  if (verboseLevel>1) {
147  G4cout << "PhysicsList::AddPhysicsList: <" << name << ">" << G4endl;
148  }
149  if (name == emName) return;
150 
152  // ELECTROMAGNETIC MODELS
154 
155  if (name == "standard_opt3") {
156  emName = name;
157  delete emPhysicsList;
158  emPhysicsList = new G4EmStandardPhysics_option3();
159  G4cout << "THE FOLLOWING ELECTROMAGNETIC PHYSICS LIST HAS BEEN ACTIVATED: G4EmStandardPhysics_option3" << G4endl;
160 
161  } else if (name == "LowE_Livermore") {
162  emName = name;
163  delete emPhysicsList;
164  emPhysicsList = new G4EmLivermorePhysics();
165  G4cout << "THE FOLLOWING ELECTROMAGNETIC PHYSICS LIST HAS BEEN ACTIVATED: G4EmLivermorePhysics" << G4endl;
166 
167  } else if (name == "LowE_Penelope") {
168  emName = name;
169  delete emPhysicsList;
170  emPhysicsList = new G4EmPenelopePhysics();
171  G4cout << "THE FOLLOWING ELECTROMAGNETIC PHYSICS LIST HAS BEEN ACTIVATED: G4EmLivermorePhysics" << G4endl;
172 
174  // HADRONIC MODELS
176  } else if (name == "elastic" && !helIsRegisted) {
177  G4cout << "THE FOLLOWING HADRONIC ELASTIC PHYSICS LIST HAS BEEN ACTIVATED: G4HadronElasticPhysics()" << G4endl;
178  hadronPhys.push_back( new G4HadronElasticPhysics());
179  helIsRegisted = true;
180 
181  } else if (name == "DElastic" && !helIsRegisted) {
182  hadronPhys.push_back( new G4HadronDElasticPhysics());
183  helIsRegisted = true;
184 
185  } else if (name == "HPElastic" && !helIsRegisted) {
186  hadronPhys.push_back( new G4HadronElasticPhysicsHP());
187  helIsRegisted = true;
188 
189  } else if (name == "binary" && !bicIsRegisted) {
190  hadronPhys.push_back(new G4HadronPhysicsQGSP_BIC_HP());
191  bicIsRegisted = true;
192  G4cout << "THE FOLLOWING HADRONIC INELASTIC PHYSICS LIST HAS BEEN ACTIVATED: HadronPhysicsQGSP_BIC_HP()" << G4endl;
193 
194  } else if (name == "binary_ion" && !biciIsRegisted) {
195  hadronPhys.push_back(new G4IonBinaryCascadePhysics());
196  biciIsRegisted = true;
197  G4cout << "THE FOLLOWING HADRONIC INELASTIC PHYSICS LIST HAS BEEN ACTIVATED: G4IonBinaryCascadePhysics()" << G4endl;
198  } else if (name == "radioactive_decay" && !radioactiveDecayIsRegisted ) {
199  hadronPhys.push_back(new G4RadioactiveDecayPhysics());
200  radioactiveDecayIsRegisted = true;
201  G4cout << "THE FOLLOWING HADRONIC INELASTIC PHYSICS LIST HAS BEEN ACTIVATED: G4RadioactiveDecayPhysics()" << G4endl;
202  } else {
203 
204  G4cout << "PhysicsList::AddPhysicsList: <" << name << ">"
205  << " is not defined"
206  << G4endl;
207  }
208 }
209 
210 
212 {
213  // Step limitation seen as a process
214 
216  particleIterator->reset();
217 
218  while ((*particleIterator)())
219 {
220  G4ParticleDefinition* particle = particleIterator->value();
221  G4ProcessManager* pmanager = particle->GetProcessManager();
222  pmanager -> AddProcess(new G4StepLimiter(), -1,-1,3);
223  }
224 
225 }
226 
228 {
229 
230  if (verboseLevel >0){
231  G4cout << "PhysicsList::SetCuts:";
232  G4cout << "CutLength : " << G4BestUnit(defaultCutValue,"Length") << G4endl;
233  }
234 
235  G4double lowLimit = 250. * eV;
236  G4double highLimit = 100. * GeV;
238 
239  // set cut values for gamma at first and for e- second and next for e+,
240  // because some processes for e+/e- need cut values for gamma
241  SetCutValue(cutForGamma, "gamma");
242  SetCutValue(cutForElectron, "e-");
243  SetCutValue(cutForPositron, "e+");
244 
246 }
247 
249 {
250  cutForGamma = cut;
251  SetParticleCuts(cutForGamma, G4Gamma::Gamma());
252 }
253 
255 {
256  cutForElectron = cut;
257  SetParticleCuts(cutForElectron, G4Electron::Electron());
258 }
259 
261 {
262  cutForPositron = cut;
263  SetParticleCuts(cutForPositron, G4Positron::Positron());
264 }
265 
const XML_Char * name
Definition: expat.h:151
void ConstructParticle()
Definition: PhysicsList.cc:117
void RegisterPhysics(G4VPhysicsConstructor *)
static G4LossTableManager * Instance()
void SetCutValue(G4double aCut, const G4String &pname)
void SetCutForGamma(G4double)
Definition: PhysicsList.cc:231
void SetEnergyRange(G4double lowedge, G4double highedge)
void AddPackage(const G4String &name)
Definition: PhysicsList.cc:95
void SetCutForPositron(G4double)
Definition: PhysicsList.cc:247
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
int G4int
Definition: G4Types.hh:78
void SetCutForElectron(G4double)
Definition: PhysicsList.cc:239
void AddPhysicsList(const G4String &name)
Definition: PhysicsList.cc:191
virtual void ConstructParticle()=0
void SetParticleCuts(G4double cut, G4ParticleDefinition *particle, G4Region *region=0)
void DumpCutValuesTable(G4int flag=1)
G4GLOB_DLL std::ostream G4cout
const G4VPhysicsConstructor * GetPhysics(G4int index) const
G4VModularPhysicsList * GetReferencePhysList(const G4String &)
static constexpr double eV
Definition: G4SIunits.hh:215
void SetVerboseLevel(G4int value)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void SetCuts()
Definition: PhysicsList.cc:219
void AddStepMax()
Definition: PhysicsList.cc:172
static G4ProductionCutsTable * GetProductionCutsTable()
static G4Positron * Positron()
Definition: G4Positron.cc:94
G4ProcessManager * GetProcessManager() const
static constexpr double GeV
Definition: G4SIunits.hh:217
virtual void ConstructProcess()=0
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void ConstructProcess()
Definition: PhysicsList.cc:170
static constexpr double micrometer
Definition: G4SIunits.hh:100