Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
examples/advanced/radioprotection/src/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 "G4ProcessManager.hh"
50 
51 #include "G4IonFluctuations.hh"
53 #include "G4EmProcessOptions.hh"
56 
57 /////////////////////////////////////////////////////////////////////////////
59 {
62  cutForGamma = defaultCutValue;
63  cutForElectron = defaultCutValue;
64  cutForPositron = defaultCutValue;
65 
66  helIsRegisted = false;
67  bicIsRegisted = false;
68  biciIsRegisted = false;
69  locIonIonInelasticIsRegistered = false;
70  radioactiveDecayIsRegisted = false;
71 
72  pMessenger = new PhysicsListMessenger(this);
73 
74  SetVerboseLevel(1);
75 
76  // EM physics
77  emPhysicsList = new G4EmStandardPhysics_option3(1);
78  emName = G4String("emstandard_opt3");
79 
80  // Decay physics and all particles
81  decPhysicsList = new G4DecayPhysics();
82 }
83 
84 /////////////////////////////////////////////////////////////////////////////
86 {
87  delete pMessenger;
88  delete emPhysicsList;
89  delete decPhysicsList;
90  for(size_t i=0; i<hadronPhys.size(); i++) {delete hadronPhys[i];}
91 }
92 
93 /////////////////////////////////////////////////////////////////////////////
95 {
96  G4PhysListFactory factory;
97  G4VModularPhysicsList* phys =factory.GetReferencePhysList(name);
98  G4int i=0;
99  const G4VPhysicsConstructor* elem= phys->GetPhysics(i);
100  G4VPhysicsConstructor* tmp = const_cast<G4VPhysicsConstructor*> (elem);
101  while (elem !=0)
102  {
103  RegisterPhysics(tmp);
104  elem= phys->GetPhysics(++i) ;
105  tmp = const_cast<G4VPhysicsConstructor*> (elem);
106  }
107 }
108 
109 /////////////////////////////////////////////////////////////////////////////
111 {
112  decPhysicsList->ConstructParticle();
113 }
114 
115 /////////////////////////////////////////////////////////////////////////////
117 {
118  // transportation
119  //
121 
122  // electromagnetic physics list
123  //
124  emPhysicsList->ConstructProcess();
125  em_config.AddModels();
126 
127  // decay physics list
128  //
129  decPhysicsList->ConstructProcess();
130 
131  // hadronic physics lists
132  for(size_t i=0; i<hadronPhys.size(); i++) {
133  hadronPhys[i]->ConstructProcess();
134  }
135 
136  // step limitation (as a full process)
137  //
138  // AddStepMax();
139 }
140 
141 /////////////////////////////////////////////////////////////////////////////
143 {
144 
145  if (verboseLevel>1) {
146  G4cout << "PhysicsList::AddPhysicsList: <" << name << ">" << G4endl;
147  }
148  if (name == emName) return;
149 
150  /////////////////////////////////////////////////////////////////////////////
151  // ELECTROMAGNETIC MODELS
152  /////////////////////////////////////////////////////////////////////////////
153 
154  if (name == "standard_opt3") {
155  emName = name;
156  delete emPhysicsList;
157  emPhysicsList = new G4EmStandardPhysics_option3();
158  G4cout << "THE FOLLOWING ELECTROMAGNETIC PHYSICS LIST HAS BEEN ACTIVATED: G4EmStandardPhysics_option3" << G4endl;
159 
160  } else if (name == "LowE_Livermore") {
161  emName = name;
162  delete emPhysicsList;
163  emPhysicsList = new G4EmLivermorePhysics();
164  G4cout << "THE FOLLOWING ELECTROMAGNETIC PHYSICS LIST HAS BEEN ACTIVATED: G4EmLivermorePhysics" << G4endl;
165 
166  } else if (name == "LowE_Penelope") {
167  emName = name;
168  delete emPhysicsList;
169  emPhysicsList = new G4EmPenelopePhysics();
170  G4cout << "THE FOLLOWING ELECTROMAGNETIC PHYSICS LIST HAS BEEN ACTIVATED: G4EmLivermorePhysics" << G4endl;
171 
172  /////////////////////////////////////////////////////////////////////////////
173  // HADRONIC MODELS
174  /////////////////////////////////////////////////////////////////////////////
175  } else if (name == "elastic" && !helIsRegisted) {
176  G4cout << "THE FOLLOWING HADRONIC ELASTIC PHYSICS LIST HAS BEEN ACTIVATED: G4HadronElasticPhysics()" << G4endl;
177  hadronPhys.push_back( new G4HadronElasticPhysics());
178  helIsRegisted = true;
179 
180  } else if (name == "DElastic" && !helIsRegisted) {
181  hadronPhys.push_back( new G4HadronDElasticPhysics());
182  helIsRegisted = true;
183 
184  } else if (name == "HPElastic" && !helIsRegisted) {
185  hadronPhys.push_back( new G4HadronElasticPhysicsHP());
186  helIsRegisted = true;
187 
188  } else if (name == "binary" && !bicIsRegisted) {
189  hadronPhys.push_back(new G4HadronPhysicsQGSP_BIC_HP());
190  bicIsRegisted = true;
191  G4cout << "THE FOLLOWING HADRONIC INELASTIC PHYSICS LIST HAS BEEN ACTIVATED: HadronPhysicsQGSP_BIC_HP()" << G4endl;
192 
193  } else if (name == "binary_ion" && !biciIsRegisted) {
194  hadronPhys.push_back(new G4IonBinaryCascadePhysics());
195  biciIsRegisted = true;
196  G4cout << "THE FOLLOWING HADRONIC INELASTIC PHYSICS LIST HAS BEEN ACTIVATED: G4IonBinaryCascadePhysics()" << G4endl;
197  } else if (name == "radioactive_decay" && !radioactiveDecayIsRegisted ) {
198  hadronPhys.push_back(new G4RadioactiveDecayPhysics());
199  radioactiveDecayIsRegisted = true;
200  G4cout << "THE FOLLOWING HADRONIC INELASTIC PHYSICS LIST HAS BEEN ACTIVATED: G4RadioactiveDecayPhysics()" << G4endl;
201  } else {
202 
203  G4cout << "PhysicsList::AddPhysicsList: <" << name << ">"
204  << " is not defined"
205  << G4endl;
206  }
207 }
208 
209 
211 {
212  // Step limitation seen as a process
213 
214  theParticleIterator->reset();
215 
216  while ((*theParticleIterator)())
217 {
218  G4ParticleDefinition* particle = theParticleIterator->value();
219  G4ProcessManager* pmanager = particle->GetProcessManager();
220  pmanager -> AddProcess(new G4StepLimiter(), -1,-1,3);
221  }
222 
223 }
224 
226 {
227 
228  if (verboseLevel >0){
229  G4cout << "PhysicsList::SetCuts:";
230  G4cout << "CutLength : " << G4BestUnit(defaultCutValue,"Length") << G4endl;
231  }
232 
233  G4double lowLimit = 250. * eV;
234  G4double highLimit = 100. * GeV;
236 
237  // set cut values for gamma at first and for e- second and next for e+,
238  // because some processes for e+/e- need cut values for gamma
239  SetCutValue(cutForGamma, "gamma");
240  SetCutValue(cutForElectron, "e-");
241  SetCutValue(cutForPositron, "e+");
242 
244 }
245 
247 {
248  cutForGamma = cut;
249  SetParticleCuts(cutForGamma, G4Gamma::Gamma());
250 }
251 
253 {
254  cutForElectron = cut;
255  SetParticleCuts(cutForElectron, G4Electron::Electron());
256 }
257 
259 {
260  cutForPositron = cut;
261  SetParticleCuts(cutForPositron, G4Positron::Positron());
262 }
263 
void RegisterPhysics(G4VPhysicsConstructor *)
static G4LossTableManager * Instance()
void SetCutValue(G4double aCut, const G4String &pname)
void SetEnergyRange(G4double lowedge, G4double highedge)
const XML_Char * name
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4ProcessManager * GetProcessManager() const
int G4int
Definition: G4Types.hh:78
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 &)
void SetVerboseLevel(G4int value)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4ProductionCutsTable * GetProductionCutsTable()
static G4Positron * Positron()
Definition: G4Positron.cc:94
virtual void ConstructProcess()=0
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
int micrometer
Definition: hepunit.py:34
#define theParticleIterator