source: trunk/examples/advanced/hadrontherapy/src/HadrontherapyPhysicsList.cc @ 1358

Last change on this file since 1358 was 1313, checked in by garnier, 14 years ago

geant4.9.4 beta rc0

File size: 13.4 KB
Line 
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// HadrontherapyPhysicsList.cc
27// See more at: http://g4advancedexamples.lngs.infn.it/Examples/hadrontherapy
28
29// This class provides all the physic models that can be activated inside Hadrontherapy;
30// Each model can be setted via macro commands;
31// Inside Hadrontherapy the models can be activate with three different complementar methods:
32//
33// 1. Use of the *Packages*.
34//    Packages (that are contained inside the
35//    Geant4 distribution at $G4INSTALL/source/physics_lists/lists) provide a full set
36//    of models (both electromagnetic and hadronic).
37//    The User can use this method simply add the line /physic/addPackage <nameOfPackage>
38//    in his/her macro file. No other action is required.
39//    For Hadrontherapy applications we suggest the use of the QGSP_BIC package
40//    for proton beams. The same can be used
41//    also for ligth ion beam.
42//    Example of use of package can be found in the packageQGSP_BIC.mac file.
43//
44// 2. Use of the *Physic Lists*.
45//    Physic lists are also already ready to use inside the Geant4 distribution
46//    ($G4INSTALL/source/physics_lists/builders). To use them the simple
47//    /physic/addPhysics <nameOfPhysicList> command must be used in the macro.
48//    In Hadrontherapy we provide physics list to activate Electromagnetic,
49//    Hadronic elastic and Hadronic inelastic models.
50//
51//    For Hadrontherapy we suggest the use of:
52//
53//    /physic/addPhysic/emstandard_option3 (electromagnetic model)
54//    /physic/addPhysic/QElastic (hadronic elastic model)
55//    /physic/addPhysic/binary (hadronic inelastic models for proton and neutrons)
56//    /physic/addPhysic/binary_ion (hadronic inelastic models for ions)
57//
58//    Example of the use of physics lists can be found in the macro files included in the
59//    'macro' folder .
60//
61// 3. Use of a *local* physics. In this case the models are implemented in local files
62//    contained in the Hadrontherapy folder. The use of local physic is recommended
63//    to more expert Users.
64//    We provide as local, only the LocalStandardICRU73EmPhysic.cc (an Elecromagnetic
65//    implementation containing the new ICRU73 data table for ions stopping powers)
66//    and the LocalIonIonInelasticPhysic.cc (physic list to use for the ion-ion interaction
67//    case)
68//    The *local* physics can be activated with the same /physic/addPhysic <nameOfPhysic> command;
69//
70//    While Packages approch must be used exclusively, Physics List and Local physics can
71//    be activated, if necessary, contemporaneously in the same simulation run.
72//
73//    AT MOMENT, IF ACCURATE RESULTS ARE NEDED, WE STRONGLY RECOMMEND THE USE OF THE MACROS:
74//    proton_therapy.mac: use of the built-in Geant4 physics list for proton beams)
75//    ion_therapy.mac   : use of mixed combination of native Geant4 physic lists
76//                        and local physic for ion-ion enelastic processes)
77
78#include "G4RunManager.hh"
79#include "G4Region.hh"
80#include "G4RegionStore.hh"
81#include "HadrontherapyPhysicsList.hh"
82#include "HadrontherapyPhysicsListMessenger.hh"
83#include "HadrontherapyStepMax.hh"
84#include "G4PhysListFactory.hh"
85#include "G4VPhysicsConstructor.hh"
86
87// Local physic directly implemented in the Hadronthrapy directory
88#include "LocalIonIonInelasticPhysic.hh"             // Physic dedicated to the ion-ion inelastic processes
89#include "LocalINCLIonIonInelasticPhysic.hh"         // Physic dedicated to the ion-ion inelastic processes using INCL/ABLA
90
91// Physic lists (contained inside the Geant4 source code, in the 'physicslists folder')
92#include "G4EmStandardPhysics_option3.hh"
93#include "G4EmLivermorePhysics.hh"
94#include "G4EmPenelopePhysics.hh"
95#include "G4DecayPhysics.hh"
96#include "G4HadronElasticPhysics.hh"
97#include "G4HadronDElasticPhysics.hh"
98#include "G4HadronHElasticPhysics.hh"
99#include "G4HadronQElasticPhysics.hh"
100#include "G4HadronInelasticQBBC.hh"
101#include "G4IonBinaryCascadePhysics.hh"
102#include "G4Decay.hh"
103
104#include "G4LossTableManager.hh"
105#include "G4UnitsTable.hh"
106#include "G4ProcessManager.hh"
107
108#include "G4IonFluctuations.hh"
109#include "G4IonParametrisedLossModel.hh"
110#include "G4EmProcessOptions.hh"
111
112#include "G4RadioactiveDecayPhysics.hh"
113
114/////////////////////////////////////////////////////////////////////////////
115HadrontherapyPhysicsList::HadrontherapyPhysicsList() : G4VModularPhysicsList()
116{
117  G4LossTableManager::Instance();
118  defaultCutValue = 1.*mm;
119  cutForGamma     = defaultCutValue;
120  cutForElectron  = defaultCutValue;
121  cutForPositron  = defaultCutValue;
122
123  helIsRegisted  = false;
124  bicIsRegisted  = false;
125  biciIsRegisted = false;
126  locIonIonInelasticIsRegistered = false;
127  radioactiveDecayIsRegisted = false;
128
129  stepMaxProcess  = 0;
130
131  pMessenger = new HadrontherapyPhysicsListMessenger(this);
132
133  SetVerboseLevel(1);
134
135  // EM physics
136  emPhysicsList = new G4EmStandardPhysics_option3(1);
137  emName = G4String("emstandard_opt3");
138
139  // Deacy physics and all particles
140  decPhysicsList = new G4DecayPhysics();
141}
142
143/////////////////////////////////////////////////////////////////////////////
144HadrontherapyPhysicsList::~HadrontherapyPhysicsList()
145{
146  delete pMessenger;
147  delete emPhysicsList;
148  delete decPhysicsList;
149  for(size_t i=0; i<hadronPhys.size(); i++) {delete hadronPhys[i];}
150}
151
152/////////////////////////////////////////////////////////////////////////////
153void HadrontherapyPhysicsList::AddPackage(const G4String& name)
154{
155  G4PhysListFactory factory;
156  G4VModularPhysicsList* phys =factory.GetReferencePhysList(name);
157  G4int i=0;
158  const G4VPhysicsConstructor* elem= phys->GetPhysics(i);
159  G4VPhysicsConstructor* tmp = const_cast<G4VPhysicsConstructor*> (elem);
160  while (elem !=0)
161        {
162          RegisterPhysics(tmp);
163          elem= phys->GetPhysics(++i) ;
164          tmp = const_cast<G4VPhysicsConstructor*> (elem);
165        }
166}
167
168/////////////////////////////////////////////////////////////////////////////
169void HadrontherapyPhysicsList::ConstructParticle()
170{
171  decPhysicsList->ConstructParticle();
172}
173
174/////////////////////////////////////////////////////////////////////////////
175void HadrontherapyPhysicsList::ConstructProcess()
176{
177  // transportation
178  AddTransportation();
179
180  // electromagnetic physics list
181  emPhysicsList->ConstructProcess();
182  em_config.AddModels();
183
184  // decay physics list
185  decPhysicsList->ConstructProcess();
186
187  // hadronic physics lists
188  for(size_t i=0; i<hadronPhys.size(); i++) {
189    hadronPhys[i]->ConstructProcess();
190  }
191
192  // step limitation (as a full process)
193  //
194  AddStepMax();
195}
196
197/////////////////////////////////////////////////////////////////////////////
198void HadrontherapyPhysicsList::AddPhysicsList(const G4String& name)
199{
200
201  if (verboseLevel>1) {
202    G4cout << "PhysicsList::AddPhysicsList: <" << name << ">" << G4endl;
203  }
204  if (name == emName) return;
205
206  /////////////////////////////////////////////////////////////////////////////
207  //   ELECTROMAGNETIC MODELS
208  /////////////////////////////////////////////////////////////////////////////
209    if (name == "standard_opt3") {
210    emName = name;
211    delete emPhysicsList;
212    emPhysicsList = new G4EmStandardPhysics_option3();
213    G4RunManager::GetRunManager() -> PhysicsHasBeenModified();
214    G4cout << "THE FOLLOWING ELECTROMAGNETIC PHYSICS LIST HAS BEEN ACTIVATED: G4EmStandardPhysics_option3" << G4endl;
215
216  } else if (name == "LowE_Livermore") {
217    emName = name;
218    delete emPhysicsList;
219    emPhysicsList = new G4EmLivermorePhysics();
220    G4RunManager::GetRunManager()-> PhysicsHasBeenModified();
221    G4cout << "THE FOLLOWING ELECTROMAGNETIC PHYSICS LIST HAS BEEN ACTIVATED: G4EmLivermorePhysics" << G4endl;
222
223  } else if (name == "LowE_Penelope") {
224    emName = name;
225    delete emPhysicsList;
226    emPhysicsList = new G4EmPenelopePhysics();
227    G4cout << "THE FOLLOWING ELECTROMAGNETIC PHYSICS LIST HAS BEEN ACTIVATED: G4EmLivermorePhysics" << G4endl;
228
229    /////////////////////////////////////////////////////////////////////////////
230    //   HADRONIC MODELS
231    /////////////////////////////////////////////////////////////////////////////
232  } else if (name == "elastic" && !helIsRegisted) {
233    G4cout << "THE FOLLOWING HADRONIC ELASTIC PHYSICS LIST HAS BEEN ACTIVATED: G4HadronElasticPhysics()" << G4endl;
234    hadronPhys.push_back( new G4HadronElasticPhysics());
235    helIsRegisted = true;
236
237  } else if (name == "DElastic" && !helIsRegisted) {
238    hadronPhys.push_back( new G4HadronDElasticPhysics());
239    helIsRegisted = true;
240
241  } else if (name == "HElastic" && !helIsRegisted) {
242    hadronPhys.push_back( new G4HadronHElasticPhysics());
243    helIsRegisted = true;
244
245  } else if (name == "QElastic" && !helIsRegisted) {
246    hadronPhys.push_back( new G4HadronQElasticPhysics());
247    helIsRegisted = true;
248
249  } else if (name == "binary" && !bicIsRegisted) {
250    hadronPhys.push_back(new G4HadronInelasticQBBC());
251    bicIsRegisted = true;
252    G4cout << "THE FOLLOWING HADRONIC INELASTIC PHYSICS LIST HAS BEEN ACTIVATED: G4HadronInelasticQBBC()" << G4endl;
253
254  } else if (name == "binary_ion" && !biciIsRegisted) {
255    hadronPhys.push_back(new G4IonBinaryCascadePhysics());
256    biciIsRegisted = true;
257
258  } else if (name == "local_ion_ion_inelastic" && !locIonIonInelasticIsRegistered) {
259    hadronPhys.push_back(new LocalIonIonInelasticPhysic());
260    locIonIonInelasticIsRegistered = true;
261
262  } else if (name == "local_incl_ion_ion_inelastic" && !locIonIonInelasticIsRegistered) {
263    hadronPhys.push_back(new LocalINCLIonIonInelasticPhysic());
264    locIonIonInelasticIsRegistered = true;
265
266  } else if (name == "radioactive_decay" && !radioactiveDecayIsRegisted ) {
267    hadronPhys.push_back(new G4RadioactiveDecayPhysics());
268    radioactiveDecayIsRegisted = true;
269
270  } else {
271
272    G4cout << "PhysicsList::AddPhysicsList: <" << name << ">"
273           << " is not defined"
274           << G4endl;
275  }
276}
277
278/////////////////////////////////////////////////////////////////////////////
279void HadrontherapyPhysicsList::AddStepMax()
280{
281  // Step limitation seen as a process
282  stepMaxProcess = new HadrontherapyStepMax();
283
284  theParticleIterator->reset();
285  while ((*theParticleIterator)()){
286    G4ParticleDefinition* particle = theParticleIterator->value();
287    G4ProcessManager* pmanager = particle->GetProcessManager();
288
289    if (stepMaxProcess->IsApplicable(*particle) && pmanager)
290      {
291        pmanager ->AddDiscreteProcess(stepMaxProcess);
292      }
293  }
294}
295
296/////////////////////////////////////////////////////////////////////////////
297void HadrontherapyPhysicsList::SetCuts()
298{
299
300  if (verboseLevel >0){
301    G4cout << "PhysicsList::SetCuts:";
302    G4cout << "CutLength : " << G4BestUnit(defaultCutValue,"Length") << G4endl;
303  }
304
305  // set cut values for gamma at first and for e- second and next for e+,
306  // because some processes for e+/e- need cut values for gamma
307  SetCutValue(cutForGamma, "gamma");
308  SetCutValue(cutForElectron, "e-");
309  SetCutValue(cutForPositron, "e+");
310
311  // Set cuts for detector
312  SetDetectorCut(defaultCutValue); 
313  if (verboseLevel>0) DumpCutValuesTable();
314}
315
316/////////////////////////////////////////////////////////////////////////////
317void HadrontherapyPhysicsList::SetCutForGamma(G4double cut)
318{
319  cutForGamma = cut;
320  SetParticleCuts(cutForGamma, G4Gamma::Gamma());
321}
322
323/////////////////////////////////////////////////////////////////////////////
324void HadrontherapyPhysicsList::SetCutForElectron(G4double cut)
325{
326  cutForElectron = cut;
327  SetParticleCuts(cutForElectron, G4Electron::Electron());
328}
329
330/////////////////////////////////////////////////////////////////////////////
331void HadrontherapyPhysicsList::SetCutForPositron(G4double cut)
332{
333  cutForPositron = cut;
334  SetParticleCuts(cutForPositron, G4Positron::Positron());
335}
336
337void HadrontherapyPhysicsList::SetDetectorCut(G4double cut)
338{
339
340  G4String regionName = "DetectorLog";
341  G4Region* region = G4RegionStore::GetInstance()->GetRegion(regionName);
342
343  G4ProductionCuts* cuts = new G4ProductionCuts ;
344  cuts -> SetProductionCut(cut,G4ProductionCuts::GetIndex("gamma"));
345  cuts -> SetProductionCut(cut,G4ProductionCuts::GetIndex("e-"));
346  cuts -> SetProductionCut(cut,G4ProductionCuts::GetIndex("e+"));
347  region -> SetProductionCuts(cuts);
348}
349
Note: See TracBrowser for help on using the repository browser.